Actual source code: test21.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 11: static char help[] = "Illustrates region filtering. "
 12:   "Based on ex5.\n"
 13:   "The command line options are:\n"
 14:   "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";

 16: #include <slepceps.h>

 18: /*
 19:    User-defined routines
 20: */
 21: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);

 23: int main(int argc,char **argv)
 24: {
 25:   Mat            A;
 26:   EPS            eps;
 27:   ST             st;
 28:   RG             rg;
 29:   PetscReal      radius,tol=PETSC_SMALL;
 30:   PetscScalar    target=0.5,kr,ki;
 31:   PetscComplex   *eigs,eval;
 32:   PetscInt       N,m=15,nev,i,nconv;
 33:   PetscBool      checkfile;
 34:   char           filename[PETSC_MAX_PATH_LEN];
 35:   PetscViewer    viewer;
 37:   char           str[50];

 39:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 40: #if defined(PETSC_HAVE_COMPLEX)
 41:   PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL);
 42:   N = m*(m+1)/2;
 43:   PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%D (m=%D)\n",N,m);
 44:   PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL);
 45:   SlepcSNPrintfScalar(str,sizeof(str),target,PETSC_FALSE);
 46:   PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str);

 48:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:      Compute the operator matrix that defines the eigensystem, Ax=kx
 50:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 52:   MatCreate(PETSC_COMM_WORLD,&A);
 53:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N);
 54:   MatSetFromOptions(A);
 55:   MatSetUp(A);
 56:   MatMarkovModel(m,A);

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:                 Create a standalone spectral transformation
 60:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 62:   STCreate(PETSC_COMM_WORLD,&st);
 63:   STSetType(st,STSINVERT);
 64:   STSetShift(st,target);

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:                     Create a region for filtering
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 70:   RGCreate(PETSC_COMM_WORLD,&rg);
 71:   RGSetType(rg,RGELLIPSE);
 72:   radius = (1.0-PetscRealPart(target))/2.0;
 73:   RGEllipseSetParameters(rg,target+radius,radius,0.1);

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:                 Create the eigensolver and set various options
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 79:   EPSCreate(PETSC_COMM_WORLD,&eps);
 80:   EPSSetST(eps,st);
 81:   EPSSetRG(eps,rg);
 82:   EPSSetOperators(eps,A,NULL);
 83:   EPSSetProblemType(eps,EPS_NHEP);
 84:   EPSSetTolerances(eps,tol,PETSC_DEFAULT);
 85:   EPSSetTarget(eps,target);
 86:   EPSSetFromOptions(eps);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:                Solve the eigensystem and display solution
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   EPSSolve(eps);
 93:   EPSGetDimensions(eps,&nev,NULL,NULL);
 94:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
 95:   EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL);

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:                    Check file containing the eigenvalues
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100:   PetscOptionsGetString(NULL,NULL,"-checkfile",filename,sizeof(filename),&checkfile);
101:   if (checkfile) {
102:     EPSGetConverged(eps,&nconv);
103:     PetscMalloc1(nconv,&eigs);
104:     PetscViewerBinaryOpen(PETSC_COMM_WORLD,filename,FILE_MODE_READ,&viewer);
105:     PetscViewerBinaryRead(viewer,eigs,nconv,NULL,PETSC_COMPLEX);
106:     PetscViewerDestroy(&viewer);
107:     for (i=0;i<nconv;i++) {
108:       EPSGetEigenpair(eps,i,&kr,&ki,NULL,NULL);
109: #if defined(PETSC_USE_COMPLEX)
110:       eval = kr;
111: #else
112:       eval = PetscCMPLX(kr,ki);
113: #endif
114:       if (eval!=eigs[i]) SETERRQ(PETSC_COMM_WORLD,1,"Eigenvalues in the file do not match");
115:     }
116:     PetscFree(eigs);
117:   }

119:   EPSDestroy(&eps);
120:   STDestroy(&st);
121:   RGDestroy(&rg);
122:   MatDestroy(&A);
123: #else
124:   SETERRQ(PETSC_COMM_WORLD,1,"This example requires C99 complex numbers");
125: #endif
126:   SlepcFinalize();
127:   return ierr;
128: }

130: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
131: {
132:   const PetscReal cst = 0.5/(PetscReal)(m-1);
133:   PetscReal       pd,pu;
134:   PetscInt        Istart,Iend,i,j,jmax,ix=0;
135:   PetscErrorCode  ierr;

138:   MatGetOwnershipRange(A,&Istart,&Iend);
139:   for (i=1;i<=m;i++) {
140:     jmax = m-i+1;
141:     for (j=1;j<=jmax;j++) {
142:       ix = ix + 1;
143:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
144:       if (j!=jmax) {
145:         pd = cst*(PetscReal)(i+j-1);
146:         /* north */
147:         if (i==1) {
148:           MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES);
149:         } else {
150:           MatSetValue(A,ix-1,ix,pd,INSERT_VALUES);
151:         }
152:         /* east */
153:         if (j==1) {
154:           MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES);
155:         } else {
156:           MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES);
157:         }
158:       }
159:       /* south */
160:       pu = 0.5 - cst*(PetscReal)(i+j-3);
161:       if (j>1) {
162:         MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES);
163:       }
164:       /* west */
165:       if (i>1) {
166:         MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES);
167:       }
168:     }
169:   }
170:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
171:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
172:   return(0);
173: }

175: /*TEST

177:    test:
178:       suffix: 1
179:       args: -eps_nev 4 -eps_ncv 20 -eps_view_values binary:myvalues.bin -checkfile myvalues.bin
180:       output_file: output/test11_1.out
181:       requires: !single c99_complex

183: TEST*/