Actual source code: ex1.c

  1: /*$Id: ex1.c,v 1.50 2001/09/07 20:12:09 bsmith Exp $*/
  2: /*
  3:        Formatted test for TS routines.

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

 10: static char help[] = "Solves 1D heat equation.nn";

 12:  #include petscda.h
 13:  #include petscsys.h
 14:  #include petscts.h

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

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

 28: extern int Monitor(TS,int,PetscReal,Vec,void *);
 29: extern int RHSFunctionHeat(TS,PetscReal,Vec,Vec,void*);
 30: extern int RHSMatrixFree(Mat,Vec,Vec);
 31: extern int Initial(Vec,void*);
 32: extern int RHSMatrixHeat(TS,PetscReal,Mat *,Mat *,MatStructure *,void *);
 33: extern int RHSJacobianHeat(TS,PetscReal,Vec,Mat*,Mat*,MatStructure *,void*);

 35: #define linear_no_matrix       0
 36: #define linear_no_time         1
 37: #define linear                 2
 38: #define nonlinear_no_jacobian  3
 39: #define nonlinear              4

 41: int main(int argc,char **argv)
 42: {
 43:   int           ierr,time_steps = 100,steps,size,m;
 44:   int           problem = linear_no_matrix;
 45:   PetscTruth    flg;
 46:   AppCtx        appctx;
 47:   PetscReal     dt,ftime;
 48:   TS            ts;
 49:   Mat           A = 0;
 50:   MatStructure  A_structure;
 51:   TSProblemType tsproblem = TS_LINEAR;
 52:   PetscDraw     draw;
 53:   PetscViewer   viewer;
 54:   char          tsinfo[120];
 55: 
 56:   PetscInitialize(&argc,&argv,(char*)0,help);
 57:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

 59:   appctx.M = 60;
 60:   PetscOptionsGetInt(PETSC_NULL,"-M",&appctx.M,PETSC_NULL);
 61:   PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
 62: 
 63:   PetscOptionsHasName(PETSC_NULL,"-nox",&flg);
 64:   if (flg) appctx.nox = 1; else appctx.nox = 0;
 65:   appctx.norm_2 = 0.0; appctx.norm_max = 0.0;

 67:   /* Set up the ghost point communication pattern */
 68:   DACreate1d(PETSC_COMM_WORLD,DA_NONPERIODIC,appctx.M,1,1,PETSC_NULL,&appctx.da);
 69:   DACreateGlobalVector(appctx.da,&appctx.global);
 70:   VecGetLocalSize(appctx.global,&m);
 71:   DACreateLocalVector(appctx.da,&appctx.local);

 73:   /* Set up display to show wave graph */

 75:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,380,400,160,&appctx.viewer1);
 76:   PetscViewerDrawGetDraw(appctx.viewer1,0,&draw);
 77:   PetscDrawSetDoubleBuffer(draw);
 78:   PetscViewerDrawOpen(PETSC_COMM_WORLD,0,"",80,0,400,160,&appctx.viewer2);
 79:   PetscViewerDrawGetDraw(appctx.viewer2,0,&draw);
 80:   PetscDrawSetDoubleBuffer(draw);


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

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

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

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

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

