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.\n\n";

 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

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

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

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

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

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


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

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

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

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

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

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

176:   TSSetFromOptions(ts);

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


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

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

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

217:   PetscFinalize();
218:   return 0;
219: }

221: /* -------------------------------------------------------------------*/
224: int Initial(Vec global,void *ctx)
225: {
226:   AppCtx *appctx = (AppCtx*) ctx;
227:   PetscScalar *localptr,h = appctx->h;
228:   int    i,mybase,myend,ierr;

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

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

244: /*
245:        Exact solution 
246: */
247: int Solution(PetscReal t,Vec solution,void *ctx)
248: {
249:   AppCtx *appctx = (AppCtx*) ctx;
250:   PetscScalar *localptr,h = appctx->h,ex1,ex2,sc1,sc2;
251:   int    i,mybase,myend,ierr;

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

256:   ex1 = exp(-36.*PETSC_PI*PETSC_PI*t);
257:   ex2 = exp(-4.*PETSC_PI*PETSC_PI*t);
258:   sc1 = PETSC_PI*6.*h;                 sc2 = PETSC_PI*2.*h;
259:   VecGetArray(solution,&localptr);
260:   for (i=mybase; i<myend; i++) {
261:     localptr[i-mybase] = PetscSinScalar(sc1*(PetscReal)i)*ex1 + 3.*PetscSinScalar(sc2*(PetscReal)i)*ex2;
262:   }
263:   VecRestoreArray(solution,&localptr);
264:   return 0;
265: }

269: int Monitor(TS ts,int step,PetscReal ltime,Vec global,void *ctx)
270: {
271:   AppCtx       *appctx = (AppCtx*) ctx;
272:   int          ierr;
273:   PetscReal    norm_2,norm_max;
274:   PetscScalar  mone = -1.0;
275:   MPI_Comm     comm;

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

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

281:   Solution(ltime,appctx->solution,ctx);
282:   VecAXPY(&mone,global,appctx->solution);
283:   VecNorm(appctx->solution,NORM_2,&norm_2);
284:   norm_2 = sqrt(appctx->h)*norm_2;
285:   VecNorm(appctx->solution,NORM_MAX,&norm_max);

287:   if (!appctx->nox) {
288:     PetscPrintf(comm,"timestep %d time %g norm of error %g %g\n",step,ltime,norm_2,norm_max);
289:   }

291:   appctx->norm_2   += norm_2;
292:   appctx->norm_max += norm_max;

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

296:   return 0;
297: }

299: /* -----------------------------------------------------------------------*/
302: int RHSMatrixFree(Mat mat,Vec x,Vec y)
303: {
304:   int  ierr;
305:   void *ctx;

307:   MatShellGetContext(mat,(void **)&ctx);
308:   RHSFunctionHeat(0,0.0,x,y,ctx);
309:   return 0;
310: }

314: int RHSFunctionHeat(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
315: {
316:   AppCtx      *appctx = (AppCtx*) ctx;
317:   DA          da = appctx->da;
318:   Vec         local = appctx->local,localwork = appctx->localwork;
319:   int         ierr,i,localsize;
320:   PetscScalar *copyptr,*localptr,sc;

322:   /*Extract local array */
323:   DAGlobalToLocalBegin(da,globalin,INSERT_VALUES,local);
324:   DAGlobalToLocalEnd(da,globalin,INSERT_VALUES,local);
325:   VecGetArray(local,&localptr);

327:   /* Extract work vector */
328:   VecGetArray(localwork,&copyptr);

330:   /* Update Locally - Make array of new values */
331:   /* Note: For the first and last entry I copy the value */
332:   /* if this is an interior node it is irrelevant */
333:   sc = 1.0/(appctx->h*appctx->h);
334:   VecGetLocalSize(local,&localsize);
335:   copyptr[0] = localptr[0];
336:   for (i=1; i<localsize-1; i++) {
337:     copyptr[i] = sc * (localptr[i+1] + localptr[i-1] - 2.0*localptr[i]);
338:   }
339:   copyptr[localsize-1] = localptr[localsize-1];
340:   VecRestoreArray(local,&localptr);
341:   VecRestoreArray(localwork,&copyptr);

343:   /* Local to Global */
344:   DALocalToGlobal(da,localwork,INSERT_VALUES,globalout);
345:   return 0;
346: }

348: /* ---------------------------------------------------------------------*/
351: int RHSMatrixHeat(TS ts,PetscReal t,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
352: {
353:   Mat    A = *AA;
354:   AppCtx *appctx = (AppCtx*) ctx;
355:   int    ierr,i,mstart,mend,rank,size,idx[3];
356:   PetscScalar v[3],stwo = -2./(appctx->h*appctx->h),sone = -.5*stwo;

358:   *str = SAME_NONZERO_PATTERN;

360:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
361:   MPI_Comm_size(PETSC_COMM_WORLD,&size);

363:   MatGetOwnershipRange(A,&mstart,&mend);
364:   if (mstart == 0) {
365:     v[0] = 1.0;
366:     MatSetValues(A,1,&mstart,1,&mstart,v,INSERT_VALUES);
367:     mstart++;
368:   }
369:   if (mend == appctx->M) {
370:     mend--;
371:     v[0] = 1.0;
372:     MatSetValues(A,1,&mend,1,&mend,v,INSERT_VALUES);
373:   }

375:   /*
376:      Construct matrice one row at a time
377:   */
378:   v[0] = sone; v[1] = stwo; v[2] = sone;
379:   for (i=mstart; i<mend; i++) {
380:     idx[0] = i-1; idx[1] = i; idx[2] = i+1;
381:     MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
382:   }

384:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
385:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
386:   return 0;
387: }

391: int RHSJacobianHeat(TS ts,PetscReal t,Vec x,Mat *AA,Mat *BB,MatStructure *str,void *ctx)
392: {
393:   return RHSMatrixHeat(ts,t,AA,BB,str,ctx);
394: }