Actual source code: ex13.c


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

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

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

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

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

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

 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;
 50:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 51:   DMCreateGlobalVector(da,&u);
 52:   VecDuplicate(u,&r);

 54:   /* Initialize user application context */
 55:   user.c = -30.0;

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

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

 70:   ftime = 1.0;
 71:   TSSetMaxTime(ts,ftime);
 72:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

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

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

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

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

102:   PetscFinalize();
103:   return ierr;
104: }
105: /* ------------------------------------------------------------------- */
106: /*
107:    RHSFunction - Evaluates nonlinear function, F(u).

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

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

128:   TSGetDM(ts,&da);
129:   DMGetLocalVector(da,&localU);
130:   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);

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

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

144:   /* Get pointers to vector data */
145:   DMDAVecGetArrayRead(da,localU,&uarray);
146:   DMDAVecGetArray(da,F,&f);

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

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

165:   /* Restore vectors */
166:   DMDAVecRestoreArrayRead(da,localU,&uarray);
167:   DMDAVecRestoreArray(da,F,&f);
168:   DMRestoreLocalVector(da,&localU);
169:   PetscLogFlops(11.0*ym*xm);
170:   return(0);
171: }

173: /* --------------------------------------------------------------------- */
174: /*
175:    RHSJacobian - User-provided routine to compute the Jacobian of
176:    the nonlinear right-hand-side function of the ODE.

178:    Input Parameters:
179:    ts - the TS context
180:    t - current time
181:    U - global input vector
182:    dummy - optional user-defined context, as set by TSetRHSJacobian()

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

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

229: /* ------------------------------------------------------------------- */
230: PetscErrorCode FormInitialSolution(DM da,Vec U,void* ptr)
231: {
232:   AppCtx         *user=(AppCtx*)ptr;
233:   PetscReal      c=user->c;
235:   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
236:   PetscScalar    **u;
237:   PetscReal      hx,hy,x,y,r;

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

242:   hx = 1.0/(PetscReal)(Mx-1);
243:   hy = 1.0/(PetscReal)(My-1);

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

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

251:   /* Compute function over the locally owned part of the grid */
252:   for (j=ys; j<ys+ym; j++) {
253:     y = j*hy;
254:     for (i=xs; i<xs+xm; i++) {
255:       x = i*hx;
256:       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
257:       if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
258:       else u[j][i] = 0.0;
259:     }
260:   }

262:   /* Restore vectors */
263:   DMDAVecRestoreArray(da,U,&u);
264:   return(0);
265: }

267: /*TEST

269:     test:
270:       args: -ts_max_steps 5 -ts_monitor

272:     test:
273:       suffix: 2
274:       args: -ts_max_steps 5 -ts_monitor

276:     test:
277:       suffix: 3
278:       args: -ts_max_steps 5 -snes_fd_color -ts_monitor

280: TEST*/