Actual source code: ex17.c

petsc-dev 2014-02-02
Report Typos and Errors
  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


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

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

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

 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*,MatStructure*,void*);
 40: static PetscErrorCode FormInitialSolution(TS,Vec,void*);

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

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

 60:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 61:      Create distributed array (DMDA) to manage parallel grid and vectors
 62:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 63:   DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-11,1,1,NULL,&da);

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

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

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

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

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

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

 93:   ftime = 1.0;
 94:   TSSetDuration(ts,maxsteps,ftime);

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


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

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

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

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:    Compute diagnostics of the solution
134:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   VecNorm(u,NORM_1,&norm);
136:   VecMax(u,NULL,&vmax);
137:   VecMin(u,NULL,&vmin);
138:   TSGetTimeStepNumber(ts,&nsteps);
139:   TSGetTime(ts,&ftime);
140:   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);

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

166:   TSGetDM(ts,&da);
167:   DMGetLocalVector(da,&localU);
168:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
169:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

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

182:   /* Get pointers to vector data */
183:   DMDAVecGetArray(da,localU,&u);
184:   DMDAVecGetArray(da,Udot,&udot);
185:   DMDAVecGetArray(da,F,&f);

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

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

202:   /* Restore vectors */
203:   DMDAVecRestoreArray(da,localU,&u);
204:   DMDAVecRestoreArray(da,Udot,&udot);
205:   DMDAVecRestoreArray(da,F,&f);
206:   DMRestoreLocalVector(da,&localU);
207:   return(0);
208: }

210: /* --------------------------------------------------------------------- */
211: /*
212:   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
213: */
216: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ctx)
217: {
219:   PetscInt       i,rstart,rend,Mx;
220:   PetscReal      hx,sx;
221:   AppCtx         *user = (AppCtx*)ctx;
222:   DM             da;
223:   MatStencil     col[3],row;
224:   PetscInt       nc;
225:   PetscScalar    vals[3];

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

252:   MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);
253:   MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);
254:   if (*J != *Jpre) {
255:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
256:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
257:   }
258:   if (user->viewJacobian) {
259:     PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n");
260:     MatView(*Jpre,PETSC_VIEWER_STDOUT_WORLD);
261:   }
262:   return(0);
263: }

265: /* ------------------------------------------------------------------- */
268: PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr)
269: {
270:   AppCtx         *user=(AppCtx*)ptr;
271:   PetscReal      c    =user->c;
272:   DM             da;
274:   PetscInt       i,xs,xm,Mx;
275:   PetscScalar    *u;
276:   PetscReal      hx,x,r;

279:   TSGetDM(ts,&da);
280:   DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
281:                      PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);

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

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

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

291:   /* Compute function over the locally owned part of the grid */
292:   for (i=xs; i<xs+xm; i++) {
293:     x = i*hx;
294:     r = PetscSqrtReal((x-.5)*(x-.5));
295:     if (r < .125) u[i] = PetscExpReal(c*r*r*r);
296:     else          u[i] = 0.0;
297:   }

299:   /* Restore vectors */
300:   DMDAVecRestoreArray(da,U,&u);
301:   return(0);
302: }