Actual source code: ex3.c

petsc-3.7.2 2016-06-05
Report Typos and Errors
  2: static char help[] ="Model Equations for Advection-Diffusion\n";

  4: /*
  5:     Page 9, Section 1.2 Model Equations for Advection-Diffusion

  7:           u_t = a u_x + d u_xx

  9:    The initial conditions used here different then in the book.

 11: */

 13: /*
 14:      Helpful runtime linear solver options:
 15:            -pc_type mg -da_refine 2 -snes_monitor -ksp_monitor -ts_view   (geometric multigrid with three levels)

 17: */

 19: /*
 20:    Include "petscts.h" so that we can use TS solvers.  Note that this file
 21:    automatically includes:
 22:      petscsys.h       - base PETSc routines   petscvec.h  - vectors
 23:      petscmat.h  - matrices
 24:      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
 25:      petscviewer.h - viewers               petscpc.h   - preconditioners
 26:      petscksp.h   - linear solvers        petscsnes.h - nonlinear solvers
 27: */

 29: #include <petscts.h>
 30: #include <petscdm.h>
 31: #include <petscdmda.h>

 33: /*
 34:    User-defined application context - contains data needed by the
 35:    application-provided call-back routines.
 36: */
 37: typedef struct {
 38:   PetscScalar a,d;   /* advection and diffusion strength */
 39:   PetscBool   upwind;
 40: } AppCtx;

 42: /*
 43:    User-defined routines
 44: */
 45: extern PetscErrorCode InitialConditions(TS,Vec,AppCtx*);
 46: extern PetscErrorCode RHSMatrixHeat(TS,PetscReal,Vec,Mat,Mat,void*);
 47: extern PetscErrorCode Solution(TS,PetscReal,Vec,AppCtx*);

 51: int main(int argc,char **argv)
 52: {
 53:   AppCtx         appctx;                 /* user-defined application context */
 54:   TS             ts;                     /* timestepping context */
 55:   Vec            U;                      /* approximate solution vector */
 57:   PetscReal      dt;
 58:   DM             da;
 59:   PetscInt       M;

 61:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 62:      Initialize program and set problem parameters
 63:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 65:   PetscInitialize(&argc,&argv,(char*)0,help);
 66:   appctx.a      = 1.0;
 67:   appctx.d      = 0.0;
 68:   PetscOptionsGetScalar(NULL,NULL,"-a",&appctx.a,NULL);
 69:   PetscOptionsGetScalar(NULL,NULL,"-d",&appctx.d,NULL);
 70:   appctx.upwind = PETSC_TRUE;
 71:   PetscOptionsGetBool(NULL,NULL,"-upwind",&appctx.upwind,NULL);

 73:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_PERIODIC, -60, 1, 1,NULL,&da);
 74:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 75:      Create vector data structures
 76:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 78:   /*
 79:      Create vector data structures for approximate and exact solutions
 80:   */
 81:   DMCreateGlobalVector(da,&U);

 83:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 84:      Create timestepping solver context
 85:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 87:   TSCreate(PETSC_COMM_WORLD,&ts);
 88:   TSSetDM(ts,da);

 90:   /*
 91:       For linear problems with a time-dependent f(U,t) in the equation
 92:      u_t = f(u,t), the user provides the discretized right-hand-side
 93:       as a time-dependent matrix.
 94:   */
 95:   TSSetRHSFunction(ts,NULL,TSComputeRHSFunctionLinear,&appctx);
 96:   TSSetRHSJacobian(ts,NULL,NULL,RHSMatrixHeat,&appctx);
 97:   TSSetSolutionFunction(ts,(PetscErrorCode (*)(TS,PetscReal,Vec,void*))Solution,&appctx);

 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100:      Customize timestepping solver:
101:        - Set timestepping duration info
102:      Then set runtime options, which can override these defaults.
103:      For example,
104:           -ts_max_steps <maxsteps> -ts_final_time <maxtime>
105:      to override the defaults set by TSSetDuration().
106:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

108:   DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);
109:   dt   = .48/(M*M);
110:   TSSetInitialTimeStep(ts,0.0,dt);
111:   TSSetDuration(ts,1000,100.0);
112:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
113:   TSSetType(ts,TSARKIMEX);
114:   TSSetFromOptions(ts);

116:   /*
117:      Evaluate initial conditions
118:   */
119:   InitialConditions(ts,U,&appctx);

