Actual source code: ex13.c

petsc-3.6.1 2015-08-06
Report Typos and Errors

  3: static char help[] = "Time-dependent PDE in 2d. Simplified from ex7.c for illustrating how to use TS on a structured domain. \n";
  4: /*
  5:    u_t = uxx + uyy
  6:    0 < x < 1, 0 < y < 1;
  7:    At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125
  8:            u(x,y) = 0.0           if r >= .125

 10:     mpiexec -n 2 ./ex13 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
 11:     mpiexec -n 1 ./ex13 -snes_fd_color -ts_monitor_draw_solution
 12:     mpiexec -n 2 ./ex13 -ts_type sundials -ts_monitor 
 13: */

 15: #include <petscdm.h>
 16: #include <petscdmda.h>
 17: #include <petscts.h>

 19: /*
 20:    User-defined data structures and routines
 21: */
 22: typedef struct {
 23:   PetscReal c;
 24: } AppCtx;

 26: extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
 27: extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
 28: extern PetscErrorCode FormInitialSolution(DM,Vec,void*);

 32: int main(int argc,char **argv)
 33: {
 34:   TS             ts;                   /* nonlinear solver */
 35:   Vec            u,r;                  /* solution, residual vector */
 36:   Mat            J;                    /* Jacobian matrix */
 37:   PetscInt       steps,maxsteps = 1000;     /* iterations for convergence */
 39:   DM             da;
 40:   PetscReal      ftime,dt;
 41:   AppCtx         user;              /* user-defined work context */

 43:   PetscInitialize(&argc,&argv,(char*)0,help);
 44:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 45:      Create distributed array (DMDA) to manage parallel grid and vectors
 46:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 47:   DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-8,-8,PETSC_DECIDE,PETSC_DECIDE,
 48:                       1,1,NULL,NULL,&da);

 50:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:      Extract global vectors from DMDA;
 52:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 53:   DMCreateGlobalVector(da,&u);
 54:   VecDuplicate(u,&r);

 56:   /* Initialize user application context */
 57:   user.c = -30.0;

 59:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 60:      Create timestepping solver context
 61:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 62:   TSCreate(PETSC_COMM_WORLD,&ts);
 63:   TSSetDM(ts,da);
 64:   TSSetType(ts,TSBEULER);
 65:   TSSetRHSFunction(ts,r,RHSFunction,&user);

 67:   /* Set Jacobian */
 68:   DMSetMatType(da,MATAIJ);
 69:   DMCreateMatrix(da,&J);
 70:   TSSetRHSJacobian(ts,J,J,RHSJacobian,NULL);

 72:   ftime = 1.0;
 73:   TSSetDuration(ts,maxsteps,ftime);

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:      Set initial conditions
 77:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 78:   FormInitialSolution(da,u,&user);
 79:   dt   = .01;
 80:   TSSetInitialTimeStep(ts,0.0,dt);

 82:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 83:      Set runtime options
 84:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 85:   TSSetFromOptions(ts);

 87:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 88:      Solve nonlinear system
 89:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 90:   TSSolve(ts,u);
 91:   TSGetSolveTime(ts,&ftime);
 92:   TSGetTimeStepNumber(ts,&steps);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:      Free work space.
 96:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   MatDestroy(&J);
 98:   VecDestroy(&u);
 99:   VecDestroy(&r);
100:   TSDestroy(&ts);
101:   DMDestroy(&da);

103:   PetscFinalize();
104:   return(0);
105: }
106: /* ------------------------------------------------------------------- */
109: /*
110:    RHSFunction - Evaluates nonlinear function, F(u).

112:    Input Parameters:
113: .  ts - the TS context
114: .  U - input vector
115: .  ptr - optional user-defined context, as set by TSSetFunction()

117:    Output Parameter:
118: .  F - function vector
119:  */
120: PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
121: {
122:   /* PETSC_UNUSED AppCtx *user=(AppCtx*)ptr; */
123:   DM             da;
125:   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
126:   PetscReal      two = 2.0,hx,hy,sx,sy;
127:   PetscScalar    u,uxx,uyy,**uarray,**f;
128:   Vec            localU;

131:   TSGetDM(ts,&da);
132:   DMGetLocalVector(da,&localU);
133:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
134:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

139:   /*
140:      Scatter ghost points to local vector,using the 2-step process
141:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
142:      By placing code between these two statements, computations can be
143:      done while messages are in transition.
144:   */
145:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
146:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);

