Actual source code: ex4.c

  1: /*$Id: ex4.c,v 1.13 2001/08/07 21:31:27 bsmith Exp $*/
  2: /*
  3:        The Problem:
  4:            Solve the convection-diffusion equation:
  5:            
  6:              u_t+a*(u_x+u_y)=epsilon*(u_xx+u_yy)
  7:              u=0 at x=0, y=0
  8:              u_x=0 at x=1
  9:              u_y=0 at y=1
 10:         
 11:        This program tests the routine of computing the Jacobian by the 
 12:        finite difference method as well as PETSc with PVODE.

 14: */

 16: static char help[] = "Solve the convection-diffusion equation. nn";

 18:  #include petscsys.h
 19:  #include petscts.h

 21: extern int Monitor(TS,int,PetscReal,Vec,void *);
 22: extern int Initial(Vec,void *);

 24: typedef struct
 25: {
 26:   int                 m;        /* the number of mesh points in x-direction */
 27:   int                 n;      /* the number of mesh points in y-direction */
 28:   PetscReal         dx;     /* the grid space in x-direction */
 29:   PetscReal     dy;     /* the grid space in y-direction */
 30:   PetscReal     a;      /* the convection coefficient    */
 31:   PetscReal     epsilon; /* the diffusion coefficient    */
 32: } Data;

 34: /* two temporal functions */
 35: extern int FormJacobian(SNES,Vec,Mat*,Mat*,MatStructure*,void*);
 36: extern int FormFunction(SNES,Vec,Vec,void*);
 37: extern int RHSFunction(TS,PetscReal,Vec,Vec,void*);
 38: extern int RHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure *,void*);

 40: /* the initial function */
 41: PetscReal f_ini(PetscReal x,PetscReal y)
 42: {
 43:   PetscReal f;
 44:   f=exp(-20.0*(pow(x-0.5,2.0)+pow(y-0.5,2.0)));
 45:   return f;
 46: }


 49: #define linear_no_matrix       0
 50: #define linear_no_time         1
 51: #define linear                 2
 52: #define nonlinear_no_jacobian  3
 53: #define nonlinear              4

 55: int main(int argc,char **argv)
 56: {
 57:   int           ierr,time_steps = 100,steps,size;
 58:   Vec           global;
 59:   PetscReal     dt,ftime;
 60:   TS            ts;
 61:   PetscViewer        viewfile;
 62:   MatStructure  A_structure;
 63:   Mat           A = 0;
 64:   TSProblemType tsproblem = TS_NONLINEAR; /* Need to be TS_NONLINEAR */
 65:   SNES          snes;
 66:   Vec                 x;
 67:   Data                data;
 68:   int                 mn;
 69: #if defined(PETSC_HAVE_PVODE) && !defined(__cplusplus)
 70:   PC                pc;
 71:   PetscViewer   viewer;
 72:   char          pcinfo[120],tsinfo[120];
 73: #endif

 75:   PetscInitialize(&argc,&argv,(char*)0,help);
 76:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
 77: 
 78:   /* set Data */
 79:   data.m = 9;
 80:   data.n = 9;
 81:   data.a = 1.0;
 82:   data.epsilon = 0.1;
 83:   data.dx = 1.0/(data.m+1.0);
 84:   data.dy = 1.0/(data.n+1.0);
 85:   mn = (data.m)*(data.n);

 87:   PetscOptionsGetInt(PETSC_NULL,"-time",&time_steps,PETSC_NULL);
 88: 
 89:   /* set initial conditions */
 90:   VecCreate(PETSC_COMM_WORLD,&global);
 91:   VecSetSizes(global,PETSC_DECIDE,mn);
 92:   VecSetFromOptions(global);
 93:   Initial(global,&data);
 94:   VecDuplicate(global,&x);
 95: 
 96:   /* make timestep context */
 97:   TSCreate(PETSC_COMM_WORLD,&ts);
 98:   TSSetProblemType(ts,tsproblem);
 99:   TSSetMonitor(ts,Monitor,PETSC_NULL,PETSC_NULL);

101:   dt = 0.1;

103:   /*
104:     The user provides the RHS and Jacobian
105:   */
106:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,mn,mn,&A);
107:   MatSetFromOptions(A);

