Actual source code: test1.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[] = "Simple 1-D nonlinear eigenproblem.\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*);

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

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

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

 54:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 55:                Prepare nonlinear eigensolver context
 56:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

 60:   /*
 61:      Create Function and Jacobian matrices; set evaluation routines
 62:   */

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

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

 78:   PetscCall(NEPSetFromOptions(nep));

 80:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 81:                       Solve the eigensystem
 82:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 84:   PetscCall(NEPSolve(nep));

 86:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 87:                     Display solution and clean up
 88:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 90:   /* show detailed info unless -terse option is given by user */
 91:   PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
 92:   if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
 93:   else {
 94:     PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
 95:     PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
 96:     PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
 97:     PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
 98:   }

100:   PetscCall(NEPDestroy(&nep));
101:   PetscCall(MatDestroy(&F));
102:   PetscCall(MatDestroy(&J));
103:   PetscCall(SlepcFinalize());
104:   return 0;
105: }

107: /* ------------------------------------------------------------------- */
108: /*
109:    FormFunction - Computes Function matrix  T(lambda)

111:    Input Parameters:
112: .  nep    - the NEP context
113: .  lambda - the scalar argument
114: .  ctx    - optional user-defined context, as set by NEPSetFunction()

116:    Output Parameters:
117: .  fun - Function matrix
118: .  B   - optionally different preconditioning matrix
119: */
120: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
121: {
122:   ApplicationCtx *user = (ApplicationCtx*)ctx;
123:   PetscScalar    A[3],c,d;
124:   PetscReal      h;
125:   PetscInt       i,n,j[3],Istart,Iend;
126:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

128:   PetscFunctionBeginUser;
129:   /*
130:      Compute Function entries and insert into matrix
131:   */
132:   PetscCall(MatGetSize(fun,&n,NULL));
133:   PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
134:   if (Istart==0) FirstBlock=PETSC_TRUE;
135:   if (Iend==n) LastBlock=PETSC_TRUE;
136:   h = user->h;
137:   c = user->kappa/(lambda-user->kappa);
138:   d = n;

140:   /*
141:      Interior grid points
142:   */
143:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
144:     j[0] = i-1; j[1] = i; j[2] = i+1;
145:     A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
146:     PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES));
147:   }

149:   /*
150:      Boundary points
151:   */
152:   if (FirstBlock) {
153:     i = 0;
154:     j[0] = 0; j[1] = 1;
155:     A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
156:     PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
157:   }

159:   if (LastBlock) {
160:     i = n-1;
161:     j[0] = n-2; j[1] = n-1;
162:     A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
163:     PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
164:   }

166:   /*
167:      Assemble matrix
168:   */
169:   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
170:   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
171:   if (fun != B) {
172:     PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
173:     PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
174:   }
175:   PetscFunctionReturn(PETSC_SUCCESS);
176: }

178: /* ------------------------------------------------------------------- */
179: /*
180:    FormJacobian - Computes Jacobian matrix  T'(lambda)

182:    Input Parameters:
183: .  nep    - the NEP context
184: .  lambda - the scalar argument
185: .  ctx    - optional user-defined context, as set by NEPSetJacobian()

187:    Output Parameters:
188: .  jac - Jacobian matrix
189: .  B   - optionally different preconditioning matrix
190: */
191: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
192: {
193:   ApplicationCtx *user = (ApplicationCtx*)ctx;
194:   PetscScalar    A[3],c;
195:   PetscReal      h;
196:   PetscInt       i,n,j[3],Istart,Iend;
197:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

199:   PetscFunctionBeginUser;
200:   /*
201:      Compute Jacobian entries and insert into matrix
202:   */
203:   PetscCall(MatGetSize(jac,&n,NULL));
204:   PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
205:   if (Istart==0) FirstBlock=PETSC_TRUE;
206:   if (Iend==n) LastBlock=PETSC_TRUE;
207:   h = user->h;
208:   c = user->kappa/(lambda-user->kappa);

210:   /*
211:      Interior grid points
212:   */
213:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
214:     j[0] = i-1; j[1] = i; j[2] = i+1;
215:     A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
216:     PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
217:   }

219:   /*
220:      Boundary points
221:   */
222:   if (FirstBlock) {
223:     i = 0;
224:     j[0] = 0; j[1] = 1;
225:     A[0] = -2.0*h/3.0; A[1] = -h/6.0;
226:     PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
227:   }

229:   if (LastBlock) {
230:     i = n-1;
231:     j[0] = n-2; j[1] = n-1;
232:     A[0] = -h/6.0; A[1] = -h/3.0-c*c;
233:     PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
234:   }

236:   /*
237:      Assemble matrix
238:   */
239:   PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
240:   PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
241:   PetscFunctionReturn(PETSC_SUCCESS);
242: }

244: /*TEST

246:    testset:
247:       args: -nep_type {{rii slp}} -nep_target 21 -terse -nep_view_vectors ::ascii_info
248:       filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/" | sed -e "s/[+-]0\.0*i//g"
249:       test:
250:          suffix: 1_real
251:          requires: !single !complex
252:       test:
253:          suffix: 1
254:          requires: !single complex

256:    test:
257:       suffix: 2_cuda
258:       args: -nep_type {{rii slp}} -nep_target 21 -mat_type aijcusparse -terse
259:       requires: cuda !single
260:       filter: sed -e "s/[+-]0\.0*i//"
261:       output_file: output/test3_1.out

263:    testset:
264:       args: -nep_type slp -nep_two_sided -nep_target 21 -terse -nep_view_vectors ::ascii_info
265:       filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/"
266:       test:
267:          suffix: 3_real
268:          requires: !single !complex
269:       test:
270:          suffix: 3
271:          requires: !single complex

273: TEST*/