Line data Source code
1 : /*
2 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3 : SLEPc - Scalable Library for Eigenvalue Problem Computations
4 : Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
5 :
6 : This file is part of SLEPc.
7 : SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9 : */
10 :
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";
15 :
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 : */
22 :
23 : #include <slepcnep.h>
24 :
25 : /*
26 : User-defined routines
27 : */
28 : PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*);
29 : PetscErrorCode FormJacobian(NEP,PetscScalar,Mat,void*);
30 :
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;
38 :
39 3 : int main(int argc,char **argv)
40 : {
41 3 : NEP nep; /* nonlinear eigensolver context */
42 3 : Mat F,J; /* Function and Jacobian matrices */
43 3 : ApplicationCtx ctx; /* user-defined context */
44 3 : PetscInt n=128;
45 3 : PetscBool terse;
46 :
47 3 : PetscFunctionBeginUser;
48 3 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
49 3 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
50 3 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
51 3 : ctx.h = 1.0/(PetscReal)n;
52 3 : ctx.kappa = 1.0;
53 :
54 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55 : Prepare nonlinear eigensolver context
56 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57 :
58 3 : PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
59 :
60 : /*
61 : Create Function and Jacobian matrices; set evaluation routines
62 : */
63 :
64 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
65 3 : PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
66 3 : PetscCall(MatSetFromOptions(F));
67 3 : PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
68 3 : PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
69 3 : PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
70 :
71 3 : PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
72 3 : PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
73 3 : PetscCall(MatSetFromOptions(J));
74 3 : PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
75 3 : PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
76 3 : PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
77 :
78 3 : PetscCall(NEPSetFromOptions(nep));
79 :
80 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81 : Solve the eigensystem
82 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83 :
84 3 : PetscCall(NEPSolve(nep));
85 :
86 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87 : Display solution and clean up
88 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89 :
90 : /* show detailed info unless -terse option is given by user */
91 3 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
92 3 : if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
93 : else {
94 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
95 0 : PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
96 0 : PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
97 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
98 : }
99 :
100 3 : PetscCall(NEPDestroy(&nep));
101 3 : PetscCall(MatDestroy(&F));
102 3 : PetscCall(MatDestroy(&J));
103 3 : PetscCall(SlepcFinalize());
104 : return 0;
105 : }
106 :
107 : /* ------------------------------------------------------------------- */
108 : /*
109 : FormFunction - Computes Function matrix T(lambda)
110 :
111 : Input Parameters:
112 : . nep - the NEP context
113 : . lambda - the scalar argument
114 : . ctx - optional user-defined context, as set by NEPSetFunction()
115 :
116 : Output Parameters:
117 : . fun - Function matrix
118 : . B - optionally different preconditioning matrix
119 : */
120 33 : PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
121 : {
122 33 : ApplicationCtx *user = (ApplicationCtx*)ctx;
123 33 : PetscScalar A[3],c,d;
124 33 : PetscReal h;
125 33 : PetscInt i,n,j[3],Istart,Iend;
126 33 : PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
127 :
128 33 : PetscFunctionBeginUser;
129 : /*
130 : Compute Function entries and insert into matrix
131 : */
132 33 : PetscCall(MatGetSize(fun,&n,NULL));
133 33 : PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
134 33 : if (Istart==0) FirstBlock=PETSC_TRUE;
135 33 : if (Iend==n) LastBlock=PETSC_TRUE;
136 33 : h = user->h;
137 33 : c = user->kappa/(lambda-user->kappa);
138 33 : d = n;
139 :
140 : /*
141 : Interior grid points
142 : */
143 4191 : for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
144 4158 : j[0] = i-1; j[1] = i; j[2] = i+1;
145 4158 : A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
146 4158 : PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES));
147 : }
148 :
149 : /*
150 : Boundary points
151 : */
152 33 : if (FirstBlock) {
153 33 : i = 0;
154 33 : j[0] = 0; j[1] = 1;
155 33 : A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
156 33 : PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
157 : }
158 :
159 33 : if (LastBlock) {
160 33 : i = n-1;
161 33 : j[0] = n-2; j[1] = n-1;
162 33 : A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
163 33 : PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
164 : }
165 :
166 : /*
167 : Assemble matrix
168 : */
169 33 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
170 33 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
171 33 : if (fun != B) {
172 0 : PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
173 0 : PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
174 : }
175 33 : PetscFunctionReturn(PETSC_SUCCESS);
176 : }
177 :
178 : /* ------------------------------------------------------------------- */
179 : /*
180 : FormJacobian - Computes Jacobian matrix T'(lambda)
181 :
182 : Input Parameters:
183 : . nep - the NEP context
184 : . lambda - the scalar argument
185 : . ctx - optional user-defined context, as set by NEPSetJacobian()
186 :
187 : Output Parameters:
188 : . jac - Jacobian matrix
189 : . B - optionally different preconditioning matrix
190 : */
191 23 : PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
192 : {
193 23 : ApplicationCtx *user = (ApplicationCtx*)ctx;
194 23 : PetscScalar A[3],c;
195 23 : PetscReal h;
196 23 : PetscInt i,n,j[3],Istart,Iend;
197 23 : PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
198 :
199 23 : PetscFunctionBeginUser;
200 : /*
201 : Compute Jacobian entries and insert into matrix
202 : */
203 23 : PetscCall(MatGetSize(jac,&n,NULL));
204 23 : PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
205 23 : if (Istart==0) FirstBlock=PETSC_TRUE;
206 23 : if (Iend==n) LastBlock=PETSC_TRUE;
207 23 : h = user->h;
208 23 : c = user->kappa/(lambda-user->kappa);
209 :
210 : /*
211 : Interior grid points
212 : */
213 2921 : for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
214 2898 : j[0] = i-1; j[1] = i; j[2] = i+1;
215 2898 : A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
216 2898 : PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
217 : }
218 :
219 : /*
220 : Boundary points
221 : */
222 23 : if (FirstBlock) {
223 23 : i = 0;
224 23 : j[0] = 0; j[1] = 1;
225 23 : A[0] = -2.0*h/3.0; A[1] = -h/6.0;
226 23 : PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
227 : }
228 :
229 23 : if (LastBlock) {
230 23 : i = n-1;
231 23 : j[0] = n-2; j[1] = n-1;
232 23 : A[0] = -h/6.0; A[1] = -h/3.0-c*c;
233 23 : PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
234 : }
235 :
236 : /*
237 : Assemble matrix
238 : */
239 23 : PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
240 23 : PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
241 23 : PetscFunctionReturn(PETSC_SUCCESS);
242 : }
243 :
244 : /*TEST
245 :
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
255 :
256 : testset:
257 : args: -nep_type {{rii slp}} -nep_target 21 -terse
258 : filter: sed -e "s/[+-]0\.0*i//"
259 : output_file: output/test3_1.out
260 : test:
261 : suffix: 2_cuda
262 : args: -mat_type aijcusparse
263 : requires: cuda !single
264 : test:
265 : suffix: 2_hip
266 : args: -mat_type aijhipsparse
267 : requires: hip !single
268 :
269 : testset:
270 : args: -nep_type slp -nep_two_sided -nep_target 21 -terse -nep_view_vectors ::ascii_info
271 : filter: sed -e "s/\(0x[0-9a-fA-F]*\)/objectid/"
272 : test:
273 : suffix: 3_real
274 : requires: !single !complex
275 : test:
276 : suffix: 3
277 : requires: !single complex
278 :
279 : TEST*/
|