121:   /*
122:      Run the timestepping solver
123:   */
124:   TSSolve(ts,U);


127:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
128:      Free work space.  All PETSc objects should be destroyed when they
129:      are no longer needed.
130:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

132:   TSDestroy(&ts);
133:   VecDestroy(&U);
134:   DMDestroy(&da);

136:   /*
137:      Always call PetscFinalize() before exiting a program.  This routine
138:        - finalizes the PETSc libraries as well as MPI
139:        - provides summary and diagnostic information if certain runtime
140:          options are chosen (e.g., -log_summary).
141:   */
142:   PetscFinalize();
143:   return 0;
144: }
145: /* --------------------------------------------------------------------- */
148: /*
149:    InitialConditions - Computes the solution at the initial time.

151:    Input Parameter:
152:    u - uninitialized solution vector (global)
153:    appctx - user-defined application context

155:    Output Parameter:
156:    u - vector with solution at initial time (global)
157: */
158: PetscErrorCode InitialConditions(TS ts,Vec U,AppCtx *appctx)
159: {
160:   PetscScalar    *u,h;
162:   PetscInt       i,mstart,mend,xm,M;
163:   DM             da;

165:   TSGetDM(ts,&da);
166:   DMDAGetCorners(da,&mstart,0,0,&xm,0,0);
167:   DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);
168:   h    = 1.0/M;
169:   mend = mstart + xm;
170:   /*
171:     Get a pointer to vector data.
172:     - For default PETSc vectors, VecGetArray() returns a pointer to
173:       the data array.  Otherwise, the routine is implementation dependent.
174:     - You MUST call VecRestoreArray() when you no longer need access to
175:       the array.
176:     - Note that the Fortran interface to VecGetArray() differs from the
177:       C version.  See the users manual for details.
178:   */
179:   DMDAVecGetArray(da,U,&u);

181:   /*
182:      We initialize the solution array by simply writing the solution
183:      directly into the array locations.  Alternatively, we could use
184:      VecSetValues() or VecSetValuesLocal().
185:   */
186:   for (i=mstart; i<mend; i++) u[i] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);

188:   /*
189:      Restore vector
190:   */
191:   DMDAVecRestoreArray(da,U,&u);
192:   return 0;
193: }
194: /* --------------------------------------------------------------------- */
197: /*
198:    Solution - Computes the exact solution at a given time.

200:    Input Parameters:
201:    t - current time
202:    solution - vector in which exact solution will be computed
203:    appctx - user-defined application context

205:    Output Parameter:
206:    solution - vector with the newly computed exact solution
207: */
208: PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *appctx)
209: {
210:   PetscScalar    *u,ex1,ex2,sc1,sc2,h;
212:   PetscInt       i,mstart,mend,xm,M;
213:   DM             da;

215:   TSGetDM(ts,&da);
216:   DMDAGetCorners(da,&mstart,0,0,&xm,0,0);
217:   DMDAGetInfo(da,PETSC_IGNORE,&M,0,0,0,0,0,0,0,0,0,0,0);
218:   h    = 1.0/M;
219:   mend = mstart + xm;
220:   /*
221:      Get a pointer to vector data.
222:   */
223:   DMDAVecGetArray(da,U,&u);

225:   /*
226:      Simply write the solution directly into the array locations.
227:      Alternatively, we culd use VecSetValues() or VecSetValuesLocal().
228:   */
229:   ex1 = PetscExpScalar(-36.*PETSC_PI*PETSC_PI*appctx->d*t);
230:   ex2 = PetscExpScalar(-4.*PETSC_PI*PETSC_PI*appctx->d*t);
231:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
232:   for (i=mstart; i<mend; i++) u[i] = PetscSinScalar(sc1*(PetscReal)i + appctx->a*PETSC_PI*6.*t)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i + appctx->a*PETSC_PI*2.*t)*ex2;

234:   /*
235:      Restore vector
236:   */
237:   DMDAVecRestoreArray(da,U,&u);
238:   return 0;
239: }

