Actual source code: ex35.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: static const char help[] = "-Laplacian u = b as a nonlinear problem.\n\n";

  3: /*T
  4:    Concepts: SNES^parallel Bratu example
  5:    Concepts: DMDA^using distributed arrays;
  6:    Concepts: IS coloirng types;
  7:    Processors: n
  8: T*/

 10: /*

 12:     The linear and nonlinear versions of these should give almost identical results on this problem

 14:     Richardson
 15:       Nonlinear:
 16:         -snes_rtol 1.e-12 -snes_monitor -snes_type nrichardson -snes_linesearch_monitor

 18:       Linear:
 19:         -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

 21:     GMRES
 22:       Nonlinear:
 23:        -snes_rtol 1.e-12 -snes_monitor  -snes_type ngmres

 25:       Linear:
 26:        -snes_rtol 1.e-12 -snes_monitor  -ksp_type gmres -ksp_monitor -ksp_rtol 1.e-12 -pc_type none

 28:     CG
 29:        Nonlinear:
 30:             -snes_rtol 1.e-12 -snes_monitor  -snes_type ncg -snes_linesearch_monitor

 32:        Linear:
 33:              -snes_rtol 1.e-12 -snes_monitor  -ksp_type cg -ksp_monitor -ksp_rtol 1.e-12 -pc_type none

 35:     Multigrid
 36:        Linear:
 37:           1 level:
 38:             -snes_rtol 1.e-12 -snes_monitor  -pc_type mg -mg_levels_ksp_type richardson -mg_levels_pc_type none -mg_levels_ksp_monitor
 39:             -mg_levels_ksp_richardson_self_scale -ksp_type richardson -ksp_monitor -ksp_rtol 1.e-12  -ksp_monitor_true_residual

 41:           n levels:
 42:             -da_refine n

 44:        Nonlinear:
 45:          1 level:
 46:            -snes_rtol 1.e-12 -snes_monitor  -snes_type fas -fas_levels_snes_monitor

 48:           n levels:
 49:             -da_refine n  -fas_coarse_snes_type newtonls -fas_coarse_pc_type lu -fas_coarse_ksp_type preonly

 51: */

 53: /*
 54:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 55:    Include "petscsnes.h" so that we can use SNES solvers.  Note that this
 56: */
 57: #include <petscdm.h>
 58: #include <petscdmda.h>
 59: #include <petscsnes.h>

 61: /*
 62:    User-defined routines
 63: */
 64: extern PetscErrorCode FormMatrix(DM,Mat);
 65: extern PetscErrorCode MyComputeFunction(SNES,Vec,Vec,void*);
 66: extern PetscErrorCode MyComputeJacobian(SNES,Vec,Mat,Mat,void*);
 67: extern PetscErrorCode NonlinearGS(SNES,Vec);

 71: int main(int argc,char **argv)
 72: {
 73:   SNES           snes;                                 /* nonlinear solver */
 74:   SNES           psnes;                                /* nonlinear Gauss-Seidel approximate solver */
 75:   Vec            x,b;                                  /* solution vector */
 76:   PetscInt       its;                                  /* iterations for convergence */
 78:   DM             da;
 79:   PetscBool      use_ngs = PETSC_FALSE;                /* use the nonlinear Gauss-Seidel approximate solver */

 81:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 82:      Initialize program
 83:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 85:   PetscInitialize(&argc,&argv,(char*)0,help);

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:      Create nonlinear solver context
 89:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 90:   SNESCreate(PETSC_COMM_WORLD,&snes);

 92:   PetscOptionsGetBool(NULL,"-use_ngs",&use_ngs,0);

 94:   if (use_ngs) {
 95:     SNESGetNPC(snes,&psnes);
 96:     SNESSetType(psnes,SNESSHELL);
 97:     SNESShellSetSolve(psnes,NonlinearGS);
 98:   }

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Create distributed array (DMDA) to manage parallel grid and vectors
102:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
104:   DMDASetUniformCoordinates(da, 0.0, 1.0, 0.0, 1.0, 0.0, 1.0);
105:   SNESSetDM(snes,da);
106:   if (use_ngs) {
107:     SNESShellSetContext(psnes,da);
108:   }
109:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110:      Extract global vectors from DMDA; then duplicate for remaining
111:      vectors that are the same types
112:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113:   DMCreateGlobalVector(da,&x);
114:   DMCreateGlobalVector(da,&b);
115:   VecSetRandom(b,NULL);

117:   SNESSetFunction(snes,NULL,MyComputeFunction,NULL);
118:   SNESSetJacobian(snes,NULL,NULL,MyComputeJacobian,NULL);

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:      Customize nonlinear solver; set runtime options
122:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123:   SNESSetFromOptions(snes);

125:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126:      Solve nonlinear system
127:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128:   SNESSolve(snes,b,x);
129:   SNESGetIterationNumber(snes,&its);

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

134:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135:      Free work space.  All PETSc objects should be destroyed when they
136:      are no longer needed.
137:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138:   VecDestroy(&x);
139:   VecDestroy(&b);
140:   SNESDestroy(&snes);
141:   DMDestroy(&da);
142:   PetscFinalize();
143:   return(0);
144: }

