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