109:   /* Create SNES context  */
110:   SNESCreate(PETSC_COMM_WORLD,&snes);
111: 

113:   /* setting the RHS function and the Jacobian's non-zero structutre */
114:   SNESSetFunction(snes,global,FormFunction,&data);
115:   SNESSetJacobian(snes,A,A,FormJacobian,&data);

117:   /* set TSPVodeRHSFunction and TSPVodeRHSJacobian, so PETSc will pick up the 
118:      RHS function from SNES and compute the Jacobian by FD */
119:   /*
120:   TSSetRHSFunction(ts,TSPVodeSetRHSFunction,snes);
121:   TSPVodeSetRHSJacobian(ts,0.0,global,&A,&A,&A_structure,snes);
122:   TSSetRHSJacobian(ts,A,A,TSPVodeSetRHSJacobian,snes);
123:   */
124: 
125:   TSSetRHSFunction(ts,RHSFunction,&data);
126:   RHSJacobian(ts,0.0,global,&A,&A,&A_structure,&data);
127:   TSSetRHSJacobian(ts,A,A,RHSJacobian,&data);

129:   /* Use PVODE */
130:   TSSetType(ts,TS_PVODE);

132:   TSSetFromOptions(ts);

134:   TSSetInitialTimeStep(ts,0.0,dt);
135:   TSSetDuration(ts,time_steps,1);
136:   TSSetSolution(ts,global);


139:   /* Pick up a Petsc preconditioner */
140:   /* one can always set method or preconditioner during the run time */
141: #if defined(PETSC_HAVE_PVODE) && !defined(__cplusplus)
142:   TSPVodeGetPC(ts,&pc);
143:   PCSetType(pc,PCJACOBI);
144:   TSPVodeSetType(ts,PVODE_BDF);
145:   /* TSPVodeSetMethodFromOptions(ts); */
146: #endif

148:   TSSetUp(ts);
149:   TSStep(ts,&steps,&ftime);

151:   TSGetSolution(ts,&global);
152:   PetscViewerASCIIOpen(PETSC_COMM_SELF,"out.m",&viewfile);
153:   PetscViewerSetFormat(viewfile,PETSC_VIEWER_ASCII_MATLAB);
154:   VecView(global,viewfile);

156: #if defined(PETSC_HAVE_PVODE) && !defined(__cplusplus)
157:   /* extracts the PC  from ts */
158:   TSPVodeGetPC(ts,&pc);
159:   PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,120,&viewer);
160:   TSView(ts,viewer);
161:   PetscViewerStringOpen(PETSC_COMM_WORLD,pcinfo,120,&viewer);
162:   PCView(pc,viewer);
163:   PetscPrintf(PETSC_COMM_WORLD,"%d Procs,%s Preconditioner,%sn",
164:                      size,tsinfo,pcinfo);
165:   PCDestroy(pc);
166: #endif

168:   /* free the memories */
169:   TSDestroy(ts);
170:   VecDestroy(global);
171:   if (A) {ierr= MatDestroy(A);}
172:   PetscFinalize();
173:   return 0;
174: }

176: /* -------------------------------------------------------------------*/
177: int Initial(Vec global,void *ctx)
178: {
179:   Data        *data = (Data*)ctx;
180:   int         m;
181:   int         row,col;
182:   PetscReal   x,y,dx,dy;
183:   PetscScalar *localptr;
184:   int         i,mybase,myend,ierr,locsize;

186:   /* make the local  copies of parameters */
187:   m = data->m;
188:   dx = data->dx;
189:   dy = data->dy;

191:   /* determine starting point of each processor */
192:   VecGetOwnershipRange(global,&mybase,&myend);
193:   VecGetLocalSize(global,&locsize);

195:   /* Initialize the array */
196:   VecGetArray(global,&localptr);

198:   for (i=0; i<locsize; i++) {
199:     row = 1+(mybase+i)-((mybase+i)/m)*m;
200:     col = (mybase+i)/m+1;
201:     x = dx*row;
202:     y = dy*col;
203:     localptr[i] = f_ini(x,y);
204:   }
205: 
206:   VecRestoreArray(global,&localptr);
207:   return 0;
208: }

