Actual source code: ex7.c

petsc-master 2017-05-26
Report Typos and Errors

  2: static char help[] = "Nonlinear, time-dependent PDE in 2d.\n";


  5: /*
  6:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
  7:    Include "petscts.h" so that we can use SNES solvers.  Note that this
  8:    file automatically includes:
  9:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 10:      petscmat.h - matrices
 11:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 12:      petscviewer.h - viewers               petscpc.h  - preconditioners
 13:      petscksp.h   - linear solvers
 14: */
 15:  #include <petscdm.h>
 16:  #include <petscdmda.h>
 17:  #include <petscts.h>


 20: /*
 21:    User-defined routines
 22: */
 23: extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec);
 24: extern PetscErrorCode MyTSMonitor(TS,PetscInt,PetscReal,Vec,void*);
 25: extern PetscErrorCode MySNESMonitor(SNES,PetscInt,PetscReal,PetscViewerAndFormat*);

 27: int main(int argc,char **argv)
 28: {
 29:   TS                   ts;                         /* time integrator */
 30:   SNES                 snes;
 31:   Vec                  x,r;                        /* solution, residual vectors */
 32:   PetscInt             maxsteps = 100;             /* iterations for convergence */
 33:   PetscErrorCode       ierr;
 34:   DM                   da;
 35:   PetscViewerAndFormat *vf;

 37:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 38:      Initialize program
 39:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 40:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 41:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 42:      Create distributed array (DMDA) to manage parallel grid and vectors
 43:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 44:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
 45:   DMSetFromOptions(da);
 46:   DMSetUp(da);

 48:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 49:      Extract global vectors from DMDA; then duplicate for remaining
 50:      vectors that are the same types
 51:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 52:   DMCreateGlobalVector(da,&x);
 53:   VecDuplicate(x,&r);

 55:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 56:      Create timestepping solver context
 57:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 58:   TSCreate(PETSC_COMM_WORLD,&ts);
 59:   TSSetProblemType(ts,TS_NONLINEAR);
 60:   TSSetRHSFunction(ts,NULL,FormFunction,da);

 62:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:      Create matrix data structure; set Jacobian evaluation routine

 65:      Set Jacobian matrix data structure and default Jacobian evaluation
 66:      routine. User can override with:
 67:      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
 68:                 (unless user explicitly sets preconditioner)
 69:      -snes_mf_operator : form preconditioning matrix as set by the user,
 70:                          but use matrix-free approx for Jacobian-vector
 71:                          products within Newton-Krylov method

 73:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 75:   TSSetDuration(ts,maxsteps,1.0);
 76:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
 77:   TSMonitorSet(ts,MyTSMonitor,PETSC_VIEWER_STDOUT_WORLD,NULL);
 78:   TSSetDM(ts,da);
 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:      Customize nonlinear solver
 81:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 82:   TSSetType(ts,TSBEULER);
 83:   TSGetSNES(ts,&snes);
 84:   PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf);
 85:   SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))MySNESMonitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:      Set initial conditions
 89:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 90:   FormInitialSolution(da,x);
 91:   TSSetInitialTimeStep(ts,0.0,.0001);
 92:   TSSetSolution(ts,x);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:      Set runtime options
 96:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   TSSetFromOptions(ts);

 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100:      Solve nonlinear system
101:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102:   TSSolve(ts,x);

104:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
105:      Free work space.  All PETSc objects should be destroyed when they
106:      are no longer needed.
107:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
108:   VecDestroy(&x);
109:   VecDestroy(&r);
110:   TSDestroy(&ts);
111:   DMDestroy(&da);

113:   PetscFinalize();
114:   return ierr;
115: }
116: /* ------------------------------------------------------------------- */
117: /*
118:    FormFunction - Evaluates nonlinear function, F(x).

120:    Input Parameters:
121: .  ts - the TS context
122: .  X - input vector
123: .  ptr - optional user-defined context, as set by SNESSetFunction()

125:    Output Parameter:
126: .  F - function vector
127:  */
128: PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
129: {
130:   DM             da;
132:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
133:   PetscReal      two = 2.0,hx,hy,sx,sy;
134:   PetscScalar    u,uxx,uyy,**x,**f;
135:   Vec            localX;

138:   TSGetDM(ts,&da);
139:   DMGetLocalVector(da,&localX);
140:   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);

