Actual source code: test7.c

slepc-3.16.3 2022-04-11
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[] = "Test the NLEIGS solver with shell matrices.\n\n"
 12:   "This is based on ex27.\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = matrix dimension.\n"
 15:   "  -split <0/1>, to select the split form in the problem definition (enabled by default).\n";

 17: /*
 18:    Solve T(lambda)x=0 using NLEIGS solver
 19:       with T(lambda) = -D+sqrt(lambda)*I
 20:       where D is the Laplacian operator in 1 dimension
 21:       and with the interpolation interval [.01,16]
 22: */

 24: #include <slepcnep.h>

 26: /* User-defined routines */
 27: PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*);
 28: PetscErrorCode ComputeSingularities(NEP,PetscInt*,PetscScalar*,void*);
 29: PetscErrorCode MatMult_A0(Mat,Vec,Vec);
 30: PetscErrorCode MatGetDiagonal_A0(Mat,Vec);
 31: PetscErrorCode MatDuplicate_A0(Mat,MatDuplicateOption,Mat*);
 32: PetscErrorCode MatMult_A1(Mat,Vec,Vec);
 33: PetscErrorCode MatGetDiagonal_A1(Mat,Vec);
 34: PetscErrorCode MatDuplicate_A1(Mat,MatDuplicateOption,Mat*);
 35: PetscErrorCode MatMult_F(Mat,Vec,Vec);
 36: PetscErrorCode MatGetDiagonal_F(Mat,Vec);
 37: PetscErrorCode MatDuplicate_F(Mat,MatDuplicateOption,Mat*);
 38: PetscErrorCode MatDestroy_F(Mat);

 40: typedef struct {
 41:   PetscScalar t;  /* square root of lambda */
 42: } MatCtx;

 44: int main(int argc,char **argv)
 45: {
 46:   NEP            nep;
 47:   KSP            *ksp;
 48:   PC             pc;
 49:   Mat            F,A[2];
 50:   NEPType        type;
 51:   PetscInt       i,n=100,nev,its,nsolve;
 52:   PetscReal      keep,tol=PETSC_SQRT_MACHINE_EPSILON/10;
 54:   RG             rg;
 55:   FN             f[2];
 56:   PetscBool      terse,flg,lock,split=PETSC_TRUE;
 57:   PetscScalar    coeffs;
 58:   MatCtx         *ctx;

 60:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 61:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 62:   PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL);
 63:   PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%D%s\n\n",n,split?" (in split form)":"");

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:      Create NEP context, configure NLEIGS with appropriate options
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 69:   NEPCreate(PETSC_COMM_WORLD,&nep);
 70:   NEPSetType(nep,NEPNLEIGS);
 71:   NEPNLEIGSSetSingularitiesFunction(nep,ComputeSingularities,NULL);
 72:   NEPGetRG(nep,&rg);
 73:   RGSetType(rg,RGINTERVAL);
 74: #if defined(PETSC_USE_COMPLEX)
 75:   RGIntervalSetEndpoints(rg,0.01,16.0,-0.001,0.001);
 76: #else
 77:   RGIntervalSetEndpoints(rg,0.01,16.0,0,0);
 78: #endif
 79:   NEPSetTarget(nep,1.1);
 80:   NEPNLEIGSGetKSPs(nep,&nsolve,&ksp);
 81:   for (i=0;i<nsolve;i++) {
 82:    KSPSetType(ksp[i],KSPBICG);
 83:    KSPGetPC(ksp[i],&pc);
 84:    PCSetType(pc,PCJACOBI);
 85:    KSPSetTolerances(ksp[i],tol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
 86:   }

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:      Define the nonlinear problem
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 92:   if (split) {
 93:     /* Create matrix A0 (tridiagonal) */
 94:     MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,NULL,&A[0]);
 95:     MatShellSetOperation(A[0],MATOP_MULT,(void(*)(void))MatMult_A0);
 96:     MatShellSetOperation(A[0],MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_A0);
 97:     MatShellSetOperation(A[0],MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_A0);
 98:     MatShellSetOperation(A[0],MATOP_DUPLICATE,(void(*)(void))MatDuplicate_A0);

100:     /* Create matrix A0 (identity) */
101:     MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,NULL,&A[1]);
102:     MatShellSetOperation(A[1],MATOP_MULT,(void(*)(void))MatMult_A1);
103:     MatShellSetOperation(A[1],MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_A1);
104:     MatShellSetOperation(A[1],MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_A1);
105:     MatShellSetOperation(A[1],MATOP_DUPLICATE,(void(*)(void))MatDuplicate_A1);

