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