Actual source code: ex25.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2014, 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[] = "Spectrum slicing on generalized symmetric eigenproblem.\n\n"
 23:   "The problem is similar to ex13.c.\n\n"
 24:   "The command line options are:\n"
 25:   "  -n <n>, where <n> = number of grid subdivisions in x dimension.\n"
 26:   "  -m <m>, where <m> = number of grid subdivisions in y dimension.\n";

 28: #include <slepceps.h>

 32: int main(int argc,char **argv)
 33: {
 34:   Mat            A,B;         /* matrices */
 35:   EPS            eps;         /* eigenproblem solver context */
 36:   ST             st;          /* spectral transformation context */
 37:   KSP            ksp;
 38:   PC             pc;
 39:   EPSType        type;
 40:   PetscInt       N,n=10,m,Istart,Iend,II,nev,i,j;
 41:   PetscReal      inta,intb;
 42:   PetscBool      flag;
 44: #if !defined(PETSC_HAVE_MUMPS)
 45:   PetscMPIInt    size;
 46: #endif

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

 50:   PetscOptionsGetInt(NULL,"-n",&n,NULL);
 51:   PetscOptionsGetInt(NULL,"-m",&m,&flag);
 52:   if (!flag) m=n;
 53:   N = n*m;
 54:   PetscPrintf(PETSC_COMM_WORLD,"\nSpectrum slicing on GHEP, N=%D (%Dx%D grid)\n\n",N,n,m);

 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 57:      Compute the matrices that define the eigensystem, Ax=kBx
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

 81:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 82:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 83:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:                 Create the eigensolver and set various options
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 95:   /*
 96:      Set operators and set problem type
 97:   */
 98:   EPSSetOperators(eps,A,B);
 99:   EPSSetProblemType(eps,EPS_GHEP);

101:   /*
102:      Set interval for spectrum slicing
103:   */
104:   inta = 0.1;
105:   intb = 0.2;
106:   EPSSetInterval(eps,inta,intb);
107:   EPSSetWhichEigenpairs(eps,EPS_ALL);

109:   /*
110:      Set shift-and-invert with Cholesky; select MUMPS if available
111:   */

113:   EPSGetST(eps,&st);
114:   STSetType(st,STSINVERT);
115:   
116:   STGetKSP(st,&ksp);
117:   KSPSetType(ksp,KSPPREONLY);
118:   KSPGetPC(ksp,&pc);
119:   PCSetType(pc,PCCHOLESKY);
120:     
121: #if defined(PETSC_HAVE_MUMPS)
122: #if defined(PETSC_USE_COMPLEX)
123:   SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Spectrum slicing with MUMPS is not available for complex scalars");
124: #endif
125:   PCFactorSetMatSolverPackage(pc,MATSOLVERMUMPS); 
126:   /*
127:      must use runtime option '-mat_mumps_icntl_13 1' (turn off scaLAPACK for
128:      matrix inertia), currently there is no better way of setting this in program
129:   */
130:   PetscOptionsInsertString("-mat_mumps_icntl_13 1"); 
131: #else
132:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
133:   if (size>1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Configure with MUMPS if you want to run this example in parallel");
134: #endif

136:   /*
137:      Set solver parameters at runtime
138:   */
139:   EPSSetFromOptions(eps);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:                       Solve the eigensystem
143:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144:   EPSSolve(eps);

146:   /*
147:      Show eigenvalues in interval and print solution
148:   */
149:   EPSGetType(eps,&type);
150:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n\n",type);
151:   EPSGetDimensions(eps,&nev,NULL,NULL);
152:   EPSGetInterval(eps,&inta,&intb);
153:   PetscPrintf(PETSC_COMM_WORLD," %D eigenvalues found in [%g, %g]\n",nev,(double)inta,(double)intb);
154:   EPSPrintSolution(eps,NULL);

156:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157:                     Clean up
158:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159:   EPSDestroy(&eps);
160:   MatDestroy(&A);
161:   MatDestroy(&B);
162:   SlepcFinalize();
163:   return 0;
164: }