Actual source code: test16.c

slepc-3.21.2 2024-09-25
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[] = "Illustrates use of NEPSetEigenvalueComparison().\n\n"
 12:   "This is a simplified version of ex20.\n"
 13:   "The command line options are:\n"
 14:   "  -n <n>, where <n> = number of grid subdivisions.\n";

 16: /*
 17:    Solve 1-D PDE
 18:             -u'' = lambda*u
 19:    on [0,1] subject to
 20:             u(0)=0, u'(1)=u(1)*lambda*kappa/(kappa-lambda)
 21: */

 23: #include <slepcnep.h>

 25: /*
 26:    User-defined routines
 27: */
 28: PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*);
 29: PetscErrorCode FormJacobian(NEP,PetscScalar,Mat,void*);
 30: PetscErrorCode MyEigenSort(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);

 32: /*
 33:    User-defined application context
 34: */
 35: typedef struct {
 36:   PetscScalar kappa;   /* ratio between stiffness of spring and attached mass */
 37:   PetscReal   h;       /* mesh spacing */
 38: } ApplicationCtx;

 40: int main(int argc,char **argv)
 41: {
 42:   NEP            nep;             /* nonlinear eigensolver context */
 43:   Mat            F,J;             /* Function and Jacobian matrices */
 44:   ApplicationCtx ctx;             /* user-defined context */
 45:   PetscScalar    target;
 46:   RG             rg;
 47:   PetscInt       n=128;
 48:   PetscBool      terse;

 50:   PetscFunctionBeginUser;
 51:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 52:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 53:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
 54:   ctx.h = 1.0/(PetscReal)n;
 55:   ctx.kappa = 1.0;

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:                Prepare nonlinear eigensolver context
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 61:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));

 63:   PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
 64:   PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
 65:   PetscCall(MatSetFromOptions(F));
 66:   PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
 67:   PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
 68:   PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));

 70:   PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
 71:   PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
 72:   PetscCall(MatSetFromOptions(J));
 73:   PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
 74:   PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
 75:   PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));

 77:   PetscCall(NEPSetType(nep,NEPNLEIGS));
 78:   PetscCall(NEPGetRG(nep,&rg));
 79:   PetscCall(RGSetType(rg,RGINTERVAL));
 80: #if defined(PETSC_USE_COMPLEX)
 81:   PetscCall(RGIntervalSetEndpoints(rg,2.0,400.0,-0.001,0.001));
 82: #else
 83:   PetscCall(RGIntervalSetEndpoints(rg,2.0,400.0,0,0));
 84: #endif
 85:   PetscCall(NEPSetTarget(nep,25.0));
 86:   PetscCall(NEPSetEigenvalueComparison(nep,MyEigenSort,&target));
 87:   PetscCall(NEPSetTolerances(nep,PETSC_SMALL,PETSC_DEFAULT));
 88:   PetscCall(NEPSetFromOptions(nep));
 89:   PetscCall(NEPGetTarget(nep,&target));

 91:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 92:               Solve the eigensystem and display the solution
 93:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 95:   PetscCall(NEPSolve(nep));

 97:   /* show detailed info unless -terse option is given by user */
 98:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
 99:   if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
100:   else {
101:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
102:     PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
103:     PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
104:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
105:   }

107:   PetscCall(NEPDestroy(&nep));
108:   PetscCall(MatDestroy(&F));
109:   PetscCall(MatDestroy(&J));
110:   PetscCall(SlepcFinalize());
111:   return 0;
112: }

114: /* ------------------------------------------------------------------- */
115: /*
116:    FormFunction - Computes Function matrix  T(lambda)

118:    Input Parameters:
119: .  nep    - the NEP context
120: .  lambda - the scalar argument
121: .  ctx    - optional user-defined context, as set by NEPSetFunction()

123:    Output Parameters:
124: .  fun - Function matrix
125: .  B   - optionally different preconditioning matrix
126: */
127: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
128: {
129:   ApplicationCtx *user = (ApplicationCtx*)ctx;
130:   PetscScalar    A[3],c,d;
131:   PetscReal      h;
132:   PetscInt       i,n,j[3],Istart,Iend;
133:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

135:   PetscFunctionBeginUser;
136:   /*
137:      Compute Function entries and insert into matrix
138:   */
139:   PetscCall(MatGetSize(fun,&n,NULL));
140:   PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
141:   if (Istart==0) FirstBlock=PETSC_TRUE;
142:   if (Iend==n) LastBlock=PETSC_TRUE;
143:   h = user->h;
144:   c = user->kappa/(lambda-user->kappa);
145:   d = n;

147:   /*
148:      Interior grid points
149:   */
150:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
151:     j[0] = i-1; j[1] = i; j[2] = i+1;
152:     A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
153:     PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES));
154:   }