132:   if (problem == linear_no_matrix) {
133:     /*
134:          The user provides the RHS as a Shell matrix.
135:     */
136:     MatCreateShell(PETSC_COMM_WORLD,m,appctx.M,appctx.M,appctx.M,&appctx,&A);
137:     MatShellSetOperation(A,MATOP_MULT,(void(*)(void))RHSMatrixFree);
138:     TSSetRHSMatrix(ts,A,A,PETSC_NULL,&appctx);
139:   } else if (problem == linear_no_time) {
140:     /*
141:          The user provides the RHS as a matrix
142:     */
143:     MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,appctx.M,appctx.M,&A);
144:     MatSetFromOptions(A);
145:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);
146:     TSSetRHSMatrix(ts,A,A,PETSC_NULL,&appctx);
147:   } else if (problem == linear) {
148:     /*
149:          The user provides the RHS as a time dependent matrix
150:     */
151:     MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,appctx.M,appctx.M,&A);
152:     MatSetFromOptions(A);
153:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);
154:     TSSetRHSMatrix(ts,A,A,RHSMatrixHeat,&appctx);
155:   } else if (problem == nonlinear_no_jacobian) {
156:     /*
157:          The user provides the RHS and a Shell Jacobian
158:     */
159:     TSSetRHSFunction(ts,RHSFunctionHeat,&appctx);
160:     MatCreateShell(PETSC_COMM_WORLD,m,appctx.M,appctx.M,appctx.M,&appctx,&A);
161:     MatShellSetOperation(A,MATOP_MULT,(void(*)(void))RHSMatrixFree);
162:     TSSetRHSJacobian(ts,A,A,PETSC_NULL,&appctx);
163:   } else if (problem == nonlinear) {
164:     /*
165:          The user provides the RHS and Jacobian
166:     */
167:     TSSetRHSFunction(ts,RHSFunctionHeat,&appctx);
168:     MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,appctx.M,appctx.M,&A);
169:     MatSetFromOptions(A);
170:     RHSMatrixHeat(ts,0.0,&A,&A,&A_structure,&appctx);
171:     TSSetRHSJacobian(ts,A,A,RHSJacobianHeat,&appctx);
172:   }

174:   TSSetFromOptions(ts);

176:   TSSetInitialTimeStep(ts,0.0,dt);
177:   TSSetDuration(ts,time_steps,100.);
178:   TSSetSolution(ts,appctx.global);


181:   TSSetUp(ts);
182:   TSStep(ts,&steps,&ftime);
183:   PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,120,&viewer);
184:   TSView(ts,viewer);

186:   PetscOptionsHasName(PETSC_NULL,"-test",&flg);
187:   if (flg) {
188:     PetscTruth iseuler;
189:     PetscTypeCompare((PetscObject)ts,"euler",&iseuler);
190:     if (iseuler) {
191:       if (!PETSC_NEAR(appctx.norm_2/steps,0.00257244,1.e-4)) {
192:         fprintf(stdout,"Error in Euler method: 2-norm %g expecting: 0.00257244n",appctx.norm_2/steps);
193:       }
194:     } else {
195:       if (!PETSC_NEAR(appctx.norm_2/steps,0.00506174,1.e-4)) {
196:         fprintf(stdout,"Error in %s method: 2-norm %g expecting: 0.00506174n",tsinfo,appctx.norm_2/steps);
197:       }
198:     }
199:   } else {
200:     PetscPrintf(PETSC_COMM_WORLD,"%d Procs Avg. error 2 norm %g max norm %g %sn",
201:                 size,appctx.norm_2/steps,appctx.norm_max/steps,tsinfo);
202:   }

204:   PetscViewerDestroy(viewer);
205:   TSDestroy(ts);
206:   PetscViewerDestroy(appctx.viewer1);
207:   PetscViewerDestroy(appctx.viewer2);
208:   VecDestroy(appctx.localwork);
209:   VecDestroy(appctx.solution);
210:   VecDestroy(appctx.local);
211:   VecDestroy(appctx.global);
212:   DADestroy(appctx.da);
213:   if (A) {ierr= MatDestroy(A);}

215:   PetscFinalize();
216:   return 0;
217: }

219: /* -------------------------------------------------------------------*/
220: int Initial(Vec global,void *ctx)
221: {
222:   AppCtx *appctx = (AppCtx*) ctx;
223:   PetscScalar *localptr,h = appctx->h;
224:   int    i,mybase,myend,ierr;

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

229:   /* Initialize the array */
230:   VecGetArray(global,&localptr);
231:   for (i=mybase; i<myend; i++) {
232:     localptr[i-mybase] = PetscSinScalar(PETSC_PI*i*6.*h) + 3.*PetscSinScalar(PETSC_PI*i*2.*h);
233:   }
234:   VecRestoreArray(global,&localptr);
235:   return 0;
236: }

238: /*
239:        Exact solution 
240: */
241: int Solution(PetscReal t,Vec solution,void *ctx)
242: {
243:   AppCtx *appctx = (AppCtx*) ctx;
244:   PetscScalar *localptr,h = appctx->h,ex1,ex2,sc1,sc2;
245:   int    i,mybase,myend,ierr;

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

250:   ex1 = exp(-36.*PETSC_PI*PETSC_PI*t);
251:   ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
252:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
253:   VecGetArray(solution,&localptr);
254:   for (i=mybase; i<myend; i++) {
255:     localptr[i-mybase] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
256:   }
257:   VecRestoreArray(solution,&localptr);
258:   return 0;
259: }

