Actual source code: ex4.c

  1: /*
  2:        The Problem:
  3:            Solve the convection-diffusion equation:
  4:            
  5:              u_t+a*(u_x+u_y)=epsilon*(u_xx+u_yy)
  6:              u=0 at x=0, y=0
  7:              u_x=0 at x=1
  8:              u_y=0 at y=1
  9:         
 10:        This program tests the routine of computing the Jacobian by the 
 11:        finite difference method as well as PETSc with PVODE.

 13: */

 15: static char help[] = "Solve the convection-diffusion equation. \n\n";

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


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

 33: /* two temporal functions */

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


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

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

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

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

104:   dt = 0.1;

106:   /*
107:     The user provides the RHS and Jacobian
108:   */
109:   MatCreate(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,mn,mn,&A);
110:   MatSetFromOptions(A);

112:   /* Create SNES context  */
113:   SNESCreate(PETSC_COMM_WORLD,&snes);
114: 

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

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

132:   /* Use PVODE */
133:   TSSetType(ts,TS_PVODE);

135:   TSSetFromOptions(ts);

137:   TSSetInitialTimeStep(ts,0.0,dt);
138:   TSSetDuration(ts,time_steps,1);
139:   TSSetSolution(ts,global);


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

151:   TSSetUp(ts);
152:   TSStep(ts,&steps,&ftime);

154:   TSGetSolution(ts,&global);
155:   PetscViewerASCIIOpen(PETSC_COMM_SELF,"out.m",&viewfile);
156:   PetscViewerSetFormat(viewfile,PETSC_VIEWER_ASCII_MATLAB);
157:   VecView(global,viewfile);

159: #if defined(PETSC_HAVE_PVODE) && !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
160:   /* extracts the PC  from ts */
161:   TSPVodeGetPC(ts,&pc);
162:   PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,120,&viewer);
163:   TSView(ts,viewer);
164:   PetscViewerStringOpen(PETSC_COMM_WORLD,pcinfo,120,&viewer);
165:   PCView(pc,viewer);
166:   PetscPrintf(PETSC_COMM_WORLD,"%d Procs,%s Preconditioner,%s\n",
167:                      size,tsinfo,pcinfo);
168:   PCDestroy(pc);
169: #endif

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

179: /* -------------------------------------------------------------------*/
182: PetscErrorCode Initial(Vec global,void *ctx)
183: {
184:   Data           *data = (Data*)ctx;
185:   PetscInt       m,row,col;
186:   PetscReal      x,y,dx,dy;
187:   PetscScalar    *localptr;
188:   PetscInt       i,mybase,myend,locsize;

191:   /* make the local  copies of parameters */
192:   m = data->m;
193:   dx = data->dx;
194:   dy = data->dy;

196:   /* determine starting point of each processor */
197:   VecGetOwnershipRange(global,&mybase,&myend);
198:   VecGetLocalSize(global,&locsize);

200:   /* Initialize the array */
201:   VecGetArray(global,&localptr);

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

217: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
218: {
219:   VecScatter     scatter;
220:   IS             from,to;
221:   PetscInt       i,n,*idx;
222:   Vec            tmp_vec;
224:   PetscScalar    *tmp;

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

229:   /* Set the index sets */
230:   PetscMalloc(n*sizeof(PetscInt),&idx);
231:   for(i=0; i<n; i++) idx[i]=i;
232: 
233:   /* Create local sequential vectors */
234:   VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);

236:   /* Create scatter context */
237:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&from);
238:   ISCreateGeneral(PETSC_COMM_SELF,n,idx,&to);
239:   VecScatterCreate(global,from,tmp_vec,to,&scatter);
240:   VecScatterBegin(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);
241:   VecScatterEnd(global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD,scatter);

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

247:   PetscFree(idx);
248:   return 0;
249: }