210: int Monitor(TS ts,int step,PetscReal time,Vec global,void *ctx)
211: {
212:   VecScatter scatter;
213:   IS from,to;
214:   int i,n,*idx;
215:   Vec tmp_vec;
216:   int      ierr;
217:   PetscScalar   *tmp;

219:   /* Get the size of the vector */
220:   VecGetSize(global,&n);

222:   /* Set the index sets */
223:   PetscMalloc(n*sizeof(int),&idx);
224:   for(i=0; i<n; i++) idx[i]=i;
225: 
226:   /* Create local sequential vectors */
227:   VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);

229:   /* Create scatter context */
230:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&from);
231:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&to);
232:   VecScatterCreate(global,from,tmp_vec,to,&scatter);
233:   VecScatterBegin(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);
234:   VecScatterEnd(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);

236:   VecGetArray(tmp_vec,&tmp);
237:   PetscPrintf(PETSC_COMM_WORLD,"At t =%14.6e u= %14.6e at the center n",time,PetscRealPart(tmp[n/2]));
238:   VecRestoreArray(tmp_vec,&tmp);

240:   PetscFree(idx);
241:   return 0;
242: }


245: int FormFunction(SNES snes,Vec globalin,Vec globalout,void *ptr)
246: {
247:   Data        *data = (Data*)ptr;
248:   int         m,n,mn;
249:   PetscReal   dx,dy;
250:   PetscReal   xc,xl,xr,yl,yr;
251:   PetscReal   a,epsilon;
252:   PetscScalar *inptr,*outptr;
253:   int         i,j,len,ierr;

255:   IS          from,to;
256:   int         *idx;
257:   VecScatter  scatter;
258:   Vec         tmp_in,tmp_out;

260:   m = data->m;
261:   n = data->n;
262:   mn = m*n;
263:   dx = data->dx;
264:   dy = data->dy;
265:   a = data->a;
266:   epsilon = data->epsilon;

268:   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
269:   xl = 0.5*a/dx+epsilon/(dx*dx);
270:   xr = -0.5*a/dx+epsilon/(dx*dx);
271:   yl = 0.5*a/dy+epsilon/(dy*dy);
272:   yr = -0.5*a/dy+epsilon/(dy*dy);
273: 
274:   /* Get the length of parallel vector */
275:   VecGetSize(globalin,&len);

277:   /* Set the index sets */
278:   PetscMalloc(len*sizeof(int),&idx);
279:   for(i=0; i<len; i++) idx[i]=i;
280: 
281:   /* Create local sequential vectors */
282:   VecCreateSeq(PETSC_COMM_SELF,len,&tmp_in);
283:   VecDuplicate(tmp_in,&tmp_out);

285:   /* Create scatter context */
286:   ISCreateGeneral(PETSC_COMM_SELF,len,idx,&from);
287:   ISCreateGeneral(PETSC_COMM_SELF,len,idx,&to);
288:   VecScatterCreate(globalin,from,tmp_in,to,&scatter);
289:   VecScatterBegin(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
290: 

292: 
293:   VecScatterEnd(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
294: 

296:   /*Extract income array */
297:   VecGetArray(tmp_in,&inptr);

299:   /* Extract outcome array*/
300:   VecGetArray(tmp_out,&outptr);

302:   outptr[0] = xc*inptr[0]+xr*inptr[1]+yr*inptr[m];
303:   outptr[m-1] = 2.0*xl*inptr[m-2]+xc*inptr[m-1]+yr*inptr[m-1+m];
304:   for (i=1; i<m-1; i++) {
305:     outptr[i] = xc*inptr[i]+xl*inptr[i-1]+xr*inptr[i+1]
306:       +yr*inptr[i+m];
307:   }

309:   for (j=1; j<n-1; j++) {
310:     outptr[j*m] = xc*inptr[j*m]+xr*inptr[j*m+1]+
311:       yl*inptr[j*m-m]+yr*inptr[j*m+m];
312:     outptr[j*m+m-1] = xc*inptr[j*m+m-1]+2.0*xl*inptr[j*m+m-1-1]+
313:       yl*inptr[j*m+m-1-m]+yr*inptr[j*m+m-1+m];
314:     for(i=1; i<m-1; i++) {
315:       outptr[j*m+i] = xc*inptr[j*m+i]+xl*inptr[j*m+i-1]+xr*inptr[j*m+i+1]
316:         +yl*inptr[j*m+i-m]+yr*inptr[j*m+i+m];
317:     }
318:   }

320:   outptr[mn-m] = xc*inptr[mn-m]+xr*inptr[mn-m+1]+2.0*yl*inptr[mn-m-m];
321:   outptr[mn-1] = 2.0*xl*inptr[mn-2]+xc*inptr[mn-1]+2.0*yl*inptr[mn-1-m];
322:   for (i=1; i<m-1; i++) {
323:     outptr[mn-m+i] = xc*inptr[mn-m+i]+xl*inptr[mn-m+i-1]+xr*inptr[mn-m+i+1]
324:       +2*yl*inptr[mn-m+i-m];
325:   }

327:   VecRestoreArray(tmp_in,&inptr);
328:   VecRestoreArray(tmp_out,&outptr);

330:   VecScatterCreate(tmp_out,from,globalout,to,&scatter);
331:   VecScatterBegin(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter
332: );
333: 
334:   VecScatterEnd(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter);
335: 
336: 
337:   /* Destroy idx aand scatter */
338:   ISDestroy(from);
339:   ISDestroy(to);
340:   VecScatterDestroy(scatter);

342:   PetscFree(idx);
343:   return 0;
344: }

346: int FormJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *flag,void *ptr)
347: {
348:   Data *data = (Data*)ptr;
349:   Mat A = *AA;
350:   PetscScalar v[1],one = 1.0;
351:   int idx[1],i,j,row,ierr;
352:   int m,n,mn;

354:   m = data->m;
355:   n = data->n;
356:   mn = (data->m)*(data->n);
357: 
358:   for(i=0; i<mn; i++) {
359:     idx[0] = i;
360:     v[0] = one;
361:     MatSetValues(A,1,&i,1,idx,v,INSERT_VALUES);
362:   }

364:   for(i=0; i<mn-m; i++) {
365:     idx[0] = i+m;
366:     v[0]= one;
367:     MatSetValues(A,1,&i,1,idx,v,INSERT_VALUES);
368:   }

370:   for(i=m; i<mn; i++) {
371:     idx[0] = i-m;
372:     v[0] = one;
373:     MatSetValues(A,1,&i,1,idx,v,INSERT_VALUES);
374:   }

376:   for(j=0; j<n; j++) {
377:     for(i=0; i<m-1; i++) {
378:       row = j*m+i;
379:       idx[0] = j*m+i+1;
380:       v[0] = one;
381:       MatSetValues(A,1,&row,1,idx,v,INSERT_VALUES);
382:     }
383:   }

385:   for(j=0; j<n; j++) {
386:     for(i=1; i<m; i++) {
387:       row = j*m+i;
388:       idx[0] = j*m+i-1;
389:       v[0] = one;
390:       MatSetValues(A,1,&row,1,idx,v,INSERT_VALUES);
391:     }
392:   }
393: 
394:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
395:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);


398:   *flag = SAME_NONZERO_PATTERN;
399:   return 0;
400: }

