Actual source code: test1.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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;

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

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

 58:   NEPCreate(PETSC_COMM_WORLD,&nep);

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

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

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

 80:   NEPSetFromOptions(nep);

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:                       Solve the eigensystem
 84:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 86:   NEPSolve(nep);

 88:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 89:                     Display solution and clean up
 90:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

103:   NEPDestroy(&nep);
104:   MatDestroy(&F);
105:   MatDestroy(&J);
106:   SlepcFinalize();
107:   return ierr;
108: }

110: /* ------------------------------------------------------------------- */
111: /*
112:    FormFunction - Computes Function matrix  T(lambda)

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

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

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

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

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

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

170:   /*
171:      Assemble matrix
172:   */
173:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
174:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
175:   if (fun != B) {
176:     MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
177:     MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
178:   }
179:   return(0);
180: }

182: /* ------------------------------------------------------------------- */
183: /*
184:    FormJacobian - Computes Jacobian matrix  T'(lambda)

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

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

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

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

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

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

241:   /*
242:      Assemble matrix
243:   */
244:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
245:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
246:   return(0);
247: }

249: /*TEST

251:    test:
252:       suffix: 1
253:       args: -nep_type {{rii slp}} -nep_target 21 -terse
254:       requires: !single

256:    test:
257:       suffix: 2_cuda
258:       args: -nep_type {{rii slp}} -nep_target 21 -mat_type aijcusparse -terse
259:       requires: cuda
260:       output_file: output/test1_1.out

262: TEST*/