107:     /* Define functions for the split form */
108:     FNCreate(PETSC_COMM_WORLD,&f[0]);
109:     FNSetType(f[0],FNRATIONAL);
110:     coeffs = 1.0;
111:     FNRationalSetNumerator(f[0],1,&coeffs);
112:     FNCreate(PETSC_COMM_WORLD,&f[1]);
113:     FNSetType(f[1],FNSQRT);
114:     NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN);
115:   } else {
116:     /* Callback form: create shell matrix for F=A0+sqrt(lambda)*A1  */
117:     PetscNew(&ctx);
118:     MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,(void*)ctx,&F);
119:     MatShellSetOperation(F,MATOP_MULT,(void(*)(void))MatMult_F);
120:     MatShellSetOperation(F,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_F);
121:     MatShellSetOperation(F,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_F);
122:     MatShellSetOperation(F,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_F);
123:     MatShellSetOperation(F,MATOP_DESTROY,(void(*)(void))MatDestroy_F);
124:     /* Set Function evaluation routine */
125:     NEPSetFunction(nep,F,F,FormFunction,NULL);
126:   }

128:   /* Set solver parameters at runtime */
129:   NEPSetFromOptions(nep);

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:                       Solve the eigensystem
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134:   NEPSolve(nep);
135:   NEPGetType(nep,&type);
136:   PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type);
137:   NEPGetDimensions(nep,&nev,NULL,NULL);
138:   PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %D\n",nev);
139:   PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
140:   if (flg) {
141:     NEPNLEIGSGetRestart(nep,&keep);
142:     NEPNLEIGSGetLocking(nep,&lock);
143:     NEPNLEIGSGetInterpolation(nep,&tol,&its);
144:     PetscPrintf(PETSC_COMM_WORLD," Restart factor is %3.2f",(double)keep);
145:     if (lock) { PetscPrintf(PETSC_COMM_WORLD," (locking activated)"); }
146:     PetscPrintf(PETSC_COMM_WORLD,"\n Divided diferences with tol=%6.2g maxit=%D\n",(double)tol,its);
147:     PetscPrintf(PETSC_COMM_WORLD,"\n");
148:   }

150:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151:                     Display solution and clean up
152:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

154:   /* show detailed info unless -terse option is given by user */
155:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
156:   if (terse) {
157:     NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
158:   } else {
159:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
160:     NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
161:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
162:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
163:   }
164:   NEPDestroy(&nep);
165:   if (split) {
166:     MatDestroy(&A[0]);
167:     MatDestroy(&A[1]);
168:     FNDestroy(&f[0]);
169:     FNDestroy(&f[1]);
170:   } else {
171:     MatDestroy(&F);
172:   }
173:   SlepcFinalize();
174:   return ierr;
175: }

177: /*
178:    FormFunction - Computes Function matrix  T(lambda)
179: */
180: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
181: {
183:   MatCtx         *ctxF;

186:   MatShellGetContext(fun,&ctxF);
187:   ctxF->t = PetscSqrtScalar(lambda);
188:   return(0);
189: }

191: /*
192:    ComputeSingularities - Computes maxnp points (at most) in the complex plane where
193:    the function T(.) is not analytic.

195:    In this case, we discretize the singularity region (-inf,0)~(-10e+6,-10e-6)
196: */
197: PetscErrorCode ComputeSingularities(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *pt)
198: {
199:   PetscReal h;
200:   PetscInt  i;

203:   h = 12.0/(*maxnp-1);
204:   xi[0] = -1e-6; xi[*maxnp-1] = -1e+6;
205:   for (i=1;i<*maxnp-1;i++) xi[i] = -PetscPowReal(10,-6+h*i);
206:   return(0);
207: }

