Actual source code: ex34.c
slepc-3.9.1 2018-05-02
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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: */
10: /*
11: This is a nonlinear eigenvalue problem. When p=2, it is reduced to a linear Laplace eigenvalue
12: problem.
14: -\nabla\cdot(|\nabla u|^{p-2} \nabla u) = k |u|^{p-2} u in (0,1)x(0,1),
16: u = 0 on the entire boundary.
18: The code is implemented based on DMPlex using Q1 FEM on a quadrilateral mesh. In this code, we consider p=3.
20: Contributed by Fande Kong fdkong.jd@gmail.com
21: */
23: static char help[] = "Nonlinear inverse iteration for A(x)*x=lambda*B(x)*x.\n\n";
26: #include <slepceps.h>
27: #include <petscdmplex.h>
28: #include <petscds.h>
30: PetscErrorCode CreateSquareMesh(MPI_Comm,DM*);
31: PetscErrorCode SetupDiscretization(DM);
32: PetscErrorCode FormJacobianA(SNES,Vec,Mat,Mat,void*);
33: PetscErrorCode FormFunctionA(SNES,Vec,Vec,void*);
34: PetscErrorCode FormJacobianB(SNES,Vec,Mat,Mat,void*);
35: PetscErrorCode FormFunctionB(SNES,Vec,Vec,void*);
36: PetscErrorCode FormFunctionAB(SNES,Vec,Vec,Vec,void*);
37: PetscErrorCode BoundaryGlobalIndex(DM,const char*,IS*);
39: typedef struct {
40: IS bdis; /* global indices for boundary DoFs */
41: } AppCtx;
43: int main(int argc,char **argv)
44: {
45: DM dm;
46: MPI_Comm comm;
47: AppCtx user;
48: EPS eps; /* eigenproblem solver context */
49: EPSType type;
50: Mat A,B;
51: PetscContainer container;
52: PetscInt nev,nconv;
53: PetscBool nonlin;
54: SNES snes;
57: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
58: comm = PETSC_COMM_WORLD;
59: /* Create a quadrilateral mesh on domain (0,1)x(0,1) */
60: CreateSquareMesh(comm,&dm);
61: /* Setup basis function */
62: SetupDiscretization(dm);
63: BoundaryGlobalIndex(dm,"marker",&user.bdis);
65: DMCreateMatrix(dm,&A);
66: MatDuplicate(A,MAT_COPY_VALUES,&B);
68: /*
69: Compose callback functions and context that will be needed by the solver
70: */
71: PetscObjectComposeFunction((PetscObject)A,"formFunction",FormFunctionA);
72: PetscObjectComposeFunction((PetscObject)A,"formFunctionAB",FormFunctionAB);
73: PetscObjectComposeFunction((PetscObject)A,"formJacobian",FormJacobianA);
74: PetscObjectComposeFunction((PetscObject)B,"formFunction",FormFunctionB);
75: PetscContainerCreate(comm,&container);
76: PetscContainerSetPointer(container,&user);
77: PetscObjectCompose((PetscObject)A,"formFunctionCtx",(PetscObject)container);
78: PetscObjectCompose((PetscObject)A,"formJacobianCtx",(PetscObject)container);
79: PetscObjectCompose((PetscObject)B,"formFunctionCtx",(PetscObject)container);
80: PetscContainerDestroy(&container);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Create the eigensolver and set various options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86: EPSCreate(comm,&eps);
87: EPSSetOperators(eps,A,B);
88: EPSSetProblemType(eps,EPS_GNHEP);
89: /*
90: Use nonlinear inverse iteration
91: */
92: EPSSetType(eps,EPSPOWER);
93: EPSPowerSetNonlinear(eps,PETSC_TRUE);
94: /*
95: Attach DM to SNES
96: */
97: EPSPowerGetSNES(eps,&snes);
98: SNESSetDM(snes,dm);
99: EPSSetFromOptions(eps);
101: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102: Solve the eigensystem
103: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105: EPSSolve(eps);
107: /*
108: Optional: Get some information from the solver and display it
109: */
110: EPSGetType(eps,&type);
111: EPSPowerGetNonlinear(eps,&nonlin);
112: PetscPrintf(comm," Solution method: %s%s\n\n",type,nonlin?" (nonlinear)":"");
113: EPSGetDimensions(eps,&nev,NULL,NULL);
114: PetscPrintf(comm," Number of requested eigenvalues: %D\n",nev);
116: /* print eigenvalue and error */
117: EPSGetConverged(eps,&nconv);
118: if (nconv>0) {
119: PetscScalar k;
120: PetscReal na,nb;
121: Vec a,b,eigen;
122: DMCreateGlobalVector(dm,&a);
123: VecDuplicate(a,&b);
124: VecDuplicate(a,&eigen);
125: EPSGetEigenpair(eps,0,&k,NULL,eigen,NULL);
126: FormFunctionA(snes,eigen,a,&user);
127: FormFunctionB(snes,eigen,b,&user);
128: VecAXPY(a,-k,b);
129: VecNorm(a,NORM_2,&na);
130: VecNorm(b,NORM_2,&nb);
131: PetscPrintf(comm,"k: %g error: %g\n",(double)PetscRealPart(k),(double)na/nb);
132: VecDestroy(&a);
133: VecDestroy(&b);
134: VecDestroy(&eigen);
135: } else {
136: PetscPrintf(comm,"Solver did not converge\n");
137: }
139: MatDestroy(&A);
140: MatDestroy(&B);
141: DMDestroy(&dm);
142: EPSDestroy(&eps);
143: ISDestroy(&user.bdis);
144: SlepcFinalize();
145: return ierr;
146: }
148: /* <|u|u, v> */
149: static void f0_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
150: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
151: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
152: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f0[])
153: {
154: PetscScalar cof = PetscAbsScalar(u[0]);
156: f0[0] = cof*u[0];
157: }
159: /* <|\nabla u| \nabla u, \nabla v> */
160: static void f1_u(PetscInt dim, PetscInt Nf, PetscInt NfAux,
161: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
162: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
163: PetscReal t, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar f1[])
164: {
165: PetscInt d;
166: PetscScalar cof = 0;
167: for (d = 0; d < dim; ++d) cof += u_x[d]*u_x[d];
169: cof = PetscSqrtScalar(cof);
171: for (d = 0; d < dim; ++d) f1[d] = u_x[d]*cof;
172: }
174: /* approximate Jacobian for <|u|u, v> */
175: static void g0_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
176: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
177: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
178: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g0[])
179: {
180: g0[0] = 1.0*PetscAbsScalar(u[0]);
181: }
183: /* approximate Jacobian for <|\nabla u| \nabla u, \nabla v> */
184: static void g3_uu(PetscInt dim, PetscInt Nf, PetscInt NfAux,
185: const PetscInt uOff[], const PetscInt uOff_x[], const PetscScalar u[], const PetscScalar u_t[], const PetscScalar u_x[],
186: const PetscInt aOff[], const PetscInt aOff_x[], const PetscScalar a[], const PetscScalar a_t[], const PetscScalar a_x[],
187: PetscReal t, PetscReal u_tShift, const PetscReal x[], PetscInt numConstants, const PetscScalar constants[], PetscScalar g3[])
188: {
189: PetscInt d;
190: for (d = 0; d < dim; ++d) g3[d*dim+d] = 1.0;
191: }
193: PetscErrorCode SetupDiscretization(DM dm)
194: {
195: PetscFE fe;
196: PetscDS prob;
197: PetscInt totDim;
201: /* Create finite element */
202: PetscFECreateDefault(dm,2,1,PETSC_FALSE,NULL,-1,&fe);
203: PetscObjectSetName((PetscObject) fe, "u");
204: DMGetDS(dm,&prob);
205: PetscDSSetDiscretization(prob,0,(PetscObject) fe);
206: DMSetDS(dm,prob);
207: PetscDSGetTotalDimension(prob,&totDim);
208: PetscFEDestroy(&fe);
209: return(0);
210: }
212: PetscErrorCode CreateSquareMesh(MPI_Comm comm,DM *dm)
213: {
214: PetscInt cells[] = {8,8};
215: PetscInt dim = 2;
216: DM pdm;
217: PetscMPIInt size;
221: DMPlexCreateBoxMesh(comm,dim,PETSC_FALSE,cells,NULL,NULL,NULL,PETSC_TRUE,dm);
222: DMSetFromOptions(*dm);
223: DMSetUp(*dm);
224: MPI_Comm_size(comm,&size);
225: if (size > 1) {
226: DMPlexDistribute(*dm,0,NULL,&pdm);
227: DMDestroy(dm);
228: *dm = pdm;
229: }
230: return(0);
231: }
233: PetscErrorCode BoundaryGlobalIndex(DM dm,const char labelname[],IS *bdis)
234: {
235: IS bdpoints;
236: PetscInt nindices,*indices,numDof,offset,npoints,i,j;
237: const PetscInt *bdpoints_indices;
238: DMLabel bdmarker;
239: PetscSection gsection;
243: DMGetDefaultGlobalSection(dm,&gsection);
244: DMGetLabel(dm,labelname,&bdmarker);
245: DMLabelGetStratumIS(bdmarker,1,&bdpoints);
246: ISGetLocalSize(bdpoints,&npoints);
247: ISGetIndices(bdpoints,&bdpoints_indices);
248: nindices = 0;
249: for (i=0;i<npoints;i++) {
250: PetscSectionGetDof(gsection,bdpoints_indices[i],&numDof);
251: if (numDof<=0) continue;
252: nindices += numDof;
253: }
254: PetscCalloc1(nindices,&indices);
255: nindices = 0;
256: for (i=0;i<npoints;i++) {
257: PetscSectionGetDof(gsection,bdpoints_indices[i],&numDof);
258: if (numDof<=0) continue;
259: PetscSectionGetOffset(gsection,bdpoints_indices[i],&offset);
260: for (j=0;j<numDof;j++) indices[nindices++] = offset+j;
261: }
262: ISRestoreIndices(bdpoints,&bdpoints_indices);
263: ISDestroy(&bdpoints);
264: ISCreateGeneral(PetscObjectComm((PetscObject)dm),nindices,indices,PETSC_OWN_POINTER,bdis);
265: return(0);
266: }
268: static PetscErrorCode FormJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx)
269: {
270: DM dm;
271: Vec Xloc;
275: SNESGetDM(snes,&dm);
276: DMGetLocalVector(dm,&Xloc);
277: VecZeroEntries(Xloc);
278: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
279: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
280: CHKMEMQ;
281: DMPlexSNESComputeJacobianFEM(dm,Xloc,A,B,ctx);
282: CHKMEMQ;
283: DMRestoreLocalVector(dm,&Xloc);
284: if (A != B) {
285: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
286: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
287: }
288: return(0);
289: }
291: PetscErrorCode FormJacobianA(SNES snes,Vec X,Mat A,Mat B,void *ctx)
292: {
294: DM dm;
295: PetscDS prob;
296: AppCtx *userctx = (AppCtx *)ctx;
299: MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
300: SNESGetDM(snes,&dm);
301: DMGetDS(dm,&prob);
302: PetscDSSetJacobian(prob,0,0,NULL,NULL,NULL,g3_uu);
303: FormJacobian(snes,X,A,B,ctx);
304: MatZeroRowsIS(B,userctx->bdis,1.0,NULL,NULL);
305: return(0);
306: }
308: PetscErrorCode FormJacobianB(SNES snes,Vec X,Mat A,Mat B,void *ctx)
309: {
311: DM dm;
312: PetscDS prob;
313: AppCtx *userctx = (AppCtx *)ctx;
316: MatSetOption(B,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE);
317: SNESGetDM(snes,&dm);
318: DMGetDS(dm,&prob);
319: PetscDSSetJacobian(prob,0,0,g0_uu,NULL,NULL,NULL);
320: FormJacobian(snes,X,A,B,ctx);
321: MatZeroRowsIS(B,userctx->bdis,0.0,NULL,NULL);
322: return(0);
323: }
325: PetscErrorCode FormFunctionAB(SNES snes,Vec x,Vec Ax,Vec Bx,void *ctx)
326: {
330: /*
331: * In real applications, users should have a generic formFunctionAB which
332: * forms Ax and Bx simultaneously for an more efficient calculation.
333: * In this example, we just call FormFunctionA+FormFunctionB to mimic how
334: * to use FormFunctionAB
335: */
336: FormFunctionA(snes,x,Ax,ctx);
337: FormFunctionB(snes,x,Bx,ctx);
338: return(0);
339: }
342: static PetscErrorCode FormFunction(SNES snes,Vec X,Vec F,void *ctx)
343: {
344: DM dm;
345: Vec Xloc,Floc;
349: SNESGetDM(snes,&dm);
350: DMGetLocalVector(dm,&Xloc);
351: DMGetLocalVector(dm,&Floc);
352: VecZeroEntries(Xloc);
353: VecZeroEntries(Floc);
354: DMGlobalToLocalBegin(dm,X,INSERT_VALUES,Xloc);
355: DMGlobalToLocalEnd(dm,X,INSERT_VALUES,Xloc);
356: CHKMEMQ;
357: DMPlexSNESComputeResidualFEM(dm,Xloc,Floc,ctx);
358: CHKMEMQ;
359: VecZeroEntries(F);
360: DMLocalToGlobalBegin(dm,Floc,ADD_VALUES,F);
361: DMLocalToGlobalEnd(dm,Floc,ADD_VALUES,F);
362: DMRestoreLocalVector(dm,&Xloc);
363: DMRestoreLocalVector(dm,&Floc);
364: return(0);
365: }
367: PetscErrorCode FormFunctionA(SNES snes,Vec X,Vec F,void *ctx)
368: {
370: DM dm;
371: PetscDS prob;
372: PetscInt nindices,iStart,iEnd,i;
373: AppCtx *userctx = (AppCtx *)ctx;
374: PetscScalar *array,value;
375: const PetscInt *indices;
376: PetscInt vecstate;
379: SNESGetDM(snes,&dm);
380: DMGetDS(dm,&prob);
381: /* hook functions */
382: PetscDSSetResidual(prob,0,NULL,f1_u);
383: FormFunction(snes,X,F,ctx);
384: /* Boundary condition */
385: VecLockGet(X,&vecstate);
386: if (vecstate>0) {
387: VecLockPop(X);
388: }
389: VecGetOwnershipRange(X,&iStart,&iEnd);
390: VecGetArray(X,&array);
391: ISGetLocalSize(userctx->bdis,&nindices);
392: ISGetIndices(userctx->bdis,&indices);
393: for (i=0;i<nindices;i++) {
394: value = array[indices[i]-iStart] - 0.0;
395: VecSetValue(F,indices[i],value,INSERT_VALUES);
396: }
397: ISRestoreIndices(userctx->bdis,&indices);
398: VecRestoreArray(X,&array);
399: if (vecstate>0) {
400: VecLockPush(X);
401: }
402: VecAssemblyBegin(F);
403: VecAssemblyEnd(F);
404: return(0);
405: }
407: PetscErrorCode FormFunctionB(SNES snes,Vec X,Vec F,void *ctx)
408: {
410: DM dm;
411: PetscDS prob;
412: PetscInt nindices,iStart,iEnd,i;
413: AppCtx *userctx = (AppCtx *)ctx;
414: PetscScalar value;
415: const PetscInt *indices;
418: SNESGetDM(snes,&dm);
419: DMGetDS(dm,&prob);
420: /* hook functions */
421: PetscDSSetResidual(prob,0,f0_u,NULL);
422: FormFunction(snes,X,F,ctx);
423: /* Boundary condition */
424: VecGetOwnershipRange(F,&iStart,&iEnd);
425: ISGetLocalSize(userctx->bdis,&nindices);
426: ISGetIndices(userctx->bdis,&indices);
427: for (i=0;i<nindices;i++) {
428: value = 0.0;
429: VecSetValue(F,indices[i],value,INSERT_VALUES);
430: }
431: ISRestoreIndices(userctx->bdis,&indices);
432: VecAssemblyBegin(F);
433: VecAssemblyEnd(F);
434: return(0);
435: }