261: int Monitor(TS ts,int step,PetscReal ltime,Vec global,void *ctx)
262: {
263:   AppCtx       *appctx = (AppCtx*) ctx;
264:   int          ierr;
265:   PetscReal    norm_2,norm_max;
266:   PetscScalar  mone = -1.0;
267:   MPI_Comm     comm;

269:   PetscObjectGetComm((PetscObject)ts,&comm);

271:   VecView(global,appctx->viewer2);

273:   Solution(ltime,appctx->solution,ctx);
274:   VecAXPY(&mone,global,appctx->solution);
275:   VecNorm(appctx->solution,NORM_2,&norm_2);
276:   norm_2 = sqrt(appctx->h)*norm_2;
277:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

279:   if (!appctx->nox) {
280:     PetscPrintf(comm,"timestep %d time %g norm of error %g %gn",step,ltime,norm_2,norm_max);
281:   }

283:   appctx->norm_2   += norm_2;
284:   appctx->norm_max += norm_max;

286:   VecView(appctx->solution,appctx->viewer1);

288:   return 0;
289: }

291: /* -----------------------------------------------------------------------*/
292: int RHSMatrixFree(Mat mat,Vec x,Vec y)
293: {
294:   int  ierr;
295:   void *ctx;

297:   MatShellGetContext(mat,(void **)&ctx);
298:   RHSFunctionHeat(0,0.0,x,y,ctx);
299:   return 0;
300: }

302: int RHSFunctionHeat(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
303: {
304:   AppCtx      *appctx = (AppCtx*) ctx;
305:   DA          da = appctx->da;
306:   Vec         local = appctx->local,localwork = appctx->localwork;
307:   int         ierr,i,localsize;
308:   PetscScalar *copyptr,*localptr,sc;

310:   /*Extract local array */
311:   DAGlobalToLocalBegin(da,globalin,INSERT_VALUES,local);
312:   DAGlobalToLocalEnd(da,globalin,INSERT_VALUES,local);
313:   VecGetArray(local,&localptr);

315:   /* Extract work vector */
316:   VecGetArray(localwork,&copyptr);

318:   /* Update Locally - Make array of new values */
319:   /* Note: For the first and last entry I copy the value */
320:   /* if this is an interior node it is irrelevant */
321:   sc = 1.0/(appctx->h*appctx->h);
322:   VecGetLocalSize(local,&localsize);
323:   copyptr[0] = localptr[0];
324:   for (i=1; i<localsize-1; i++) {
325:     copyptr[i] = sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);
326:   }
327:   copyptr[localsize-1] = localptr[localsize-1];
328:   VecRestoreArray(local,&localptr);
329:   VecRestoreArray(localwork,&copyptr);

331:   /* Local to Global */
332:   DALocalToGlobal(da,localwork,INSERT_VALUES,globalout);
333:   return 0;
334: }

336: /* ---------------------------------------------------------------------*/
337: int RHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
338: {
339:   Mat    A = *AA;
340:   AppCtx *appctx = (AppCtx*) ctx;
341:   int    ierr,i,mstart,mend,rank,size,idx[3];
342:   PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;

344:   *str = SAME_NONZERO_PATTERN;

346:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
347:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

349:   MatGetOwnershipRange(A,&mstart,&mend);
350:   if (mstart == 0) {
351:     v[0] = 1.0;
352:     MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
353:     mstart++;
354:   }
355:   if (mend == appctx->M) {
356:     mend--;
357:     v[0] = 1.0;
358:     MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
359:   }

361:   /*
362:      Construct matrice one row at a time
363:   */
364:   v[0] = sone; v[1] = stwo; v[2] = sone;
365:   for (i=mstart; i<mend; i++) {
366:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
367:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
368:   }

370:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
371:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
372:   return 0;
373: }

375: int RHSJacobianHeat(TS ts,PetscReal t,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
376: {
377:   return RHSMatrixHeat(ts,t,AA,BB,str,ctx);
378: }