142:   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
143:   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);

145:   /*
146:      Scatter ghost points to local vector,using the 2-step process
147:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
148:      By placing code between these two statements, computations can be
149:      done while messages are in transition.
150:   */
151:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX);
152:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX);

154:   /*
155:      Get pointers to vector data
156:   */
157:   DMDAVecGetArrayRead(da,localX,&x);
158:   DMDAVecGetArray(da,F,&f);

160:   /*
161:      Get local grid boundaries
162:   */
163:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

165:   /*
166:      Compute function over the locally owned part of the grid
167:   */
168:   for (j=ys; j<ys+ym; j++) {
169:     for (i=xs; i<xs+xm; i++) {
170:       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
171:         f[j][i] = x[j][i];
172:         continue;
173:       }
174:       u   = x[j][i];
175:       uxx = (two*u - x[j][i-1] - x[j][i+1])*sx;
176:       uyy = (two*u - x[j-1][i] - x[j+1][i])*sy;
177:       /*      f[j][i] = -(uxx + uyy); */
178:       f[j][i] = -u*(uxx + uyy) - (4.0 - 1.0)*((x[j][i+1] - x[j][i-1])*(x[j][i+1] - x[j][i-1])*.25*sx +
179:                                               (x[j+1][i] - x[j-1][i])*(x[j+1][i] - x[j-1][i])*.25*sy);
180:     }
181:   }

183:   /*
184:      Restore vectors
185:   */
186:   DMDAVecRestoreArrayRead(da,localX,&x);
187:   DMDAVecRestoreArray(da,F,&f);
188:   DMRestoreLocalVector(da,&localX);
189:   PetscLogFlops(11.0*ym*xm);
190:   return(0);
191: }

193: /* ------------------------------------------------------------------- */
194: PetscErrorCode FormInitialSolution(DM da,Vec U)
195: {
197:   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
198:   PetscScalar    **u;
199:   PetscReal      hx,hy,x,y,r;

202:   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);

204:   hx = 1.0/(PetscReal)(Mx-1);
205:   hy = 1.0/(PetscReal)(My-1);

207:   /*
208:      Get pointers to vector data
209:   */
210:   DMDAVecGetArray(da,U,&u);

212:   /*
213:      Get local grid boundaries
214:   */
215:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

217:   /*
218:      Compute function over the locally owned part of the grid
219:   */
220:   for (j=ys; j<ys+ym; j++) {
221:     y = j*hy;
222:     for (i=xs; i<xs+xm; i++) {
223:       x = i*hx;
224:       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
225:       if (r < .125) u[j][i] = PetscExpReal(-30.0*r*r*r);
226:       else          u[j][i] = 0.0;
227:     }
228:   }

230:   /*
231:      Restore vectors
232:   */
233:   DMDAVecRestoreArray(da,U,&u);
234:   return(0);
235: }

237: PetscErrorCode MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
238: {
240:   PetscReal      norm;
241:   MPI_Comm       comm;

244:   if (step < 0) return(0); /* step of -1 indicates an interpolated solution */
245:   VecNorm(v,NORM_2,&norm);
246:   PetscObjectGetComm((PetscObject)ts,&comm);
247:   PetscPrintf(comm,"timestep %D time %g norm %g\n",step,(double)ptime,(double)norm);
248:   return(0);
249: }

251: /*
252:    MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
253:    Input Parameters:
254:      snes - the SNES context
255:      its - iteration number
256:      fnorm - 2-norm function value (may be estimated)
257:      ctx - optional user-defined context for private data for the
258:          monitor routine, as set by SNESMonitorSet()
259:  */
260: PetscErrorCode MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *vf)
261: {

265:   SNESMonitorDefaultShort(snes,its,fnorm,vf);
266:   return(0);
267: }