xref: /petsc/src/snes/tutorials/ex35.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   DM             da;
76c4762a1bSJed Brown   PetscBool      use_ngs_as_npc = PETSC_FALSE;                /* use the nonlinear Gauss-Seidel approximate solver */
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79c4762a1bSJed Brown      Initialize program
80c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81c4762a1bSJed Brown 
82*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85c4762a1bSJed Brown      Create nonlinear solver context
86c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
875f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESCreate(PETSC_COMM_WORLD,&snes));
88c4762a1bSJed Brown 
895f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-use_ngs_as_npc",&use_ngs_as_npc,0));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   if (use_ngs_as_npc) {
925f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESGetNPC(snes,&psnes));
935f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESSetType(psnes,SNESSHELL));
945f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESShellSetSolve(psnes,NonlinearGS));
95c4762a1bSJed Brown   }
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
99c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1005f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetDM(snes,da));
105c4762a1bSJed Brown   if (use_ngs_as_npc) {
1065f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESShellSetContext(psnes,da));
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
110c4762a1bSJed Brown      vectors that are the same types
111c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&b));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(b,1.0));
115c4762a1bSJed Brown 
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFunction(snes,NULL,MyComputeFunction,NULL));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetJacobian(snes,NULL,NULL,MyComputeJacobian,NULL));
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120c4762a1bSJed Brown      Customize nonlinear solver; set runtime options
121c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSetFromOptions(snes));
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125c4762a1bSJed Brown      Solve nonlinear system
126c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESSolve(snes,b,x));
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetIterationNumber(snes,&its));
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
135c4762a1bSJed Brown      are no longer needed.
136c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&b));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESDestroy(&snes));
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
141*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
142*b122ec5aSJacob Faibussowitsch   return 0;
143c4762a1bSJed Brown }
144c4762a1bSJed Brown 
145c4762a1bSJed Brown /* ------------------------------------------------------------------- */
146c4762a1bSJed Brown PetscErrorCode MyComputeFunction(SNES snes,Vec x,Vec F,void *ctx)
147c4762a1bSJed Brown {
148c4762a1bSJed Brown   Mat            J;
149c4762a1bSJed Brown   DM             dm;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   PetscFunctionBeginUser;
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetApplicationContext(dm,&J));
154c4762a1bSJed Brown   if (!J) {
1555f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetMatType(dm,MATAIJ));
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(DMCreateMatrix(dm,&J));
1575f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetDM(J, NULL));
1585f80ce2aSJacob Faibussowitsch     CHKERRQ(FormMatrix(dm,J));
1595f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetApplicationContext(dm,J));
1605f80ce2aSJacob Faibussowitsch     CHKERRQ(DMSetApplicationContextDestroy(dm,(PetscErrorCode (*)(void**))MatDestroy));
161c4762a1bSJed Brown   }
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(MatMult(J,x,F));
163c4762a1bSJed Brown   PetscFunctionReturn(0);
164c4762a1bSJed Brown }
165c4762a1bSJed Brown 
166c4762a1bSJed Brown PetscErrorCode MyComputeJacobian(SNES snes,Vec x,Mat J,Mat Jp,void *ctx)
167c4762a1bSJed Brown {
168c4762a1bSJed Brown   DM             dm;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   PetscFunctionBeginUser;
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetDM(snes,&dm));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(FormMatrix(dm,Jp));
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown PetscErrorCode FormMatrix(DM da,Mat jac)
177c4762a1bSJed Brown {
178c4762a1bSJed Brown   PetscInt       i,j,nrows = 0;
179c4762a1bSJed Brown   MatStencil     col[5],row,*rows;
180c4762a1bSJed Brown   PetscScalar    v[5],hx,hy,hxdhy,hydhx;
181c4762a1bSJed Brown   DMDALocalInfo  info;
182c4762a1bSJed Brown 
183c4762a1bSJed Brown   PetscFunctionBeginUser;
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetLocalInfo(da,&info));
185c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(info.mx-1);
186c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(info.my-1);
187c4762a1bSJed Brown   hxdhy = hx/hy;
188c4762a1bSJed Brown   hydhx = hy/hx;
189c4762a1bSJed Brown 
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscMalloc1(info.ym*info.xm,&rows));
191c4762a1bSJed Brown   /*
192c4762a1bSJed Brown      Compute entries for the locally owned part of the Jacobian.
193c4762a1bSJed Brown       - Currently, all PETSc parallel matrix formats are partitioned by
194c4762a1bSJed Brown         contiguous chunks of rows across the processors.
195c4762a1bSJed Brown       - Each processor needs to insert only elements that it owns
196c4762a1bSJed Brown         locally (but any non-local elements will be sent to the
197c4762a1bSJed Brown         appropriate processor during matrix assembly).
198c4762a1bSJed Brown       - Here, we set all entries for a particular row at once.
199c4762a1bSJed Brown       - We can set matrix entries either using either
200c4762a1bSJed Brown         MatSetValuesLocal() or MatSetValues(), as discussed above.
201c4762a1bSJed Brown   */
202c4762a1bSJed Brown   for (j=info.ys; j<info.ys+info.ym; j++) {
203c4762a1bSJed Brown     for (i=info.xs; i<info.xs+info.xm; i++) {
204c4762a1bSJed Brown       row.j = j; row.i = i;
205c4762a1bSJed Brown       /* boundary points */
206c4762a1bSJed Brown       if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
207c4762a1bSJed Brown         v[0]            = 2.0*(hydhx + hxdhy);
2085f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES));
209c4762a1bSJed Brown         rows[nrows].i   = i;
210c4762a1bSJed Brown         rows[nrows++].j = j;
211c4762a1bSJed Brown       } else {
212c4762a1bSJed Brown         /* interior grid points */
213c4762a1bSJed Brown         v[0] = -hxdhy;                                           col[0].j = j - 1; col[0].i = i;
214c4762a1bSJed Brown         v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
215c4762a1bSJed Brown         v[2] = 2.0*(hydhx + hxdhy);                              col[2].j = row.j; col[2].i = row.i;
216c4762a1bSJed Brown         v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
217c4762a1bSJed Brown         v[4] = -hxdhy;                                           col[4].j = j + 1; col[4].i = i;
2185f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES));
219c4762a1bSJed Brown       }
220c4762a1bSJed Brown     }
221c4762a1bSJed Brown   }
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /*
224c4762a1bSJed Brown      Assemble matrix, using the 2-step process:
225c4762a1bSJed Brown        MatAssemblyBegin(), MatAssemblyEnd().
226c4762a1bSJed Brown   */
2275f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroRowsColumnsStencil(jac,nrows,rows,2.0*(hydhx + hxdhy),NULL,NULL));
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(rows));
231c4762a1bSJed Brown   /*
232c4762a1bSJed Brown      Tell the matrix we will never add a new nonzero location to the
233c4762a1bSJed Brown      matrix. If we do, it will generate an error.
234c4762a1bSJed Brown   */
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE));
236c4762a1bSJed Brown   PetscFunctionReturn(0);
237c4762a1bSJed Brown }
238c4762a1bSJed Brown 
239c4762a1bSJed Brown /* ------------------------------------------------------------------- */
240c4762a1bSJed Brown /*
241c4762a1bSJed Brown       Applies some sweeps on nonlinear Gauss-Seidel on each process
242c4762a1bSJed Brown 
243c4762a1bSJed Brown  */
244c4762a1bSJed Brown PetscErrorCode NonlinearGS(SNES snes,Vec X)
245c4762a1bSJed Brown {
246c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym,its,l;
247c4762a1bSJed Brown   PetscReal      hx,hy,hxdhy,hydhx;
248c4762a1bSJed Brown   PetscScalar    **x,F,J,u,uxx,uyy;
249c4762a1bSJed Brown   DM             da;
250c4762a1bSJed Brown   Vec            localX;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   PetscFunctionBeginUser;
2535f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetTolerances(snes,NULL,NULL,NULL,&its,NULL));
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESShellGetContext(snes,&da));
255c4762a1bSJed Brown 
2565f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   hx    = 1.0/(PetscReal)(Mx-1);
259c4762a1bSJed Brown   hy    = 1.0/(PetscReal)(My-1);
260c4762a1bSJed Brown   hxdhy = hx/hy;
261c4762a1bSJed Brown   hydhx = hy/hx;
262c4762a1bSJed Brown 
2635f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   for (l=0; l<its; l++) {
266c4762a1bSJed Brown 
2675f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
2685f80ce2aSJacob Faibussowitsch     CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
269c4762a1bSJed Brown     /*
270c4762a1bSJed Brown      Get a pointer to vector data.
271c4762a1bSJed Brown      - For default PETSc vectors, VecGetArray() returns a pointer to
272c4762a1bSJed Brown      the data array.  Otherwise, the routine is implementation dependent.
273c4762a1bSJed Brown      - You MUST call VecRestoreArray() when you no longer need access to
274c4762a1bSJed Brown      the array.
275c4762a1bSJed Brown      */
2765f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecGetArray(da,localX,&x));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown     /*
279c4762a1bSJed Brown      Get local grid boundaries (for 2-dimensional DMDA):
280c4762a1bSJed Brown      xs, ys   - starting grid indices (no ghost points)
281c4762a1bSJed Brown      xm, ym   - widths of local grid (no ghost points)
282c4762a1bSJed Brown 
283c4762a1bSJed Brown      */
2845f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
285c4762a1bSJed Brown 
286c4762a1bSJed Brown     for (j=ys; j<ys+ym; j++) {
287c4762a1bSJed Brown       for (i=xs; i<xs+xm; i++) {
288c4762a1bSJed Brown         if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
289c4762a1bSJed Brown           /* boundary conditions are all zero Dirichlet */
290c4762a1bSJed Brown           x[j][i] = 0.0;
291c4762a1bSJed Brown         } else {
292c4762a1bSJed Brown           u   = x[j][i];
293c4762a1bSJed Brown           uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx;
294c4762a1bSJed Brown           uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy;
295c4762a1bSJed Brown           F   = uxx + uyy;
296c4762a1bSJed Brown           J   = 2.0*(hydhx + hxdhy);
297c4762a1bSJed Brown           u   = u - F/J;
298c4762a1bSJed Brown 
299c4762a1bSJed Brown           x[j][i] = u;
300c4762a1bSJed Brown         }
301c4762a1bSJed Brown       }
302c4762a1bSJed Brown     }
303c4762a1bSJed Brown 
304c4762a1bSJed Brown     /*
305c4762a1bSJed Brown      Restore vector
306c4762a1bSJed Brown      */
3075f80ce2aSJacob Faibussowitsch     CHKERRQ(DMDAVecRestoreArray(da,localX,&x));
3085f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X));
3095f80ce2aSJacob Faibussowitsch     CHKERRQ(DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X));
310c4762a1bSJed Brown   }
3115f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
312c4762a1bSJed Brown   PetscFunctionReturn(0);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown 
315c4762a1bSJed Brown /*TEST
316c4762a1bSJed Brown 
317c4762a1bSJed Brown    test:
318c4762a1bSJed Brown       args: -snes_monitor_short -snes_type nrichardson
319c4762a1bSJed Brown       requires: !single
320c4762a1bSJed Brown 
321c4762a1bSJed Brown    test:
322c4762a1bSJed Brown       suffix: 2
323c4762a1bSJed Brown       args: -snes_monitor_short -ksp_monitor_short -ksp_type richardson -pc_type none -ksp_richardson_self_scale
324c4762a1bSJed Brown       requires: !single
325c4762a1bSJed Brown 
326c4762a1bSJed Brown    test:
327c4762a1bSJed Brown       suffix: 3
328c4762a1bSJed Brown       args: -snes_monitor_short -snes_type ngmres
329c4762a1bSJed Brown 
330c4762a1bSJed Brown    test:
331c4762a1bSJed Brown       suffix: 4
332c4762a1bSJed Brown       args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none
333c4762a1bSJed Brown 
334c4762a1bSJed Brown    test:
335c4762a1bSJed Brown       suffix: 5
336c4762a1bSJed Brown       args: -snes_monitor_short -snes_type ncg
337c4762a1bSJed Brown 
338c4762a1bSJed Brown    test:
339c4762a1bSJed Brown       suffix: 6
340c4762a1bSJed Brown       args: -snes_monitor_short -ksp_type cg -ksp_monitor_short -pc_type none
341c4762a1bSJed Brown 
342c4762a1bSJed Brown    test:
343c4762a1bSJed Brown       suffix: 7
344c4762a1bSJed 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
345c4762a1bSJed Brown       requires: !single
346c4762a1bSJed Brown 
347c4762a1bSJed Brown    test:
348c4762a1bSJed Brown       suffix: 8
349c4762a1bSJed 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
350c4762a1bSJed Brown 
35141ba4c6cSHeeho Park    test:
35241ba4c6cSHeeho Park       suffix: 9
35341ba4c6cSHeeho Park       args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc
35441ba4c6cSHeeho Park 
35541ba4c6cSHeeho Park    test:
35641ba4c6cSHeeho Park       suffix: 10
35741ba4c6cSHeeho Park       args: -snes_monitor_short -ksp_type gmres -ksp_monitor_short -pc_type none -snes_type newtontrdc -snes_trdc_use_cauchy false
35841ba4c6cSHeeho Park 
359c4762a1bSJed Brown TEST*/
360