xref: /petsc/src/snes/tutorials/ex35.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static const char help[] = "-Laplacian u = b as a nonlinear problem.\n\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*T
4*c4762a1bSJed Brown    Concepts: SNES^parallel Bratu example
5*c4762a1bSJed Brown    Concepts: DMDA^using distributed arrays;
6*c4762a1bSJed Brown    Concepts: IS coloirng types;
7*c4762a1bSJed Brown    Processors: n
8*c4762a1bSJed Brown T*/
9*c4762a1bSJed Brown 
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown /*
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown     The linear and nonlinear versions of these should give almost identical results on this problem
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown     Richardson
17*c4762a1bSJed Brown       Nonlinear:
18*c4762a1bSJed Brown         -snes_rtol 1.e-12 -snes_monitor -snes_type nrichardson -snes_linesearch_monitor
19*c4762a1bSJed Brown 
20*c4762a1bSJed Brown       Linear:
21*c4762a1bSJed Brown         -snes_rtol 1.e-12 -snes_monitor -ksp_rtol 1.e-12  -ksp_monitor -ksp_type richardson -pc_type none -ksp_richardson_self_scale -info
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown     GMRES
24*c4762a1bSJed Brown       Nonlinear:
25*c4762a1bSJed Brown        -snes_rtol 1.e-12 -snes_monitor  -snes_type ngmres
26*c4762a1bSJed Brown 
27*c4762a1bSJed Brown       Linear:
28*c4762a1bSJed Brown        -snes_rtol 1.e-12 -snes_monitor  -ksp_type gmres -ksp_monitor -ksp_rtol 1.e-12 -pc_type none
29*c4762a1bSJed Brown 
30*c4762a1bSJed Brown     CG
31*c4762a1bSJed Brown        Nonlinear:
32*c4762a1bSJed Brown             -snes_rtol 1.e-12 -snes_monitor  -snes_type ncg -snes_linesearch_monitor
33*c4762a1bSJed Brown 
34*c4762a1bSJed Brown        Linear:
35*c4762a1bSJed Brown              -snes_rtol 1.e-12 -snes_monitor  -ksp_type cg -ksp_monitor -ksp_rtol 1.e-12 -pc_type none
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown     Multigrid
38*c4762a1bSJed Brown        Linear:
39*c4762a1bSJed Brown           1 level:
40*c4762a1bSJed Brown             -snes_rtol 1.e-12 -snes_monitor  -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type none -mg_levels_ksp_monitor
41*c4762a1bSJed Brown             -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor -ksp_rtol 1.e-12  -ksp_monitor_true_residual
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown           n levels:
44*c4762a1bSJed Brown             -da_refine n
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown        Nonlinear:
47*c4762a1bSJed Brown          1 level:
48*c4762a1bSJed Brown            -snes_rtol 1.e-12 -snes_monitor  -snes_type fas -fas_levels_snes_monitor
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown           n levels:
51*c4762a1bSJed Brown             -da_refine n  -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown */
54*c4762a1bSJed Brown 
55*c4762a1bSJed Brown /*
56*c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
57*c4762a1bSJed Brown    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
58*c4762a1bSJed Brown */
59*c4762a1bSJed Brown #include <petscdm.h>
60*c4762a1bSJed Brown #include <petscdmda.h>
61*c4762a1bSJed Brown #include <petscsnes.h>
62*c4762a1bSJed Brown 
63*c4762a1bSJed Brown /*
64*c4762a1bSJed Brown    User-defined routines
65*c4762a1bSJed Brown */
66*c4762a1bSJed Brown extern PetscErrorCode FormMatrix(DM,Mat);
67*c4762a1bSJed Brown extern PetscErrorCode MyComputeFunction(SNES,Vec,Vec,void*);
68*c4762a1bSJed Brown extern PetscErrorCode MyComputeJacobian(SNES,Vec,Mat,Mat,void*);
69*c4762a1bSJed Brown extern PetscErrorCode NonlinearGS(SNES,Vec);
70*c4762a1bSJed Brown 
71*c4762a1bSJed Brown int main(int argc,char **argv)
72*c4762a1bSJed Brown {
73*c4762a1bSJed Brown   SNES           snes;                                 /* nonlinear solver */
74*c4762a1bSJed Brown   SNES           psnes;                                /* nonlinear Gauss-Seidel approximate solver */
75*c4762a1bSJed Brown   Vec            x,b;                                  /* solution vector */
76*c4762a1bSJed Brown   PetscInt       its;                                  /* iterations for convergence */
77*c4762a1bSJed Brown   PetscErrorCode ierr;
78*c4762a1bSJed Brown   DM             da;
79*c4762a1bSJed Brown   PetscBool      use_ngs_as_npc = PETSC_FALSE;                /* use the nonlinear Gauss-Seidel approximate solver */
80*c4762a1bSJed Brown 
81*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82*c4762a1bSJed Brown      Initialize program
83*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
86*c4762a1bSJed Brown 
87*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88*c4762a1bSJed Brown      Create nonlinear solver context
89*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90*c4762a1bSJed Brown   ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-use_ngs_as_npc",&use_ngs_as_npc,0);CHKERRQ(ierr);
93*c4762a1bSJed Brown 
94*c4762a1bSJed Brown   if (use_ngs_as_npc) {
95*c4762a1bSJed Brown     ierr = SNESGetNPC(snes,&psnes);CHKERRQ(ierr);
96*c4762a1bSJed Brown     ierr = SNESSetType(psnes,SNESSHELL);CHKERRQ(ierr);
97*c4762a1bSJed Brown     ierr = SNESShellSetSolve(psnes,NonlinearGS);CHKERRQ(ierr);
98*c4762a1bSJed Brown   }
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101*c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
102*c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103*c4762a1bSJed Brown   ierr = DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = SNESSetDM(snes,da);CHKERRQ(ierr);
108*c4762a1bSJed Brown   if (use_ngs_as_npc) {
109*c4762a1bSJed Brown     ierr = SNESShellSetContext(psnes,da);CHKERRQ(ierr);
110*c4762a1bSJed Brown   }
111*c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112*c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
113*c4762a1bSJed Brown      vectors that are the same types
114*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&x);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&b);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = VecSet(b,1.0);CHKERRQ(ierr);
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown   ierr = SNESSetFunction(snes,NULL,MyComputeFunction,NULL);CHKERRQ(ierr);
120*c4762a1bSJed Brown   ierr = SNESSetJacobian(snes,NULL,NULL,MyComputeJacobian,NULL);CHKERRQ(ierr);
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123*c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
124*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125*c4762a1bSJed Brown   ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
126*c4762a1bSJed Brown 
127*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128*c4762a1bSJed Brown      Solve nonlinear system
129*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130*c4762a1bSJed Brown   ierr = SNESSolve(snes,b,x);CHKERRQ(ierr);
131*c4762a1bSJed Brown   ierr = SNESGetIterationNumber(snes,&its);CHKERRQ(ierr);
132*c4762a1bSJed Brown 
133*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135*c4762a1bSJed Brown 
136*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
137*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
138*c4762a1bSJed Brown      are no longer needed.
139*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
140*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
141*c4762a1bSJed Brown   ierr = VecDestroy(&b);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = SNESDestroy(&snes);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = PetscFinalize();
145*c4762a1bSJed Brown   return ierr;
146*c4762a1bSJed Brown }
147*c4762a1bSJed Brown 
148*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
149*c4762a1bSJed Brown PetscErrorCode MyComputeFunction(SNES snes,Vec x,Vec F,void *ctx)
150*c4762a1bSJed Brown {
151*c4762a1bSJed Brown   PetscErrorCode ierr;
152*c4762a1bSJed Brown   Mat            J;
153*c4762a1bSJed Brown   DM             dm;
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown   PetscFunctionBeginUser;
156*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
157*c4762a1bSJed Brown   ierr = DMGetApplicationContext(dm,&J);CHKERRQ(ierr);
158*c4762a1bSJed Brown   if (!J) {
159*c4762a1bSJed Brown     ierr = DMSetMatType(dm,MATAIJ);CHKERRQ(ierr);
160*c4762a1bSJed Brown     ierr = DMCreateMatrix(dm,&J);CHKERRQ(ierr);
161*c4762a1bSJed Brown     ierr = MatSetDM(J, NULL);CHKERRQ(ierr);
162*c4762a1bSJed Brown     ierr = FormMatrix(dm,J);CHKERRQ(ierr);
163*c4762a1bSJed Brown     ierr = DMSetApplicationContext(dm,J);CHKERRQ(ierr);
164*c4762a1bSJed Brown     ierr = DMSetApplicationContextDestroy(dm,(PetscErrorCode (*)(void**))MatDestroy);CHKERRQ(ierr);
165*c4762a1bSJed Brown   }
166*c4762a1bSJed Brown   ierr = MatMult(J,x,F);CHKERRQ(ierr);
167*c4762a1bSJed Brown   PetscFunctionReturn(0);
168*c4762a1bSJed Brown }
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown PetscErrorCode MyComputeJacobian(SNES snes,Vec x,Mat J,Mat Jp,void *ctx)
171*c4762a1bSJed Brown {
172*c4762a1bSJed Brown   PetscErrorCode ierr;
173*c4762a1bSJed Brown   DM             dm;
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown   PetscFunctionBeginUser;
176*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
177*c4762a1bSJed Brown   ierr = FormMatrix(dm,Jp);CHKERRQ(ierr);
178*c4762a1bSJed Brown   PetscFunctionReturn(0);
179*c4762a1bSJed Brown }
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown PetscErrorCode FormMatrix(DM da,Mat jac)
182*c4762a1bSJed Brown {
183*c4762a1bSJed Brown   PetscErrorCode ierr;
184*c4762a1bSJed Brown   PetscInt       i,j,nrows = 0;
185*c4762a1bSJed Brown   MatStencil     col[5],row,*rows;
186*c4762a1bSJed Brown   PetscScalar    v[5],hx,hy,hxdhy,hydhx;
187*c4762a1bSJed Brown   DMDALocalInfo  info;
188*c4762a1bSJed Brown 
189*c4762a1bSJed Brown   PetscFunctionBeginUser;
190*c4762a1bSJed Brown   ierr  = DMDAGetLocalInfo(da,&info);CHKERRQ(ierr);
191*c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(info.mx-1);
192*c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(info.my-1);
193*c4762a1bSJed Brown   hxdhy = hx/hy;
194*c4762a1bSJed Brown   hydhx = hy/hx;
195*c4762a1bSJed Brown 
196*c4762a1bSJed Brown   ierr = PetscMalloc1(info.ym*info.xm,&rows);CHKERRQ(ierr);
197*c4762a1bSJed Brown   /*
198*c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
199*c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
200*c4762a1bSJed Brown         contiguous chunks of rows across the processors.
201*c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
202*c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
203*c4762a1bSJed Brown         appropriate processor during matrix assembly).
204*c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
205*c4762a1bSJed Brown       - We can set matrix entries either using either
206*c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues(), as discussed above.
207*c4762a1bSJed Brown   */
208*c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
209*c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
210*c4762a1bSJed Brown       row.j = j; row.i = i;
211*c4762a1bSJed Brown       /* boundary points */
212*c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
213*c4762a1bSJed Brown         v[0]            = 2.0*(hydhx + hxdhy);
214*c4762a1bSJed Brown         ierr            = MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);CHKERRQ(ierr);
215*c4762a1bSJed Brown         rows[nrows].i   = i;
216*c4762a1bSJed Brown         rows[nrows++].j = j;
217*c4762a1bSJed Brown       } else {
218*c4762a1bSJed Brown         /* interior grid points */
219*c4762a1bSJed Brown         v[0] = -hxdhy;                                           col[0].j = j - 1; col[0].i = i;
220*c4762a1bSJed Brown         v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
221*c4762a1bSJed Brown         v[2] = 2.0*(hydhx + hxdhy);                              col[2].j = row.j; col[2].i = row.i;
222*c4762a1bSJed Brown         v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
223*c4762a1bSJed Brown         v[4] = -hxdhy;                                           col[4].j = j + 1; col[4].i = i;
224*c4762a1bSJed Brown         ierr = MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);CHKERRQ(ierr);
225*c4762a1bSJed Brown       }
226*c4762a1bSJed Brown     }
227*c4762a1bSJed Brown   }
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   /*
230*c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
231*c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
232*c4762a1bSJed Brown   */
233*c4762a1bSJed Brown   ierr = MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = MatZeroRowsColumnsStencil(jac,nrows,rows,2.0*(hydhx + hxdhy),NULL,NULL);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = PetscFree(rows);CHKERRQ(ierr);
237*c4762a1bSJed Brown   /*
238*c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
239*c4762a1bSJed Brown      matrix. If we do, it will generate an error.
240*c4762a1bSJed Brown   */
241*c4762a1bSJed Brown   ierr = MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
242*c4762a1bSJed Brown   PetscFunctionReturn(0);
243*c4762a1bSJed Brown }
244*c4762a1bSJed Brown 
245*c4762a1bSJed Brown 
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown /* ------------------------------------------------------------------- */
248*c4762a1bSJed Brown /*
249*c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown  */
252*c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X)
253*c4762a1bSJed Brown {
254*c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym,its,l;
255*c4762a1bSJed Brown   PetscErrorCode ierr;
256*c4762a1bSJed Brown   PetscReal      hx,hy,hxdhy,hydhx;
257*c4762a1bSJed Brown   PetscScalar    **x,F,J,u,uxx,uyy;
258*c4762a1bSJed Brown   DM             da;
259*c4762a1bSJed Brown   Vec            localX;
260*c4762a1bSJed Brown 
261*c4762a1bSJed Brown   PetscFunctionBeginUser;
262*c4762a1bSJed Brown   ierr = SNESGetTolerances(snes,NULL,NULL,NULL,&its,NULL);CHKERRQ(ierr);
263*c4762a1bSJed Brown   ierr = SNESShellGetContext(snes,(void**)&da);CHKERRQ(ierr);
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown   ierr = DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);CHKERRQ(ierr);
266*c4762a1bSJed Brown 
267*c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(Mx-1);
268*c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(My-1);
269*c4762a1bSJed Brown   hxdhy = hx/hy;
270*c4762a1bSJed Brown   hydhx = hy/hx;
271*c4762a1bSJed Brown 
272*c4762a1bSJed Brown 
273*c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localX);CHKERRQ(ierr);
274*c4762a1bSJed Brown 
275*c4762a1bSJed Brown   for (l=0; l<its; l++) {
276*c4762a1bSJed Brown 
277*c4762a1bSJed Brown     ierr = DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
278*c4762a1bSJed Brown     ierr = DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);CHKERRQ(ierr);
279*c4762a1bSJed Brown     /*
280*c4762a1bSJed Brown      Get a pointer to vector data.
281*c4762a1bSJed Brown      - For default PETSc vectors, VecGetArray() returns a pointer to
282*c4762a1bSJed Brown      the data array.  Otherwise, the routine is implementation dependent.
283*c4762a1bSJed Brown      - You MUST call VecRestoreArray() when you no longer need access to
284*c4762a1bSJed Brown      the array.
285*c4762a1bSJed Brown      */
286*c4762a1bSJed Brown     ierr = DMDAVecGetArray(da,localX,&x);CHKERRQ(ierr);
287*c4762a1bSJed Brown 
288*c4762a1bSJed Brown     /*
289*c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
290*c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
291*c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
292*c4762a1bSJed Brown 
293*c4762a1bSJed Brown      */
294*c4762a1bSJed Brown     ierr = DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);CHKERRQ(ierr);
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
297*c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
298*c4762a1bSJed Brown         if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
299*c4762a1bSJed Brown           /* boundary conditions are all zero Dirichlet */
300*c4762a1bSJed Brown           x[j][i] = 0.0;
301*c4762a1bSJed Brown         } else {
302*c4762a1bSJed Brown           u   = x[j][i];
303*c4762a1bSJed Brown           uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx;
304*c4762a1bSJed Brown           uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy;
305*c4762a1bSJed Brown           F   = uxx + uyy;
306*c4762a1bSJed Brown           J   = 2.0*(hydhx + hxdhy);
307*c4762a1bSJed Brown           u   = u - F/J;
308*c4762a1bSJed Brown 
309*c4762a1bSJed Brown           x[j][i] = u;
310*c4762a1bSJed Brown         }
311*c4762a1bSJed Brown       }
312*c4762a1bSJed Brown     }
313*c4762a1bSJed Brown 
314*c4762a1bSJed Brown     /*
315*c4762a1bSJed Brown      Restore vector
316*c4762a1bSJed Brown      */
317*c4762a1bSJed Brown     ierr = DMDAVecRestoreArray(da,localX,&x);CHKERRQ(ierr);
318*c4762a1bSJed Brown     ierr = DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);CHKERRQ(ierr);
319*c4762a1bSJed Brown     ierr = DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);CHKERRQ(ierr);
320*c4762a1bSJed Brown   }
321*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localX);CHKERRQ(ierr);
322*c4762a1bSJed Brown   PetscFunctionReturn(0);
323*c4762a1bSJed Brown }
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown 
326*c4762a1bSJed Brown /*TEST
327*c4762a1bSJed Brown 
328*c4762a1bSJed Brown    test:
329*c4762a1bSJed Brown       args: -snes_monitor_short -snes_type nrichardson
330*c4762a1bSJed Brown       requires: !single
331*c4762a1bSJed Brown 
332*c4762a1bSJed Brown    test:
333*c4762a1bSJed Brown       suffix: 2
334*c4762a1bSJed Brown       args: -snes_monitor_short -ksp_monitor_short -ksp_type richardson -pc_type none -ksp_richardson_self_scale
335*c4762a1bSJed Brown       requires: !single
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown    test:
338*c4762a1bSJed Brown       suffix: 3
339*c4762a1bSJed Brown       args: -snes_monitor_short -snes_type ngmres
340*c4762a1bSJed Brown 
341*c4762a1bSJed Brown    test:
342*c4762a1bSJed Brown       suffix: 4
343*c4762a1bSJed Brown       args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none
344*c4762a1bSJed Brown 
345*c4762a1bSJed Brown    test:
346*c4762a1bSJed Brown       suffix: 5
347*c4762a1bSJed Brown       args: -snes_monitor_short -snes_type ncg
348*c4762a1bSJed Brown 
349*c4762a1bSJed Brown    test:
350*c4762a1bSJed Brown       suffix: 6
351*c4762a1bSJed Brown       args: -snes_monitor_short -ksp_type cg -ksp_monitor_short -pc_type none
352*c4762a1bSJed Brown 
353*c4762a1bSJed Brown    test:
354*c4762a1bSJed Brown       suffix: 7
355*c4762a1bSJed Brown       args: -da_refine 2 -snes_monitor_short -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type none -mg_levels_ksp_monitor_short -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor_short
356*c4762a1bSJed Brown       requires: !single
357*c4762a1bSJed Brown 
358*c4762a1bSJed Brown    test:
359*c4762a1bSJed Brown       suffix: 8
360*c4762a1bSJed Brown       args: -da_refine 2 -snes_monitor_short -snes_type fas -fas_levels_snes_monitor_short -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly -snes_type fas -snes_rtol 1.e-5
361*c4762a1bSJed Brown 
362*c4762a1bSJed Brown TEST*/
363