Actual source code: test11.c

slepc-3.22.1 2024-10-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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[] = "Solves the same problem as in ex5, but with a user-defined sorting criterion."
 12:   "It is a standard nonsymmetric eigenproblem with real eigenvalues and the rightmost eigenvalue is known to be 1.\n"
 13:   "This example illustrates how the user can set a custom spectrum selection.\n\n"
 14:   "The command line options are:\n"
 15:   "  -m <m>, where <m> = number of grid subdivisions in each dimension.\n\n";

 17: #include <slepceps.h>

 19: /*
 20:    User-defined routines
 21: */

 23: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx);
 24: PetscErrorCode MatMarkovModel(PetscInt m,Mat A);

 26: int main(int argc,char **argv)
 27: {
 28:   Vec            v0;              /* initial vector */
 29:   Mat            A;               /* operator matrix */
 30:   EPS            eps;             /* eigenproblem solver context */
 31:   ST             st;              /* spectral transformation associated */
 32:   PetscReal      tol=PETSC_SMALL;
 33:   PetscScalar    target=0.5;
 34:   PetscInt       N,m=15,nev;
 35:   char           str[50];

 37:   PetscFunctionBeginUser;
 38:   PetscCall(SlepcInitialize(&argc,&argv,NULL,help));

 40:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
 41:   N = m*(m+1)/2;
 42:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nMarkov Model, N=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n",N,m));
 43:   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-target",&target,NULL));
 44:   PetscCall(SlepcSNPrintfScalar(str,sizeof(str),target,PETSC_FALSE));
 45:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Searching closest eigenvalues to the right of %s.\n\n",str));

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

 51:   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
 52:   PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,N,N));
 53:   PetscCall(MatSetFromOptions(A));
 54:   PetscCall(MatMarkovModel(m,A));

 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 57:                 Create the eigensolver and set various options
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 60:   /*
 61:      Create eigensolver context
 62:   */
 63:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));

 65:   /*
 66:      Set operators. In this case, it is a standard eigenvalue problem
 67:   */
 68:   PetscCall(EPSSetOperators(eps,A,NULL));
 69:   PetscCall(EPSSetProblemType(eps,EPS_NHEP));
 70:   PetscCall(EPSSetTolerances(eps,tol,PETSC_CURRENT));

 72:   /*
 73:      Set the custom comparing routine in order to obtain the eigenvalues
 74:      closest to the target on the right only
 75:   */
 76:   PetscCall(EPSSetEigenvalueComparison(eps,MyEigenSort,&target));

 78:   /*
 79:      Set solver parameters at runtime
 80:   */
 81:   PetscCall(EPSSetFromOptions(eps));

 83:   /*
 84:      Set the preconditioner based on A - target * I
 85:   */
 86:   PetscCall(EPSGetST(eps,&st));
 87:   PetscCall(STSetShift(st,target));

 89:   /*
 90:      Set the initial vector. This is optional, if not done the initial
 91:      vector is set to random values
 92:   */
 93:   PetscCall(MatCreateVecs(A,&v0,NULL));
 94:   PetscCall(VecSet(v0,1.0));
 95:   PetscCall(EPSSetInitialSpace(eps,1,&v0));

 97:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 98:                       Solve the eigensystem
 99:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

101:   PetscCall(EPSSolve(eps));
102:   PetscCall(EPSGetDimensions(eps,&nev,NULL,NULL));
103:   PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));

105:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106:                     Display solution and clean up
107:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

109:   PetscCall(EPSErrorView(eps,EPS_ERROR_RELATIVE,NULL));
110:   PetscCall(EPSDestroy(&eps));
111:   PetscCall(MatDestroy(&A));
112:   PetscCall(VecDestroy(&v0));
113:   PetscCall(SlepcFinalize());
114:   return 0;
115: }

117: PetscErrorCode MatMarkovModel(PetscInt m,Mat A)
118: {
119:   const PetscReal cst = 0.5/(PetscReal)(m-1);
120:   PetscReal       pd,pu;
121:   PetscInt        Istart,Iend,i,j,jmax,ix=0;

123:   PetscFunctionBeginUser;
124:   PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
125:   for (i=1;i<=m;i++) {
126:     jmax = m-i+1;
127:     for (j=1;j<=jmax;j++) {
128:       ix = ix + 1;
129:       if (ix-1<Istart || ix>Iend) continue;  /* compute only owned rows */
130:       if (j!=jmax) {
131:         pd = cst*(PetscReal)(i+j-1);
132:         /* north */
133:         if (i==1) PetscCall(MatSetValue(A,ix-1,ix,2*pd,INSERT_VALUES));
134:         else PetscCall(MatSetValue(A,ix-1,ix,pd,INSERT_VALUES));
135:         /* east */
136:         if (j==1) PetscCall(MatSetValue(A,ix-1,ix+jmax-1,2*pd,INSERT_VALUES));
137:         else PetscCall(MatSetValue(A,ix-1,ix+jmax-1,pd,INSERT_VALUES));
138:       }
139:       /* south */
140:       pu = 0.5 - cst*(PetscReal)(i+j-3);
141:       if (j>1) PetscCall(MatSetValue(A,ix-1,ix-2,pu,INSERT_VALUES));
142:       /* west */
143:       if (i>1) PetscCall(MatSetValue(A,ix-1,ix-jmax-2,pu,INSERT_VALUES));
144:     }
145:   }
146:   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
147:   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
148:   PetscFunctionReturn(PETSC_SUCCESS);
149: }

151: /*
152:     Function for user-defined eigenvalue ordering criterion.

154:     Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
155:     one of them as the preferred one according to the criterion.
156:     In this example, the preferred value is the one closest to the target,
157:     but on the right side.
158: */
159: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
160: {
161:   PetscScalar target = *(PetscScalar*)ctx;
162:   PetscReal   da,db;
163:   PetscBool   aisright,bisright;

165:   PetscFunctionBeginUser;
166:   if (PetscRealPart(target) < PetscRealPart(ar)) aisright = PETSC_TRUE;
167:   else aisright = PETSC_FALSE;
168:   if (PetscRealPart(target) < PetscRealPart(br)) bisright = PETSC_TRUE;
169:   else bisright = PETSC_FALSE;
170:   if (aisright == bisright) {
171:     /* both are on the same side of the target */
172:     da = SlepcAbsEigenvalue(ar-target,ai);
173:     db = SlepcAbsEigenvalue(br-target,bi);
174:     if (da < db) *r = -1;
175:     else if (da > db) *r = 1;
176:     else *r = 0;
177:   } else if (aisright && !bisright) *r = -1; /* 'a' is on the right */
178:   else *r = 1;  /* 'b' is on the right */
179:   PetscFunctionReturn(PETSC_SUCCESS);
180: }

182: /*TEST

184:    testset:
185:       args: -eps_nev 4
186:       requires: !single
187:       output_file: output/test11_1.out
188:       test:
189:          suffix: 1
190:          args: -eps_type {{krylovschur arnoldi lapack}} -st_type sinvert
191:       test:
192:          suffix: 1_ks_cayley
193:          args: -st_type cayley -st_cayley_antishift 1

195:    test:
196:       suffix: 2
197:       args: -target 0.77 -eps_type gd -eps_nev 4 -eps_tol 1e-7 -eps_gd_krylov_start -eps_gd_blocksize 3
198:       requires: double
199:       filter: sed -e "s/[+-]0\.0*i//g"

201: TEST*/