402: int RHSJacobian(TS ts,PetscReal t,Vec x,Mat *AA,Mat *BB,MatStructure *flag,void *ptr)
403: {
404:   Data        *data = (Data*)ptr;
405:   Mat         A = *AA;
406:   PetscScalar v[5];
407:   int         idx[5],i,j,row,ierr;
408:   int         m,n,mn;
409:   PetscReal   dx,dy,a,epsilon,xc,xl,xr,yl,yr;

411:   m = data->m;
412:   n = data->n;
413:   mn = m*n;
414:   dx = data->dx;
415:   dy = data->dy;
416:   a = data->a;
417:   epsilon = data->epsilon;

419:   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
420:   xl = 0.5*a/dx+epsilon/(dx*dx);
421:   xr = -0.5*a/dx+epsilon/(dx*dx);
422:   yl = 0.5*a/dy+epsilon/(dy*dy);
423:   yr = -0.5*a/dy+epsilon/(dy*dy);

425:   row=0;
426:   v[0] = xc; v[1]=xr; v[2]=yr;
427:   idx[0]=0; idx[1]=2; idx[2]=m;
428:   MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);

430:   row=m-1;
431:   v[0]=2.0*xl; v[1]=xc; v[2]=yr;
432:   idx[0]=m-2; idx[1]=m-1; idx[2]=m-1+m;
433:   MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);

