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