Subversion Repositories slepc-dev

Rev

Details | Compare with Previous | Last modification | View Log | RSS feed

Rev Author Line No. Line
1921 jroman 1
/*
2
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3
   SLEPc - Scalable Library for Eigenvalue Problem Computations
2575 eromero 4
   Copyright (c) 2002-2011, Universitat Politecnica de Valencia, Spain
1921 jroman 5
 
6
   This file is part of SLEPc.
7
 
8
   SLEPc is free software: you can redistribute it and/or modify it under  the
9
   terms of version 3 of the GNU Lesser General Public License as published by
10
   the Free Software Foundation.
11
 
12
   SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
13
   WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
14
   FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
15
   more details.
16
 
17
   You  should have received a copy of the GNU Lesser General  Public  License
18
   along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
19
   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
20
*/
21
 
22
static char help[] = "Solves a quadratic eigenproblem (l^2*M + l*C + K)*x = 0 with matrices loaded from a file.\n\n"
23
  "The command line options are:\n"
24
  "  -M <filename>, where <filename> = matrix (M) file in PETSc binary form.\n"
25
  "  -C <filename>, where <filename> = matrix (C) file in PETSc binary form.\n"
26
  "  -K <filename>, where <filename> = matrix (K) file in PETSc binary form.\n\n";
27
 
2283 jroman 28
#include <slepcqep.h>
1921 jroman 29
 
30
#undef __FUNCT__
31
#define __FUNCT__ "main"
2331 jroman 32
int main(int argc,char **argv)
1921 jroman 33
{
2331 jroman 34
  Mat            M,C,K;           /* problem matrices */
2316 jroman 35
  QEP            qep;             /* quadratic eigenproblem solver context */
1921 jroman 36
  const QEPType  type;
2513 jroman 37
  PetscReal      tol;
38
  PetscInt       nev,maxit,its;
2374 jroman 39
  char           filename[PETSC_MAX_PATH_LEN];
2316 jroman 40
  PetscViewer    viewer;
41
  PetscBool      flg;
2317 jroman 42
  PetscErrorCode ierr;
1921 jroman 43
 
44
  SlepcInitialize(&argc,&argv,(char*)0,help);
45
 
46
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47
        Load the matrices that define the quadratic eigenproblem
48
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49
 
50
  ierr = PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic eigenproblem stored in file.\n\n");CHKERRQ(ierr);
51
#if defined(PETSC_USE_COMPLEX)
52
  ierr = PetscPrintf(PETSC_COMM_WORLD," Reading COMPLEX matrices from binary files...\n");CHKERRQ(ierr);
53
#else
54
  ierr = PetscPrintf(PETSC_COMM_WORLD," Reading REAL matrices from binary files...\n");CHKERRQ(ierr);
55
#endif
56
 
2374 jroman 57
  ierr = PetscOptionsGetString(PETSC_NULL,"-M",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1921 jroman 58
  if (!flg) {
2214 jroman 59
    SETERRQ(PETSC_COMM_WORLD,1,"Must indicate a file name for matrix M with the -M option.");
1921 jroman 60
  }
61
  ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
2246 eromero 62
  ierr = MatCreate(PETSC_COMM_WORLD,&M);CHKERRQ(ierr);
63
  ierr = MatSetFromOptions(M);CHKERRQ(ierr);
64
  ierr = MatLoad(M,viewer);CHKERRQ(ierr);
2305 jroman 65
  ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1921 jroman 66
 
2374 jroman 67
  ierr = PetscOptionsGetString(PETSC_NULL,"-C",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1921 jroman 68
  if (!flg) {
2214 jroman 69
    SETERRQ(PETSC_COMM_WORLD,1,"Must indicate a file name for matrix C with the -C option.");
1921 jroman 70
  }
71
  ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
2246 eromero 72
  ierr = MatCreate(PETSC_COMM_WORLD,&C);CHKERRQ(ierr);
73
  ierr = MatSetFromOptions(C);CHKERRQ(ierr);
74
  ierr = MatLoad(C,viewer);CHKERRQ(ierr);
2305 jroman 75
  ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1921 jroman 76
 
2374 jroman 77
  ierr = PetscOptionsGetString(PETSC_NULL,"-K",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr);
1921 jroman 78
  if (!flg) {
2214 jroman 79
    SETERRQ(PETSC_COMM_WORLD,1,"Must indicate a file name for matrix K with the -K option.");
1921 jroman 80
  }
81
  ierr = PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);CHKERRQ(ierr);
2246 eromero 82
  ierr = MatCreate(PETSC_COMM_WORLD,&K);CHKERRQ(ierr);
83
  ierr = MatSetFromOptions(K);CHKERRQ(ierr);
84
  ierr = MatLoad(K,viewer);CHKERRQ(ierr);
2305 jroman 85
  ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr);
1921 jroman 86
 
87
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88
                Create the eigensolver and set various options
89
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90
 
91
  /*
92
     Create eigensolver context
93
  */
94
  ierr = QEPCreate(PETSC_COMM_WORLD,&qep);CHKERRQ(ierr);
95
 
96
  /*
97
     Set matrices
98
  */
99
  ierr = QEPSetOperators(qep,M,C,K);CHKERRQ(ierr);
100
 
101
  /*
102
     Set solver parameters at runtime
103
  */
104
  ierr = QEPSetFromOptions(qep);CHKERRQ(ierr);
105
 
106
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107
                      Solve the eigensystem
108
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109
 
110
  ierr = QEPSolve(qep);CHKERRQ(ierr);
2331 jroman 111
  ierr = QEPGetIterationNumber(qep,&its);CHKERRQ(ierr);
2513 jroman 112
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of iterations of the method: %D\n",its);CHKERRQ(ierr);
1921 jroman 113
 
114
  /*
115
     Optional: Get some information from the solver and display it
116
  */
117
  ierr = QEPGetType(qep,&type);CHKERRQ(ierr);
118
  ierr = PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);CHKERRQ(ierr);
119
  ierr = QEPGetDimensions(qep,&nev,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);
2513 jroman 120
  ierr = PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);CHKERRQ(ierr);
1921 jroman 121
  ierr = QEPGetTolerances(qep,&tol,&maxit);CHKERRQ(ierr);
2513 jroman 122
  ierr = PetscPrintf(PETSC_COMM_WORLD," Stopping condition: tol=%.4G, maxit=%D\n",tol,maxit);CHKERRQ(ierr);
1921 jroman 123
 
124
  /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125
                    Display solution and clean up
126
     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127
 
2513 jroman 128
  ierr = QEPPrintSolution(qep,PETSC_NULL);CHKERRQ(ierr);
2312 jroman 129
  ierr = QEPDestroy(&qep);CHKERRQ(ierr);
2305 jroman 130
  ierr = MatDestroy(&M);CHKERRQ(ierr);
131
  ierr = MatDestroy(&C);CHKERRQ(ierr);
132
  ierr = MatDestroy(&K);CHKERRQ(ierr);
1921 jroman 133
  ierr = SlepcFinalize();CHKERRQ(ierr);
134
  return 0;
135
}
136