254: PetscErrorCode FormFunction(SNES snes,Vec globalin,Vec globalout,void *ptr)
255: {
256:   Data           *data = (Data*)ptr;
257:   PetscInt       m,n,mn;
258:   PetscReal      dx,dy;
259:   PetscReal      xc,xl,xr,yl,yr;
260:   PetscReal      a,epsilon;
261:   PetscScalar    *inptr,*outptr;
262:   PetscInt       i,j,len;
264:   IS             from,to;
265:   PetscInt       *idx;
266:   VecScatter     scatter;
267:   Vec            tmp_in,tmp_out;

269:   m = data->m;
270:   n = data->n;
271:   mn = m*n;
272:   dx = data->dx;
273:   dy = data->dy;
274:   a = data->a;
275:   epsilon = data->epsilon;

277:   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
278:   xl = 0.5*a/dx+epsilon/(dx*dx);
279:   xr = -0.5*a/dx+epsilon/(dx*dx);
280:   yl = 0.5*a/dy+epsilon/(dy*dy);
281:   yr = -0.5*a/dy+epsilon/(dy*dy);
282: 
283:   /* Get the length of parallel vector */
284:   VecGetSize(globalin,&len);

286:   /* Set the index sets */
287:   PetscMalloc(len*sizeof(PetscInt),&idx);
288:   for(i=0; i<len; i++) idx[i]=i;
289: 
290:   /* Create local sequential vectors */
291:   VecCreateSeq(PETSC_COMM_SELF,len,&tmp_in);
292:   VecDuplicate(tmp_in,&tmp_out);

294:   /* Create scatter context */
295:   ISCreateGeneral(PETSC_COMM_SELF,len,idx,&from);
296:   ISCreateGeneral(PETSC_COMM_SELF,len,idx,&to);
297:   VecScatterCreate(globalin,from,tmp_in,to,&scatter);
298:   VecScatterBegin(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
299: 

301: 
302:   VecScatterEnd(globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD,scatter);
303: 

305:   /*Extract income array */
306:   VecGetArray(tmp_in,&inptr);

308:   /* Extract outcome array*/
309:   VecGetArray(tmp_out,&outptr);

311:   outptr[0] = xc*inptr[0]+xr*inptr[1]+yr*inptr[m];
312:   outptr[m-1] = 2.0*xl*inptr[m-2]+xc*inptr[m-1]+yr*inptr[m-1+m];
313:   for (i=1; i<m-1; i++) {
314:     outptr[i] = xc*inptr[i]+xl*inptr[i-1]+xr*inptr[i+1]
315:       +yr*inptr[i+m];
316:   }

318:   for (j=1; j<n-1; j++) {
319:     outptr[j*m] = xc*inptr[j*m]+xr*inptr[j*m+1]+
320:       yl*inptr[j*m-m]+yr*inptr[j*m+m];
321:     outptr[j*m+m-1] = xc*inptr[j*m+m-1]+2.0*xl*inptr[j*m+m-1-1]+
322:       yl*inptr[j*m+m-1-m]+yr*inptr[j*m+m-1+m];
323:     for(i=1; i<m-1; i++) {
324:       outptr[j*m+i] = xc*inptr[j*m+i]+xl*inptr[j*m+i-1]+xr*inptr[j*m+i+1]
325:         +yl*inptr[j*m+i-m]+yr*inptr[j*m+i+m];
326:     }
327:   }

329:   outptr[mn-m] = xc*inptr[mn-m]+xr*inptr[mn-m+1]+2.0*yl*inptr[mn-m-m];
330:   outptr[mn-1] = 2.0*xl*inptr[mn-2]+xc*inptr[mn-1]+2.0*yl*inptr[mn-1-m];
331:   for (i=1; i<m-1; i++) {
332:     outptr[mn-m+i] = xc*inptr[mn-m+i]+xl*inptr[mn-m+i-1]+xr*inptr[mn-m+i+1]
333:       +2*yl*inptr[mn-m+i-m];
334:   }

336:   VecRestoreArray(tmp_in,&inptr);
337:   VecRestoreArray(tmp_out,&outptr);

339:   VecScatterCreate(tmp_out,from,globalout,to,&scatter);
340:   VecScatterBegin(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter
341: );
342: 
343:   VecScatterEnd(tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD,scatter);
344: 
345: 
346:   /* Destroy idx aand scatter */
347:   ISDestroy(from);
348:   ISDestroy(to);
349:   VecScatterDestroy(scatter);

351:   PetscFree(idx);
352:   return 0;
353: }

357: PetscErrorCode FormJacobian(SNES snes,Vec x,Mat *AA,Mat *BB,MatStructure *flag,void *ptr)
358: {
359:   Data           *data = (Data*)ptr;
360:   Mat            A = *AA;
361:   PetscScalar    v[1],one = 1.0;
362:   PetscInt       idx[1],i,j,row;
364:   PetscInt       m,n,mn;

366:   m = data->m;
367:   n = data->n;
368:   mn = (data->m)*(data->n);
369: 
370:   for(i=0; i<mn; i++) {
371:     idx[0] = i;
372:     v[0] = one;
373:     MatSetValues(A,1,&i,1,idx,v,INSERT_VALUES);
374:   }

376:   for(i=0; i<mn-m; i++) {
377:     idx[0] = i+m;
378:     v[0]= one;
379:     MatSetValues(A,1,&i,1,idx,v,INSERT_VALUES);
380:   }

382:   for(i=m; i<mn; i++) {
383:     idx[0] = i-m;
384:     v[0] = one;
385:     MatSetValues(A,1,&i,1,idx,v,INSERT_VALUES);
386:   }

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

397:   for(j=0; j<n; j++) {
398:     for(i=1; i<m; i++) {
399:       row = j*m+i;
400:       idx[0] = j*m+i-1;
401:       v[0] = one;
402:       MatSetValues(A,1,&row,1,idx,v,INSERT_VALUES);
403:     }
404:   }
405: 
406:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
407:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);


410:   *flag = SAME_NONZERO_PATTERN;
411:   return 0;
412: }

