xref: /petsc/src/snes/tutorials/ex5.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Bratu nonlinear PDE in 2d.\n\
3*c4762a1bSJed Brown We solve the  Bratu (SFI - solid fuel ignition) problem in a 2D rectangular\n\
4*c4762a1bSJed Brown domain, using distributed arrays (DMDAs) to partition the parallel grid.\n\
5*c4762a1bSJed Brown The command line options include:\n\
6*c4762a1bSJed Brown   -par <parameter>, where <parameter> indicates the problem's nonlinearity\n\
7*c4762a1bSJed Brown      problem SFI:  <parameter> = Bratu parameter (0 <= par <= 6.81)\n\n\
8*c4762a1bSJed Brown   -m_par/n_par <parameter>, where <parameter> indicates an integer\n \
9*c4762a1bSJed Brown       that MMS3 will be evaluated with 2^m_par, 2^n_par";
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown /*T
12*c4762a1bSJed Brown    Concepts: SNES^parallel Bratu example
13*c4762a1bSJed Brown    Concepts: DMDA^using distributed arrays;
14*c4762a1bSJed Brown    Concepts: IS coloirng types;
15*c4762a1bSJed Brown    Processors: n
16*c4762a1bSJed Brown T*/
17*c4762a1bSJed Brown 
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown /* ------------------------------------------------------------------------
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown     Solid Fuel Ignition (SFI) problem.  This problem is modeled by
23*c4762a1bSJed Brown     the partial differential equation
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown             -Laplacian u - lambda*exp(u) = 0,  0 < x,y < 1,
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown     with boundary conditions
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown              u = 0  for  x = 0, x = 1, y = 0, y = 1.
30*c4762a1bSJed Brown 
31*c4762a1bSJed Brown     A finite difference approximation with the usual 5-point stencil
32*c4762a1bSJed Brown     is used to discretize the boundary value problem to obtain a nonlinear
33*c4762a1bSJed Brown     system of equations.
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown       This example shows how geometric multigrid can be run transparently with a nonlinear solver so long
37*c4762a1bSJed Brown       as SNESSetDM() is provided. Example usage
38*c4762a1bSJed Brown 
39*c4762a1bSJed Brown       ./ex5 -pc_type mg -ksp_monitor  -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17
40*c4762a1bSJed Brown              -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown       or to run with grid sequencing on the nonlinear problem (note that you do not need to provide the number of
43*c4762a1bSJed Brown          multigrid levels, it will be determined automatically based on the number of refinements done)
44*c4762a1bSJed Brown 
45*c4762a1bSJed Brown       ./ex5 -pc_type mg -ksp_monitor  -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3
46*c4762a1bSJed Brown              -mg_levels_ksp_monitor -snes_monitor -mg_levels_pc_type sor -pc_mg_type full
47*c4762a1bSJed Brown 
48*c4762a1bSJed Brown   ------------------------------------------------------------------------- */
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown /*
51*c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
52*c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
53*c4762a1bSJed Brown */
54*c4762a1bSJed Brown #include <petscdm.h>
55*c4762a1bSJed Brown #include <petscdmda.h>
56*c4762a1bSJed Brown #include <petscsnes.h>
57*c4762a1bSJed Brown #include <petscmatlab.h>
58*c4762a1bSJed Brown #include <petsc/private/snesimpl.h> /* For SNES_Solve event */
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown /*
61*c4762a1bSJed Brown    User-defined application context - contains data needed by the
62*c4762a1bSJed Brown    application-provided call-back routines, FormJacobianLocal() and
63*c4762a1bSJed Brown    FormFunctionLocal().
64*c4762a1bSJed Brown */
65*c4762a1bSJed Brown typedef struct AppCtx AppCtx;
66*c4762a1bSJed Brown struct AppCtx {
67*c4762a1bSJed Brown   PetscReal param;          /* test problem parameter */
68*c4762a1bSJed Brown   PetscInt  m,n;            /* MMS3 parameters */
69*c4762a1bSJed Brown   PetscErrorCode (*mms_solution)(AppCtx*,const DMDACoor2d*,PetscScalar*);
70*c4762a1bSJed Brown   PetscErrorCode (*mms_forcing)(AppCtx*,const DMDACoor2d*,PetscScalar*);
71*c4762a1bSJed Brown };
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown /*
74*c4762a1bSJed Brown    User-defined routines
75*c4762a1bSJed Brown */
76*c4762a1bSJed Brown extern PetscErrorCode FormInitialGuess(DM,AppCtx*,Vec);
77*c4762a1bSJed Brown extern PetscErrorCode FormFunctionLocal(DMDALocalInfo*,PetscScalar**,PetscScalar**,AppCtx*);
78*c4762a1bSJed Brown extern PetscErrorCode FormExactSolution(DM,AppCtx*,Vec);
79*c4762a1bSJed Brown extern PetscErrorCode ZeroBCSolution(AppCtx*,const DMDACoor2d*,PetscScalar*);
80*c4762a1bSJed Brown extern PetscErrorCode MMSSolution1(AppCtx*,const DMDACoor2d*,PetscScalar*);
81*c4762a1bSJed Brown extern PetscErrorCode MMSForcing1(AppCtx*,const DMDACoor2d*,PetscScalar*);
82*c4762a1bSJed Brown extern PetscErrorCode MMSSolution2(AppCtx*,const DMDACoor2d*,PetscScalar*);
83*c4762a1bSJed Brown extern PetscErrorCode MMSForcing2(AppCtx*,const DMDACoor2d*,PetscScalar*);
84*c4762a1bSJed Brown extern PetscErrorCode MMSSolution3(AppCtx*,const DMDACoor2d*,PetscScalar*);
85*c4762a1bSJed Brown extern PetscErrorCode MMSForcing3(AppCtx*,const DMDACoor2d*,PetscScalar*);
86*c4762a1bSJed Brown extern PetscErrorCode MMSSolution4(AppCtx*,const DMDACoor2d*,PetscScalar*);
87*c4762a1bSJed Brown extern PetscErrorCode MMSForcing4(AppCtx*,const DMDACoor2d*,PetscScalar*);
88*c4762a1bSJed Brown extern PetscErrorCode FormJacobianLocal(DMDALocalInfo*,PetscScalar**,Mat,Mat,AppCtx*);
89*c4762a1bSJed Brown extern PetscErrorCode FormObjectiveLocal(DMDALocalInfo*,PetscScalar**,PetscReal*,AppCtx*);
90*c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
91*c4762a1bSJed Brown extern PetscErrorCode FormFunctionMatlab(SNES,Vec,Vec,void*);
92*c4762a1bSJed Brown #endif
93*c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES,Vec,Vec,void*);
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown int main(int argc,char **argv)
96*c4762a1bSJed Brown {
97*c4762a1bSJed Brown   SNES           snes;                         /* nonlinear solver */
98*c4762a1bSJed Brown   Vec            x;                            /* solution vector */
99*c4762a1bSJed Brown   AppCtx         user;                         /* user-defined work context */
100*c4762a1bSJed Brown   PetscInt       its;                          /* iterations for convergence */
101*c4762a1bSJed Brown   PetscErrorCode ierr;
102*c4762a1bSJed Brown   PetscReal      bratu_lambda_max = 6.81;
103*c4762a1bSJed Brown   PetscReal      bratu_lambda_min = 0.;
104*c4762a1bSJed Brown   PetscInt       MMS              = 0;
105*c4762a1bSJed Brown   PetscBool      flg              = PETSC_FALSE;
106*c4762a1bSJed Brown   DM             da;
107*c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
108*c4762a1bSJed Brown   Vec            r               = NULL;
109*c4762a1bSJed Brown   PetscBool      matlab_function = PETSC_FALSE;
110*c4762a1bSJed Brown #endif
111*c4762a1bSJed Brown   KSP            ksp;
112*c4762a1bSJed Brown   PetscInt       lits,slits;
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115*c4762a1bSJed Brown      Initialize program
116*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121*c4762a1bSJed Brown      Initialize problem parameters
122*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123*c4762a1bSJed Brown   user.param = 6.0;
124*c4762a1bSJed Brown   ierr       = PetscOptionsGetReal(NULL,NULL,"-par",&user.param,NULL);CHKERRQ(ierr);
125*c4762a1bSJed Brown   if (user.param > bratu_lambda_max || user.param < bratu_lambda_min) SETERRQ3(PETSC_COMM_SELF,1,"Lambda, %g, is out of range, [%g, %g]", user.param, bratu_lambda_min, bratu_lambda_max);
126*c4762a1bSJed Brown   ierr       = PetscOptionsGetInt(NULL,NULL,"-mms",&MMS,NULL);CHKERRQ(ierr);
127*c4762a1bSJed Brown   if (MMS == 3) {
128*c4762a1bSJed Brown     PetscInt mPar = 2, nPar = 1;
129*c4762a1bSJed Brown     ierr = PetscOptionsGetInt(NULL,NULL,"-m_par",&mPar,NULL);CHKERRQ(ierr);
130*c4762a1bSJed Brown     ierr = PetscOptionsGetInt(NULL,NULL,"-n_par",&nPar,NULL);CHKERRQ(ierr);
131*c4762a1bSJed Brown     user.m = PetscPowInt(2,mPar);
132*c4762a1bSJed Brown     user.n = PetscPowInt(2,nPar);
133*c4762a1bSJed Brown   }
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136*c4762a1bSJed Brown      Create nonlinear solver context
137*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = SNESSetCountersReset(snes,PETSC_FALSE);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = SNESSetNGS(snes, NonlinearGS, NULL);CHKERRQ(ierr);
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
144*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);
149*c4762a1bSJed Brown   ierr = DMSetApplicationContext(da,&user);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = SNESSetDM(snes,da);CHKERRQ(ierr);
151*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152*c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
153*c4762a1bSJed Brown      vectors that are the same types
154*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158*c4762a1bSJed Brown      Set local function evaluation routine
159*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160*c4762a1bSJed Brown   user.mms_solution = ZeroBCSolution;
161*c4762a1bSJed Brown   switch (MMS) {
162*c4762a1bSJed Brown   case 0: user.mms_solution = NULL; user.mms_forcing = NULL; CHKERRQ(ierr);
163*c4762a1bSJed Brown   case 1: user.mms_solution = MMSSolution1; user.mms_forcing = MMSForcing1; break;
164*c4762a1bSJed Brown   case 2: user.mms_solution = MMSSolution2; user.mms_forcing = MMSForcing2; break;
165*c4762a1bSJed Brown   case 3: user.mms_solution = MMSSolution3; user.mms_forcing = MMSForcing3; break;
166*c4762a1bSJed Brown   case 4: user.mms_solution = MMSSolution4; user.mms_forcing = MMSForcing4; break;
167*c4762a1bSJed Brown   default: SETERRQ1(PETSC_COMM_WORLD,PETSC_ERR_USER,"Unknown MMS type %d",MMS);
168*c4762a1bSJed Brown   }
169*c4762a1bSJed Brown   ierr = DMDASNESSetFunctionLocal(da,INSERT_VALUES,(DMDASNESFunction)FormFunctionLocal,&user);CHKERRQ(ierr);
170*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-fd",&flg,NULL);CHKERRQ(ierr);
171*c4762a1bSJed Brown   if (!flg) {
172*c4762a1bSJed Brown     ierr = DMDASNESSetJacobianLocal(da,(DMDASNESJacobian)FormJacobianLocal,&user);CHKERRQ(ierr);
173*c4762a1bSJed Brown   }
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-obj",&flg,NULL);CHKERRQ(ierr);
176*c4762a1bSJed Brown   if (flg) {
177*c4762a1bSJed Brown     ierr = DMDASNESSetObjectiveLocal(da,(DMDASNESObjective)FormObjectiveLocal,&user);CHKERRQ(ierr);
178*c4762a1bSJed Brown   }
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
181*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-matlab_function",&matlab_function,0);CHKERRQ(ierr);
182*c4762a1bSJed Brown   if (matlab_function) {
183*c4762a1bSJed Brown     ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
184*c4762a1bSJed Brown     ierr = SNESSetFunction(snes,r,FormFunctionMatlab,&user);CHKERRQ(ierr);
185*c4762a1bSJed Brown   }
186*c4762a1bSJed Brown #endif
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189*c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
190*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
192*c4762a1bSJed Brown 
193*c4762a1bSJed Brown   ierr = FormInitialGuess(da,&user,x);CHKERRQ(ierr);
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196*c4762a1bSJed Brown      Solve nonlinear system
197*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198*c4762a1bSJed Brown   ierr = SNESSolve(snes,NULL,x);CHKERRQ(ierr);
199*c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
200*c4762a1bSJed Brown 
201*c4762a1bSJed Brown   ierr = SNESGetLinearSolveIterations(snes,&slits);CHKERRQ(ierr);
202*c4762a1bSJed Brown   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
203*c4762a1bSJed Brown   ierr = KSPGetTotalIterations(ksp,&lits);CHKERRQ(ierr);
204*c4762a1bSJed Brown   if (lits != slits) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_PLIB,"Number of total linear iterations reported by SNES %D does not match reported by KSP %D",slits,lits);
205*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207*c4762a1bSJed Brown 
208*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209*c4762a1bSJed Brown      If using MMS, check the l_2 error
210*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211*c4762a1bSJed Brown   if (MMS) {
212*c4762a1bSJed Brown     Vec       e;
213*c4762a1bSJed Brown     PetscReal errorl2, errorinf;
214*c4762a1bSJed Brown     PetscInt  N;
215*c4762a1bSJed Brown 
216*c4762a1bSJed Brown     ierr = VecDuplicate(x, &e);CHKERRQ(ierr);
217*c4762a1bSJed Brown     ierr = PetscObjectViewFromOptions((PetscObject) x, NULL, "-sol_view");CHKERRQ(ierr);
218*c4762a1bSJed Brown     ierr = FormExactSolution(da, &user, e);CHKERRQ(ierr);
219*c4762a1bSJed Brown     ierr = PetscObjectViewFromOptions((PetscObject) e, NULL, "-exact_view");CHKERRQ(ierr);
220*c4762a1bSJed Brown     ierr = VecAXPY(e, -1.0, x);CHKERRQ(ierr);
221*c4762a1bSJed Brown     ierr = PetscObjectViewFromOptions((PetscObject) e, NULL, "-error_view");CHKERRQ(ierr);
222*c4762a1bSJed Brown     ierr = VecNorm(e, NORM_2, &errorl2);CHKERRQ(ierr);
223*c4762a1bSJed Brown     ierr = VecNorm(e, NORM_INFINITY, &errorinf);CHKERRQ(ierr);
224*c4762a1bSJed Brown     ierr = VecGetSize(e, &N);CHKERRQ(ierr);
225*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "N: %D error L2 %g inf %g\n", N, (double) errorl2/PetscSqrtReal(N), (double) errorinf);CHKERRQ(ierr);
226*c4762a1bSJed Brown     ierr = VecDestroy(&e);CHKERRQ(ierr);
227*c4762a1bSJed Brown     ierr = PetscLogEventSetDof(SNES_Solve, 0, N);CHKERRQ(ierr);
228*c4762a1bSJed Brown     ierr = PetscLogEventSetError(SNES_Solve, 0, errorl2/PetscSqrtReal(N));CHKERRQ(ierr);
229*c4762a1bSJed Brown   }
230*c4762a1bSJed Brown 
231*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
233*c4762a1bSJed Brown      are no longer needed.
234*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
235*c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
236*c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
237*c4762a1bSJed Brown #endif
238*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
239*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = PetscFinalize();
242*c4762a1bSJed Brown   return ierr;
243*c4762a1bSJed Brown }
244*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
245*c4762a1bSJed Brown /*
246*c4762a1bSJed Brown    FormInitialGuess - Forms initial approximation.
247*c4762a1bSJed Brown 
248*c4762a1bSJed Brown    Input Parameters:
249*c4762a1bSJed Brown    da - The DM
250*c4762a1bSJed Brown    user - user-defined application context
251*c4762a1bSJed Brown 
252*c4762a1bSJed Brown    Output Parameter:
253*c4762a1bSJed Brown    X - vector
254*c4762a1bSJed Brown  */
255*c4762a1bSJed Brown PetscErrorCode FormInitialGuess(DM da,AppCtx *user,Vec X)
256*c4762a1bSJed Brown {
257*c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
258*c4762a1bSJed Brown   PetscErrorCode ierr;
259*c4762a1bSJed Brown   PetscReal      lambda,temp1,temp,hx,hy;
260*c4762a1bSJed Brown   PetscScalar    **x;
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   PetscFunctionBeginUser;
263*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown   lambda = user->param;
266*c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(Mx-1);
267*c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(My-1);
268*c4762a1bSJed Brown   temp1  = lambda/(lambda + 1.0);
269*c4762a1bSJed Brown 
270*c4762a1bSJed Brown   /*
271*c4762a1bSJed Brown      Get a pointer to vector data.
272*c4762a1bSJed Brown        - For default PETSc vectors, VecGetArray() returns a pointer to
273*c4762a1bSJed Brown          the data array.  Otherwise, the routine is implementation dependent.
274*c4762a1bSJed Brown        - You MUST call VecRestoreArray() when you no longer need access to
275*c4762a1bSJed Brown          the array.
276*c4762a1bSJed Brown   */
277*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da,X,&x);CHKERRQ(ierr);
278*c4762a1bSJed Brown 
279*c4762a1bSJed Brown   /*
280*c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
281*c4762a1bSJed Brown        xs, ys   - starting grid indices (no ghost points)
282*c4762a1bSJed Brown        xm, ym   - widths of local grid (no ghost points)
283*c4762a1bSJed Brown 
284*c4762a1bSJed Brown   */
285*c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
286*c4762a1bSJed Brown 
287*c4762a1bSJed Brown   /*
288*c4762a1bSJed Brown      Compute initial guess over the locally owned part of the grid
289*c4762a1bSJed Brown   */
290*c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
291*c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j,My-j-1))*hy;
292*c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
293*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
294*c4762a1bSJed Brown         /* boundary conditions are all zero Dirichlet */
295*c4762a1bSJed Brown         x[j][i] = 0.0;
296*c4762a1bSJed Brown       } else {
297*c4762a1bSJed Brown         x[j][i] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,Mx-i-1))*hx,temp));
298*c4762a1bSJed Brown       }
299*c4762a1bSJed Brown     }
300*c4762a1bSJed Brown   }
301*c4762a1bSJed Brown 
302*c4762a1bSJed Brown   /*
303*c4762a1bSJed Brown      Restore vector
304*c4762a1bSJed Brown   */
305*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da,X,&x);CHKERRQ(ierr);
306*c4762a1bSJed Brown   PetscFunctionReturn(0);
307*c4762a1bSJed Brown }
308*c4762a1bSJed Brown 
309*c4762a1bSJed Brown /*
310*c4762a1bSJed Brown   FormExactSolution - Forms MMS solution
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown   Input Parameters:
313*c4762a1bSJed Brown   da - The DM
314*c4762a1bSJed Brown   user - user-defined application context
315*c4762a1bSJed Brown 
316*c4762a1bSJed Brown   Output Parameter:
317*c4762a1bSJed Brown   X - vector
318*c4762a1bSJed Brown  */
319*c4762a1bSJed Brown PetscErrorCode FormExactSolution(DM da, AppCtx *user, Vec U)
320*c4762a1bSJed Brown {
321*c4762a1bSJed Brown   DM             coordDA;
322*c4762a1bSJed Brown   Vec            coordinates;
323*c4762a1bSJed Brown   DMDACoor2d   **coords;
324*c4762a1bSJed Brown   PetscScalar  **u;
325*c4762a1bSJed Brown   PetscInt       xs, ys, xm, ym, i, j;
326*c4762a1bSJed Brown   PetscErrorCode ierr;
327*c4762a1bSJed Brown 
328*c4762a1bSJed Brown   PetscFunctionBeginUser;
329*c4762a1bSJed Brown   ierr = DMDAGetCorners(da, &xs, &ys, NULL, &xm, &ym, NULL);CHKERRQ(ierr);
330*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(da, &coordDA);CHKERRQ(ierr);
331*c4762a1bSJed Brown   ierr = DMGetCoordinates(da, &coordinates);CHKERRQ(ierr);
332*c4762a1bSJed Brown   ierr = DMDAVecGetArray(coordDA, coordinates, &coords);CHKERRQ(ierr);
333*c4762a1bSJed Brown   ierr = DMDAVecGetArray(da, U, &u);CHKERRQ(ierr);
334*c4762a1bSJed Brown   for (j = ys; j < ys+ym; ++j) {
335*c4762a1bSJed Brown     for (i = xs; i < xs+xm; ++i) {
336*c4762a1bSJed Brown       user->mms_solution(user,&coords[j][i],&u[j][i]);
337*c4762a1bSJed Brown     }
338*c4762a1bSJed Brown   }
339*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(da, U, &u);CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(coordDA, coordinates, &coords);CHKERRQ(ierr);
341*c4762a1bSJed Brown   PetscFunctionReturn(0);
342*c4762a1bSJed Brown }
343*c4762a1bSJed Brown 
344*c4762a1bSJed Brown PetscErrorCode ZeroBCSolution(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
345*c4762a1bSJed Brown {
346*c4762a1bSJed Brown   u[0] = 0.;
347*c4762a1bSJed Brown   return 0;
348*c4762a1bSJed Brown }
349*c4762a1bSJed Brown 
350*c4762a1bSJed Brown /* The functions below evaluate the MMS solution u(x,y) and associated forcing
351*c4762a1bSJed Brown 
352*c4762a1bSJed Brown      f(x,y) = -u_xx - u_yy - lambda exp(u)
353*c4762a1bSJed Brown 
354*c4762a1bSJed Brown   such that u(x,y) is an exact solution with f(x,y) as the right hand side forcing term.
355*c4762a1bSJed Brown  */
356*c4762a1bSJed Brown PetscErrorCode MMSSolution1(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
357*c4762a1bSJed Brown {
358*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
359*c4762a1bSJed Brown   u[0] = x*(1 - x)*y*(1 - y);
360*c4762a1bSJed Brown   PetscLogFlops(5);
361*c4762a1bSJed Brown   return 0;
362*c4762a1bSJed Brown }
363*c4762a1bSJed Brown PetscErrorCode MMSForcing1(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
364*c4762a1bSJed Brown {
365*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
366*c4762a1bSJed Brown   f[0] = 2*x*(1 - x) + 2*y*(1 - y) - user->param*PetscExpReal(x*(1 - x)*y*(1 - y));
367*c4762a1bSJed Brown   return 0;
368*c4762a1bSJed Brown }
369*c4762a1bSJed Brown 
370*c4762a1bSJed Brown PetscErrorCode MMSSolution2(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
371*c4762a1bSJed Brown {
372*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
373*c4762a1bSJed Brown   u[0] = PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y);
374*c4762a1bSJed Brown   PetscLogFlops(5);
375*c4762a1bSJed Brown   return 0;
376*c4762a1bSJed Brown }
377*c4762a1bSJed Brown PetscErrorCode MMSForcing2(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
378*c4762a1bSJed Brown {
379*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
380*c4762a1bSJed Brown   f[0] = 2*PetscSqr(PETSC_PI)*PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y) - user->param*PetscExpReal(PetscSinReal(PETSC_PI*x)*PetscSinReal(PETSC_PI*y));
381*c4762a1bSJed Brown   return 0;
382*c4762a1bSJed Brown }
383*c4762a1bSJed Brown 
384*c4762a1bSJed Brown PetscErrorCode MMSSolution3(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
385*c4762a1bSJed Brown {
386*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
387*c4762a1bSJed Brown   u[0] = PetscSinReal(user->m*PETSC_PI*x*(1-y))*PetscSinReal(user->n*PETSC_PI*y*(1-x));
388*c4762a1bSJed Brown   PetscLogFlops(5);
389*c4762a1bSJed Brown   return 0;
390*c4762a1bSJed Brown }
391*c4762a1bSJed Brown PetscErrorCode MMSForcing3(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
392*c4762a1bSJed Brown {
393*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
394*c4762a1bSJed Brown   PetscReal m = user->m, n = user->n, lambda = user->param;
395*c4762a1bSJed Brown   f[0] = (-(PetscExpReal(PetscSinReal(m*PETSC_PI*x*(1 - y))*PetscSinReal(n*PETSC_PI*(1 - x)*y))*lambda)
396*c4762a1bSJed Brown           + PetscSqr(PETSC_PI)*(-2*m*n*((-1 + x)*x + (-1 + y)*y)*PetscCosReal(m*PETSC_PI*x*(-1 + y))*PetscCosReal(n*PETSC_PI*(-1 + x)*y)
397*c4762a1bSJed Brown                                 + (PetscSqr(m)*(PetscSqr(x) + PetscSqr(-1 + y)) + PetscSqr(n)*(PetscSqr(-1 + x) + PetscSqr(y)))
398*c4762a1bSJed Brown                                 *PetscSinReal(m*PETSC_PI*x*(-1 + y))*PetscSinReal(n*PETSC_PI*(-1 + x)*y)));
399*c4762a1bSJed Brown   return 0;
400*c4762a1bSJed Brown }
401*c4762a1bSJed Brown 
402*c4762a1bSJed Brown PetscErrorCode MMSSolution4(AppCtx *user,const DMDACoor2d *c,PetscScalar *u)
403*c4762a1bSJed Brown {
404*c4762a1bSJed Brown   const PetscReal Lx = 1.,Ly = 1.;
405*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
406*c4762a1bSJed Brown   u[0] = (PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y));
407*c4762a1bSJed Brown   PetscLogFlops(9);
408*c4762a1bSJed Brown   return 0;
409*c4762a1bSJed Brown }
410*c4762a1bSJed Brown PetscErrorCode MMSForcing4(AppCtx *user,const DMDACoor2d *c,PetscScalar *f)
411*c4762a1bSJed Brown {
412*c4762a1bSJed Brown   const PetscReal Lx = 1.,Ly = 1.;
413*c4762a1bSJed Brown   PetscReal x = PetscRealPart(c->x), y = PetscRealPart(c->y);
414*c4762a1bSJed Brown   f[0] = (2*PetscSqr(x)*(PetscSqr(x)-PetscSqr(Lx))*(PetscSqr(Ly)-6*PetscSqr(y))
415*c4762a1bSJed Brown           + 2*PetscSqr(y)*(PetscSqr(Lx)-6*PetscSqr(x))*(PetscSqr(y)-PetscSqr(Ly))
416*c4762a1bSJed Brown           - user->param*PetscExpReal((PetscPowReal(x,4)-PetscSqr(Lx)*PetscSqr(x))*(PetscPowReal(y,4)-PetscSqr(Ly)*PetscSqr(y))));
417*c4762a1bSJed Brown   return 0;
418*c4762a1bSJed Brown }
419*c4762a1bSJed Brown 
420*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
421*c4762a1bSJed Brown /*
422*c4762a1bSJed Brown    FormFunctionLocal - Evaluates nonlinear function, F(x) on local process patch
423*c4762a1bSJed Brown 
424*c4762a1bSJed Brown 
425*c4762a1bSJed Brown  */
426*c4762a1bSJed Brown PetscErrorCode FormFunctionLocal(DMDALocalInfo *info,PetscScalar **x,PetscScalar **f,AppCtx *user)
427*c4762a1bSJed Brown {
428*c4762a1bSJed Brown   PetscErrorCode ierr;
429*c4762a1bSJed Brown   PetscInt       i,j;
430*c4762a1bSJed Brown   PetscReal      lambda,hx,hy,hxdhy,hydhx;
431*c4762a1bSJed Brown   PetscScalar    u,ue,uw,un,us,uxx,uyy,mms_solution,mms_forcing;
432*c4762a1bSJed Brown   DMDACoor2d     c;
433*c4762a1bSJed Brown 
434*c4762a1bSJed Brown   PetscFunctionBeginUser;
435*c4762a1bSJed Brown   lambda = user->param;
436*c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info->mx-1);
437*c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info->my-1);
438*c4762a1bSJed Brown   hxdhy  = hx/hy;
439*c4762a1bSJed Brown   hydhx  = hy/hx;
440*c4762a1bSJed Brown   /*
441*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
442*c4762a1bSJed Brown   */
443*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
444*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
445*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
446*c4762a1bSJed Brown         c.x = i*hx; c.y = j*hy;
447*c4762a1bSJed Brown         ierr = user->mms_solution(user,&c,&mms_solution);CHKERRQ(ierr);
448*c4762a1bSJed Brown         f[j][i] = 2.0*(hydhx+hxdhy)*(x[j][i] - mms_solution);
449*c4762a1bSJed Brown       } else {
450*c4762a1bSJed Brown         u  = x[j][i];
451*c4762a1bSJed Brown         uw = x[j][i-1];
452*c4762a1bSJed Brown         ue = x[j][i+1];
453*c4762a1bSJed Brown         un = x[j-1][i];
454*c4762a1bSJed Brown         us = x[j+1][i];
455*c4762a1bSJed Brown 
456*c4762a1bSJed Brown         /* Enforce boundary conditions at neighboring points -- setting these values causes the Jacobian to be symmetric. */
457*c4762a1bSJed Brown         if (i-1 == 0) {c.x = (i-1)*hx; c.y = j*hy; ierr = user->mms_solution(user,&c,&uw);CHKERRQ(ierr);}
458*c4762a1bSJed Brown         if (i+1 == info->mx-1) {c.x = (i+1)*hx; c.y = j*hy; ierr = user->mms_solution(user,&c,&ue);CHKERRQ(ierr);}
459*c4762a1bSJed Brown         if (j-1 == 0) {c.x = i*hx; c.y = (j-1)*hy; ierr = user->mms_solution(user,&c,&un);CHKERRQ(ierr);}
460*c4762a1bSJed Brown         if (j+1 == info->my-1) {c.x = i*hx; c.y = (j+1)*hy; ierr = user->mms_solution(user,&c,&us);CHKERRQ(ierr);}
461*c4762a1bSJed Brown 
462*c4762a1bSJed Brown         uxx     = (2.0*u - uw - ue)*hydhx;
463*c4762a1bSJed Brown         uyy     = (2.0*u - un - us)*hxdhy;
464*c4762a1bSJed Brown         mms_forcing = 0;
465*c4762a1bSJed Brown         c.x = i*hx; c.y = j*hy;
466*c4762a1bSJed Brown         if (user->mms_forcing) {ierr = user->mms_forcing(user,&c,&mms_forcing);CHKERRQ(ierr);}
467*c4762a1bSJed Brown         f[j][i] = uxx + uyy - hx*hy*(lambda*PetscExpScalar(u) + mms_forcing);
468*c4762a1bSJed Brown       }
469*c4762a1bSJed Brown     }
470*c4762a1bSJed Brown   }
471*c4762a1bSJed Brown   ierr = PetscLogFlops(11.0*info->ym*info->xm);CHKERRQ(ierr);
472*c4762a1bSJed Brown   PetscFunctionReturn(0);
473*c4762a1bSJed Brown }
474*c4762a1bSJed Brown 
475*c4762a1bSJed Brown /* FormObjectiveLocal - Evaluates nonlinear function, F(x) on local process patch */
476*c4762a1bSJed Brown PetscErrorCode FormObjectiveLocal(DMDALocalInfo *info,PetscScalar **x,PetscReal *obj,AppCtx *user)
477*c4762a1bSJed Brown {
478*c4762a1bSJed Brown   PetscErrorCode ierr;
479*c4762a1bSJed Brown   PetscInt       i,j;
480*c4762a1bSJed Brown   PetscReal      lambda,hx,hy,hxdhy,hydhx,sc,lobj=0;
481*c4762a1bSJed Brown   PetscScalar    u,ue,uw,un,us,uxux,uyuy;
482*c4762a1bSJed Brown   MPI_Comm       comm;
483*c4762a1bSJed Brown 
484*c4762a1bSJed Brown   PetscFunctionBeginUser;
485*c4762a1bSJed Brown   *obj   = 0;
486*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)info->da,&comm);CHKERRQ(ierr);
487*c4762a1bSJed Brown   lambda = user->param;
488*c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(info->mx-1);
489*c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(info->my-1);
490*c4762a1bSJed Brown   sc     = hx*hy*lambda;
491*c4762a1bSJed Brown   hxdhy  = hx/hy;
492*c4762a1bSJed Brown   hydhx  = hy/hx;
493*c4762a1bSJed Brown   /*
494*c4762a1bSJed Brown      Compute function over the locally owned part of the grid
495*c4762a1bSJed Brown   */
496*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
497*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
498*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
499*c4762a1bSJed Brown         lobj += PetscRealPart((hydhx + hxdhy)*x[j][i]*x[j][i]);
500*c4762a1bSJed Brown       } else {
501*c4762a1bSJed Brown         u  = x[j][i];
502*c4762a1bSJed Brown         uw = x[j][i-1];
503*c4762a1bSJed Brown         ue = x[j][i+1];
504*c4762a1bSJed Brown         un = x[j-1][i];
505*c4762a1bSJed Brown         us = x[j+1][i];
506*c4762a1bSJed Brown 
507*c4762a1bSJed Brown         if (i-1 == 0) uw = 0.;
508*c4762a1bSJed Brown         if (i+1 == info->mx-1) ue = 0.;
509*c4762a1bSJed Brown         if (j-1 == 0) un = 0.;
510*c4762a1bSJed Brown         if (j+1 == info->my-1) us = 0.;
511*c4762a1bSJed Brown 
512*c4762a1bSJed Brown         /* F[u] = 1/2\int_{\omega}\nabla^2u(x)*u(x)*dx */
513*c4762a1bSJed Brown 
514*c4762a1bSJed Brown         uxux = u*(2.*u - ue - uw)*hydhx;
515*c4762a1bSJed Brown         uyuy = u*(2.*u - un - us)*hxdhy;
516*c4762a1bSJed Brown 
517*c4762a1bSJed Brown         lobj += PetscRealPart(0.5*(uxux + uyuy) - sc*PetscExpScalar(u));
518*c4762a1bSJed Brown       }
519*c4762a1bSJed Brown     }
520*c4762a1bSJed Brown   }
521*c4762a1bSJed Brown   ierr = PetscLogFlops(12.0*info->ym*info->xm);CHKERRQ(ierr);
522*c4762a1bSJed Brown   ierr = MPI_Allreduce(&lobj,obj,1,MPIU_REAL,MPIU_SUM,comm);CHKERRQ(ierr);
523*c4762a1bSJed Brown   PetscFunctionReturn(0);
524*c4762a1bSJed Brown }
525*c4762a1bSJed Brown 
526*c4762a1bSJed Brown /*
527*c4762a1bSJed Brown    FormJacobianLocal - Evaluates Jacobian matrix on local process patch
528*c4762a1bSJed Brown */
529*c4762a1bSJed Brown PetscErrorCode FormJacobianLocal(DMDALocalInfo *info,PetscScalar **x,Mat jac,Mat jacpre,AppCtx *user)
530*c4762a1bSJed Brown {
531*c4762a1bSJed Brown   PetscErrorCode ierr;
532*c4762a1bSJed Brown   PetscInt       i,j,k;
533*c4762a1bSJed Brown   MatStencil     col[5],row;
534*c4762a1bSJed Brown   PetscScalar    lambda,v[5],hx,hy,hxdhy,hydhx,sc;
535*c4762a1bSJed Brown   DM             coordDA;
536*c4762a1bSJed Brown   Vec            coordinates;
537*c4762a1bSJed Brown   DMDACoor2d   **coords;
538*c4762a1bSJed Brown 
539*c4762a1bSJed Brown   PetscFunctionBeginUser;
540*c4762a1bSJed Brown   lambda = user->param;
541*c4762a1bSJed Brown   /* Extract coordinates */
542*c4762a1bSJed Brown   ierr = DMGetCoordinateDM(info->da, &coordDA);CHKERRQ(ierr);
543*c4762a1bSJed Brown   ierr = DMGetCoordinates(info->da, &coordinates);CHKERRQ(ierr);
544*c4762a1bSJed Brown   ierr = DMDAVecGetArray(coordDA, coordinates, &coords);CHKERRQ(ierr);
545*c4762a1bSJed Brown   hx     = info->xm > 1 ? PetscRealPart(coords[info->ys][info->xs+1].x) - PetscRealPart(coords[info->ys][info->xs].x) : 1.0;
546*c4762a1bSJed Brown   hy     = info->ym > 1 ? PetscRealPart(coords[info->ys+1][info->xs].y) - PetscRealPart(coords[info->ys][info->xs].y) : 1.0;
547*c4762a1bSJed Brown   ierr = DMDAVecRestoreArray(coordDA, coordinates, &coords);CHKERRQ(ierr);
548*c4762a1bSJed Brown   hxdhy  = hx/hy;
549*c4762a1bSJed Brown   hydhx  = hy/hx;
550*c4762a1bSJed Brown   sc     = hx*hy*lambda;
551*c4762a1bSJed Brown 
552*c4762a1bSJed Brown 
553*c4762a1bSJed Brown   /*
554*c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
555*c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
556*c4762a1bSJed Brown         contiguous chunks of rows across the processors.
557*c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
558*c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
559*c4762a1bSJed Brown         appropriate processor during matrix assembly).
560*c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
561*c4762a1bSJed Brown       - We can set matrix entries either using either
562*c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues(), as discussed above.
563*c4762a1bSJed Brown   */
564*c4762a1bSJed Brown   for (j=info->ys; j<info->ys+info->ym; j++) {
565*c4762a1bSJed Brown     for (i=info->xs; i<info->xs+info->xm; i++) {
566*c4762a1bSJed Brown       row.j = j; row.i = i;
567*c4762a1bSJed Brown       /* boundary points */
568*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info->mx-1 || j == info->my-1) {
569*c4762a1bSJed Brown         v[0] =  2.0*(hydhx + hxdhy);
570*c4762a1bSJed Brown         ierr = MatSetValuesStencil(jacpre,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
571*c4762a1bSJed Brown       } else {
572*c4762a1bSJed Brown         k = 0;
573*c4762a1bSJed Brown         /* interior grid points */
574*c4762a1bSJed Brown         if (j-1 != 0) {
575*c4762a1bSJed Brown           v[k]     = -hxdhy;
576*c4762a1bSJed Brown           col[k].j = j - 1; col[k].i = i;
577*c4762a1bSJed Brown           k++;
578*c4762a1bSJed Brown         }
579*c4762a1bSJed Brown         if (i-1 != 0) {
580*c4762a1bSJed Brown           v[k]     = -hydhx;
581*c4762a1bSJed Brown           col[k].j = j;     col[k].i = i-1;
582*c4762a1bSJed Brown           k++;
583*c4762a1bSJed Brown         }
584*c4762a1bSJed Brown 
585*c4762a1bSJed Brown         v[k] = 2.0*(hydhx + hxdhy) - sc*PetscExpScalar(x[j][i]); col[k].j = row.j; col[k].i = row.i; k++;
586*c4762a1bSJed Brown 
587*c4762a1bSJed Brown         if (i+1 != info->mx-1) {
588*c4762a1bSJed Brown           v[k]     = -hydhx;
589*c4762a1bSJed Brown           col[k].j = j;     col[k].i = i+1;
590*c4762a1bSJed Brown           k++;
591*c4762a1bSJed Brown         }
592*c4762a1bSJed Brown         if (j+1 != info->mx-1) {
593*c4762a1bSJed Brown           v[k]     = -hxdhy;
594*c4762a1bSJed Brown           col[k].j = j + 1; col[k].i = i;
595*c4762a1bSJed Brown           k++;
596*c4762a1bSJed Brown         }
597*c4762a1bSJed Brown         ierr = MatSetValuesStencil(jacpre,1,&row,k,col,v,INSERT_VALUES);CHKERRQ(ierr);
598*c4762a1bSJed Brown       }
599*c4762a1bSJed Brown     }
600*c4762a1bSJed Brown   }
601*c4762a1bSJed Brown 
602*c4762a1bSJed Brown   /*
603*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
604*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
605*c4762a1bSJed Brown   */
606*c4762a1bSJed Brown   ierr = MatAssemblyBegin(jacpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
607*c4762a1bSJed Brown   ierr = MatAssemblyEnd(jacpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
608*c4762a1bSJed Brown   /*
609*c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
610*c4762a1bSJed Brown      matrix. If we do, it will generate an error.
611*c4762a1bSJed Brown   */
612*c4762a1bSJed Brown   ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
613*c4762a1bSJed Brown   PetscFunctionReturn(0);
614*c4762a1bSJed Brown }
615*c4762a1bSJed Brown 
616*c4762a1bSJed Brown #if defined(PETSC_HAVE_MATLAB_ENGINE)
617*c4762a1bSJed Brown PetscErrorCode FormFunctionMatlab(SNES snes,Vec X,Vec F,void *ptr)
618*c4762a1bSJed Brown {
619*c4762a1bSJed Brown   AppCtx         *user = (AppCtx*)ptr;
620*c4762a1bSJed Brown   PetscErrorCode ierr;
621*c4762a1bSJed Brown   PetscInt       Mx,My;
622*c4762a1bSJed Brown   PetscReal      lambda,hx,hy;
623*c4762a1bSJed Brown   Vec            localX,localF;
624*c4762a1bSJed Brown   MPI_Comm       comm;
625*c4762a1bSJed Brown   DM             da;
626*c4762a1bSJed Brown 
627*c4762a1bSJed Brown   PetscFunctionBeginUser;
628*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&da);CHKERRQ(ierr);
629*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
630*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localF);CHKERRQ(ierr);
631*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)localX,"localX");CHKERRQ(ierr);
632*c4762a1bSJed Brown   ierr = PetscObjectSetName((PetscObject)localF,"localF");CHKERRQ(ierr);
633*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
634*c4762a1bSJed Brown 
635*c4762a1bSJed Brown   lambda = user->param;
636*c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(Mx-1);
637*c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(My-1);
638*c4762a1bSJed Brown 
639*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
640*c4762a1bSJed Brown   /*
641*c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
642*c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
643*c4762a1bSJed Brown      By placing code between these two statements, computations can be
644*c4762a1bSJed Brown      done while messages are in transition.
645*c4762a1bSJed Brown   */
646*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
647*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
648*c4762a1bSJed Brown   ierr = PetscMatlabEnginePut(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localX);CHKERRQ(ierr);
649*c4762a1bSJed Brown   ierr = PetscMatlabEngineEvaluate(PETSC_MATLAB_ENGINE_(comm),"localF=ex5m(localX,%18.16e,%18.16e,%18.16e)",hx,hy,lambda);CHKERRQ(ierr);
650*c4762a1bSJed Brown   ierr = PetscMatlabEngineGet(PETSC_MATLAB_ENGINE_(comm),(PetscObject)localF);CHKERRQ(ierr);
651*c4762a1bSJed Brown 
652*c4762a1bSJed Brown   /*
653*c4762a1bSJed Brown      Insert values into global vector
654*c4762a1bSJed Brown   */
655*c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(da,localF,INSERT_VALUES,F);CHKERRQ(ierr);
656*c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(da,localF,INSERT_VALUES,F);CHKERRQ(ierr);
657*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
658*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localF);CHKERRQ(ierr);
659*c4762a1bSJed Brown   PetscFunctionReturn(0);
660*c4762a1bSJed Brown }
661*c4762a1bSJed Brown #endif
662*c4762a1bSJed Brown 
663*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
664*c4762a1bSJed Brown /*
665*c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
666*c4762a1bSJed Brown 
667*c4762a1bSJed Brown  */
668*c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X, Vec B, void *ctx)
669*c4762a1bSJed Brown {
670*c4762a1bSJed Brown   PetscInt       i,j,k,Mx,My,xs,ys,xm,ym,its,tot_its,sweeps,l;
671*c4762a1bSJed Brown   PetscErrorCode ierr;
672*c4762a1bSJed Brown   PetscReal      lambda,hx,hy,hxdhy,hydhx,sc;
673*c4762a1bSJed Brown   PetscScalar    **x,**b,bij,F,F0=0,J,u,un,us,ue,eu,uw,uxx,uyy,y;
674*c4762a1bSJed Brown   PetscReal      atol,rtol,stol;
675*c4762a1bSJed Brown   DM             da;
676*c4762a1bSJed Brown   AppCtx         *user;
677*c4762a1bSJed Brown   Vec            localX,localB;
678*c4762a1bSJed Brown 
679*c4762a1bSJed Brown   PetscFunctionBeginUser;
680*c4762a1bSJed Brown   tot_its = 0;
681*c4762a1bSJed Brown   ierr    = SNESNGSGetSweeps(snes,&sweeps);CHKERRQ(ierr);
682*c4762a1bSJed Brown   ierr    = SNESNGSGetTolerances(snes,&atol,&rtol,&stol,&its);CHKERRQ(ierr);
683*c4762a1bSJed Brown   ierr    = SNESGetDM(snes,&da);CHKERRQ(ierr);
684*c4762a1bSJed Brown   ierr    = DMGetApplicationContext(da,(void**)&user);CHKERRQ(ierr);
685*c4762a1bSJed Brown 
686*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
687*c4762a1bSJed Brown 
688*c4762a1bSJed Brown   lambda = user->param;
689*c4762a1bSJed Brown   hx     = 1.0/(PetscReal)(Mx-1);
690*c4762a1bSJed Brown   hy     = 1.0/(PetscReal)(My-1);
691*c4762a1bSJed Brown   sc     = hx*hy*lambda;
692*c4762a1bSJed Brown   hxdhy  = hx/hy;
693*c4762a1bSJed Brown   hydhx  = hy/hx;
694*c4762a1bSJed Brown 
695*c4762a1bSJed Brown 
696*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
697*c4762a1bSJed Brown   if (B) {
698*c4762a1bSJed Brown     ierr = DMGetLocalVector(da,&localB);CHKERRQ(ierr);
699*c4762a1bSJed Brown   }
700*c4762a1bSJed Brown   for (l=0; l<sweeps; l++) {
701*c4762a1bSJed Brown 
702*c4762a1bSJed Brown     ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
703*c4762a1bSJed Brown     ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
704*c4762a1bSJed Brown     if (B) {
705*c4762a1bSJed Brown       ierr = DMGlobalToLocalBegin(da,B,INSERT_VALUES,localB);CHKERRQ(ierr);
706*c4762a1bSJed Brown       ierr = DMGlobalToLocalEnd(da,B,INSERT_VALUES,localB);CHKERRQ(ierr);
707*c4762a1bSJed Brown     }
708*c4762a1bSJed Brown     /*
709*c4762a1bSJed Brown      Get a pointer to vector data.
710*c4762a1bSJed Brown      - For default PETSc vectors, VecGetArray() returns a pointer to
711*c4762a1bSJed Brown      the data array.  Otherwise, the routine is implementation dependent.
712*c4762a1bSJed Brown      - You MUST call VecRestoreArray() when you no longer need access to
713*c4762a1bSJed Brown      the array.
714*c4762a1bSJed Brown      */
715*c4762a1bSJed Brown     ierr = DMDAVecGetArray(da,localX,&x);CHKERRQ(ierr);
716*c4762a1bSJed Brown     if (B) ierr = DMDAVecGetArray(da,localB,&b);CHKERRQ(ierr);
717*c4762a1bSJed Brown     /*
718*c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
719*c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
720*c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
721*c4762a1bSJed Brown      */
722*c4762a1bSJed Brown     ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
723*c4762a1bSJed Brown 
724*c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
725*c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
726*c4762a1bSJed Brown         if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
727*c4762a1bSJed Brown           /* boundary conditions are all zero Dirichlet */
728*c4762a1bSJed Brown           x[j][i] = 0.0;
729*c4762a1bSJed Brown         } else {
730*c4762a1bSJed Brown           if (B) bij = b[j][i];
731*c4762a1bSJed Brown           else   bij = 0.;
732*c4762a1bSJed Brown 
733*c4762a1bSJed Brown           u  = x[j][i];
734*c4762a1bSJed Brown           un = x[j-1][i];
735*c4762a1bSJed Brown           us = x[j+1][i];
736*c4762a1bSJed Brown           ue = x[j][i-1];
737*c4762a1bSJed Brown           uw = x[j][i+1];
738*c4762a1bSJed Brown 
739*c4762a1bSJed Brown           for (k=0; k<its; k++) {
740*c4762a1bSJed Brown             eu  = PetscExpScalar(u);
741*c4762a1bSJed Brown             uxx = (2.0*u - ue - uw)*hydhx;
742*c4762a1bSJed Brown             uyy = (2.0*u - un - us)*hxdhy;
743*c4762a1bSJed Brown             F   = uxx + uyy - sc*eu - bij;
744*c4762a1bSJed Brown             if (k == 0) F0 = F;
745*c4762a1bSJed Brown             J  = 2.0*(hydhx + hxdhy) - sc*eu;
746*c4762a1bSJed Brown             y  = F/J;
747*c4762a1bSJed Brown             u -= y;
748*c4762a1bSJed Brown             tot_its++;
749*c4762a1bSJed Brown 
750*c4762a1bSJed Brown             if (atol > PetscAbsReal(PetscRealPart(F)) ||
751*c4762a1bSJed Brown                 rtol*PetscAbsReal(PetscRealPart(F0)) > PetscAbsReal(PetscRealPart(F)) ||
752*c4762a1bSJed Brown                 stol*PetscAbsReal(PetscRealPart(u)) > PetscAbsReal(PetscRealPart(y))) {
753*c4762a1bSJed Brown               break;
754*c4762a1bSJed Brown             }
755*c4762a1bSJed Brown           }
756*c4762a1bSJed Brown           x[j][i] = u;
757*c4762a1bSJed Brown         }
758*c4762a1bSJed Brown       }
759*c4762a1bSJed Brown     }
760*c4762a1bSJed Brown     /*
761*c4762a1bSJed Brown      Restore vector
762*c4762a1bSJed Brown      */
763*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(da,localX,&x);CHKERRQ(ierr);
764*c4762a1bSJed Brown     ierr = DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);CHKERRQ(ierr);
765*c4762a1bSJed Brown     ierr = DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);CHKERRQ(ierr);
766*c4762a1bSJed Brown   }
767*c4762a1bSJed Brown   ierr = PetscLogFlops(tot_its*(21.0));CHKERRQ(ierr);
768*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
769*c4762a1bSJed Brown   if (B) {
770*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(da,localB,&b);CHKERRQ(ierr);
771*c4762a1bSJed Brown     ierr = DMRestoreLocalVector(da,&localB);CHKERRQ(ierr);
772*c4762a1bSJed Brown   }
773*c4762a1bSJed Brown   PetscFunctionReturn(0);
774*c4762a1bSJed Brown }
775*c4762a1bSJed Brown 
776*c4762a1bSJed Brown /*TEST
777*c4762a1bSJed Brown 
778*c4762a1bSJed Brown    test:
779*c4762a1bSJed Brown      suffix: asm_0
780*c4762a1bSJed Brown      requires: !single
781*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu
782*c4762a1bSJed Brown 
783*c4762a1bSJed Brown    test:
784*c4762a1bSJed Brown      suffix: msm_0
785*c4762a1bSJed Brown      requires: !single
786*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu
787*c4762a1bSJed Brown 
788*c4762a1bSJed Brown    test:
789*c4762a1bSJed Brown      suffix: asm_1
790*c4762a1bSJed Brown      requires: !single
791*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
792*c4762a1bSJed Brown 
793*c4762a1bSJed Brown    test:
794*c4762a1bSJed Brown      suffix: msm_1
795*c4762a1bSJed Brown      requires: !single
796*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
797*c4762a1bSJed Brown 
798*c4762a1bSJed Brown    test:
799*c4762a1bSJed Brown      suffix: asm_2
800*c4762a1bSJed Brown      requires: !single
801*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
802*c4762a1bSJed Brown 
803*c4762a1bSJed Brown    test:
804*c4762a1bSJed Brown      suffix: msm_2
805*c4762a1bSJed Brown      requires: !single
806*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 3 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
807*c4762a1bSJed Brown 
808*c4762a1bSJed Brown    test:
809*c4762a1bSJed Brown      suffix: asm_3
810*c4762a1bSJed Brown      requires: !single
811*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
812*c4762a1bSJed Brown 
813*c4762a1bSJed Brown    test:
814*c4762a1bSJed Brown      suffix: msm_3
815*c4762a1bSJed Brown      requires: !single
816*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
817*c4762a1bSJed Brown 
818*c4762a1bSJed Brown    test:
819*c4762a1bSJed Brown      suffix: asm_4
820*c4762a1bSJed Brown      requires: !single
821*c4762a1bSJed Brown      nsize: 2
822*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
823*c4762a1bSJed Brown 
824*c4762a1bSJed Brown    test:
825*c4762a1bSJed Brown      suffix: msm_4
826*c4762a1bSJed Brown      requires: !single
827*c4762a1bSJed Brown      nsize: 2
828*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 2 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
829*c4762a1bSJed Brown 
830*c4762a1bSJed Brown    test:
831*c4762a1bSJed Brown      suffix: asm_5
832*c4762a1bSJed Brown      requires: !single
833*c4762a1bSJed Brown      nsize: 2
834*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type additive -sub_pc_type lu -da_grid_x 8
835*c4762a1bSJed Brown 
836*c4762a1bSJed Brown    test:
837*c4762a1bSJed Brown      suffix: msm_5
838*c4762a1bSJed Brown      requires: !single
839*c4762a1bSJed Brown      nsize: 2
840*c4762a1bSJed Brown      args: -mms 1 -par 0.0 -snes_monitor_short -snes_converged_reason -snes_view -ksp_rtol 1.0e-9 -ksp_monitor_short -ksp_type richardson -pc_type asm -pc_asm_blocks 4 -pc_asm_overlap 0 -pc_asm_local_type multiplicative -sub_pc_type lu -da_grid_x 8
841*c4762a1bSJed Brown 
842*c4762a1bSJed Brown    test:
843*c4762a1bSJed Brown      args: -snes_rtol 1.e-5 -pc_type mg -ksp_monitor_short -snes_view -pc_mg_levels 3 -pc_mg_galerkin pmat -da_grid_x 17 -da_grid_y 17 -mg_levels_ksp_monitor_short -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full
844*c4762a1bSJed Brown      requires: !single
845*c4762a1bSJed Brown 
846*c4762a1bSJed Brown    test:
847*c4762a1bSJed Brown      suffix: 2
848*c4762a1bSJed Brown      args: -pc_type mg -ksp_converged_reason -snes_view -pc_mg_galerkin pmat -snes_grid_sequence 3 -mg_levels_ksp_norm_type unpreconditioned -snes_monitor_short -mg_levels_ksp_chebyshev_esteig 0.5,1.1 -mg_levels_pc_type sor -pc_mg_type full -ksp_atol -1.
849*c4762a1bSJed Brown 
850*c4762a1bSJed Brown    test:
851*c4762a1bSJed Brown      suffix: 3
852*c4762a1bSJed Brown      nsize: 2
853*c4762a1bSJed Brown      args: -snes_grid_sequence 2 -snes_mf_operator -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -snes_rtol 1.e-2
854*c4762a1bSJed Brown      filter: grep -v "otal number of function evaluations"
855*c4762a1bSJed Brown 
856*c4762a1bSJed Brown    test:
857*c4762a1bSJed Brown      suffix: 4
858*c4762a1bSJed Brown      nsize: 2
859*c4762a1bSJed Brown      args: -snes_grid_sequence 2 -snes_monitor_short -ksp_converged_reason -snes_converged_reason -snes_view -pc_type mg -snes_atol -1 -ksp_atol -1
860*c4762a1bSJed Brown 
861*c4762a1bSJed Brown    test:
862*c4762a1bSJed Brown      suffix: 5_anderson
863*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type anderson
864*c4762a1bSJed Brown 
865*c4762a1bSJed Brown    test:
866*c4762a1bSJed Brown      suffix: 5_aspin
867*c4762a1bSJed Brown      nsize: 4
868*c4762a1bSJed Brown      args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type aspin -snes_view
869*c4762a1bSJed Brown 
870*c4762a1bSJed Brown    test:
871*c4762a1bSJed Brown      suffix: 5_broyden
872*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_qn_type broyden -snes_qn_m 10
873*c4762a1bSJed Brown 
874*c4762a1bSJed Brown    test:
875*c4762a1bSJed Brown      suffix: 5_fas
876*c4762a1bSJed Brown      args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6
877*c4762a1bSJed Brown      requires: !single
878*c4762a1bSJed Brown 
879*c4762a1bSJed Brown    test:
880*c4762a1bSJed Brown      suffix: 5_fas_additive
881*c4762a1bSJed Brown      args: -fas_coarse_snes_max_it 1 -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_monitor_short -snes_type fas -fas_coarse_ksp_type richardson -da_refine 6 -snes_fas_type additive -snes_max_it 50
882*c4762a1bSJed Brown 
883*c4762a1bSJed Brown    test:
884*c4762a1bSJed Brown      suffix: 5_fas_monitor
885*c4762a1bSJed Brown      args: -da_refine 1 -snes_type fas -snes_fas_monitor
886*c4762a1bSJed Brown      requires: !single
887*c4762a1bSJed Brown 
888*c4762a1bSJed Brown    test:
889*c4762a1bSJed Brown      suffix: 5_ls
890*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls
891*c4762a1bSJed Brown 
892*c4762a1bSJed Brown    test:
893*c4762a1bSJed Brown      suffix: 5_ls_sell_sor
894*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type newtonls -dm_mat_type sell -pc_type sor
895*c4762a1bSJed Brown      output_file: output/ex5_5_ls.out
896*c4762a1bSJed Brown 
897*c4762a1bSJed Brown    test:
898*c4762a1bSJed Brown      suffix: 5_nasm
899*c4762a1bSJed Brown      nsize: 4
900*c4762a1bSJed Brown      args: -snes_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type nasm -snes_nasm_type restrict -snes_max_it 10
901*c4762a1bSJed Brown 
902*c4762a1bSJed Brown    test:
903*c4762a1bSJed Brown      suffix: 5_ncg
904*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ncg -snes_ncg_type fr
905*c4762a1bSJed Brown 
906*c4762a1bSJed Brown    test:
907*c4762a1bSJed Brown      suffix: 5_newton_asm_dmda
908*c4762a1bSJed Brown      nsize: 4
909*c4762a1bSJed Brown      args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type asm -pc_asm_dm_subdomains -malloc_dump
910*c4762a1bSJed Brown      requires: !single
911*c4762a1bSJed Brown 
912*c4762a1bSJed Brown    test:
913*c4762a1bSJed Brown      suffix: 5_newton_gasm_dmda
914*c4762a1bSJed Brown      nsize: 4
915*c4762a1bSJed Brown      args: -snes_monitor_short -ksp_monitor_short -snes_converged_reason -da_refine 4 -da_overlap 3 -snes_type newtonls -pc_type gasm -malloc_dump
916*c4762a1bSJed Brown      requires: !single
917*c4762a1bSJed Brown 
918*c4762a1bSJed Brown    test:
919*c4762a1bSJed Brown      suffix: 5_ngmres
920*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10
921*c4762a1bSJed Brown 
922*c4762a1bSJed Brown    test:
923*c4762a1bSJed Brown      suffix: 5_ngmres_fas
924*c4762a1bSJed Brown      args: -snes_rtol 1.e-4 -snes_type ngmres -npc_fas_coarse_snes_max_it 1 -npc_fas_coarse_snes_type newtonls -npc_fas_coarse_pc_type lu -npc_fas_coarse_ksp_type preonly -snes_ngmres_m 10 -snes_monitor_short -npc_snes_max_it 1 -npc_snes_type fas -npc_fas_coarse_ksp_type richardson -da_refine 6
925*c4762a1bSJed Brown 
926*c4762a1bSJed Brown    test:
927*c4762a1bSJed Brown      suffix: 5_ngmres_ngs
928*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -npc_snes_type ngs -npc_snes_max_it 1
929*c4762a1bSJed Brown 
930*c4762a1bSJed Brown    test:
931*c4762a1bSJed Brown      suffix: 5_ngmres_nrichardson
932*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type ngmres -snes_ngmres_m 10 -npc_snes_type nrichardson -npc_snes_max_it 3
933*c4762a1bSJed Brown      output_file: output/ex5_5_ngmres_richardson.out
934*c4762a1bSJed Brown 
935*c4762a1bSJed Brown    test:
936*c4762a1bSJed Brown      suffix: 5_nrichardson
937*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type nrichardson
938*c4762a1bSJed Brown 
939*c4762a1bSJed Brown    test:
940*c4762a1bSJed Brown      suffix: 5_qn
941*c4762a1bSJed Brown      args: -da_grid_x 81 -da_grid_y 81 -snes_monitor_short -snes_max_it 50 -par 6.0 -snes_type qn -snes_linesearch_type cp -snes_qn_m 10
942*c4762a1bSJed Brown 
943*c4762a1bSJed Brown    test:
944*c4762a1bSJed Brown      suffix: 6
945*c4762a1bSJed Brown      nsize: 4
946*c4762a1bSJed Brown      args: -snes_converged_reason -ksp_converged_reason -da_grid_x 129 -da_grid_y 129 -pc_type mg -pc_mg_levels 8 -mg_levels_ksp_type chebyshev -mg_levels_ksp_chebyshev_esteig 0,0.5,0,1.1 -mg_levels_ksp_max_it 2
947*c4762a1bSJed Brown 
948*c4762a1bSJed Brown    test:
949*c4762a1bSJed Brown      requires: complex !single
950*c4762a1bSJed Brown      suffix: complex
951*c4762a1bSJed Brown      args: -snes_mf_operator -mat_mffd_complex -snes_monitor
952*c4762a1bSJed Brown 
953*c4762a1bSJed Brown TEST*/
954