156:   /*
157:      Boundary points
158:   */
159:   if (FirstBlock) {
160:     i = 0;
161:     j[0] = 0; j[1] = 1;
162:     A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
163:     PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
164:   }

166:   if (LastBlock) {
167:     i = n-1;
168:     j[0] = n-2; j[1] = n-1;
169:     A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
170:     PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
171:   }

173:   /*
174:      Assemble matrix
175:   */
176:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
177:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
178:   if (fun != B) {
179:     PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
180:     PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
181:   }
182:   PetscFunctionReturn(PETSC_SUCCESS);
183: }

185: /* ------------------------------------------------------------------- */
186: /*
187:    FormJacobian - Computes Jacobian matrix  T'(lambda)

189:    Input Parameters:
190: .  nep    - the NEP context
191: .  lambda - the scalar argument
192: .  ctx    - optional user-defined context, as set by NEPSetJacobian()

194:    Output Parameters:
195: .  jac - Jacobian matrix
196: .  B   - optionally different preconditioning matrix
197: */
198: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
199: {
200:   ApplicationCtx *user = (ApplicationCtx*)ctx;
201:   PetscScalar    A[3],c;
202:   PetscReal      h;
203:   PetscInt       i,n,j[3],Istart,Iend;
204:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

206:   PetscFunctionBeginUser;
207:   /*
208:      Compute Jacobian entries and insert into matrix
209:   */
210:   PetscCall(MatGetSize(jac,&n,NULL));
211:   PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
212:   if (Istart==0) FirstBlock=PETSC_TRUE;
213:   if (Iend==n) LastBlock=PETSC_TRUE;
214:   h = user->h;
215:   c = user->kappa/(lambda-user->kappa);

217:   /*
218:      Interior grid points
219:   */
220:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
221:     j[0] = i-1; j[1] = i; j[2] = i+1;
222:     A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
223:     PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
224:   }

226:   /*
227:      Boundary points
228:   */
229:   if (FirstBlock) {
230:     i = 0;
231:     j[0] = 0; j[1] = 1;
232:     A[0] = -2.0*h/3.0; A[1] = -h/6.0;
233:     PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
234:   }

236:   if (LastBlock) {
237:     i = n-1;
238:     j[0] = n-2; j[1] = n-1;
239:     A[0] = -h/6.0; A[1] = -h/3.0-c*c;
240:     PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
241:   }

243:   /*
244:      Assemble matrix
245:   */
246:   PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
247:   PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
248:   PetscFunctionReturn(PETSC_SUCCESS);
249: }

251: /*
252:     Function for user-defined eigenvalue ordering criterion.

254:     Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
255:     one of them as the preferred one according to the criterion.
256:     In this example, eigenvalues are sorted with respect to the target,
257:     but those on the right of the target are preferred.
258: */
259: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
260: {
261:   PetscReal   a,b;
262:   PetscScalar target = *(PetscScalar*)ctx;

264:   PetscFunctionBeginUser;
265:   if (PetscRealPart(ar-target)<0.0 && PetscRealPart(br-target)>0.0) *r = 1;
266:   else {
267:     a = SlepcAbsEigenvalue(ar-target,ai);
268:     b = SlepcAbsEigenvalue(br-target,bi);
269:     if (a>b) *r = 1;
270:     else if (a<b) *r = -1;
271:     else *r = 0;
272:   }
273:   PetscFunctionReturn(PETSC_SUCCESS);
274: }

276: /*TEST

278:    test:
279:       suffix: 1
280:       args: -nep_nev 4 -nep_ncv 8 -terse
281:       requires: double !complex

283: TEST*/