Actual source code: ex17.c

  1: static const char help[] = "Time-dependent PDE in 1d. Simplified from ex15.c for illustrating how to solve DAEs. \n";
  2: /*
  3:    u_t = uxx
  4:    0 < x < 1;
  5:    At t=0: u(x) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5)) < .125
  6:            u(x) = 0.0           if r >= .125

  8:    Boundary conditions:
  9:    Dirichlet BC:
 10:    At x=0, x=1, u = 0.0

 12:    Neumann BC:
 13:    At x=0, x=1: du(x,t)/dx = 0

 15:    mpiexec -n 2 ./ex17 -da_grid_x 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
 16:          ./ex17 -da_grid_x 40 -monitor_solution
 17:          ./ex17 -da_grid_x 100  -ts_type theta -ts_theta_theta 0.5 # Midpoint is not L-stable
 18:          ./ex17 -jac_type fd_coloring  -da_grid_x 500 -boundary 1
 19:          ./ex17 -da_grid_x 100  -ts_type gl -ts_adapt_type none -ts_max_steps 2
 20: */

 22: #include <petscdm.h>
 23: #include <petscdmda.h>
 24: #include <petscts.h>

 26: typedef enum {JACOBIAN_ANALYTIC,JACOBIAN_FD_COLORING,JACOBIAN_FD_FULL} JacobianType;
 27: static const char *const JacobianTypes[] = {"analytic","fd_coloring","fd_full","JacobianType","fd_",0};

 29: /*
 30:    User-defined data structures and routines
 31: */
 32: typedef struct {
 33:   PetscReal c;
 34:   PetscInt  boundary;            /* Type of boundary condition */
 35:   PetscBool viewJacobian;
 36: } AppCtx;

 38: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 39: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
 40: static PetscErrorCode FormInitialSolution(TS,Vec,void*);

 42: int main(int argc,char **argv)
 43: {
 44:   TS             ts;                   /* nonlinear solver */
 45:   Vec            u;                    /* solution, residual vectors */
 46:   Mat            J;                    /* Jacobian matrix */
 47:   PetscInt       nsteps;
 48:   PetscReal      vmin,vmax,norm;
 50:   DM             da;
 51:   PetscReal      ftime,dt;
 52:   AppCtx         user;              /* user-defined work context */
 53:   JacobianType   jacType;

 55:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:      Create distributed array (DMDA) to manage parallel grid and vectors
 59:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,1,1,NULL,&da);
 61:   DMSetFromOptions(da);
 62:   DMSetUp(da);

 64:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 65:      Extract global vectors from DMDA;
 66:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 67:   DMCreateGlobalVector(da,&u);

 69:   /* Initialize user application context */
 70:   user.c            = -30.0;
 71:   user.boundary     = 0;  /* 0: Dirichlet BC; 1: Neumann BC */
 72:   user.viewJacobian = PETSC_FALSE;

 74:   PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL);
 75:   PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian);

 77:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 78:      Create timestepping solver context
 79:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 80:   TSCreate(PETSC_COMM_WORLD,&ts);
 81:   TSSetProblemType(ts,TS_NONLINEAR);
 82:   TSSetType(ts,TSTHETA);
 83:   TSThetaSetTheta(ts,1.0); /* Make the Theta method behave like backward Euler */
 84:   TSSetIFunction(ts,NULL,FormIFunction,&user);

 86:   DMSetMatType(da,MATAIJ);
 87:   DMCreateMatrix(da,&J);
 88:   jacType = JACOBIAN_ANALYTIC; /* use user-provide Jacobian */

 90:   TSSetDM(ts,da); /* Use TSGetDM() to access. Setting here allows easy use of geometric multigrid. */

 92:   ftime = 1.0;
 93:   TSSetMaxTime(ts,ftime);
 94:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

 96:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 97:      Set initial conditions
 98:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 99:   FormInitialSolution(ts,u,&user);
100:   TSSetSolution(ts,u);
101:   dt   = .01;
102:   TSSetTimeStep(ts,dt);

104:   /* Use slow fd Jacobian or fast fd Jacobian with colorings.
105:      Note: this requirs snes which is not created until TSSetUp()/TSSetFromOptions() is called */
106:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for Jacobian evaluation",NULL);
107:   PetscOptionsEnum("-jac_type","Type of Jacobian","",JacobianTypes,(PetscEnum)jacType,(PetscEnum*)&jacType,0);
108:   PetscOptionsEnd();
109:   if (jacType == JACOBIAN_ANALYTIC) {
110:     TSSetIJacobian(ts,J,J,FormIJacobian,&user);
111:   } else if (jacType == JACOBIAN_FD_COLORING) {
112:     SNES snes;
113:     TSGetSNES(ts,&snes);
114:     SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0);
115:   } else if (jacType == JACOBIAN_FD_FULL) {
116:     SNES snes;
117:     TSGetSNES(ts,&snes);
118:     SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,&user);
119:   }

121:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122:      Set runtime options
123:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124:   TSSetFromOptions(ts);

126:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127:      Integrate ODE system
128:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129:   TSSolve(ts,u);

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:    Compute diagnostics of the solution
133:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134:   VecNorm(u,NORM_1,&norm);
135:   VecMax(u,NULL,&vmax);
136:   VecMin(u,NULL,&vmin);
137:   TSGetStepNumber(ts,&nsteps);
138:   TSGetTime(ts,&ftime);
139:   PetscPrintf(PETSC_COMM_WORLD,"timestep %D: time %g, solution norm %g, max %g, min %g\n",nsteps,(double)ftime,(double)norm,(double)vmax,(double)vmin);

141:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142:      Free work space.
143:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144:   MatDestroy(&J);
145:   VecDestroy(&u);
146:   TSDestroy(&ts);
147:   DMDestroy(&da);
148:   PetscFinalize();
149:   return ierr;
150: }
151: /* ------------------------------------------------------------------- */
152: static PetscErrorCode FormIFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
153: {
154:   AppCtx         *user=(AppCtx*)ptr;
155:   DM             da;
157:   PetscInt       i,Mx,xs,xm;
158:   PetscReal      hx,sx;
159:   PetscScalar    *u,*udot,*f;
160:   Vec            localU;

163:   TSGetDM(ts,&da);
164:   DMGetLocalVector(da,&localU);
165:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

167:   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);

169:   /*
170:      Scatter ghost points to local vector,using the 2-step process
171:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
172:      By placing code between these two statements, computations can be
173:      done while messages are in transition.
174:   */
175:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
176:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);

178:   /* Get pointers to vector data */
179:   DMDAVecGetArrayRead(da,localU,&u);
180:   DMDAVecGetArrayRead(da,Udot,&udot);
181:   DMDAVecGetArray(da,F,&f);

183:   /* Get local grid boundaries */
184:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

186:   /* Compute function over the locally owned part of the grid */
187:   for (i=xs; i<xs+xm; i++) {
188:     if (user->boundary == 0) { /* Dirichlet BC */
189:       if (i == 0 || i == Mx-1) f[i] = u[i]; /* F = U */
190:       else                     f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
191:     } else { /* Neumann BC */
192:       if (i == 0)         f[i] = u[0] - u[1];
193:       else if (i == Mx-1) f[i] = u[i] - u[i-1];
194:       else                f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
195:     }
196:   }

198:   /* Restore vectors */
199:   DMDAVecRestoreArrayRead(da,localU,&u);
200:   DMDAVecRestoreArrayRead(da,Udot,&udot);
201:   DMDAVecRestoreArray(da,F,&f);
202:   DMRestoreLocalVector(da,&localU);
203:   return(0);
204: }

206: /* --------------------------------------------------------------------- */
207: /*
208:   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
209: */
210: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
211: {
213:   PetscInt       i,rstart,rend,Mx;
214:   PetscReal      hx,sx;
215:   AppCtx         *user = (AppCtx*)ctx;
216:   DM             da;
217:   MatStencil     col[3],row;
218:   PetscInt       nc;
219:   PetscScalar    vals[3];

222:   TSGetDM(ts,&da);
223:   MatGetOwnershipRange(Jpre,&rstart,&rend);
224:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
225:   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
226:   for (i=rstart; i<rend; i++) {
227:     nc    = 0;
228:     row.i = i;
229:     if (user->boundary == 0 && (i == 0 || i == Mx-1)) {
230:       col[nc].i = i; vals[nc++] = 1.0;
231:     } else if (user->boundary > 0 && i == 0) { /* Left Neumann */
232:       col[nc].i = i;   vals[nc++] = 1.0;
233:       col[nc].i = i+1; vals[nc++] = -1.0;
234:     } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
235:       col[nc].i = i-1; vals[nc++] = -1.0;
236:       col[nc].i = i;   vals[nc++] = 1.0;
237:     } else {                    /* Interior */
238:       col[nc].i = i-1; vals[nc++] = -1.0*sx;
239:       col[nc].i = i;   vals[nc++] = 2.0*sx + a;
240:       col[nc].i = i+1; vals[nc++] = -1.0*sx;
241:     }
242:     MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES);
243:   }

245:   MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
246:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
247:   if (J != Jpre) {
248:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
249:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
250:   }
251:   if (user->viewJacobian) {
252:     PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n");
253:     MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);
254:   }
255:   return(0);
256: }

258: /* ------------------------------------------------------------------- */
259: PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr)
260: {
261:   AppCtx         *user=(AppCtx*)ptr;
262:   PetscReal      c    =user->c;
263:   DM             da;
265:   PetscInt       i,xs,xm,Mx;
266:   PetscScalar    *u;
267:   PetscReal      hx,x,r;

270:   TSGetDM(ts,&da);
271:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

273:   hx = 1.0/(PetscReal)(Mx-1);

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

278:   /* Get local grid boundaries */
279:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

281:   /* Compute function over the locally owned part of the grid */
282:   for (i=xs; i<xs+xm; i++) {
283:     x = i*hx;
284:     r = PetscSqrtReal((x-.5)*(x-.5));
285:     if (r < .125) u[i] = PetscExpReal(c*r*r*r);
286:     else          u[i] = 0.0;
287:   }

289:   /* Restore vectors */
290:   DMDAVecRestoreArray(da,U,&u);
291:   return(0);
292: }

294: /*TEST

296:     test:
297:       requires: !single
298:       args: -da_grid_x 40 -ts_max_steps 2 -snes_monitor_short -ksp_monitor_short -ts_monitor

300:     test:
301:       suffix: 2
302:       requires: !single
303:       args: -da_grid_x 100 -ts_type theta -ts_theta_theta 0.5

305: TEST*/