Actual source code: ex1.c

  1: /*
  2:        Formatted test for TS routines.

  4:           Solves U_t = U_xx 
  5:      F(t,u) = (u_i+1 - 2u_i + u_i-1)/h^2
  6:        using several different schemes. 
  7: */

  9: /* Usage: 
 10:    ./ex1 -nox -ts_type beuler -ts_view 
 11:    ./ex1 -nox -linear_constant_matrix -ts_type beuler -pc_type lu
 12:    ./ex1 -nox -linear_variable_matrix -ts_type beuler
 13: */

 15: static char help[] = "Solves 1D heat equation.\n\n";

 17: #include <petscdmda.h>
 18: #include <petscts.h>

 20: #define PETSC_NEAR(a,b,c) (!(PetscAbsReal((a)-(b)) > (c)*PetscMax(PetscAbsReal(a),PetscAbsReal(b))))

 22: typedef struct {
 23:   Vec         global,local,localwork,solution;    /* location for local work (with ghost points) vector */
 24:   DM          da;                    /* manages ghost point communication */
 25:   PetscViewer viewer1,viewer2;
 26:   PetscInt    M;                     /* total number of grid points */
 27:   PetscReal   h;                     /* mesh width h = 1/(M-1) */
 28:   PetscReal   norm_2,norm_max;
 29:   PetscBool   nox;                   /* indicates problem is to be run without graphics */
 30: } AppCtx;


 40: #define linear_no_matrix        0
 41: #define linear_constant_matrix  1
 42: #define linear_variable_matrix  2
 43: #define nonlinear_no_jacobian   3
 44: #define nonlinear_jacobian      4

 48: int main(int argc,char **argv)
 49: {
 51:   PetscInt       maxsteps = 100,steps,m;
 52:   PetscMPIInt    size;
 53:   PetscInt       problem = linear_no_matrix;
 54:   PetscBool      flg;
 55:   AppCtx         appctx;
 56:   PetscReal      dt,ftime,maxtime=100.;
 57:   TS             ts;
 58:   Mat            A=0,Alhs=0;
 59:   MatStructure   A_structure;
 60:   TSProblemType  tsproblem = TS_LINEAR;
 61:   PetscDraw      draw;
 62:   char           tsinfo[120];
 63: 
 64:   PetscInitialize(&argc,&argv,(char*)0,help);
 65:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 67:   appctx.M = 60;
 68:   PetscOptionsGetInt(PETSC_NULL,"-M",&appctx.M,PETSC_NULL);
 69:   PetscOptionsGetInt(PETSC_NULL,"-time",&maxsteps,PETSC_NULL);
 70: 
 71:   PetscOptionsHasName(PETSC_NULL,"-nox",&appctx.nox);
 72:   appctx.norm_2 = 0.0; appctx.norm_max = 0.0;

 74:   /* Set up the ghost point communication pattern */
 75:   DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,appctx.M,1,1,PETSC_NULL,&appctx.da);
 76:   DMCreateGlobalVector(appctx.da,&appctx.global);
 77:   VecGetLocalSize(appctx.global,&m);
 78:   DMCreateLocalVector(appctx.da,&appctx.local);

 80:   /* Set up display to show wave graph */
 81:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,380,400,160,&appctx.viewer1);
 82:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
 83:   PetscDrawSetDoubleBuffer(draw);
 84:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,0,400,160,&appctx.viewer2);
 85:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
 86:   PetscDrawSetDoubleBuffer(draw);

 88:   /* make work array for evaluating right hand side function */
 89:   VecDuplicate(appctx.local,&appctx.localwork);

 91:   /* make work array for storing exact solution */
 92:   VecDuplicate(appctx.global,&appctx.solution);

 94:   appctx.h = 1.0/(appctx.M-1.0);

 96:   /* set initial conditions */
 97:   Initial(appctx.global,&appctx);
 98: 
 99:   /*
100:      This example is written to allow one to easily test parts 
101:     of TS, we do not expect users to generally need to use more
102:     then a single TSProblemType
103:   */
104:   tsproblem = TS_NONLINEAR;
105:   problem   = nonlinear_no_jacobian;
106:   PetscOptionsHasName(PETSC_NULL,"-linear_no_matrix",&flg);
107:   if (flg) {
108:     tsproblem = TS_LINEAR;
109:     problem   = linear_no_matrix;
110:   }
111:   PetscOptionsHasName(PETSC_NULL,"-linear_constant_matrix",&flg);
112:   if (flg) {
113:     tsproblem = TS_LINEAR;
114:     problem   = linear_constant_matrix;
115:   }
116:   PetscOptionsHasName(PETSC_NULL,"-linear_variable_matrix",&flg);
117:   if (flg) {
118:     tsproblem = TS_LINEAR;
119:     problem   = linear_variable_matrix;
120:   }
121:   PetscOptionsHasName(PETSC_NULL,"-nonlinear_no_jacobian",&flg);
122:   if (flg) {
123:     tsproblem = TS_NONLINEAR;
124:     problem   = nonlinear_jacobian;
125:   }
126:   PetscOptionsHasName(PETSC_NULL,"-nonlinear_jacobian",&flg);
127:   if (flg) {
128:     tsproblem = TS_NONLINEAR;
129:     problem   = nonlinear_jacobian;
130:   }
131: 
132:   /* make timestep context */
133:   TSCreate(PETSC_COMM_WORLD,&ts);
134:   TSSetProblemType(ts,tsproblem);
135:   PetscOptionsHasName(PETSC_NULL,"-monitor",&flg);
136:   if (flg) {
137:     TSMonitorSet(ts,Monitor,&appctx,PETSC_NULL);
138:   }

