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