435:   for (i=1; i<m-1; i++) {
436:     row=i;
437:     v[0]=xl; v[1]=xc; v[2]=xr; v[3]=yr;
438:     idx[0]=i-1; idx[1]=i; idx[2]=i+1; idx[3]=i+m;
439:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
440:   }

442:   for (j=1; j<n-1; j++) {
443:     row=j*m;
444:     v[0]=xc; v[1]=xr; v[2]=yl; v[3]=yr;
445:     idx[0]=j*m; idx[1]=j*m; idx[2]=j*m-m; idx[3]=j*m+m;
446:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
447: 
448:     row=j*m+m-1;
449:     v[0]=xc; v[1]=2.0*xl; v[2]=yl; v[3]=yr;
450:     idx[0]=j*m+m-1; idx[1]=j*m+m-1-1; idx[2]=j*m+m-1-m; idx[3]=j*m+m-1+m;
451:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);

453:     for(i=1; i<m-1; i++) {
454:       row=j*m+i;
455:       v[0]=xc; v[1]=xl; v[2]=xr; v[3]=yl; v[4]=yr;
456:       idx[0]=j*m+i; idx[1]=j*m+i-1; idx[2]=j*m+i+1; idx[3]=j*m+i-m;
457:       idx[4]=j*m+i+m;
458:       MatSetValues(A,1,&row,5,idx,v,INSERT_VALUES);
459:     }
460:   }

462:   row=mn-m;
463:   v[0] = xc; v[1]=xr; v[2]=2.0*yl;
464:   idx[0]=mn-m; idx[1]=mn-m+1; idx[2]=mn-m-m;
465:   MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
466: 
467:   row=mn-1;
468:   v[0] = xc; v[1]=2.0*xl; v[2]=2.0*yl;
469:   idx[0]=mn-1; idx[1]=mn-2; idx[2]=mn-1-m;
470:   MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);

472:   for (i=1; i<m-1; i++) {
473:     row=mn-m+i;
474:     v[0]=xl; v[1]=xc; v[2]=xr; v[3]=2.0*yl;
475:     idx[0]=mn-m+i-1; idx[1]=mn-m+i; idx[2]=mn-m+i+1; idx[3]=mn-m+i-m;
476:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
477:   }


480:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
481:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);


484:   *flag = SAME_NONZERO_PATTERN;
485:   return 0;
486: }

488: int RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
489: {
491:   SNES snes = PETSC_NULL;

493:   FormFunction(snes,globalin,globalout,ctx);


496:   return 0;
497: }