Actual source code: test1.c
slepc-3.20.1 2023-11-27
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(MatSetUp(F));
70: PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
72: PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
73: PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
74: PetscCall(MatSetFromOptions(J));
75: PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
76: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
77: PetscCall(MatSetUp(J));
78: PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
80: PetscCall(NEPSetFromOptions(nep));
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Solve the eigensystem
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: PetscCall(NEPSolve(nep));
88: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89: Display solution and clean up
90: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92: /* show detailed info unless -terse option is given by user */
93: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
94: if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
95: else {
96: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
97: PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
98: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
99: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
100: }
102: PetscCall(NEPDestroy(&nep));
103: PetscCall(MatDestroy(&F));
104: PetscCall(MatDestroy(&J));
105: PetscCall(SlepcFinalize());
106: return 0;
107: }
109: /* ------------------------------------------------------------------- */
110: /*
111: FormFunction - Computes Function matrix T(lambda)
113: Input Parameters:
114: . nep - the NEP context
115: . lambda - the scalar argument
116: . ctx - optional user-defined context, as set by NEPSetFunction()
118: Output Parameters:
119: . fun - Function matrix
120: . B - optionally different preconditioning matrix
121: */
122: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
123: {
124: ApplicationCtx *user = (ApplicationCtx*)ctx;
125: PetscScalar A[3],c,d;
126: PetscReal h;
127: PetscInt i,n,j[3],Istart,Iend;
128: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
130: PetscFunctionBeginUser;
131: /*
132: Compute Function entries and insert into matrix
133: */
134: PetscCall(MatGetSize(fun,&n,NULL));
135: PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
136: if (Istart==0) FirstBlock=PETSC_TRUE;
137: if (Iend==n) LastBlock=PETSC_TRUE;
138: h = user->h;
139: c = user->kappa/(lambda-user->kappa);
140: d = n;
142: /*
143: Interior grid points
144: */
145: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
146: j[0] = i-1; j[1] = i; j[2] = i+1;
147: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
148: PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES));
149: }
151: /*
152: Boundary points
153: */
154: if (FirstBlock) {
155: i = 0;
156: j[0] = 0; j[1] = 1;
157: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
158: PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
159: }
161: if (LastBlock) {
162: i = n-1;
163: j[0] = n-2; j[1] = n-1;
164: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
165: PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
166: }
168: /*
169: Assemble matrix
170: */
171: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
172: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
173: if (fun != B) {
174: PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
175: PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
176: }
177: PetscFunctionReturn(PETSC_SUCCESS);
178: }
180: /* ------------------------------------------------------------------- */
181: /*
182: FormJacobian - Computes Jacobian matrix T'(lambda)
184: Input Parameters:
185: . nep - the NEP context
186: . lambda - the scalar argument
187: . ctx - optional user-defined context, as set by NEPSetJacobian()
189: Output Parameters:
190: . jac - Jacobian matrix
191: . B - optionally different preconditioning matrix
192: */
193: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
194: {
195: ApplicationCtx *user = (ApplicationCtx*)ctx;
196: PetscScalar A[3],c;
197: PetscReal h;
198: PetscInt i,n,j[3],Istart,Iend;
199: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
201: PetscFunctionBeginUser;
202: /*
203: Compute Jacobian entries and insert into matrix
204: */
205: PetscCall(MatGetSize(jac,&n,NULL));
206: PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
207: if (Istart==0) FirstBlock=PETSC_TRUE;
208: if (Iend==n) LastBlock=PETSC_TRUE;
209: h = user->h;
210: c = user->kappa/(lambda-user->kappa);
212: /*
213: Interior grid points
214: */
215: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
216: j[0] = i-1; j[1] = i; j[2] = i+1;
217: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
218: PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
219: }
221: /*
222: Boundary points
223: */
224: if (FirstBlock) {
225: i = 0;
226: j[0] = 0; j[1] = 1;
227: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
228: PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
229: }
231: if (LastBlock) {
232: i = n-1;
233: j[0] = n-2; j[1] = n-1;
234: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
235: PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
236: }
238: /*
239: Assemble matrix
240: */
241: PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
242: PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
243: PetscFunctionReturn(PETSC_SUCCESS);
244: }
246: /*TEST
248: testset:
249: args: -nep_type {{rii slp}} -nep_target 21 -terse -nep_view_vectors ::ascii_info
250: filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/" | sed -e "s/[+-]0\.0*i//g"
251: test:
252: suffix: 1_real
253: requires: !single !complex
254: test:
255: suffix: 1
256: requires: !single complex
258: test:
259: suffix: 2_cuda
260: args: -nep_type {{rii slp}} -nep_target 21 -mat_type aijcusparse -terse
261: requires: cuda !single
262: filter: sed -e "s/[+-]0\.0*i//"
263: output_file: output/test3_1.out
265: testset:
266: args: -nep_type slp -nep_two_sided -nep_target 21 -terse -nep_view_vectors ::ascii_info
267: filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/"
268: test:
269: suffix: 3_real
270: requires: !single !complex
271: test:
272: suffix: 3
273: requires: !single complex
275: TEST*/