140:   dt = appctx.h*appctx.h/2.01;

142:   if (problem == linear_no_matrix) {
143:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This case needs rewritten");
144:     /*
145:          The user provides the RHS as a Shell matrix.
146:     */
147:     /*
148:     MatCreateShell(PETSC_COMM_WORLD,m,m,PETSC_DECIDE,PETSC_DECIDE,&appctx,&A);
149:     MatShellSetOperation(A,MATOP_MULT,(void(*)(void))RHSMatrixFree);
150:     TSSetMatrices(ts,A,PETSC_NULL,PETSC_NULL,PETSC_NULL,DIFFERENT_NONZERO_PATTERN,&appctx);
151:      */
152:   } else if (problem == linear_constant_matrix) {
153:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This case needs rewritten");
154:     /*
155:          The user provides the RHS as a constant matrix
156:     */
157:     /*
158:     MatCreate(PETSC_COMM_WORLD,&A);
159:     MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);
160:     MatSetFromOptions(A);
161:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx); 

163:     MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Alhs); 
164:     MatZeroEntries(Alhs);
165:     MatShift(Alhs,1.0);
166:     TSSetMatrices(ts,A,PETSC_NULL,Alhs,PETSC_NULL,DIFFERENT_NONZERO_PATTERN,&appctx);
167:      */
168:   } else if (problem == linear_variable_matrix) {
169:     SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"This case needs rewritten");
170:     /*
171:          The user provides the RHS as a time dependent matrix
172:     */
173:     /*
174:     MatCreate(PETSC_COMM_WORLD,&A);
175:     MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);
176:     MatSetFromOptions(A);
177:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);

179:     MatDuplicate(A,MAT_DO_NOT_COPY_VALUES,&Alhs);    
180:     MatZeroEntries(Alhs);
181:     MatShift(Alhs,1.0);
182:     LHSMatrixHeat(ts,0.0,&Alhs,&Alhs,&A_structure,&appctx);
183:     TSSetMatrices(ts,A,RHSMatrixHeat,Alhs,LHSMatrixHeat,DIFFERENT_NONZERO_PATTERN,&appctx);
184:      */
185:   } else if (problem == nonlinear_no_jacobian) {
186:     /*
187:          The user provides the RHS and a Shell Jacobian
188:     */
189:     TSSetRHSFunction(ts,PETSC_NULL,RHSFunctionHeat,&appctx);
190:     MatCreateShell(PETSC_COMM_WORLD,m,m,PETSC_DECIDE,PETSC_DECIDE,&appctx,&A);
191:     MatShellSetOperation(A,MATOP_MULT,(void(*)(void))RHSMatrixFree);
192:     TSSetRHSJacobian(ts,A,A,PETSC_NULL,&appctx);
193:   } else if (problem == nonlinear_jacobian) {
194:     /*
195:          The user provides the RHS and Jacobian
196:     */
197:     TSSetRHSFunction(ts,PETSC_NULL,RHSFunctionHeat,&appctx);
198:     MatCreate(PETSC_COMM_WORLD,&A);
199:     MatSetSizes(A,m,m,PETSC_DECIDE,PETSC_DECIDE);
200:     MatSetFromOptions(A);
201:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);
202:     TSSetRHSJacobian(ts,A,A,RHSJacobianHeat,&appctx);
203:   }

