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[] = "Test the SLP solver with a user-provided EPS.\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 2 : int main(int argc,char **argv)
40 : {
41 2 : NEP nep;
42 2 : EPS eps;
43 2 : KSP ksp;
44 2 : PC pc;
45 2 : Mat F,J;
46 2 : ApplicationCtx ctx;
47 2 : PetscInt n=128;
48 2 : PetscReal deftol;
49 2 : PetscBool terse,flag,ts;
50 :
51 2 : PetscFunctionBeginUser;
52 2 : PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
53 2 : PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
54 2 : PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
55 2 : ctx.h = 1.0/(PetscReal)n;
56 2 : ctx.kappa = 1.0;
57 :
58 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59 : Create a standalone EPS with appropriate settings
60 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61 :
62 2 : PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
63 2 : PetscCall(EPSSetType(eps,EPSGD));
64 2 : PetscCall(EPSSetFromOptions(eps));
65 :
66 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67 : Create a standalone KSP with appropriate settings
68 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69 :
70 2 : PetscCall(KSPCreate(PETSC_COMM_WORLD,&ksp));
71 2 : PetscCall(KSPSetType(ksp,KSPBCGS));
72 2 : PetscCall(KSPGetPC(ksp,&pc));
73 2 : PetscCall(PCSetType(pc,PCBJACOBI));
74 2 : PetscCall(KSPSetFromOptions(ksp));
75 :
76 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77 : Prepare nonlinear eigensolver context
78 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79 :
80 2 : PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
81 :
82 : /* Create Function and Jacobian matrices; set evaluation routines */
83 2 : PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
84 2 : PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
85 2 : PetscCall(MatSetFromOptions(F));
86 2 : PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
87 2 : PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
88 2 : PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
89 :
90 2 : PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
91 2 : PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
92 2 : PetscCall(MatSetFromOptions(J));
93 2 : PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
94 2 : PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
95 2 : PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
96 :
97 : /* Set options */
98 2 : PetscCall(NEPSetType(nep,NEPSLP));
99 2 : PetscCall(NEPSLPSetEPS(nep,eps));
100 2 : PetscCall(NEPSLPSetKSP(nep,ksp));
101 2 : PetscCall(NEPSetFromOptions(nep));
102 :
103 : /* Print some options */
104 2 : PetscCall(PetscObjectTypeCompare((PetscObject)nep,NEPSLP,&flag));
105 2 : if (flag) {
106 2 : PetscCall(NEPGetTwoSided(nep,&ts));
107 2 : if (ts) {
108 1 : PetscCall(NEPSLPGetDeflationThreshold(nep,&deftol));
109 1 : PetscCall(PetscPrintf(PETSC_COMM_WORLD," Two-sided solve with deflation threshold=%g\n",(double)deftol));
110 : }
111 : }
112 :
113 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114 : Solve the eigensystem
115 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
116 :
117 2 : PetscCall(NEPSolve(nep));
118 :
119 : /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120 : Display solution and clean up
121 : - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122 :
123 : /* show detailed info unless -terse option is given by user */
124 2 : PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
125 2 : if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
126 : else {
127 0 : PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
128 0 : PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
129 0 : PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
130 0 : PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
131 : }
132 :
133 2 : PetscCall(NEPDestroy(&nep));
134 2 : PetscCall(EPSDestroy(&eps));
135 2 : PetscCall(KSPDestroy(&ksp));
136 2 : PetscCall(MatDestroy(&F));
137 2 : PetscCall(MatDestroy(&J));
138 2 : PetscCall(SlepcFinalize());
139 : return 0;
140 : }
141 :
142 : /* ------------------------------------------------------------------- */
143 : /*
144 : FormFunction - Computes Function matrix T(lambda)
145 :
146 : Input Parameters:
147 : . nep - the NEP context
148 : . lambda - the scalar argument
149 : . ctx - optional user-defined context, as set by NEPSetFunction()
150 :
151 : Output Parameters:
152 : . fun - Function matrix
153 : . B - optionally different preconditioning matrix
154 : */
155 11 : PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
156 : {
157 11 : ApplicationCtx *user = (ApplicationCtx*)ctx;
158 11 : PetscScalar A[3],c,d;
159 11 : PetscReal h;
160 11 : PetscInt i,n,j[3],Istart,Iend;
161 11 : PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
162 :
163 11 : PetscFunctionBeginUser;
164 : /*
165 : Compute Function entries and insert into matrix
166 : */
167 11 : PetscCall(MatGetSize(fun,&n,NULL));
168 11 : PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
169 11 : if (Istart==0) FirstBlock=PETSC_TRUE;
170 11 : if (Iend==n) LastBlock=PETSC_TRUE;
171 11 : h = user->h;
172 11 : c = user->kappa/(lambda-user->kappa);
173 11 : d = n;
174 :
175 : /*
176 : Interior grid points
177 : */
178 1397 : for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
179 1386 : j[0] = i-1; j[1] = i; j[2] = i+1;
180 1386 : A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
181 1386 : PetscCall(MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES));
182 : }
183 :
184 : /*
185 : Boundary points
186 : */
187 11 : if (FirstBlock) {
188 11 : i = 0;
189 11 : j[0] = 0; j[1] = 1;
190 11 : A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
191 11 : PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
192 : }
193 :
194 11 : if (LastBlock) {
195 11 : i = n-1;
196 11 : j[0] = n-2; j[1] = n-1;
197 11 : A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
198 11 : PetscCall(MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES));
199 : }
200 :
201 : /*
202 : Assemble matrix
203 : */
204 11 : PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
205 11 : PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
206 11 : if (fun != B) {
207 0 : PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
208 0 : PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
209 : }
210 11 : PetscFunctionReturn(PETSC_SUCCESS);
211 : }
212 :
213 : /* ------------------------------------------------------------------- */
214 : /*
215 : FormJacobian - Computes Jacobian matrix T'(lambda)
216 :
217 : Input Parameters:
218 : . nep - the NEP context
219 : . lambda - the scalar argument
220 : . ctx - optional user-defined context, as set by NEPSetJacobian()
221 :
222 : Output Parameters:
223 : . jac - Jacobian matrix
224 : . B - optionally different preconditioning matrix
225 : */
226 4 : PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
227 : {
228 4 : ApplicationCtx *user = (ApplicationCtx*)ctx;
229 4 : PetscScalar A[3],c;
230 4 : PetscReal h;
231 4 : PetscInt i,n,j[3],Istart,Iend;
232 4 : PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
233 :
234 4 : PetscFunctionBeginUser;
235 : /*
236 : Compute Jacobian entries and insert into matrix
237 : */
238 4 : PetscCall(MatGetSize(jac,&n,NULL));
239 4 : PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
240 4 : if (Istart==0) FirstBlock=PETSC_TRUE;
241 4 : if (Iend==n) LastBlock=PETSC_TRUE;
242 4 : h = user->h;
243 4 : c = user->kappa/(lambda-user->kappa);
244 :
245 : /*
246 : Interior grid points
247 : */
248 508 : for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
249 504 : j[0] = i-1; j[1] = i; j[2] = i+1;
250 504 : A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
251 504 : PetscCall(MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES));
252 : }
253 :
254 : /*
255 : Boundary points
256 : */
257 4 : if (FirstBlock) {
258 4 : i = 0;
259 4 : j[0] = 0; j[1] = 1;
260 4 : A[0] = -2.0*h/3.0; A[1] = -h/6.0;
261 4 : PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
262 : }
263 :
264 4 : if (LastBlock) {
265 4 : i = n-1;
266 4 : j[0] = n-2; j[1] = n-1;
267 4 : A[0] = -h/6.0; A[1] = -h/3.0-c*c;
268 4 : PetscCall(MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES));
269 : }
270 :
271 : /*
272 : Assemble matrix
273 : */
274 4 : PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
275 4 : PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
276 4 : PetscFunctionReturn(PETSC_SUCCESS);
277 : }
278 :
279 : /*TEST
280 :
281 : test:
282 : args: -nep_target 21 -terse
283 : requires: !single
284 : test:
285 : suffix: 1
286 : test:
287 : suffix: 1_ts
288 : args: -nep_two_sided
289 :
290 : TEST*/
|