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)
 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,&A);
110:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,mn,mn);
111:   MatSetFromOptions(A);

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

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

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

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

136:   TSSetFromOptions(ts);

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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


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


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

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

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


514:   return 0;
515: }