xref: /petsc/src/snes/tests/ex69.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Tests recovery from domain errors in MatMult() and PCApply()\n\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown       See src/ksp/ksp/tutorials/ex19.c from which this was copied
5c4762a1bSJed Brown */
6c4762a1bSJed Brown 
7c4762a1bSJed Brown #include <petscsnes.h>
8c4762a1bSJed Brown #include <petscdm.h>
9c4762a1bSJed Brown #include <petscdmda.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown /*
12c4762a1bSJed Brown    User-defined routines and data structures
13c4762a1bSJed Brown */
14c4762a1bSJed Brown typedef struct {
15c4762a1bSJed Brown   PetscScalar u, v, omega, temp;
16c4762a1bSJed Brown } Field;
17c4762a1bSJed Brown 
18c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *, Field **, Field **, void *);
19c4762a1bSJed Brown 
20c4762a1bSJed Brown typedef struct {
21c4762a1bSJed Brown   PetscReal lidvelocity, prandtl, grashof; /* physical parameters */
22c4762a1bSJed Brown   PetscBool draw_contours;                 /* flag - 1 indicates drawing contours */
23c4762a1bSJed Brown   PetscBool errorindomain;
24c4762a1bSJed Brown   PetscBool errorindomainmf;
25c4762a1bSJed Brown   SNES      snes;
26c4762a1bSJed Brown } AppCtx;
27c4762a1bSJed Brown 
28c4762a1bSJed Brown typedef struct {
29c4762a1bSJed Brown   Mat Jmf;
30c4762a1bSJed Brown } MatShellCtx;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(AppCtx *, DM, Vec);
33c4762a1bSJed Brown extern PetscErrorCode MatMult_MyShell(Mat, Vec, Vec);
34c4762a1bSJed Brown extern PetscErrorCode MatAssemblyEnd_MyShell(Mat, MatAssemblyType);
35c4762a1bSJed Brown extern PetscErrorCode PCApply_MyShell(PC, Vec, Vec);
36c4762a1bSJed Brown extern PetscErrorCode SNESComputeJacobian_MyShell(SNES, Vec, Mat, Mat, void *);
37c4762a1bSJed Brown 
main(int argc,char ** argv)38d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
39d71ae5a4SJacob Faibussowitsch {
40c4762a1bSJed Brown   AppCtx      user; /* user-defined work context */
41c4762a1bSJed Brown   PetscInt    mx, my;
42c4762a1bSJed Brown   MPI_Comm    comm;
43c4762a1bSJed Brown   DM          da;
44c4762a1bSJed Brown   Vec         x;
45c4762a1bSJed Brown   Mat         J = NULL, Jmf = NULL;
46c4762a1bSJed Brown   MatShellCtx matshellctx;
47c4762a1bSJed Brown   PetscInt    mlocal, nlocal;
48c4762a1bSJed Brown   PC          pc;
49c4762a1bSJed Brown   KSP         ksp;
50c4762a1bSJed Brown   PetscBool   errorinmatmult = PETSC_FALSE, errorinpcapply = PETSC_FALSE, errorinpcsetup = PETSC_FALSE;
51c4762a1bSJed Brown 
52327415f7SBarry Smith   PetscFunctionBeginUser;
53c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
549566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_matmult", &errorinmatmult, NULL));
559566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_pcapply", &errorinpcapply, NULL));
569566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_pcsetup", &errorinpcsetup, NULL));
57c4762a1bSJed Brown   user.errorindomain = PETSC_FALSE;
589566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_domain", &user.errorindomain, NULL));
59c4762a1bSJed Brown   user.errorindomainmf = PETSC_FALSE;
609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-error_in_domainmf", &user.errorindomainmf, NULL));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   comm = PETSC_COMM_WORLD;
639566063dSJacob Faibussowitsch   PetscCall(SNESCreate(comm, &user.snes));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /*
66c4762a1bSJed Brown       Create distributed array object to manage parallel grid and vectors
67c4762a1bSJed Brown       for principal unknowns (x) and governing residuals (f)
68c4762a1bSJed Brown   */
699566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE, DMDA_STENCIL_STAR, 4, 4, PETSC_DECIDE, PETSC_DECIDE, 4, 1, 0, 0, &da));
709566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
719566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
729566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(user.snes, da));
73c4762a1bSJed Brown 
749371c9d4SSatish Balay   PetscCall(DMDAGetInfo(da, 0, &mx, &my, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
75c4762a1bSJed Brown   /*
76c4762a1bSJed Brown      Problem parameters (velocity of lid, prandtl, and grashof numbers)
77c4762a1bSJed Brown   */
78c4762a1bSJed Brown   user.lidvelocity = 1.0 / (mx * my);
79c4762a1bSJed Brown   user.prandtl     = 1.0;
80c4762a1bSJed Brown   user.grashof     = 1.0;
81c4762a1bSJed Brown 
829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-lidvelocity", &user.lidvelocity, NULL));
839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-prandtl", &user.prandtl, NULL));
849566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-grashof", &user.grashof, NULL));
859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsHasName(NULL, NULL, "-contours", &user.draw_contours));
86c4762a1bSJed Brown 
879566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "x_velocity"));
889566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "y_velocity"));
899566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 2, "Omega"));
909566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 3, "temperature"));
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93c4762a1bSJed Brown      Create user context, set problem data, create vector data structures.
94c4762a1bSJed Brown      Also, compute the initial guess.
95c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown      Create nonlinear solver context
99c4762a1bSJed Brown 
100c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1019566063dSJacob Faibussowitsch   PetscCall(DMSetApplicationContext(da, &user));
1029566063dSJacob Faibussowitsch   PetscCall(DMDASNESSetFunctionLocal(da, INSERT_VALUES, (PetscErrorCode (*)(DMDALocalInfo *, void *, void *, void *))FormFunctionLocal, &user));
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   if (errorinmatmult) {
1059566063dSJacob Faibussowitsch     PetscCall(MatCreateSNESMF(user.snes, &Jmf));
1069566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(Jmf));
1079566063dSJacob Faibussowitsch     PetscCall(MatGetLocalSize(Jmf, &mlocal, &nlocal));
108c4762a1bSJed Brown     matshellctx.Jmf = Jmf;
1099566063dSJacob Faibussowitsch     PetscCall(MatCreateShell(PetscObjectComm((PetscObject)Jmf), mlocal, nlocal, PETSC_DECIDE, PETSC_DECIDE, &matshellctx, &J));
11057d50842SBarry Smith     PetscCall(MatShellSetOperation(J, MATOP_MULT, (PetscErrorCodeFn *)MatMult_MyShell));
11157d50842SBarry Smith     PetscCall(MatShellSetOperation(J, MATOP_ASSEMBLY_END, (PetscErrorCodeFn *)MatAssemblyEnd_MyShell));
1129566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(user.snes, J, J, MatMFFDComputeJacobian, NULL));
113c4762a1bSJed Brown   }
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(user.snes));
1169566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(comm, "lid velocity = %g, prandtl # = %g, grashof # = %g\n", (double)user.lidvelocity, (double)user.prandtl, (double)user.grashof));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   if (errorinpcapply) {
1199566063dSJacob Faibussowitsch     PetscCall(SNESGetKSP(user.snes, &ksp));
1209566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
1219566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCSHELL));
1229566063dSJacob Faibussowitsch     PetscCall(PCShellSetApply(pc, PCApply_MyShell));
123c4762a1bSJed Brown   }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126c4762a1bSJed Brown      Solve the nonlinear system
127c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1289566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &x));
1299566063dSJacob Faibussowitsch   PetscCall(FormInitialGuess(&user, da, x));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   if (errorinpcsetup) {
1329566063dSJacob Faibussowitsch     PetscCall(SNESSetUp(user.snes));
1339566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(user.snes, NULL, NULL, SNESComputeJacobian_MyShell, NULL));
134c4762a1bSJed Brown   }
1359566063dSJacob Faibussowitsch   PetscCall(SNESSolve(user.snes, NULL, x));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
139c4762a1bSJed Brown      are no longer needed.
140c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1419566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1429566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jmf));
1439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1449566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1459566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&user.snes));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
147b122ec5aSJacob Faibussowitsch   return 0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown /*
151c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    Input Parameters:
154c4762a1bSJed Brown    user - user-defined application context
155c4762a1bSJed Brown    X - vector
156c4762a1bSJed Brown 
157c4762a1bSJed Brown    Output Parameter:
158c4762a1bSJed Brown    X - vector
159c4762a1bSJed Brown */
FormInitialGuess(AppCtx * user,DM da,Vec X)160d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(AppCtx *user, DM da, Vec X)
161d71ae5a4SJacob Faibussowitsch {
162c4762a1bSJed Brown   PetscInt  i, j, mx, xs, ys, xm, ym;
163c4762a1bSJed Brown   PetscReal grashof, dx;
164c4762a1bSJed Brown   Field   **x;
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   PetscFunctionBeginUser;
167c4762a1bSJed Brown   grashof = user->grashof;
168c4762a1bSJed Brown 
1699566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, 0, &mx, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0));
170c4762a1bSJed Brown   dx = 1.0 / (mx - 1);
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /*
173c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
174c4762a1bSJed Brown        xs, ys   - starting grid indices (no ghost points)
175c4762a1bSJed Brown        xm, ym   - widths of local grid (no ghost points)
176c4762a1bSJed Brown   */
1779566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /*
180c4762a1bSJed Brown      Get a pointer to vector data.
181c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
182c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
183c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
184c4762a1bSJed Brown          the array.
185c4762a1bSJed Brown   */
1869566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /*
189c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
190c4762a1bSJed Brown      Initial condition is motionless fluid and equilibrium temperature
191c4762a1bSJed Brown   */
192c4762a1bSJed Brown   for (j = ys; j < ys + ym; j++) {
193c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
194c4762a1bSJed Brown       x[j][i].u     = 0.0;
195c4762a1bSJed Brown       x[j][i].v     = 0.0;
196c4762a1bSJed Brown       x[j][i].omega = 0.0;
197c4762a1bSJed Brown       x[j][i].temp  = (grashof > 0) * i * dx;
198c4762a1bSJed Brown     }
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /*
202c4762a1bSJed Brown      Restore vector
203c4762a1bSJed Brown   */
2049566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
2053ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
206c4762a1bSJed Brown }
207c4762a1bSJed Brown 
FormFunctionLocal(DMDALocalInfo * info,Field ** x,Field ** f,void * ptr)208d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunctionLocal(DMDALocalInfo *info, Field **x, Field **f, void *ptr)
209d71ae5a4SJacob Faibussowitsch {
210c4762a1bSJed Brown   AppCtx         *user = (AppCtx *)ptr;
211c4762a1bSJed Brown   PetscInt        xints, xinte, yints, yinte, i, j;
212c4762a1bSJed Brown   PetscReal       hx, hy, dhx, dhy, hxdhy, hydhx;
213c4762a1bSJed Brown   PetscReal       grashof, prandtl, lid;
214c4762a1bSJed Brown   PetscScalar     u, uxx, uyy, vx, vy, avx, avy, vxp, vxm, vyp, vym;
215c4762a1bSJed Brown   static PetscInt fail = 0;
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   PetscFunctionBeginUser;
218c4762a1bSJed Brown   if ((fail++ > 7 && user->errorindomainmf) || (fail++ > 36 && user->errorindomain)) {
219c4762a1bSJed Brown     PetscMPIInt rank;
2209566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)user->snes), &rank));
22148a46eb9SPierre Jolivet     if (rank == 0) PetscCall(SNESSetFunctionDomainError(user->snes));
222c4762a1bSJed Brown   }
223c4762a1bSJed Brown   grashof = user->grashof;
224c4762a1bSJed Brown   prandtl = user->prandtl;
225c4762a1bSJed Brown   lid     = user->lidvelocity;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /*
228c4762a1bSJed Brown      Define mesh intervals ratios for uniform grid.
229c4762a1bSJed Brown 
230c4762a1bSJed Brown      Note: FD formulae below are normalized by multiplying through by
231c4762a1bSJed Brown      local volume element (i.e. hx*hy) to obtain coefficients O(1) in two dimensions.
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   */
2349371c9d4SSatish Balay   dhx   = (PetscReal)(info->mx - 1);
2359371c9d4SSatish Balay   dhy   = (PetscReal)(info->my - 1);
2369371c9d4SSatish Balay   hx    = 1.0 / dhx;
2379371c9d4SSatish Balay   hy    = 1.0 / dhy;
2389371c9d4SSatish Balay   hxdhy = hx * dhy;
2399371c9d4SSatish Balay   hydhx = hy * dhx;
240c4762a1bSJed Brown 
2419371c9d4SSatish Balay   xints = info->xs;
2429371c9d4SSatish Balay   xinte = info->xs + info->xm;
2439371c9d4SSatish Balay   yints = info->ys;
2449371c9d4SSatish Balay   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        convective coefficients for upwinding
303c4762a1bSJed Brown       */
3049371c9d4SSatish Balay       vx  = x[j][i].u;
3059371c9d4SSatish Balay       avx = PetscAbsScalar(vx);
3069371c9d4SSatish Balay       vxp = .5 * (vx + avx);
3079371c9d4SSatish Balay       vxm = .5 * (vx - avx);
3089371c9d4SSatish Balay       vy  = x[j][i].v;
3099371c9d4SSatish Balay       avy = PetscAbsScalar(vy);
3109371c9d4SSatish Balay       vyp = .5 * (vy + avy);
3119371c9d4SSatish Balay       vym = .5 * (vy - avy);
312c4762a1bSJed Brown 
313c4762a1bSJed Brown       /* U velocity */
314c4762a1bSJed Brown       u         = x[j][i].u;
315c4762a1bSJed Brown       uxx       = (2.0 * u - x[j][i - 1].u - x[j][i + 1].u) * hydhx;
316c4762a1bSJed Brown       uyy       = (2.0 * u - x[j - 1][i].u - x[j + 1][i].u) * hxdhy;
317c4762a1bSJed Brown       f[j][i].u = uxx + uyy - .5 * (x[j + 1][i].omega - x[j - 1][i].omega) * hx;
318c4762a1bSJed Brown 
319c4762a1bSJed Brown       /* V velocity */
320c4762a1bSJed Brown       u         = x[j][i].v;
321c4762a1bSJed Brown       uxx       = (2.0 * u - x[j][i - 1].v - x[j][i + 1].v) * hydhx;
322c4762a1bSJed Brown       uyy       = (2.0 * u - x[j - 1][i].v - x[j + 1][i].v) * hxdhy;
323c4762a1bSJed Brown       f[j][i].v = uxx + uyy + .5 * (x[j][i + 1].omega - x[j][i - 1].omega) * hy;
324c4762a1bSJed Brown 
325c4762a1bSJed Brown       /* Omega */
326c4762a1bSJed Brown       u             = x[j][i].omega;
327c4762a1bSJed Brown       uxx           = (2.0 * u - x[j][i - 1].omega - x[j][i + 1].omega) * hydhx;
328c4762a1bSJed Brown       uyy           = (2.0 * u - x[j - 1][i].omega - x[j + 1][i].omega) * hxdhy;
3299371c9d4SSatish Balay       f[j][i].omega = uxx + uyy + (vxp * (u - x[j][i - 1].omega) + vxm * (x[j][i + 1].omega - u)) * hy + (vyp * (u - x[j - 1][i].omega) + vym * (x[j + 1][i].omega - u)) * hx - .5 * grashof * (x[j][i + 1].temp - x[j][i - 1].temp) * hy;
330c4762a1bSJed Brown 
331c4762a1bSJed Brown       /* Temperature */
332c4762a1bSJed Brown       u            = x[j][i].temp;
333c4762a1bSJed Brown       uxx          = (2.0 * u - x[j][i - 1].temp - x[j][i + 1].temp) * hydhx;
334c4762a1bSJed Brown       uyy          = (2.0 * u - x[j - 1][i].temp - x[j + 1][i].temp) * hxdhy;
3359371c9d4SSatish Balay       f[j][i].temp = uxx + uyy + prandtl * ((vxp * (u - x[j][i - 1].temp) + vxm * (x[j][i + 1].temp - u)) * hy + (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   */
3429566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(84.0 * info->ym * info->xm));
3433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
344c4762a1bSJed Brown }
345c4762a1bSJed Brown 
MatMult_MyShell(Mat A,Vec x,Vec y)346d71ae5a4SJacob Faibussowitsch PetscErrorCode MatMult_MyShell(Mat A, Vec x, Vec y)
347d71ae5a4SJacob Faibussowitsch {
348c4762a1bSJed Brown   MatShellCtx    *matshellctx;
349c4762a1bSJed Brown   static PetscInt fail = 0;
350f480ea8aSBarry Smith   PetscMPIInt     rank;
351c4762a1bSJed Brown 
352c4762a1bSJed Brown   PetscFunctionBegin;
3539566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &matshellctx));
3549566063dSJacob Faibussowitsch   PetscCall(MatMult(matshellctx->Jmf, x, y));
3559566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)A), &rank));
356f480ea8aSBarry Smith   if (fail++ > 5) {
357f480ea8aSBarry Smith     PetscCall(VecFlag(y, rank == 0));
3583872ee93SStefano Zampini     PetscCall(VecAssemblyBegin(y));
3593872ee93SStefano Zampini     PetscCall(VecAssemblyEnd(y));
360c4762a1bSJed Brown   }
3613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
362c4762a1bSJed Brown }
363c4762a1bSJed Brown 
MatAssemblyEnd_MyShell(Mat A,MatAssemblyType tp)364d71ae5a4SJacob Faibussowitsch PetscErrorCode MatAssemblyEnd_MyShell(Mat A, MatAssemblyType tp)
365d71ae5a4SJacob Faibussowitsch {
366c4762a1bSJed Brown   MatShellCtx *matshellctx;
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   PetscFunctionBegin;
3699566063dSJacob Faibussowitsch   PetscCall(MatShellGetContext(A, &matshellctx));
3709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(matshellctx->Jmf, tp));
3713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
372c4762a1bSJed Brown }
373c4762a1bSJed Brown 
PCApply_MyShell(PC pc,Vec x,Vec y)374d71ae5a4SJacob Faibussowitsch PetscErrorCode PCApply_MyShell(PC pc, Vec x, Vec y)
375d71ae5a4SJacob Faibussowitsch {
376c4762a1bSJed Brown   static PetscInt fail = 0;
377f480ea8aSBarry Smith   PetscMPIInt     rank;
378c4762a1bSJed Brown 
379c4762a1bSJed Brown   PetscFunctionBegin;
3809566063dSJacob Faibussowitsch   PetscCall(VecCopy(x, y));
3819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
382f480ea8aSBarry Smith   if (fail++ > 3) {
383f480ea8aSBarry Smith     PetscCall(VecFlag(y, rank == 0));
3843872ee93SStefano Zampini     PetscCall(VecAssemblyBegin(y));
3853872ee93SStefano Zampini     PetscCall(VecAssemblyEnd(y));
386c4762a1bSJed Brown   }
3873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
388c4762a1bSJed Brown }
389c4762a1bSJed Brown 
390c4762a1bSJed Brown PETSC_EXTERN PetscErrorCode SNESComputeJacobian_DMDA(SNES, Vec, Mat, Mat, void *);
391c4762a1bSJed Brown 
SNESComputeJacobian_MyShell(SNES snes,Vec X,Mat A,Mat B,PetscCtx ctx)392*2a8381b2SBarry Smith PetscErrorCode SNESComputeJacobian_MyShell(SNES snes, Vec X, Mat A, Mat B, PetscCtx ctx)
393d71ae5a4SJacob Faibussowitsch {
394c4762a1bSJed Brown   static PetscInt fail = 0;
395c4762a1bSJed Brown 
396c4762a1bSJed Brown   PetscFunctionBegin;
3979566063dSJacob Faibussowitsch   PetscCall(SNESComputeJacobian_DMDA(snes, X, A, B, ctx));
39848a46eb9SPierre Jolivet   if (fail++ > 0) PetscCall(MatZeroEntries(A));
3993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
400c4762a1bSJed Brown }
401c4762a1bSJed Brown 
402c4762a1bSJed Brown /*TEST
403c4762a1bSJed Brown 
404c4762a1bSJed Brown    test:
405c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason
406c4762a1bSJed Brown 
407c4762a1bSJed Brown    test:
408c4762a1bSJed Brown       suffix: 2
409308041e4SStefano Zampini       args: -snes_converged_reason -ksp_converged_reason -error_in_matmult -fp_trap 0
410c4762a1bSJed Brown 
411c4762a1bSJed Brown    test:
412c4762a1bSJed Brown       suffix: 3
413308041e4SStefano Zampini       args: -snes_converged_reason -ksp_converged_reason -error_in_pcapply -fp_trap 0
414c4762a1bSJed Brown 
415c4762a1bSJed Brown    test:
416c4762a1bSJed Brown       suffix: 4
417308041e4SStefano Zampini       args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -fp_trap 0
418c4762a1bSJed Brown 
419c4762a1bSJed Brown    test:
420c4762a1bSJed Brown       suffix: 5
421308041e4SStefano Zampini       args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -pc_type bjacobi -fp_trap 0
422c4762a1bSJed Brown 
423c4762a1bSJed Brown    test:
424c4762a1bSJed Brown       suffix: 5_fieldsplit
425308041e4SStefano Zampini       args: -snes_converged_reason -ksp_converged_reason -error_in_pcsetup -pc_type fieldsplit -fp_trap 0
426c4762a1bSJed Brown       output_file: output/ex69_5.out
427c4762a1bSJed Brown 
428c4762a1bSJed Brown    test:
429c4762a1bSJed Brown       suffix: 6
430308041e4SStefano Zampini       args: -snes_converged_reason -ksp_converged_reason -error_in_domainmf -snes_mf -pc_type none -fp_trap 0
431c4762a1bSJed Brown 
432c4762a1bSJed Brown    test:
433c4762a1bSJed Brown       suffix: 7
434308041e4SStefano Zampini       args: -snes_converged_reason -ksp_converged_reason -error_in_domain -fp_trap 0
435c4762a1bSJed Brown 
436c4762a1bSJed Brown    test:
437c4762a1bSJed Brown       suffix: 8
438c4762a1bSJed Brown       args: -snes_converged_reason -ksp_converged_reason -error_in_domain -snes_mf -pc_type none
439c4762a1bSJed Brown 
440c4762a1bSJed Brown TEST*/
441