241: /* --------------------------------------------------------------------- */
244: /*
245:    RHSMatrixHeat - User-provided routine to compute the right-hand-side
246:    matrix for the heat equation.

248:    Input Parameters:
249:    ts - the TS context
250:    t - current time
251:    global_in - global input vector
252:    dummy - optional user-defined context, as set by TSetRHSJacobian()

254:    Output Parameters:
255:    AA - Jacobian matrix
256:    BB - optionally different preconditioning matrix
257:    str - flag indicating matrix structure

259:    Notes:
260:    Recall that MatSetValues() uses 0-based row and column numbers
261:    in Fortran as well as in C.
262: */
263: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Vec U,Mat AA,Mat BB,void *ctx)
264: {
265:   Mat            A       = AA;                /* Jacobian matrix */
266:   AppCtx         *appctx = (AppCtx*)ctx;     /* user-defined application context */
267:   PetscInt       mstart, mend;
269:   PetscInt       i,idx[3],M,xm;
270:   PetscScalar    v[3],h;
271:   DM             da;

273:   TSGetDM(ts,&da);
274:   DMDAGetInfo(da,0,&M,0,0,0,0,0,0,0,0,0,0,0);
275:   DMDAGetCorners(da,&mstart,0,0,&xm,0,0);
276:   h    = 1.0/M;
277:   mend = mstart + xm;
278:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
279:      Compute entries for the locally owned part of the matrix
280:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
281:   /*
282:      Set matrix rows corresponding to boundary data
283:   */

285:   /* diffusion */
286:   v[0] = appctx->d/(h*h);
287:   v[1] = -2.0*appctx->d/(h*h);
288:   v[2] = appctx->d/(h*h);
289:   if (!mstart) {
290:     idx[0] = M-1; idx[1] = 0; idx[2] = 1;
291:     MatSetValues(A,1,&mstart,3,idx,v,INSERT_VALUES);
292:     mstart++;
293:   }

295:   if (mend == M) {
296:     mend--;
297:     idx[0] = M-2; idx[1] = M-1; idx[2] = 0;
298:     MatSetValues(A,1,&mend,3,idx,v,INSERT_VALUES);
299:   }

301:   /*
302:      Set matrix rows corresponding to interior data.  We construct the
303:      matrix one row at a time.
304:   */
305:   for (i=mstart; i<mend; i++) {
306:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
307:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
308:   }
309:   MatAssemblyBegin(A,MAT_FLUSH_ASSEMBLY);
310:   MatAssemblyEnd(A,MAT_FLUSH_ASSEMBLY);

312:   DMDAGetCorners(da,&mstart,0,0,&xm,0,0);
313:   mend = mstart + xm;
314:   if (!appctx->upwind) {
315:     /* advection -- centered differencing */
316:     v[0] = -.5*appctx->a/(h);
317:     v[1] = .5*appctx->a/(h);
318:     if (!mstart) {
319:       idx[0] = M-1; idx[1] = 1;
320:       MatSetValues(A,1,&mstart,2,idx,v,ADD_VALUES);
321:       mstart++;
322:     }

324:     if (mend == M) {
325:       mend--;
326:       idx[0] = M-2; idx[1] = 0;
327:       MatSetValues(A,1,&mend,2,idx,v,ADD_VALUES);
328:     }

330:     for (i=mstart; i<mend; i++) {
331:       idx[0] = i-1; idx[1] = i+1;
332:       MatSetValues(A,1,&i,2,idx,v,ADD_VALUES);
333:     }
334:   } else {
335:     /* advection -- upwinding */
336:     v[0] = -appctx->a/(h);
337:     v[1] = appctx->a/(h);
338:     if (!mstart) {
339:       idx[0] = 0; idx[1] = 1;
340:       MatSetValues(A,1,&mstart,2,idx,v,ADD_VALUES);
341:       mstart++;
342:     }

344:     if (mend == M) {
345:       mend--;
346:       idx[0] = M-1; idx[1] = 0;
347:       MatSetValues(A,1,&mend,2,idx,v,ADD_VALUES);
348:     }

350:     for (i=mstart; i<mend; i++) {
351:       idx[0] = i; idx[1] = i+1;
352:       MatSetValues(A,1,&i,2,idx,v,ADD_VALUES);
353:     }
354:   }


357:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
358:      Complete the matrix assembly process and set some options
359:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
360:   /*
361:      Assemble matrix, using the 2-step process:
362:        MatAssemblyBegin(), MatAssemblyEnd()
363:      Computations can be done while messages are in transition
364:      by placing code between these two statements.
365:   */
366:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
367:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);

369:   /*
370:      Set and option to indicate that we will never add a new nonzero location
371:      to the matrix. If we do, it will generate an error.
372:   */
373:   MatSetOption(A,MAT_NEW_NONZERO_LOCATION_ERR,PETSC_TRUE);
374:   return 0;
375: }