Actual source code: ex13.c

slepc-3.7.1 2016-05-27
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2016, 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[] = "Generalized Symmetric eigenproblem.\n\n"
 23:   "The problem is Ax = lambda Bx, with:\n"
 24:   "   A = Laplacian operator in 2-D\n"
 25:   "   B = diagonal matrix with all values equal to 4 except nulldim zeros\n\n"
 26:   "The command line options are:\n"
 27:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 28:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n"
 29:   "  -nulldim <k>, where <k> = dimension of the nullspace of B.\n\n";

 31: #include <slepceps.h>

 35: int main(int argc,char **argv)
 36: {
 37:   Mat            A,B;         /* matrices */
 38:   EPS            eps;         /* eigenproblem solver context */
 39:   ST             st;          /* spectral transformation context */
 40:   EPSType        type;
 41:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j,nulldim=0;
 42:   PetscBool      flag,terse;

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

 47:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 48:   PetscOptionsGetInt(NULL,NULL,"-m",&m,&flag);
 49:   if (!flag) m=n;
 50:   N = n*m;
 51:   PetscOptionsGetInt(NULL,NULL,"-nulldim",&nulldim,NULL);
 52:   PetscPrintf(PETSC_COMM_WORLD,"\nGeneralized Symmetric Eigenproblem, N=%D (%Dx%D grid), null(B)=%D\n\n",N,n,m,nulldim);

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:      Compute the matrices that define the eigensystem, Ax=kBx
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 58:   MatCreate(PETSC_COMM_WORLD,&A);
 59:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 60:   MatSetFromOptions(A);
 61:   MatSetUp(A);

 63:   MatCreate(PETSC_COMM_WORLD,&B);
 64:   MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,N,N);
 65:   MatSetFromOptions(B);
 66:   MatSetUp(B);

 68:   MatGetOwnershipRange(A,&Istart,&Iend);
 69:   for (II=Istart;II<Iend;II++) {
 70:     i = II/n; j = II-i*n;
 71:     if (i>0) { MatSetValue(A,II,II-n,-1.0,INSERT_VALUES); }
 72:     if (i<m-1) { MatSetValue(A,II,II+n,-1.0,INSERT_VALUES); }
 73:     if (j>0) { MatSetValue(A,II,II-1,-1.0,INSERT_VALUES); }
 74:     if (j<n-1) { MatSetValue(A,II,II+1,-1.0,INSERT_VALUES); }
 75:     MatSetValue(A,II,II,4.0,INSERT_VALUES);
 76:     if (II>=nulldim) { MatSetValue(B,II,II,4.0,INSERT_VALUES); }
 77:   }

 79:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 81:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 82:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:                 Create the eigensolver and set various options
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 88:   /*
 89:      Create eigensolver context
 90:   */
 91:   EPSCreate(PETSC_COMM_WORLD,&eps);

 93:   /*
 94:      Set operators. In this case, it is a generalized eigenvalue problem
 95:   */
 96:   EPSSetOperators(eps,A,B);
 97:   EPSSetProblemType(eps,EPS_GHEP);

 99:   /*
100:      Set solver parameters at runtime
101:   */
102:   EPSSetFromOptions(eps);

104:   PetscObjectTypeCompareAny((PetscObject)eps,&flag,EPSBLOPEX,EPSRQCG,"");
105:   if (flag) {
106:     EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);
107:   } else {
108:     /*
109:        Select portion of spectrum
110:     */
111:     EPSSetTarget(eps,0.0);
112:     EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
113:     /*
114:        Use shift-and-invert to avoid solving linear systems with a singular B
115:        in case nulldim>0
116:     */
117:     PetscObjectTypeCompareAny((PetscObject)eps,&flag,EPSGD,EPSJD,"");
118:     if (!flag) {
119:       EPSGetST(eps,&st);
120:       STSetType(st,STSINVERT);
121:     }
122:   }

124:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125:                       Solve the eigensystem
126:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

128:   EPSSolve(eps);

130:   /*
131:      Optional: Get some information from the solver and display it
132:   */
133:   EPSGetType(eps,&type);
134:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
135:   EPSGetDimensions(eps,&nev,NULL,NULL);
136:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);

138:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139:                     Display solution and clean up
140:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

142:   /* show detailed info unless -terse option is given by user */
143:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
144:   if (terse) {
145:     EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);
146:   } else {
147:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
148:     EPSReasonView(eps,PETSC_VIEWER_STDOUT_WORLD);
149:     EPSErrorView(eps,EPS_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
150:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
151:   }
152:   EPSDestroy(&eps);
153:   MatDestroy(&A);
154:   MatDestroy(&B);
155:   SlepcFinalize();
156:   return ierr;
157: }