Actual source code: test1.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  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.

 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.

 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: */

 22: static char help[] = "Test the solution of a QEP without calling QEPSetFromOptions (based on ex16.c).\n\n"
 23:   "The command line options are:\n"
 24:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 25:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
 26:   "  -type <qep_type> = qep type to test.\n"
 27:   "  -epstype <eps_type> = eps type to test (for linear).\n\n";

 29: #include <slepcqep.h>

 33: int main(int argc,char **argv)
 34: {
 35:   Mat            M,C,K;           /* problem matrices */
 36:   QEP            qep;             /* quadratic eigenproblem solver context */
 37:   QEPType        type;
 38:   PetscInt       N,n=10,m,Istart,Iend,II,nev,maxit,i,j;
 39:   PetscBool      flag,isgd2;
 40:   char           qeptype[30] = "linear",epstype[30] = "";
 41:   EPS            eps;
 42:   ST             st;
 43:   KSP            ksp;
 44:   PC             pc;

 47:   SlepcInitialize(&argc,&argv,(char*)0,help);

 49:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 50:   PetscOptionsGetInt(NULL,"-m",&m,&flag);
 51:   if (!flag) m=n;
 52:   N = n*m;
 53:   PetscOptionsGetString(NULL,"-type",qeptype,30,NULL);
 54:   PetscOptionsGetString(NULL,"-epstype",epstype,30,&flag);
 55:   PetscPrintf(PETSC_COMM_WORLD,"\nQuadratic Eigenproblem, N=%D (%Dx%D grid)",N,n,m);
 56:   PetscPrintf(PETSC_COMM_WORLD,"\nQEP type: %s",qeptype);
 57:   if (flag) {
 58:     PetscPrintf(PETSC_COMM_WORLD,"\nEPS type: %s",epstype);
 59:   }
 60:   PetscPrintf(PETSC_COMM_WORLD,"\n\n");

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:      Compute the matrices that define the eigensystem, (k^2*M+k*C+K)x=0
 64:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 66:   /* K is the 2-D Laplacian */
 67:   MatCreate(PETSC_COMM_WORLD,&K);
 68:   MatSetSizes(K,PETSC_DECIDE,PETSC_DECIDE,N,N);
 69:   MatSetFromOptions(K);
 70:   MatSetUp(K);

 72:   MatGetOwnershipRange(K,&Istart,&Iend);
 73:   for (II=Istart;II<Iend;II++) {
 74:     i = II/n; j = II-i*n;
 75:     if (i>0) { MatSetValue(K,II,II-n,-1.0,INSERT_VALUES); }
 76:     if (i<m-1) { MatSetValue(K,II,II+n,-1.0,INSERT_VALUES); }
 77:     if (j>0) { MatSetValue(K,II,II-1,-1.0,INSERT_VALUES); }
 78:     if (j<n-1) { MatSetValue(K,II,II+1,-1.0,INSERT_VALUES); }
 79:     MatSetValue(K,II,II,4.0,INSERT_VALUES);
 80:   }

 82:   MatAssemblyBegin(K,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyEnd(K,MAT_FINAL_ASSEMBLY);

 85:   /* C is the zero matrix */
 86:   MatCreate(PETSC_COMM_WORLD,&C);
 87:   MatSetSizes(C,PETSC_DECIDE,PETSC_DECIDE,N,N);
 88:   MatSetFromOptions(C);
 89:   MatSetUp(C);
 90:   MatAssemblyBegin(C,MAT_FINAL_ASSEMBLY);
 91:   MatAssemblyEnd(C,MAT_FINAL_ASSEMBLY);

 93:   /* M is the identity matrix */
 94:   MatCreate(PETSC_COMM_WORLD,&M);
 95:   MatSetSizes(M,PETSC_DECIDE,PETSC_DECIDE,N,N);
 96:   MatSetFromOptions(M);
 97:   MatSetUp(M);
 98:   MatAssemblyBegin(M,MAT_FINAL_ASSEMBLY);
 99:   MatAssemblyEnd(M,MAT_FINAL_ASSEMBLY);
100:   MatShift(M,1.0);

102:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
103:                 Create the eigensolver and set various options
104:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

106:   /*
107:      Create eigensolver context
108:   */
109:   QEPCreate(PETSC_COMM_WORLD,&qep);

111:   /*
112:      Set matrices and problem type
113:   */
114:   QEPSetOperators(qep,M,C,K);
115:   QEPSetProblemType(qep,QEP_GENERAL);
116:   QEPSetDimensions(qep,4,20,PETSC_DEFAULT);
117:   QEPSetTolerances(qep,PETSC_SMALL,PETSC_DEFAULT);

119:   /*
120:      Set solver type at runtime
121:   */
122:   QEPSetType(qep,qeptype);
123:   if (flag) {
124:     PetscObjectTypeCompare((PetscObject)qep,QEPLINEAR,&flag);
125:     if (flag) {
126:       QEPLinearGetEPS(qep,&eps);
127:       PetscStrcmp(epstype,"gd2",&isgd2);
128:       if (isgd2) {
129:         EPSSetType(eps,EPSGD);
130:         EPSGDSetDoubleExpansion(eps,PETSC_TRUE);
131:       } else {
132:         EPSSetType(eps,epstype);
133:       }
134:       EPSGetST(eps,&st);
135:       STGetKSP(st,&ksp);
136:       KSPGetPC(ksp,&pc);
137:       PCSetType(pc,PCJACOBI);
138:       PetscObjectTypeCompare((PetscObject)eps,EPSGD,&flag);
139:     }
140:   }

142:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143:                       Solve the eigensystem
144:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

146:   QEPSolve(qep);

148:   /*
149:      Optional: Get some information from the solver and display it
150:   */
151:   QEPGetType(qep,&type);
152:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
153:   QEPGetDimensions(qep,&nev,NULL,NULL);
154:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
155:   QEPGetTolerances(qep,NULL,&maxit);
156:   PetscPrintf(PETSC_COMM_WORLD," Stopping condition: maxit=%D\n",maxit);

158:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159:                     Display solution and clean up
160:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

162:   QEPPrintSolution(qep,NULL);

164:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
165:      Free work space
166:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

168:   QEPDestroy(&qep);
169:   MatDestroy(&M);
170:   MatDestroy(&C);
171:   MatDestroy(&K);
172:   SlepcFinalize();
173:   return 0;
174: }