146: /* ------------------------------------------------------------------- */
149: PetscErrorCode MyComputeFunction(SNES snes,Vec x,Vec F,void *ctx)
150: {
152:   Mat            J;
153:   DM             dm;

156:   SNESGetDM(snes,&dm);
157:   DMGetApplicationContext(dm,&J);
158:   if (!J) {
159:     DMSetMatType(dm,MATAIJ);
160:     DMCreateMatrix(dm,&J);
161:     MatSetDM(J, NULL);
162:     FormMatrix(dm,J);
163:     DMSetApplicationContext(dm,J);
164:     DMSetApplicationContextDestroy(dm,(PetscErrorCode (*)(void**))MatDestroy);
165:   }
166:   MatMult(J,x,F);
167:   return(0);
168: }

172: PetscErrorCode MyComputeJacobian(SNES snes,Vec x,Mat J,Mat Jp,void *ctx)
173: {
175:   DM             dm;

178:   SNESGetDM(snes,&dm);
179:   FormMatrix(dm,Jp);
180:   return(0);
181: }

185: PetscErrorCode FormMatrix(DM da,Mat jac)
186: {
188:   PetscInt       i,j,nrows = 0;
189:   MatStencil     col[5],row,*rows;
190:   PetscScalar    v[5],hx,hy,hxdhy,hydhx;
191:   DMDALocalInfo  info;

194:   DMDAGetLocalInfo(da,&info);
195:   hx    = 1.0/(PetscReal)(info.mx-1);
196:   hy    = 1.0/(PetscReal)(info.my-1);
197:   hxdhy = hx/hy;
198:   hydhx = hy/hx;

200:   PetscMalloc1(info.ym*info.xm,&rows);
201:   /*
202:      Compute entries for the locally owned part of the Jacobian.
203:       - Currently, all PETSc parallel matrix formats are partitioned by
204:         contiguous chunks of rows across the processors.
205:       - Each processor needs to insert only elements that it owns
206:         locally (but any non-local elements will be sent to the
207:         appropriate processor during matrix assembly).
208:       - Here, we set all entries for a particular row at once.
209:       - We can set matrix entries either using either
210:         MatSetValuesLocal() or MatSetValues(), as discussed above.
211:   */
212:   for (j=info.ys; j<info.ys+info.ym; j++) {
213:     for (i=info.xs; i<info.xs+info.xm; i++) {
214:       row.j = j; row.i = i;
215:       /* boundary points */
216:       if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
217:         v[0]            = 2.0*(hydhx + hxdhy);
218:         MatSetValuesStencil(jac,1,&row,1,&row,v,INSERT_VALUES);
219:         rows[nrows].i   = i;
220:         rows[nrows++].j = j;
221:       } else {
222:         /* interior grid points */
223:         v[0] = -hxdhy;                                           col[0].j = j - 1; col[0].i = i;
224:         v[1] = -hydhx;                                           col[1].j = j;     col[1].i = i-1;
225:         v[2] = 2.0*(hydhx + hxdhy);                              col[2].j = row.j; col[2].i = row.i;
226:         v[3] = -hydhx;                                           col[3].j = j;     col[3].i = i+1;
227:         v[4] = -hxdhy;                                           col[4].j = j + 1; col[4].i = i;
228:         MatSetValuesStencil(jac,1,&row,5,col,v,INSERT_VALUES);
229:       }
230:     }
231:   }

