Actual source code: ex69.c
petsc-3.6.4 2016-04-12
2: static char help[] = "Tests recovery from domain errors in MatMult() and PCApply()\n\n";
4: /*
5: See src/ksp/ksp/examples/tutorials/ex19.c from which this was copied
6: */
8: #include <petscsnes.h>
9: #include <petscdm.h>
10: #include <petscdmda.h>
12: /*
13: User-defined routines and data structures
14: */
15: typedef struct {
16: PetscScalar u,v,omega,temp;
17: } Field;
19: PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field**,Field**,void*);
21: typedef struct {
22: PassiveReal lidvelocity,prandtl,grashof; /* physical parameters */
23: PetscBool draw_contours; /* flag - 1 indicates drawing contours */
24: PetscBool errorindomain;
25: PetscBool errorindomainmf;
26: SNES snes;
27: } AppCtx;
29: typedef struct {
30: Mat Jmf;
31: } MatShellCtx;
33: extern PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
34: extern PetscErrorCode MatMult_MyShell(Mat,Vec,Vec);
35: extern PetscErrorCode MatAssemblyEnd_MyShell(Mat,MatAssemblyType);
36: extern PetscErrorCode PCApply_MyShell(PC,Vec,Vec);
37: extern PetscErrorCode SNESComputeJacobian_MyShell(SNES,Vec,Mat,Mat,void*);
41: int main(int argc,char **argv)
42: {
43: AppCtx user; /* user-defined work context */
44: PetscInt mx,my;
46: MPI_Comm comm;
47: DM da;
48: Vec x;
49: Mat J = NULL,Jmf = NULL;
50: MatShellCtx matshellctx;
51: PetscInt mlocal,nlocal;
52: PC pc;
53: KSP ksp;
54: PetscBool errorinmatmult = PETSC_FALSE,errorinpcapply = PETSC_FALSE,errorinpcsetup = PETSC_FALSE;
56: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return(1);
59: PetscOptionsGetBool(NULL,"-error_in_matmult",&errorinmatmult,NULL);
60: PetscOptionsGetBool(NULL,"-error_in_pcapply",&errorinpcapply,NULL);
61: PetscOptionsGetBool(NULL,"-error_in_pcsetup",&errorinpcsetup,NULL);
62: PetscOptionsGetBool(NULL,"-error_in_domain",&user.errorindomain,NULL);
63: PetscOptionsGetBool(NULL,"-error_in_domainmf",&user.errorindomainmf,NULL);
65: comm = PETSC_COMM_WORLD;
66: SNESCreate(comm,&user.snes);
68: /*
69: Create distributed array object to manage parallel grid and vectors
70: for principal unknowns (x) and governing residuals (f)
71: */
72: DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da);
73: SNESSetDM(user.snes,da);
75: DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
76: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
77: /*
78: Problem parameters (velocity of lid, prandtl, and grashof numbers)
79: */
80: user.lidvelocity = 1.0/(mx*my);
81: user.prandtl = 1.0;
82: user.grashof = 1.0;
84: PetscOptionsGetReal(NULL,"-lidvelocity",&user.lidvelocity,NULL);
85: PetscOptionsGetReal(NULL,"-prandtl",&user.prandtl,NULL);
86: PetscOptionsGetReal(NULL,"-grashof",&user.grashof,NULL);
87: PetscOptionsHasName(NULL,"-contours",&user.draw_contours);
89: DMDASetFieldName(da,0,"x_velocity");
90: DMDASetFieldName(da,1,"y_velocity");
91: DMDASetFieldName(da,2,"Omega");
92: DMDASetFieldName(da,3,"temperature");
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Create user context, set problem data, create vector data structures.
96: Also, compute the initial guess.
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Create nonlinear solver context
102: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: DMSetApplicationContext(da,&user);
104: DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user);
106: if (errorinmatmult) {
107: MatCreateSNESMF(user.snes,&Jmf);
108: MatSetFromOptions(Jmf);
109: MatGetLocalSize(Jmf,&mlocal,&nlocal);
110: matshellctx.Jmf = Jmf;
111: MatCreateShell(PetscObjectComm((PetscObject)Jmf),mlocal,nlocal,PETSC_DECIDE,PETSC_DECIDE,&matshellctx,&J);
112: MatShellSetOperation(J,MATOP_MULT,(void (*)(void))MatMult_MyShell);
113: MatShellSetOperation(J,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MyShell);
114: SNESSetJacobian(user.snes,J,J,MatMFFDComputeJacobian,NULL);
115: }
117: SNESSetFromOptions(user.snes);
118: PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %g\n",(double)user.lidvelocity,(double)user.prandtl,(double)user.grashof);
120: if (errorinpcapply) {
121: SNESGetKSP(user.snes,&ksp);
122: KSPGetPC(ksp,&pc);
123: PCSetType(pc,PCSHELL);
124: PCShellSetApply(pc,PCApply_MyShell);
125: }
127: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128: Solve the nonlinear system
129: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: DMCreateGlobalVector(da,&x);
131: FormInitialGuess(&user,da,x);
133: if (errorinpcsetup) {
134: SNESSetUp(user.snes);
135: SNESSetJacobian(user.snes,NULL,NULL,SNESComputeJacobian_MyShell,NULL);
136: }
137: SNESSolve(user.snes,NULL,x);
140: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141: Free work space. All PETSc objects should be destroyed when they
142: are no longer needed.
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: MatDestroy(&J);
145: MatDestroy(&Jmf);
146: VecDestroy(&x);
147: DMDestroy(&da);
148: SNESDestroy(&user.snes);
149: PetscFinalize();
150: return 0;
151: }
153: /* ------------------------------------------------------------------- */
158: /*
159: FormInitialGuess - Forms initial approximation.
161: Input Parameters:
162: user - user-defined application context
163: X - vector
165: Output Parameter:
166: X - vector
167: */
168: PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
169: {
170: PetscInt i,j,mx,xs,ys,xm,ym;
172: PetscReal grashof,dx;
173: Field **x;
176: grashof = user->grashof;
178: DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0);
179: dx = 1.0/(mx-1);
181: /*
182: Get local grid boundaries (for 2-dimensional DMDA):
183: xs, ys - starting grid indices (no ghost points)
184: xm, ym - widths of local grid (no ghost points)
185: */
186: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
188: /*
189: Get a pointer to vector data.
190: - For default PETSc vectors, VecGetArray() returns a pointer to
191: the data array. Otherwise, the routine is implementation dependent.
192: - You MUST call VecRestoreArray() when you no longer need access to
193: the array.
194: */
195: DMDAVecGetArray(da,X,&x);
197: /*
198: Compute initial guess over the locally owned part of the grid
199: Initial condition is motionless fluid and equilibrium temperature
200: */
201: for (j=ys; j<ys+ym; j++) {
202: for (i=xs; i<xs+xm; i++) {
203: x[j][i].u = 0.0;
204: x[j][i].v = 0.0;
205: x[j][i].omega = 0.0;
206: x[j][i].temp = (grashof>0)*i*dx;
207: }
208: }
210: /*
211: Restore vector
212: */
213: DMDAVecRestoreArray(da,X,&x);
214: return(0);
215: }
219: PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
220: {
221: AppCtx *user = (AppCtx*)ptr;
222: PetscErrorCode ierr;
223: PetscInt xints,xinte,yints,yinte,i,j;
224: PetscReal hx,hy,dhx,dhy,hxdhy,hydhx;
225: PetscReal grashof,prandtl,lid;
226: PetscScalar u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
227: static PetscInt fail = 0;
230: if ((fail++ > 7 && user->errorindomainmf) || (fail++ > 36 && user->errorindomain)){
231: PetscMPIInt rank;
232: MPI_Comm_rank(PetscObjectComm((PetscObject)user->snes),&rank);
233: if (!rank) {
234: SNESSetFunctionDomainError(user->snes);
235: }
236: }
237: grashof = user->grashof;
238: prandtl = user->prandtl;
239: lid = user->lidvelocity;
241: /*
242: Define mesh intervals ratios for uniform grid.
244: Note: FD formulae below are normalized by multiplying through by
245: local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions.
248: */
249: dhx = (PetscReal)(info->mx-1); dhy = (PetscReal)(info->my-1);
250: hx = 1.0/dhx; hy = 1.0/dhy;
251: hxdhy = hx*dhy; hydhx = hy*dhx;
253: xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym;
255: /* Test whether we are on the bottom edge of the global array */
256: if (yints == 0) {
257: j = 0;
258: yints = yints + 1;
259: /* bottom edge */
260: for (i=info->xs; i<info->xs+info->xm; i++) {
261: f[j][i].u = x[j][i].u;
262: f[j][i].v = x[j][i].v;
263: f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
264: f[j][i].temp = x[j][i].temp-x[j+1][i].temp;
265: }
266: }
268: /* Test whether we are on the top edge of the global array */
269: if (yinte == info->my) {
270: j = info->my - 1;
271: yinte = yinte - 1;
272: /* top edge */
273: for (i=info->xs; i<info->xs+info->xm; i++) {
274: f[j][i].u = x[j][i].u - lid;
275: f[j][i].v = x[j][i].v;
276: f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
277: f[j][i].temp = x[j][i].temp-x[j-1][i].temp;
278: }
279: }
281: /* Test whether we are on the left edge of the global array */
282: if (xints == 0) {
283: i = 0;
284: xints = xints + 1;
285: /* left edge */
286: for (j=info->ys; j<info->ys+info->ym; j++) {
287: f[j][i].u = x[j][i].u;
288: f[j][i].v = x[j][i].v;
289: f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
290: f[j][i].temp = x[j][i].temp;
291: }
292: }
294: /* Test whether we are on the right edge of the global array */
295: if (xinte == info->mx) {
296: i = info->mx - 1;
297: xinte = xinte - 1;
298: /* right edge */
299: for (j=info->ys; j<info->ys+info->ym; j++) {
300: f[j][i].u = x[j][i].u;
301: f[j][i].v = x[j][i].v;
302: f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
303: f[j][i].temp = x[j][i].temp - (PetscReal)(grashof>0);
304: }
305: }
307: /* Compute over the interior points */
308: for (j=yints; j<yinte; j++) {
309: for (i=xints; i<xinte; i++) {
311: /*
312: convective coefficients for upwinding
313: */
314: vx = x[j][i].u; avx = PetscAbsScalar(vx);
315: vxp = .5*(vx+avx); vxm = .5*(vx-avx);
316: vy = x[j][i].v; avy = PetscAbsScalar(vy);
317: vyp = .5*(vy+avy); vym = .5*(vy-avy);
319: /* U velocity */
320: u = x[j][i].u;
321: uxx = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
322: uyy = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
323: f[j][i].u = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;
325: /* V velocity */
326: u = x[j][i].v;
327: uxx = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
328: uyy = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
329: f[j][i].v = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
331: /* Omega */
332: u = x[j][i].omega;
333: uxx = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
334: uyy = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
335: f[j][i].omega = uxx + uyy + (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
336: (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
337: .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy;
339: /* Temperature */
340: u = x[j][i].temp;
341: uxx = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
342: uyy = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
343: f[j][i].temp = uxx + uyy + prandtl*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy +
344: (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx);
345: }
346: }
348: /*
349: Flop count (multiply-adds are counted as 2 operations)
350: */
351: PetscLogFlops(84.0*info->ym*info->xm);
352: return(0);
353: }
357: PetscErrorCode MatMult_MyShell(Mat A,Vec x,Vec y)
358: {
359: PetscErrorCode ierr;
360: MatShellCtx *matshellctx;
361: static PetscInt fail = 0;
364: MatShellGetContext(A,&matshellctx);
365: MatMult(matshellctx->Jmf,x,y);
366: if (fail++ > 5) {
367: PetscMPIInt rank;
368: MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank);
369: if (!rank) {VecSetInf(y);}
370: }
371: return(0);
372: }
376: PetscErrorCode MatAssemblyEnd_MyShell(Mat A,MatAssemblyType tp)
377: {
379: MatShellCtx *matshellctx;
382: MatShellGetContext(A,&matshellctx);
383: MatAssemblyEnd(matshellctx->Jmf,tp);
384: return(0);
385: }
389: PetscErrorCode PCApply_MyShell(PC pc,Vec x,Vec y)
390: {
392: static PetscInt fail = 0;
395: VecCopy(x,y);
396: if (fail++ > 3) {
397: PetscMPIInt rank;
398: MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);
399: if (!rank) {VecSetInf(y);}
400: }
401: return(0);
402: }
404: extern PetscErrorCode SNESComputeJacobian_DMDA(SNES,Vec,Mat,Mat,void*);
408: PetscErrorCode SNESComputeJacobian_MyShell(SNES snes,Vec X,Mat A,Mat B,void *ctx)
409: {
410: static PetscInt fail = 0;
411: PetscErrorCode ierr;
414: SNESComputeJacobian_DMDA(snes,X,A,B,ctx);
415: if (fail++ > 0) {
416: MatZeroEntries(A);
417: }
418: return(0);
419: }