Actual source code: ex17.c

petsc-3.10.5 2019-03-28
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 <petscdm.h>
 24:  #include <petscdmda.h>
 25:  #include <petscts.h>

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

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

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

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

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

 58:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 59:      Create distributed array (DMDA) to manage parallel grid and vectors
 60:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 61:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,1,1,NULL,&da);
 62:   DMSetFromOptions(da);
 63:   DMSetUp(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,NULL,"-boundary",&user.boundary,NULL);
 76:   PetscOptionsHasName(NULL,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:   TSSetMaxTime(ts,ftime);
 95:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);

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


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

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

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

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

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

165:   TSGetDM(ts,&da);
166:   DMGetLocalVector(da,&localU);
167:   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);

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

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

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

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

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

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

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

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

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

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

272:   TSGetDM(ts,&da);
273:   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);

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

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

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

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

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

296: /*TEST

298:     test:
299:       requires: !single
300:       args: -da_grid_x 40 -ts_max_steps 2 -snes_monitor_short -ksp_monitor_short -ts_monitor

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


308: TEST*/