416: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat *AA,Mat *BB,MatStructure *flag,void *ptr)
417: {
418:   Data           *data = (Data*)ptr;
419:   Mat            A = *AA;
420:   PetscScalar    v[5];
421:   PetscInt       idx[5],i,j,row;
423:   PetscInt       m,n,mn;
424:   PetscReal      dx,dy,a,epsilon,xc,xl,xr,yl,yr;

426:   m = data->m;
427:   n = data->n;
428:   mn = m*n;
429:   dx = data->dx;
430:   dy = data->dy;
431:   a = data->a;
432:   epsilon = data->epsilon;

434:   xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
435:   xl = 0.5*a/dx+epsilon/(dx*dx);
436:   xr = -0.5*a/dx+epsilon/(dx*dx);
437:   yl = 0.5*a/dy+epsilon/(dy*dy);
438:   yr = -0.5*a/dy+epsilon/(dy*dy);

440:   row=0;
441:   v[0] = xc; v[1]=xr; v[2]=yr;
442:   idx[0]=0; idx[1]=2; idx[2]=m;
443:   MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);

445:   row=m-1;
446:   v[0]=2.0*xl; v[1]=xc; v[2]=yr;
447:   idx[0]=m-2; idx[1]=m-1; idx[2]=m-1+m;
448:   MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);

450:   for (i=1; i<m-1; i++) {
451:     row=i;
452:     v[0]=xl; v[1]=xc; v[2]=xr; v[3]=yr;
453:     idx[0]=i-1; idx[1]=i; idx[2]=i+1; idx[3]=i+m;
454:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
455:   }

457:   for (j=1; j<n-1; j++) {
458:     row=j*m;
459:     v[0]=xc; v[1]=xr; v[2]=yl; v[3]=yr;
460:     idx[0]=j*m; idx[1]=j*m; idx[2]=j*m-m; idx[3]=j*m+m;
461:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
462: 
463:     row=j*m+m-1;
464:     v[0]=xc; v[1]=2.0*xl; v[2]=yl; v[3]=yr;
465:     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;
466:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);

468:     for(i=1; i<m-1; i++) {
469:       row=j*m+i;
470:       v[0]=xc; v[1]=xl; v[2]=xr; v[3]=yl; v[4]=yr;
471:       idx[0]=j*m+i; idx[1]=j*m+i-1; idx[2]=j*m+i+1; idx[3]=j*m+i-m;
472:       idx[4]=j*m+i+m;
473:       MatSetValues(A,1,&row,5,idx,v,INSERT_VALUES);
474:     }
475:   }

477:   row=mn-m;
478:   v[0] = xc; v[1]=xr; v[2]=2.0*yl;
479:   idx[0]=mn-m; idx[1]=mn-m+1; idx[2]=mn-m-m;
480:   MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
481: 
482:   row=mn-1;
483:   v[0] = xc; v[1]=2.0*xl; v[2]=2.0*yl;
484:   idx[0]=mn-1; idx[1]=mn-2; idx[2]=mn-1-m;
485:   MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);

487:   for (i=1; i<m-1; i++) {
488:     row=mn-m+i;
489:     v[0]=xl; v[1]=xc; v[2]=xr; v[3]=2.0*yl;
490:     idx[0]=mn-m+i-1; idx[1]=mn-m+i; idx[2]=mn-m+i+1; idx[3]=mn-m+i-m;
491:     MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
492:   }


495:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
496:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);


499:   *flag = SAME_NONZERO_PATTERN;
500:   return 0;
501: }

505: PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
506: {
508:   SNES snes = PETSC_NULL;

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


513:   return 0;
514: }