xref: /petsc/src/snes/tests/ex69.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Tests recovery from domain errors in MatMult() and PCApply()\n\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown       See src/ksp/ksp/tutorials/ex19.c from which this was copied
6c4762a1bSJed Brown */
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include <petscsnes.h>
9c4762a1bSJed Brown #include <petscdm.h>
10c4762a1bSJed Brown #include <petscdmda.h>
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown    User-defined routines and data structures
14c4762a1bSJed Brown */
15c4762a1bSJed Brown typedef struct {
16c4762a1bSJed Brown   PetscScalar u,v,omega,temp;
17c4762a1bSJed Brown } Field;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo*,Field**,Field**,void*);
20c4762a1bSJed Brown 
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscReal   lidvelocity,prandtl,grashof;  /* physical parameters */
23c4762a1bSJed Brown   PetscBool   draw_contours;                /* flag - 1 indicates drawing contours */
24c4762a1bSJed Brown   PetscBool   errorindomain;
25c4762a1bSJed Brown   PetscBool   errorindomainmf;
26c4762a1bSJed Brown   SNES        snes;
27c4762a1bSJed Brown } AppCtx;
28c4762a1bSJed Brown 
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   Mat Jmf;
31c4762a1bSJed Brown } MatShellCtx;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx*,DM,Vec);
34c4762a1bSJed Brown extern PetscErrorCode MatMult_MyShell(Mat,Vec,Vec);
35c4762a1bSJed Brown extern PetscErrorCode MatAssemblyEnd_MyShell(Mat,MatAssemblyType);
36c4762a1bSJed Brown extern PetscErrorCode PCApply_MyShell(PC,Vec,Vec);
37c4762a1bSJed Brown extern PetscErrorCode SNESComputeJacobian_MyShell(SNES,Vec,Mat,Mat,void*);
38c4762a1bSJed Brown 
39c4762a1bSJed Brown int main(int argc,char **argv)
40c4762a1bSJed Brown {
41c4762a1bSJed Brown   AppCtx         user;                /* user-defined work context */
42c4762a1bSJed Brown   PetscInt       mx,my;
43c4762a1bSJed Brown   PetscErrorCode ierr;
44c4762a1bSJed Brown   MPI_Comm       comm;
45c4762a1bSJed Brown   DM             da;
46c4762a1bSJed Brown   Vec            x;
47c4762a1bSJed Brown   Mat            J = NULL,Jmf = NULL;
48c4762a1bSJed Brown   MatShellCtx    matshellctx;
49c4762a1bSJed Brown   PetscInt       mlocal,nlocal;
50c4762a1bSJed Brown   PC             pc;
51c4762a1bSJed Brown   KSP            ksp;
52c4762a1bSJed Brown   PetscBool      errorinmatmult = PETSC_FALSE,errorinpcapply = PETSC_FALSE,errorinpcsetup = PETSC_FALSE;
53c4762a1bSJed Brown 
54*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
555f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-error_in_matmult",&errorinmatmult,NULL));
565f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-error_in_pcapply",&errorinpcapply,NULL));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-error_in_pcsetup",&errorinpcsetup,NULL));
58c4762a1bSJed Brown   user.errorindomain = PETSC_FALSE;
595f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-error_in_domain",&user.errorindomain,NULL));
60c4762a1bSJed Brown   user.errorindomainmf = PETSC_FALSE;
615f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-error_in_domainmf",&user.errorindomainmf,NULL));
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
645f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(comm,&user.snes));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /*
67c4762a1bSJed Brown       Create distributed array object to manage parallel grid and vectors
68c4762a1bSJed Brown       for principal unknowns (x) and governing residuals (f)
69c4762a1bSJed Brown   */
705f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,4,1,0,0,&da));
715f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
735f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(user.snes,da));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   ierr = DMDAGetInfo(da,0,&mx,&my,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
76c4762a1bSJed Brown                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
77c4762a1bSJed Brown   /*
78c4762a1bSJed Brown      Problem parameters (velocity of lid, prandtl, and grashof numbers)
79c4762a1bSJed Brown   */
80c4762a1bSJed Brown   user.lidvelocity = 1.0/(mx*my);
81c4762a1bSJed Brown   user.prandtl     = 1.0;
82c4762a1bSJed Brown   user.grashof     = 1.0;
83c4762a1bSJed Brown 
845f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-lidvelocity",&user.lidvelocity,NULL));
855f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-prandtl",&user.prandtl,NULL));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-grashof",&user.grashof,NULL));
875f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsHasName(NULL,NULL,"-contours",&user.draw_contours));
88c4762a1bSJed Brown 
895f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,0,"x_velocity"));
905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,1,"y_velocity"));
915f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,2,"Omega"));
925f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,3,"temperature"));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95c4762a1bSJed Brown      Create user context, set problem data, create vector data structures.
96c4762a1bSJed Brown      Also, compute the initial guess.
97c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown      Create nonlinear solver context
101c4762a1bSJed Brown 
102c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetApplicationContext(da,&user));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASNESSetFunctionLocal(da,INSERT_VALUES,(PetscErrorCode (*)(DMDALocalInfo*,void*,void*,void*))FormFunctionLocal,&user));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   if (errorinmatmult) {
1075f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateSNESMF(user.snes,&Jmf));
1085f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetFromOptions(Jmf));
1095f80ce2aSJacob Faibussowitsch     CHKERRQ(MatGetLocalSize(Jmf,&mlocal,&nlocal));
110c4762a1bSJed Brown     matshellctx.Jmf = Jmf;
1115f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateShell(PetscObjectComm((PetscObject)Jmf),mlocal,nlocal,PETSC_DECIDE,PETSC_DECIDE,&matshellctx,&J));
1125f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShellSetOperation(J,MATOP_MULT,(void (*)(void))MatMult_MyShell));
1135f80ce2aSJacob Faibussowitsch     CHKERRQ(MatShellSetOperation(J,MATOP_ASSEMBLY_END,(void (*)(void))MatAssemblyEnd_MyShell));
1145f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetJacobian(user.snes,J,J,MatMFFDComputeJacobian,NULL));
115c4762a1bSJed Brown   }
116c4762a1bSJed Brown 
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(user.snes));
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"lid velocity = %g, prandtl # = %g, grashof # = %g\n",(double)user.lidvelocity,(double)user.prandtl,(double)user.grashof));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   if (errorinpcapply) {
1215f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetKSP(user.snes,&ksp));
1225f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
1235f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCSHELL));
1245f80ce2aSJacob Faibussowitsch     CHKERRQ(PCShellSetApply(pc,PCApply_MyShell));
125c4762a1bSJed Brown   }
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128c4762a1bSJed Brown      Solve the nonlinear system
129c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(&user,da,x));
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   if (errorinpcsetup) {
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetUp(user.snes));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetJacobian(user.snes,NULL,NULL,SNESComputeJacobian_MyShell,NULL));
136c4762a1bSJed Brown   }
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(user.snes,NULL,x));
138c4762a1bSJed Brown 
139c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
141c4762a1bSJed Brown      are no longer needed.
142c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Jmf));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&user.snes));
148*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
149*b122ec5aSJacob Faibussowitsch   return 0;
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /* ------------------------------------------------------------------- */
153c4762a1bSJed Brown 
154c4762a1bSJed Brown /*
155c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    Input Parameters:
158c4762a1bSJed Brown    user - user-defined application context
159c4762a1bSJed Brown    X - vector
160c4762a1bSJed Brown 
161c4762a1bSJed Brown    Output Parameter:
162c4762a1bSJed Brown    X - vector
163c4762a1bSJed Brown */
164c4762a1bSJed Brown PetscErrorCode FormInitialGuess(AppCtx *user,DM da,Vec X)
165c4762a1bSJed Brown {
166c4762a1bSJed Brown   PetscInt       i,j,mx,xs,ys,xm,ym;
167c4762a1bSJed Brown   PetscReal      grashof,dx;
168c4762a1bSJed Brown   Field          **x;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   PetscFunctionBeginUser;
171c4762a1bSJed Brown   grashof = user->grashof;
172c4762a1bSJed Brown 
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,0,&mx,0,0,0,0,0,0,0,0,0,0,0));
174c4762a1bSJed Brown   dx   = 1.0/(mx-1);
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /*
177c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
178c4762a1bSJed Brown        xs, ys   - starting grid indices (no ghost points)
179c4762a1bSJed Brown        xm, ym   - widths of local grid (no ghost points)
180c4762a1bSJed Brown   */
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   /*
184c4762a1bSJed Brown      Get a pointer to vector data.
185c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
186c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
187c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
188c4762a1bSJed Brown          the array.
189c4762a1bSJed Brown   */
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,X,&x));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /*
193c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
194c4762a1bSJed Brown      Initial condition is motionless fluid and equilibrium temperature
195c4762a1bSJed Brown   */
196c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
197c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
198c4762a1bSJed Brown       x[j][i].u     = 0.0;
199c4762a1bSJed Brown       x[j][i].v     = 0.0;
200c4762a1bSJed Brown       x[j][i].omega = 0.0;
201c4762a1bSJed Brown       x[j][i].temp  = (grashof>0)*i*dx;
202c4762a1bSJed Brown     }
203c4762a1bSJed Brown   }
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /*
206c4762a1bSJed Brown      Restore vector
207c4762a1bSJed Brown   */
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,X,&x));
209c4762a1bSJed Brown   PetscFunctionReturn(0);
210c4762a1bSJed Brown }
211c4762a1bSJed Brown 
212c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,Field **x,Field **f,void *ptr)
213c4762a1bSJed Brown {
214c4762a1bSJed Brown   AppCtx          *user = (AppCtx*)ptr;
215c4762a1bSJed Brown   PetscInt        xints,xinte,yints,yinte,i,j;
216c4762a1bSJed Brown   PetscReal       hx,hy,dhx,dhy,hxdhy,hydhx;
217c4762a1bSJed Brown   PetscReal       grashof,prandtl,lid;
218c4762a1bSJed Brown   PetscScalar     u,uxx,uyy,vx,vy,avx,avy,vxp,vxm,vyp,vym;
219c4762a1bSJed Brown   static PetscInt fail = 0;
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   PetscFunctionBeginUser;
222c4762a1bSJed Brown   if ((fail++ > 7 && user->errorindomainmf) || (fail++ > 36 && user->errorindomain)) {
223c4762a1bSJed Brown     PetscMPIInt rank;
2245f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)user->snes),&rank));
225dd400576SPatrick Sanan     if (rank == 0) {
2265f80ce2aSJacob Faibussowitsch       CHKERRQ(SNESSetFunctionDomainError(user->snes));
227c4762a1bSJed Brown     }
228c4762a1bSJed Brown   }
229c4762a1bSJed Brown   grashof = user->grashof;
230c4762a1bSJed Brown   prandtl = user->prandtl;
231c4762a1bSJed Brown   lid     = user->lidvelocity;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /*
234c4762a1bSJed Brown      Define mesh intervals ratios for uniform grid.
235c4762a1bSJed Brown 
236c4762a1bSJed Brown      Note: FD formulae below are normalized by multiplying through by
237c4762a1bSJed Brown      local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions.
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   */
240c4762a1bSJed Brown   dhx   = (PetscReal)(info->mx-1);  dhy = (PetscReal)(info->my-1);
241c4762a1bSJed Brown   hx    = 1.0/dhx;                   hy = 1.0/dhy;
242c4762a1bSJed Brown   hxdhy = hx*dhy;                 hydhx = hy*dhx;
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   xints = info->xs; xinte = info->xs+info->xm; yints = info->ys; yinte = info->ys+info->ym;
245c4762a1bSJed Brown 
246c4762a1bSJed Brown   /* Test whether we are on the bottom edge of the global array */
247c4762a1bSJed Brown   if (yints == 0) {
248c4762a1bSJed Brown     j     = 0;
249c4762a1bSJed Brown     yints = yints + 1;
250c4762a1bSJed Brown     /* bottom edge */
251c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
252c4762a1bSJed Brown       f[j][i].u     = x[j][i].u;
253c4762a1bSJed Brown       f[j][i].v     = x[j][i].v;
254c4762a1bSJed Brown       f[j][i].omega = x[j][i].omega + (x[j+1][i].u - x[j][i].u)*dhy;
255c4762a1bSJed Brown       f[j][i].temp  = x[j][i].temp-x[j+1][i].temp;
256c4762a1bSJed Brown     }
257c4762a1bSJed Brown   }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* Test whether we are on the top edge of the global array */
260c4762a1bSJed Brown   if (yinte == info->my) {
261c4762a1bSJed Brown     j     = info->my - 1;
262c4762a1bSJed Brown     yinte = yinte - 1;
263c4762a1bSJed Brown     /* top edge */
264c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
265c4762a1bSJed Brown       f[j][i].u     = x[j][i].u - lid;
266c4762a1bSJed Brown       f[j][i].v     = x[j][i].v;
267c4762a1bSJed Brown       f[j][i].omega = x[j][i].omega + (x[j][i].u - x[j-1][i].u)*dhy;
268c4762a1bSJed Brown       f[j][i].temp  = x[j][i].temp-x[j-1][i].temp;
269c4762a1bSJed Brown     }
270c4762a1bSJed Brown   }
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   /* Test whether we are on the left edge of the global array */
273c4762a1bSJed Brown   if (xints == 0) {
274c4762a1bSJed Brown     i     = 0;
275c4762a1bSJed Brown     xints = xints + 1;
276c4762a1bSJed Brown     /* left edge */
277c4762a1bSJed Brown     for (j=info->ys; j<info->ys+info->ym; j++) {
278c4762a1bSJed Brown       f[j][i].u     = x[j][i].u;
279c4762a1bSJed Brown       f[j][i].v     = x[j][i].v;
280c4762a1bSJed Brown       f[j][i].omega = x[j][i].omega - (x[j][i+1].v - x[j][i].v)*dhx;
281c4762a1bSJed Brown       f[j][i].temp  = x[j][i].temp;
282c4762a1bSJed Brown     }
283c4762a1bSJed Brown   }
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   /* Test whether we are on the right edge of the global array */
286c4762a1bSJed Brown   if (xinte == info->mx) {
287c4762a1bSJed Brown     i     = info->mx - 1;
288c4762a1bSJed Brown     xinte = xinte - 1;
289c4762a1bSJed Brown     /* right edge */
290c4762a1bSJed Brown     for (j=info->ys; j<info->ys+info->ym; j++) {
291c4762a1bSJed Brown       f[j][i].u     = x[j][i].u;
292c4762a1bSJed Brown       f[j][i].v     = x[j][i].v;
293c4762a1bSJed Brown       f[j][i].omega = x[j][i].omega - (x[j][i].v - x[j][i-1].v)*dhx;
294c4762a1bSJed Brown       f[j][i].temp  = x[j][i].temp - (PetscReal)(grashof>0);
295c4762a1bSJed Brown     }
296c4762a1bSJed Brown   }
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   /* Compute over the interior points */
299c4762a1bSJed Brown   for (j=yints; j<yinte; j++) {
300c4762a1bSJed Brown     for (i=xints; i<xinte; i++) {
301c4762a1bSJed Brown 
302c4762a1bSJed Brown       /*
303c4762a1bSJed Brown        convective coefficients for upwinding
304c4762a1bSJed Brown       */
305c4762a1bSJed Brown       vx  = x[j][i].u; avx = PetscAbsScalar(vx);
306c4762a1bSJed Brown       vxp = .5*(vx+avx); vxm = .5*(vx-avx);
307c4762a1bSJed Brown       vy  = x[j][i].v; avy = PetscAbsScalar(vy);
308c4762a1bSJed Brown       vyp = .5*(vy+avy); vym = .5*(vy-avy);
309c4762a1bSJed Brown 
310c4762a1bSJed Brown       /* U velocity */
311c4762a1bSJed Brown       u         = x[j][i].u;
312c4762a1bSJed Brown       uxx       = (2.0*u - x[j][i-1].u - x[j][i+1].u)*hydhx;
313c4762a1bSJed Brown       uyy       = (2.0*u - x[j-1][i].u - x[j+1][i].u)*hxdhy;
314c4762a1bSJed Brown       f[j][i].u = uxx + uyy - .5*(x[j+1][i].omega-x[j-1][i].omega)*hx;
315c4762a1bSJed Brown 
316c4762a1bSJed Brown       /* V velocity */
317c4762a1bSJed Brown       u         = x[j][i].v;
318c4762a1bSJed Brown       uxx       = (2.0*u - x[j][i-1].v - x[j][i+1].v)*hydhx;
319c4762a1bSJed Brown       uyy       = (2.0*u - x[j-1][i].v - x[j+1][i].v)*hxdhy;
320c4762a1bSJed Brown       f[j][i].v = uxx + uyy + .5*(x[j][i+1].omega-x[j][i-1].omega)*hy;
321c4762a1bSJed Brown 
322c4762a1bSJed Brown       /* Omega */
323c4762a1bSJed Brown       u             = x[j][i].omega;
324c4762a1bSJed Brown       uxx           = (2.0*u - x[j][i-1].omega - x[j][i+1].omega)*hydhx;
325c4762a1bSJed Brown       uyy           = (2.0*u - x[j-1][i].omega - x[j+1][i].omega)*hxdhy;
326c4762a1bSJed Brown       f[j][i].omega = uxx + uyy + (vxp*(u - x[j][i-1].omega) + vxm*(x[j][i+1].omega - u))*hy +
327c4762a1bSJed Brown                       (vyp*(u - x[j-1][i].omega) + vym*(x[j+1][i].omega - u))*hx -
328c4762a1bSJed Brown                       .5*grashof*(x[j][i+1].temp - x[j][i-1].temp)*hy;
329c4762a1bSJed Brown 
330c4762a1bSJed Brown       /* Temperature */
331c4762a1bSJed Brown       u            = x[j][i].temp;
332c4762a1bSJed Brown       uxx          = (2.0*u - x[j][i-1].temp - x[j][i+1].temp)*hydhx;
333c4762a1bSJed Brown       uyy          = (2.0*u - x[j-1][i].temp - x[j+1][i].temp)*hxdhy;
334c4762a1bSJed Brown       f[j][i].temp =  uxx + uyy  + prandtl*((vxp*(u - x[j][i-1].temp) + vxm*(x[j][i+1].temp - u))*hy +
335c4762a1bSJed Brown                                             (vyp*(u - x[j-1][i].temp) + vym*(x[j+1][i].temp - u))*hx);
336c4762a1bSJed Brown     }
337c4762a1bSJed Brown   }
338c4762a1bSJed Brown 
339c4762a1bSJed Brown   /*
340c4762a1bSJed Brown      Flop count (multiply-adds are counted as 2 operations)
341c4762a1bSJed Brown   */
3425f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(84.0*info->ym*info->xm));
343c4762a1bSJed Brown   PetscFunctionReturn(0);
344c4762a1bSJed Brown }
345c4762a1bSJed Brown 
346c4762a1bSJed Brown PetscErrorCode MatMult_MyShell(Mat A,Vec x,Vec y)
347c4762a1bSJed Brown {
348c4762a1bSJed Brown   MatShellCtx     *matshellctx;
349c4762a1bSJed Brown   static PetscInt fail = 0;
350c4762a1bSJed Brown 
351c4762a1bSJed Brown   PetscFunctionBegin;
3525f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&matshellctx));
3535f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(matshellctx->Jmf,x,y));
354c4762a1bSJed Brown   if (fail++ > 5) {
355c4762a1bSJed Brown     PetscMPIInt rank;
3565f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A),&rank));
3575f80ce2aSJacob Faibussowitsch     if (rank == 0) CHKERRQ(VecSetInf(y));
358c4762a1bSJed Brown   }
359c4762a1bSJed Brown   PetscFunctionReturn(0);
360c4762a1bSJed Brown }
361c4762a1bSJed Brown 
362c4762a1bSJed Brown PetscErrorCode MatAssemblyEnd_MyShell(Mat A,MatAssemblyType tp)
363c4762a1bSJed Brown {
364c4762a1bSJed Brown   MatShellCtx    *matshellctx;
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   PetscFunctionBegin;
3675f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShellGetContext(A,&matshellctx));
3685f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(matshellctx->Jmf,tp));
369c4762a1bSJed Brown   PetscFunctionReturn(0);
370c4762a1bSJed Brown }
371c4762a1bSJed Brown 
372c4762a1bSJed Brown PetscErrorCode PCApply_MyShell(PC pc,Vec x,Vec y)
373c4762a1bSJed Brown {
374c4762a1bSJed Brown   static PetscInt fail = 0;
375c4762a1bSJed Brown 
376c4762a1bSJed Brown   PetscFunctionBegin;
3775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(x,y));
378c4762a1bSJed Brown   if (fail++ > 3) {
379c4762a1bSJed Brown     PetscMPIInt rank;
3805f80ce2aSJacob Faibussowitsch     CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank));
3815f80ce2aSJacob Faibussowitsch     if (rank == 0) CHKERRQ(VecSetInf(y));
382c4762a1bSJed Brown   }
383c4762a1bSJed Brown   PetscFunctionReturn(0);
384c4762a1bSJed Brown }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES,Vec,Mat,Mat,void*);
387c4762a1bSJed Brown 
388c4762a1bSJed Brown PetscErrorCode SNESComputeJacobian_MyShell(SNES snes,Vec X,Mat A,Mat B,void *ctx)
389c4762a1bSJed Brown {
390c4762a1bSJed Brown   static PetscInt fail = 0;
391c4762a1bSJed Brown 
392c4762a1bSJed Brown   PetscFunctionBegin;
3935f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESComputeJacobian_DMDA(snes,X,A,B,ctx));
394c4762a1bSJed Brown   if (fail++ > 0) {
3955f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(A));
396c4762a1bSJed Brown   }
397c4762a1bSJed Brown   PetscFunctionReturn(0);
398c4762a1bSJed Brown }
399c4762a1bSJed Brown 
400c4762a1bSJed Brown /*TEST
401c4762a1bSJed Brown 
402c4762a1bSJed Brown    test:
403c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason
404c4762a1bSJed Brown 
405c4762a1bSJed Brown    test:
406c4762a1bSJed Brown       suffix: 2
407c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_matmult
408c4762a1bSJed Brown 
409c4762a1bSJed Brown    test:
410c4762a1bSJed Brown       suffix: 3
411c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_pcapply
412c4762a1bSJed Brown 
413c4762a1bSJed Brown    test:
414c4762a1bSJed Brown       suffix: 4
415c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup
416c4762a1bSJed Brown 
417c4762a1bSJed Brown    test:
418c4762a1bSJed Brown       suffix: 5
419c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -pc_type bjacobi
420c4762a1bSJed Brown 
421c4762a1bSJed Brown    test:
422c4762a1bSJed Brown       suffix: 5_fieldsplit
423c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -pc_type fieldsplit
424c4762a1bSJed Brown       output_file: output/ex69_5.out
425c4762a1bSJed Brown 
426c4762a1bSJed Brown    test:
427c4762a1bSJed Brown       suffix: 6
428c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_domainmf -snes_mf -pc_type none
429c4762a1bSJed Brown 
430c4762a1bSJed Brown    test:
431c4762a1bSJed Brown       suffix: 7
432c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_domain
433c4762a1bSJed Brown 
434c4762a1bSJed Brown    test:
435c4762a1bSJed Brown       suffix: 8
436c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_domain -snes_mf -pc_type none
437c4762a1bSJed Brown 
438c4762a1bSJed Brown TEST*/
439