205:   TSSetInitialTimeStep(ts,0.0,dt);
206:   TSSetDuration(ts,maxsteps,maxtime);
207:   TSSetSolution(ts,appctx.global);

209:   TSSetFromOptions(ts);

211:   TSSolve(ts,appctx.global,&ftime);

213:   PetscOptionsHasName(PETSC_NULL,"-testinfo",&flg);
214:   if (flg) {
215:     PetscBool  iseuler;
216:     TSGetTimeStepNumber(ts,&steps);

218:     PetscTypeCompare((PetscObject)ts,"euler",&iseuler);
219:     if (iseuler) {
220:       if (!PETSC_NEAR(appctx.norm_2/steps,0.00257244,1.e-4)) {
221:         fprintf(stdout,"Error in Euler method: 2-norm %G expecting: 0.00257244\n",appctx.norm_2/steps);
222:       }
223:     } else {
224:       PetscPrintf(PETSC_COMM_WORLD,"%D Procs; Avg. error 2-norm %G; max-norm %G; %s\n",
225:                 size,appctx.norm_2/steps,appctx.norm_max/steps,tsinfo);
226:     }
227:   }

229:   TSDestroy(&ts);
230:   PetscViewerDestroy(&appctx.viewer1);
231:   PetscViewerDestroy(&appctx.viewer2);
232:   VecDestroy(&appctx.localwork);
233:   VecDestroy(&appctx.solution);
234:   VecDestroy(&appctx.local);
235:   VecDestroy(&appctx.global);
236:   DMDestroy(&appctx.da);
237:   if (A) {ierr= MatDestroy(&A);}
238:   if (Alhs) {ierr= MatDestroy(&Alhs);}

240:   PetscFinalize();
241:   return 0;
242: }

244: /* -------------------------------------------------------------------*/
245: /*
246:   Set initial condition: u(t=0) = sin(6*pi*x) + 3*sin(2*pi*x)
247: */
250: PetscErrorCode Initial(Vec global,void *ctx)
251: {
252:   AppCtx         *appctx = (AppCtx*) ctx;
253:   PetscScalar    *localptr,h = appctx->h;
254:   PetscInt       i,mybase,myend;

258:   /* determine starting point of each processor */
259:   VecGetOwnershipRange(global,&mybase,&myend);

261:   /* Initialize the array */
262:   VecGetArray(global,&localptr);
263:   for (i=mybase; i<myend; i++) {
264:     localptr[i-mybase] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
265:   }
266:   VecRestoreArray(global,&localptr);
267:   return(0);
268: }

272: /*
273:    Exact solution: 
274:      u = sin(6*pi*x)*exp(-36*pi*pi*t) + 3*sin(2*pi*x)*exp(-4*pi*pi*t)
275: */
276: PetscErrorCode Solution(PetscReal t,Vec solution,void *ctx)
277: {
278:   AppCtx *       appctx = (AppCtx*) ctx;
279:   PetscScalar    *localptr,h = appctx->h,ex1,ex2,sc1,sc2;
280:   PetscInt       i,mybase,myend;

284:   /* determine starting point of each processor */
285:   VecGetOwnershipRange(solution,&mybase,&myend);

287:   ex1 = exp(-36.*PETSC_PI*PETSC_PI*t);
288:   ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
289:   sc1 = PETSC_PI*6.*h;
290:   sc2 = PETSC_PI*2.*h;
291:   VecGetArray(solution,&localptr);
292:   for (i=mybase; i<myend; i++) {
293:     localptr[i-mybase] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
294:   }
295:   VecRestoreArray(solution,&localptr);
296:   return(0);
297: }

