Actual source code: test1.c
slepc-3.17.1 2022-04-11
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: SlepcInitialize(&argc,&argv,(char*)0,help);
48: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
49: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n);
50: ctx.h = 1.0/(PetscReal)n;
51: ctx.kappa = 1.0;
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Prepare nonlinear eigensolver context
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: NEPCreate(PETSC_COMM_WORLD,&nep);
59: /*
60: Create Function and Jacobian matrices; set evaluation routines
61: */
63: MatCreate(PETSC_COMM_WORLD,&F);
64: MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
65: MatSetFromOptions(F);
66: MatSeqAIJSetPreallocation(F,3,NULL);
67: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
68: MatSetUp(F);
69: NEPSetFunction(nep,F,F,FormFunction,&ctx);
71: MatCreate(PETSC_COMM_WORLD,&J);
72: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
73: MatSetFromOptions(J);
74: MatSeqAIJSetPreallocation(J,3,NULL);
75: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
76: MatSetUp(J);
77: NEPSetJacobian(nep,J,FormJacobian,&ctx);
79: NEPSetFromOptions(nep);
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Solve the eigensystem
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: NEPSolve(nep);
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Display solution and clean up
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91: /* show detailed info unless -terse option is given by user */
92: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
93: if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
94: else {
95: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
96: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
97: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
98: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
99: }
101: NEPDestroy(&nep);
102: MatDestroy(&F);
103: MatDestroy(&J);
104: SlepcFinalize();
105: return 0;
106: }
108: /* ------------------------------------------------------------------- */
109: /*
110: FormFunction - Computes Function matrix T(lambda)
112: Input Parameters:
113: . nep - the NEP context
114: . lambda - the scalar argument
115: . ctx - optional user-defined context, as set by NEPSetFunction()
117: Output Parameters:
118: . fun - Function matrix
119: . B - optionally different preconditioning matrix
120: */
121: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
122: {
123: ApplicationCtx *user = (ApplicationCtx*)ctx;
124: PetscScalar A[3],c,d;
125: PetscReal h;
126: PetscInt i,n,j[3],Istart,Iend;
127: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
130: /*
131: Compute Function entries and insert into matrix
132: */
133: MatGetSize(fun,&n,NULL);
134: MatGetOwnershipRange(fun,&Istart,&Iend);
135: if (Istart==0) FirstBlock=PETSC_TRUE;
136: if (Iend==n) LastBlock=PETSC_TRUE;
137: h = user->h;
138: c = user->kappa/(lambda-user->kappa);
139: d = n;
141: /*
142: Interior grid points
143: */
144: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
145: j[0] = i-1; j[1] = i; j[2] = i+1;
146: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
147: MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
148: }
150: /*
151: Boundary points
152: */
153: if (FirstBlock) {
154: i = 0;
155: j[0] = 0; j[1] = 1;
156: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
157: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
158: }
160: if (LastBlock) {
161: i = n-1;
162: j[0] = n-2; j[1] = n-1;
163: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
164: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
165: }
167: /*
168: Assemble matrix
169: */
170: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
171: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
172: if (fun != B) {
173: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
174: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
175: }
176: PetscFunctionReturn(0);
177: }
179: /* ------------------------------------------------------------------- */
180: /*
181: FormJacobian - Computes Jacobian matrix T'(lambda)
183: Input Parameters:
184: . nep - the NEP context
185: . lambda - the scalar argument
186: . ctx - optional user-defined context, as set by NEPSetJacobian()
188: Output Parameters:
189: . jac - Jacobian matrix
190: . B - optionally different preconditioning matrix
191: */
192: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
193: {
194: ApplicationCtx *user = (ApplicationCtx*)ctx;
195: PetscScalar A[3],c;
196: PetscReal h;
197: PetscInt i,n,j[3],Istart,Iend;
198: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
201: /*
202: Compute Jacobian entries and insert into matrix
203: */
204: MatGetSize(jac,&n,NULL);
205: MatGetOwnershipRange(jac,&Istart,&Iend);
206: if (Istart==0) FirstBlock=PETSC_TRUE;
207: if (Iend==n) LastBlock=PETSC_TRUE;
208: h = user->h;
209: c = user->kappa/(lambda-user->kappa);
211: /*
212: Interior grid points
213: */
214: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
215: j[0] = i-1; j[1] = i; j[2] = i+1;
216: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
217: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
218: }
220: /*
221: Boundary points
222: */
223: if (FirstBlock) {
224: i = 0;
225: j[0] = 0; j[1] = 1;
226: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
227: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
228: }
230: if (LastBlock) {
231: i = n-1;
232: j[0] = n-2; j[1] = n-1;
233: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
234: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
235: }
237: /*
238: Assemble matrix
239: */
240: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
241: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
242: PetscFunctionReturn(0);
243: }
245: /*TEST
247: testset:
248: args: -nep_type {{rii slp}} -nep_target 21 -terse -nep_view_vectors ::ascii_info
249: filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/" | sed -e "s/[+-]0\.0*i//g"
250: test:
251: suffix: 1_real
252: requires: !single !complex
253: test:
254: suffix: 1
255: requires: !single complex
257: test:
258: suffix: 2_cuda
259: args: -nep_type {{rii slp}} -nep_target 21 -mat_type aijcusparse -terse
260: requires: cuda !single
261: filter: sed -e "s/[+-]0\.0*i//"
262: output_file: output/test3_1.out
264: testset:
265: args: -nep_type slp -nep_two_sided -nep_target 21 -terse -nep_view_vectors ::ascii_info
266: filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/"
267: test:
268: suffix: 3_real
269: requires: !single !complex
270: test:
271: suffix: 3
272: requires: !single complex
274: TEST*/