233:   /*
234:      Assemble matrix, using the 2-step process:
235:        MatAssemblyBegin(), MatAssemblyEnd().
236:   */
237:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
238:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
239:   MatZeroRowsColumnsStencil(jac,nrows,rows,2.0*(hydhx + hxdhy),NULL,NULL);
240:   PetscFree(rows);
241:   /*
242:      Tell the matrix we will never add a new nonzero location to the
243:      matrix. If we do, it will generate an error.
244:   */
245:   MatSetOption(jac,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
246:   return(0);
247: }



251: /* ------------------------------------------------------------------- */
254: /*
255:       Applies some sweeps on nonlinear Gauss-Seidel on each process

257:  */
258: PetscErrorCode NonlinearGS(SNES snes,Vec X)
259: {
260:   PetscInt       i,j,Mx,My,xs,ys,xm,ym,its,l;
262:   PetscReal      hx,hy,hxdhy,hydhx;
263:   PetscScalar    **x,F,J,u,uxx,uyy;
264:   DM             da;
265:   Vec            localX;

268:   SNESGetTolerances(snes,NULL,NULL,NULL,&its,NULL);
269:   SNESShellGetContext(snes,(void**)&da);

271:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
272:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

274:   hx    = 1.0/(PetscReal)(Mx-1);
275:   hy    = 1.0/(PetscReal)(My-1);
276:   hxdhy = hx/hy;
277:   hydhx = hy/hx;


280:   DMGetLocalVector(da,&localX);

282:   for (l=0; l<its; l++) {

284:     DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
285:     DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);
286:     /*
287:      Get a pointer to vector data.
288:      - For default PETSc vectors, VecGetArray() returns a pointer to
289:      the data array.  Otherwise, the routine is implementation dependent.
290:      - You MUST call VecRestoreArray() when you no longer need access to
291:      the array.
292:      */
293:     DMDAVecGetArray(da,localX,&x);

295:     /*
296:      Get local grid boundaries (for 2-dimensional DMDA):
297:      xs, ys   - starting grid indices (no ghost points)
298:      xm, ym   - widths of local grid (no ghost points)

300:      */
301:     DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

303:     for (j=ys; j<ys+ym; j++) {
304:       for (i=xs; i<xs+xm; i++) {
305:         if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
306:           /* boundary conditions are all zero Dirichlet */
307:           x[j][i] = 0.0;
308:         } else {
309:           u   = x[j][i];
310:           uxx = (2.0*u - x[j][i-1] - x[j][i+1])*hydhx;
311:           uyy = (2.0*u - x[j-1][i] - x[j+1][i])*hxdhy;
312:           F   = uxx + uyy;
313:           J   = 2.0*(hydhx + hxdhy);
314:           u   = u - F/J;

316:           x[j][i] = u;
317:         }
318:       }
319:     }

321:     /*
322:      Restore vector
323:      */
324:     DMDAVecRestoreArray(da,localX,&x);
325:     DMLocalToGlobalBegin(da,localX,INSERT_VALUES,X);
326:     DMLocalToGlobalEnd(da,localX,INSERT_VALUES,X);
327:   }
328:   DMRestoreLocalVector(da,&localX);
329:   return(0);
330: }