148:   /* Get pointers to vector data */
149:   DMDAVecGetArrayRead(da,localU,&uarray);
150:   DMDAVecGetArray(da,F,&f);

152:   /* Get local grid boundaries */
153:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

155:   /* Compute function over the locally owned part of the grid */
156:   for (j=ys; j<ys+ym; j++) {
157:     for (i=xs; i<xs+xm; i++) {
158:       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
159:         f[j][i] = uarray[j][i];
160:         continue;
161:       }
162:       u       = uarray[j][i];
163:       uxx     = (-two*u + uarray[j][i-1] + uarray[j][i+1])*sx;
164:       uyy     = (-two*u + uarray[j-1][i] + uarray[j+1][i])*sy;
165:       f[j][i] = uxx + uyy;
166:     }
167:   }

169:   /* Restore vectors */
170:   DMDAVecRestoreArrayRead(da,localU,&uarray);
171:   DMDAVecRestoreArray(da,F,&f);
172:   DMRestoreLocalVector(da,&localU);
173:   PetscLogFlops(11.0*ym*xm);
174:   return(0);
175: }

177: /* --------------------------------------------------------------------- */
180: /*
181:    RHSJacobian - User-provided routine to compute the Jacobian of
182:    the nonlinear right-hand-side function of the ODE.

184:    Input Parameters:
185:    ts - the TS context
186:    t - current time
187:    U - global input vector
188:    dummy - optional user-defined context, as set by TSetRHSJacobian()

190:    Output Parameters:
191:    J - Jacobian matrix
192:    Jpre - optionally different preconditioning matrix
193:    str - flag indicating matrix structure
194: */
195: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat J,Mat Jpre,void *ctx)
196: {
198:   DM             da;
199:   DMDALocalInfo  info;
200:   PetscInt       i,j;
201:   PetscReal      hx,hy,sx,sy;

204:   TSGetDM(ts,&da);
205:   DMDAGetLocalInfo(da,&info);
206:   hx   = 1.0/(PetscReal)(info.mx-1); sx = 1.0/(hx*hx);
207:   hy   = 1.0/(PetscReal)(info.my-1); sy = 1.0/(hy*hy);
208:   for (j=info.ys; j<info.ys+info.ym; j++) {
209:     for (i=info.xs; i<info.xs+info.xm; i++) {
210:       PetscInt    nc = 0;
211:       MatStencil  row,col[5];
212:       PetscScalar val[5];
213:       row.i = i; row.j = j;
214:       if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
215:         col[nc].i = i; col[nc].j = j; val[nc++] = 1.0;
216:       } else {
217:         col[nc].i = i-1; col[nc].j = j;   val[nc++] = sx;
218:         col[nc].i = i+1; col[nc].j = j;   val[nc++] = sx;
219:         col[nc].i = i;   col[nc].j = j-1; val[nc++] = sy;
220:         col[nc].i = i;   col[nc].j = j+1; val[nc++] = sy;
221:         col[nc].i = i;   col[nc].j = j;   val[nc++] = -2*sx - 2*sy;
222:       }
223:       MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);
224:     }
225:   }
226:   MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
227:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
228:   if (J != Jpre) {
229:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
230:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
231:   }
232:   return(0);
233: }

235: /* ------------------------------------------------------------------- */
238: PetscErrorCode FormInitialSolution(DM da,Vec U,void* ptr)
239: {
240:   AppCtx         *user=(AppCtx*)ptr;
241:   PetscReal      c=user->c;
243:   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
244:   PetscScalar    **u;
245:   PetscReal      hx,hy,x,y,r;

248:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
249:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

251:   hx = 1.0/(PetscReal)(Mx-1);
252:   hy = 1.0/(PetscReal)(My-1);

254:   /* Get pointers to vector data */
255:   DMDAVecGetArray(da,U,&u);

257:   /* Get local grid boundaries */
258:   DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);

260:   /* Compute function over the locally owned part of the grid */
261:   for (j=ys; j<ys+ym; j++) {
262:     y = j*hy;
263:     for (i=xs; i<xs+xm; i++) {
264:       x = i*hx;
265:       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
266:       if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
267:       else u[j][i] = 0.0;
268:     }
269:   }

271:   /* Restore vectors */
272:   DMDAVecRestoreArray(da,U,&u);
273:   return(0);
274: }