299: /*
300:   step   - iteration number
301:   ltime  - current time
302:   global - current iterate
303:  */
306: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal ltime,Vec global,void *ctx)
307: {
308:   AppCtx         *appctx = (AppCtx*) ctx;
310:   PetscReal      norm_2,norm_max;
311:   MPI_Comm       comm;

314:   PetscObjectGetComm((PetscObject)ts,&comm);
315:   if (!appctx->nox) {
316:     VecView(global,appctx->viewer2); /* show wave graph */
317:   }
318:   Solution(ltime,appctx->solution,ctx); /* get true solution at current time */
319:   VecAXPY(appctx->solution,-1.0,global);
320:   VecNorm(appctx->solution,NORM_2,&norm_2);
321:   norm_2 = sqrt(appctx->h)*norm_2;
322:   VecNorm(appctx->solution,NORM_MAX,&norm_max);
323:   PetscPrintf(comm,"timestep %D time %G norm of error %.5f %.5f\n",step,ltime,norm_2,norm_max);

325:   appctx->norm_2   += norm_2;
326:   appctx->norm_max += norm_max;
327:   if (!appctx->nox) {
328:     VecView(appctx->solution,appctx->viewer1);
329:   }
330:   return(0);
331: }

333: /* -----------------------------------------------------------------------*/
336: PetscErrorCode RHSMatrixFree(Mat mat,Vec x,Vec y)
337: {
338:   PetscErrorCode  ierr;
339:   void            *ctx;

342:   MatShellGetContext(mat,(void **)&ctx);
343:   RHSFunctionHeat(0,0.0,x,y,ctx);
344:   return(0);
345: }

349: PetscErrorCode RHSFunctionHeat(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
350: {
351:   AppCtx         *appctx = (AppCtx*) ctx;
352:   DM             da = appctx->da;
353:   Vec            local = appctx->local,localwork = appctx->localwork;
355:   PetscInt       i,localsize;
356:   PetscScalar    *copyptr,*localptr,sc;

359:   /*Extract local array */
360:   DMGlobalToLocalBegin(da,globalin,INSERT_VALUES,local);
361:   DMGlobalToLocalEnd(da,globalin,INSERT_VALUES,local);
362:   VecGetArray(local,&localptr);

364:   /* Extract work vector */
365:   VecGetArray(localwork,&copyptr);

367:   /* Update Locally - Make array of new values */
368:   /* Note: For the first and last entry I copy the value */
369:   /* if this is an interior node it is irrelevant */
370:   sc = 1.0/(appctx->h*appctx->h);
371:   VecGetLocalSize(local,&localsize);
372:   copyptr[0] = localptr[0];
373:   for (i=1; i<localsize-1; i++) {
374:     copyptr[i] = sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);
375:   }
376:   copyptr[localsize-1] = localptr[localsize-1];
377:   VecRestoreArray(local,&localptr);
378:   VecRestoreArray(localwork,&copyptr);

380:   /* Local to Global */
381:   DMLocalToGlobalBegin(da,localwork,INSERT_VALUES,globalout);
382:   DMLocalToGlobalEnd(da,localwork,INSERT_VALUES,globalout);
383:   return(0);
384: }

386: /* ---------------------------------------------------------------------*/
389: PetscErrorCode RHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
390: {
391:   Mat            A = *AA;
392:   AppCtx         *appctx = (AppCtx*) ctx;
394:   PetscInt       i,mstart,mend,idx[3];
395:   PetscMPIInt    size,rank;
396:   PetscScalar    v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;

399:   *str = SAME_NONZERO_PATTERN;
400: 
401:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
402:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

404:   MatGetOwnershipRange(A,&mstart,&mend);
405:   if (mstart == 0) {
406:     v[0] = 1.0;
407:     MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
408:     mstart++;
409:   }
410:   if (mend == appctx->M) {
411:     mend--;
412:     v[0] = 1.0;
413:     MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
414:   }

416:   /*
417:      Construct matrice one row at a time
418:   */
419:   v[0] = sone; v[1] = stwo; v[2] = sone;
420:   for (i=mstart; i<mend; i++) {
421:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
422:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
423:   }

425:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
426:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
427:   return(0);
428: }

432: PetscErrorCode RHSJacobianHeat(TS ts,PetscReal t,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
433: {
435:   RHSMatrixHeat(ts,t,AA,BB,str,ctx);
436:   return(0);
437: }

439: /* A = indentity matrix */
442: PetscErrorCode LHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
443: {
445:   Mat            A = *AA;

448:   *str = SAME_NONZERO_PATTERN;
449:   MatZeroEntries(A);
450:   MatShift(A,1.0);
451:   return(0);
452: }