209: /* -------------------------------- A0 ----------------------------------- */

211: PetscErrorCode MatMult_A0(Mat A,Vec x,Vec y)
212: {
213:   PetscErrorCode    ierr;
214:   PetscInt          i,n;
215:   PetscMPIInt       rank,size,next,prev;
216:   const PetscScalar *px;
217:   PetscScalar       *py,upper=0.0,lower=0.0;
218:   MPI_Comm          comm;

221:   PetscObjectGetComm((PetscObject)A,&comm);
222:   MPI_Comm_size(comm,&size);
223:   MPI_Comm_rank(comm,&rank);
224:   next = rank==size-1? MPI_PROC_NULL: rank+1;
225:   prev = rank==0? MPI_PROC_NULL: rank-1;

227:   VecGetArrayRead(x,&px);
228:   VecGetArray(y,&py);
229:   VecGetLocalSize(x,&n);

231:   MPI_Sendrecv(px,1,MPIU_SCALAR,prev,0,&lower,1,MPIU_SCALAR,next,0,comm,MPI_STATUS_IGNORE);
232:   MPI_Sendrecv(px+n-1,1,MPIU_SCALAR,next,0,&upper,1,MPIU_SCALAR,prev,0,comm,MPI_STATUS_IGNORE);

234:   py[0] = upper-2.0*px[0]+px[1];
235:   for (i=1;i<n-1;i++) py[i] = px[i-1]-2.0*px[i]+px[i+1];
236:   py[n-1] = px[n-2]-2.0*px[n-1]+lower;
237:   VecRestoreArrayRead(x,&px);
238:   VecRestoreArray(y,&py);
239:   return(0);
240: }

242: PetscErrorCode MatGetDiagonal_A0(Mat A,Vec diag)
243: {

247:   VecSet(diag,-2.0);
248:   return(0);
249: }

251: PetscErrorCode MatDuplicate_A0(Mat A,MatDuplicateOption op,Mat *B)
252: {
253:   PetscInt       m,n,M,N;
254:   MPI_Comm       comm;

258:   MatGetSize(A,&M,&N);
259:   MatGetLocalSize(A,&m,&n);
260:   PetscObjectGetComm((PetscObject)A,&comm);
261:   MatCreateShell(comm,m,n,M,N,NULL,B);
262:   MatShellSetOperation(*B,MATOP_MULT,(void(*)(void))MatMult_A0);
263:   MatShellSetOperation(*B,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_A0);
264:   MatShellSetOperation(*B,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_A0);
265:   MatShellSetOperation(*B,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_A0);
266:   return(0);
267: }

269: /* -------------------------------- A1 ----------------------------------- */

271: PetscErrorCode MatMult_A1(Mat A,Vec x,Vec y)
272: {

276:   VecCopy(x,y);
277:   return(0);
278: }

280: PetscErrorCode MatGetDiagonal_A1(Mat A,Vec diag)
281: {

285:   VecSet(diag,1.0);
286:   return(0);
287: }

289: PetscErrorCode MatDuplicate_A1(Mat A,MatDuplicateOption op,Mat *B)
290: {
291:   PetscInt       m,n,M,N;
292:   MPI_Comm       comm;

296:   MatGetSize(A,&M,&N);
297:   MatGetLocalSize(A,&m,&n);
298:   PetscObjectGetComm((PetscObject)A,&comm);
299:   MatCreateShell(comm,m,n,M,N,NULL,B);
300:   MatShellSetOperation(*B,MATOP_MULT,(void(*)(void))MatMult_A1);
301:   MatShellSetOperation(*B,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_A1);
302:   MatShellSetOperation(*B,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_A1);
303:   MatShellSetOperation(*B,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_A1);
304:   return(0);
305: }

307: /* -------------------------------- F ----------------------------------- */

309: PetscErrorCode MatMult_F(Mat A,Vec x,Vec y)
310: {
311:   PetscErrorCode    ierr;
312:   PetscInt          i,n;
313:   PetscMPIInt       rank,size,next,prev;
314:   const PetscScalar *px;
315:   PetscScalar       *py,d,upper=0.0,lower=0.0;
316:   MatCtx            *ctx;
317:   MPI_Comm          comm;

320:   PetscObjectGetComm((PetscObject)A,&comm);
321:   MPI_Comm_size(comm,&size);
322:   MPI_Comm_rank(comm,&rank);
323:   next = rank==size-1? MPI_PROC_NULL: rank+1;
324:   prev = rank==0? MPI_PROC_NULL: rank-1;

326:   MatShellGetContext(A,&ctx);
327:   VecGetArrayRead(x,&px);
328:   VecGetArray(y,&py);
329:   VecGetLocalSize(x,&n);

331:   MPI_Sendrecv(px,1,MPIU_SCALAR,prev,0,&lower,1,MPIU_SCALAR,next,0,comm,MPI_STATUS_IGNORE);
332:   MPI_Sendrecv(px+n-1,1,MPIU_SCALAR,next,0,&upper,1,MPIU_SCALAR,prev,0,comm,MPI_STATUS_IGNORE);

334:   d = -2.0+ctx->t;
335:   py[0] = upper+d*px[0]+px[1];
336:   for (i=1;i<n-1;i++) py[i] = px[i-1]+d*px[i]+px[i+1];
337:   py[n-1] = px[n-2]+d*px[n-1]+lower;
338:   VecRestoreArrayRead(x,&px);
339:   VecRestoreArray(y,&py);
340:   return(0);
341: }

343: PetscErrorCode MatGetDiagonal_F(Mat A,Vec diag)
344: {
346:   MatCtx         *ctx;

349:   MatShellGetContext(A,&ctx);
350:   VecSet(diag,-2.0+ctx->t);
351:   return(0);
352: }

354: PetscErrorCode MatDuplicate_F(Mat A,MatDuplicateOption op,Mat *B)
355: {
356:   MatCtx         *actx,*bctx;
357:   PetscInt       m,n,M,N;
358:   MPI_Comm       comm;

362:   MatShellGetContext(A,&actx);
363:   MatGetSize(A,&M,&N);
364:   MatGetLocalSize(A,&m,&n);
365:   PetscNew(&bctx);
366:   bctx->t = actx->t;
367:   PetscObjectGetComm((PetscObject)A,&comm);
368:   MatCreateShell(comm,m,n,M,N,(void*)bctx,B);
369:   MatShellSetOperation(*B,MATOP_MULT,(void(*)(void))MatMult_F);
370:   MatShellSetOperation(*B,MATOP_MULT_TRANSPOSE,(void(*)(void))MatMult_F);
371:   MatShellSetOperation(*B,MATOP_GET_DIAGONAL,(void(*)(void))MatGetDiagonal_F);
372:   MatShellSetOperation(*B,MATOP_DUPLICATE,(void(*)(void))MatDuplicate_F);
373:   MatShellSetOperation(*B,MATOP_DESTROY,(void(*)(void))MatDestroy_F);
374:   return(0);
375: }

377: PetscErrorCode MatDestroy_F(Mat A)
378: {
379:   MatCtx         *ctx;

383:   MatShellGetContext(A,&ctx);
384:   PetscFree(ctx);
385:   return(0);
386: }

388: /*TEST

390:    testset:
391:       nsize: {{1 2}}
392:       args: -nep_nev 3 -nep_tol 1e-8 -terse
393:       filter: sed -e "s/[+-]0\.0*i//g"
394:       requires: !single
395:       test:
396:          suffix: 1
397:          args: -nep_nleigs_locking 0 -nep_nleigs_interpolation_degree 90 -nep_nleigs_interpolation_tol 1e-8 -nep_nleigs_restart 0.4
398:       test:
399:          suffix: 2
400:          args: -split 0

402: TEST*/