Actual source code: ex4.c
petsc-3.7.3 2016-08-01
1: /*
2: The Problem:
3: Solve the convection-diffusion equation:
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: u = exp(-20.0*(pow(x-0.5,2.0)+pow(y-0.5,2.0))) at t=0
11: This program tests the routine of computing the Jacobian by the
12: finite difference method as well as PETSc with SUNDIALS.
14: */
16: static char help[] = "Solve the convection-diffusion equation. \n\n";
18: #include <petscts.h>
20: typedef struct
21: {
22: PetscInt m; /* the number of mesh points in x-direction */
23: PetscInt n; /* the number of mesh points in y-direction */
24: PetscReal dx; /* the grid space in x-direction */
25: PetscReal dy; /* the grid space in y-direction */
26: PetscReal a; /* the convection coefficient */
27: PetscReal epsilon; /* the diffusion coefficient */
28: PetscReal tfinal;
29: } Data;
31: extern PetscErrorCode Monitor(TS,PetscInt,PetscReal,Vec,void*);
32: extern PetscErrorCode Initial(Vec,void*);
33: extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
34: extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat,Mat,void*);
35: extern PetscErrorCode PostStep(TS);
39: int main(int argc,char **argv)
40: {
42: PetscInt time_steps=100,iout,NOUT=1;
43: PetscMPIInt size;
44: Vec global;
45: PetscReal dt,ftime,ftime_original;
46: TS ts;
47: PetscViewer viewfile;
48: Mat J = 0;
49: Vec x;
50: Data data;
51: PetscInt mn;
52: PetscBool flg;
53: MatColoring mc;
54: ISColoring iscoloring;
55: MatFDColoring matfdcoloring = 0;
56: PetscBool fd_jacobian_coloring = PETSC_FALSE;
57: SNES snes;
58: KSP ksp;
59: PC pc;
60: PetscViewer viewer;
61: char pcinfo[120],tsinfo[120];
62: TSType tstype;
63: PetscBool sundials;
65: PetscInitialize(&argc,&argv,(char*)0,help);
66: MPI_Comm_size(PETSC_COMM_WORLD,&size);
68: /* set data */
69: data.m = 9;
70: data.n = 9;
71: data.a = 1.0;
72: data.epsilon = 0.1;
73: data.dx = 1.0/(data.m+1.0);
74: data.dy = 1.0/(data.n+1.0);
75: mn = (data.m)*(data.n);
76: PetscOptionsGetInt(NULL,NULL,"-time",&time_steps,NULL);
78: /* set initial conditions */
79: VecCreate(PETSC_COMM_WORLD,&global);
80: VecSetSizes(global,PETSC_DECIDE,mn);
81: VecSetFromOptions(global);
82: Initial(global,&data);
83: VecDuplicate(global,&x);
85: /* create timestep context */
86: TSCreate(PETSC_COMM_WORLD,&ts);
87: TSMonitorSet(ts,Monitor,&data,NULL);
88: #if defined(PETSC_HAVE_SUNDIALS)
89: TSSetType(ts,TSSUNDIALS);
90: #else
91: TSSetType(ts,TSEULER);
92: #endif
93: dt = 0.1;
94: ftime_original = data.tfinal = 1.0;
96: TSSetInitialTimeStep(ts,0.0,dt);
97: TSSetDuration(ts,time_steps,ftime_original);
98: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
99: TSSetSolution(ts,global);
101: /* set user provided RHSFunction and RHSJacobian */
102: TSSetRHSFunction(ts,NULL,RHSFunction,&data);
103: MatCreate(PETSC_COMM_WORLD,&J);
104: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,mn,mn);
105: MatSetFromOptions(J);
106: MatSeqAIJSetPreallocation(J,5,NULL);
107: MatMPIAIJSetPreallocation(J,5,NULL,5,NULL);
109: PetscOptionsHasName(NULL,NULL,"-ts_fd",&flg);
110: if (!flg) {
111: TSSetRHSJacobian(ts,J,J,RHSJacobian,&data);
112: } else {
113: TSGetSNES(ts,&snes);
114: PetscOptionsHasName(NULL,NULL,"-fd_color",&fd_jacobian_coloring);
115: if (fd_jacobian_coloring) { /* Use finite differences with coloring */
116: /* Get data structure of J */
117: PetscBool pc_diagonal;
118: PetscOptionsHasName(NULL,NULL,"-pc_diagonal",&pc_diagonal);
119: if (pc_diagonal) { /* the preconditioner of J is a diagonal matrix */
120: PetscInt rstart,rend,i;
121: PetscScalar zero=0.0;
122: MatGetOwnershipRange(J,&rstart,&rend);
123: for (i=rstart; i<rend; i++) {
124: MatSetValues(J,1,&i,1,&i,&zero,INSERT_VALUES);
125: }
126: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
127: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
128: } else {
129: /* Fill the structure using the expensive SNESComputeJacobianDefault. Temporarily set up the TS so we can call this function */
130: TSSetType(ts,TSBEULER);
131: TSSetUp(ts);
132: SNESComputeJacobianDefault(snes,x,J,J,ts);
133: }
135: /* create coloring context */
136: MatColoringCreate(J,&mc);
137: MatColoringSetType(mc,MATCOLORINGSL);
138: MatColoringSetFromOptions(mc);
139: MatColoringApply(mc,&iscoloring);
140: MatColoringDestroy(&mc);
141: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
142: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
143: MatFDColoringSetFromOptions(matfdcoloring);
144: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
145: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
146: ISColoringDestroy(&iscoloring);
147: } else { /* Use finite differences (slow) */
148: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);
149: }
150: }
152: /* Pick up a Petsc preconditioner */
153: /* one can always set method or preconditioner during the run time */
154: TSGetSNES(ts,&snes);
155: SNESGetKSP(snes,&ksp);
156: KSPGetPC(ksp,&pc);
157: PCSetType(pc,PCJACOBI);
158: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
160: TSSetFromOptions(ts);
161: TSSetUp(ts);
163: /* Test TSSetPostStep() */
164: PetscOptionsHasName(NULL,NULL,"-test_PostStep",&flg);
165: if (flg) {
166: TSSetPostStep(ts,PostStep);
167: }
169: PetscOptionsGetInt(NULL,NULL,"-NOUT",&NOUT,NULL);
170: for (iout=1; iout<=NOUT; iout++) {
171: TSSetDuration(ts,time_steps,iout*ftime_original/NOUT);
172: TSSolve(ts,global);
173: TSGetSolveTime(ts,&ftime);
174: TSSetInitialTimeStep(ts,ftime,dt);
175: }
176: /* Interpolate solution at tfinal */
177: TSGetSolution(ts,&global);
178: TSInterpolate(ts,ftime_original,global);
180: PetscOptionsHasName(NULL,NULL,"-matlab_view",&flg);
181: if (flg) { /* print solution into a MATLAB file */
182: PetscViewerASCIIOpen(PETSC_COMM_WORLD,"out.m",&viewfile);
183: PetscViewerPushFormat(viewfile,PETSC_VIEWER_ASCII_MATLAB);
184: VecView(global,viewfile);
185: PetscViewerPopFormat(viewfile);
186: PetscViewerDestroy(&viewfile);
187: }
189: /* display solver info for Sundials */
190: TSGetType(ts,&tstype);
191: PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&sundials);
192: if (sundials) {
193: PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,120,&viewer);
194: TSView(ts,viewer);
195: PetscViewerDestroy(&viewer);
196: PetscViewerStringOpen(PETSC_COMM_WORLD,pcinfo,120,&viewer);
197: PCView(pc,viewer);
198: PetscPrintf(PETSC_COMM_WORLD,"%d Procs,%s TSType, %s Preconditioner\n",size,tsinfo,pcinfo);
199: PetscViewerDestroy(&viewer);
200: }
202: /* free the memories */
203: TSDestroy(&ts);
204: VecDestroy(&global);
205: VecDestroy(&x);
206: MatDestroy(&J);
207: if (fd_jacobian_coloring) {MatFDColoringDestroy(&matfdcoloring);}
208: PetscFinalize();
209: return 0;
210: }
212: /* -------------------------------------------------------------------*/
213: /* the initial function */
214: PetscReal f_ini(PetscReal x,PetscReal y)
215: {
216: PetscReal f;
218: f=PetscExpReal(-20.0*(PetscPowRealInt(x-0.5,2)+PetscPowRealInt(y-0.5,2)));
219: return f;
220: }
224: PetscErrorCode Initial(Vec global,void *ctx)
225: {
226: Data *data = (Data*)ctx;
227: PetscInt m,row,col;
228: PetscReal x,y,dx,dy;
229: PetscScalar *localptr;
230: PetscInt i,mybase,myend,locsize;
234: /* make the local copies of parameters */
235: m = data->m;
236: dx = data->dx;
237: dy = data->dy;
239: /* determine starting point of each processor */
240: VecGetOwnershipRange(global,&mybase,&myend);
241: VecGetLocalSize(global,&locsize);
243: /* Initialize the array */
244: VecGetArray(global,&localptr);
246: for (i=0; i<locsize; i++) {
247: row = 1+(mybase+i)-((mybase+i)/m)*m;
248: col = (mybase+i)/m+1;
249: x = dx*row;
250: y = dy*col;
251: localptr[i] = f_ini(x,y);
252: }
254: VecRestoreArray(global,&localptr);
255: return(0);
256: }
260: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
261: {
262: VecScatter scatter;
263: IS from,to;
264: PetscInt i,n,*idx,nsteps,maxsteps;
265: Vec tmp_vec;
266: PetscErrorCode ierr;
267: const PetscScalar *tmp;
268: PetscReal maxtime;
271: TSGetTimeStepNumber(ts,&nsteps);
272: /* display output at selected time steps */
273: TSGetDuration(ts, &maxsteps, &maxtime);
274: if (nsteps % 10 != 0) return(0);
276: /* Get the size of the vector */
277: VecGetSize(global,&n);
279: /* Set the index sets */
280: PetscMalloc1(n,&idx);
281: for (i=0; i<n; i++) idx[i]=i;
283: /* Create local sequential vectors */
284: VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);
286: /* Create scatter context */
287: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);
288: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);
289: VecScatterCreate(global,from,tmp_vec,to,&scatter);
290: VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);
291: VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);
293: VecGetArrayRead(tmp_vec,&tmp);
294: PetscPrintf(PETSC_COMM_WORLD,"At t[%D] =%14.2e u= %14.2e at the center \n",nsteps,(double)time,(double)PetscRealPart(tmp[n/2]));
295: VecRestoreArrayRead(tmp_vec,&tmp);
297: PetscFree(idx);
298: ISDestroy(&from);
299: ISDestroy(&to);
300: VecScatterDestroy(&scatter);
301: VecDestroy(&tmp_vec);
302: return(0);
303: }
307: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ptr)
308: {
309: Data *data = (Data*)ptr;
310: PetscScalar v[5];
311: PetscInt idx[5],i,j,row;
313: PetscInt m,n,mn;
314: PetscReal dx,dy,a,epsilon,xc,xl,xr,yl,yr;
317: m = data->m;
318: n = data->n;
319: mn = m*n;
320: dx = data->dx;
321: dy = data->dy;
322: a = data->a;
323: epsilon = data->epsilon;
325: xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
326: xl = 0.5*a/dx+epsilon/(dx*dx);
327: xr = -0.5*a/dx+epsilon/(dx*dx);
328: yl = 0.5*a/dy+epsilon/(dy*dy);
329: yr = -0.5*a/dy+epsilon/(dy*dy);
331: row = 0;
332: v[0] = xc; v[1] = xr; v[2] = yr;
333: idx[0] = 0; idx[1] = 2; idx[2] = m;
334: MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
336: row = m-1;
337: v[0] = 2.0*xl; v[1] = xc; v[2] = yr;
338: idx[0] = m-2; idx[1] = m-1; idx[2] = m-1+m;
339: MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
341: for (i=1; i<m-1; i++) {
342: row = i;
343: v[0] = xl; v[1] = xc; v[2] = xr; v[3] = yr;
344: idx[0] = i-1; idx[1] = i; idx[2] = i+1; idx[3] = i+m;
345: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
346: }
348: for (j=1; j<n-1; j++) {
349: row = j*m;
350: v[0] = xc; v[1] = xr; v[2] = yl; v[3] = yr;
351: idx[0] = j*m; idx[1] = j*m; idx[2] = j*m-m; idx[3] = j*m+m;
352: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
354: row = j*m+m-1;
355: v[0] = xc; v[1] = 2.0*xl; v[2] = yl; v[3] = yr;
356: 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;
357: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
359: for (i=1; i<m-1; i++) {
360: row = j*m+i;
361: v[0] = xc; v[1] = xl; v[2] = xr; v[3] = yl; v[4]=yr;
362: idx[0] = j*m+i; idx[1] = j*m+i-1; idx[2] = j*m+i+1; idx[3] = j*m+i-m;
363: idx[4] = j*m+i+m;
364: MatSetValues(A,1,&row,5,idx,v,INSERT_VALUES);
365: }
366: }
368: row = mn-m;
369: v[0] = xc; v[1] = xr; v[2] = 2.0*yl;
370: idx[0] = mn-m; idx[1] = mn-m+1; idx[2] = mn-m-m;
371: MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
373: row = mn-1;
374: v[0] = xc; v[1] = 2.0*xl; v[2] = 2.0*yl;
375: idx[0] = mn-1; idx[1] = mn-2; idx[2] = mn-1-m;
376: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
378: for (i=1; i<m-1; i++) {
379: row = mn-m+i;
380: v[0] = xl; v[1] = xc; v[2] = xr; v[3] = 2.0*yl;
381: idx[0] = mn-m+i-1; idx[1] = mn-m+i; idx[2] = mn-m+i+1; idx[3] = mn-m+i-m;
382: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
383: }
385: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
386: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
388: return(0);
389: }
391: /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
394: PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
395: {
396: Data *data = (Data*)ctx;
397: PetscInt m,n,mn;
398: PetscReal dx,dy;
399: PetscReal xc,xl,xr,yl,yr;
400: PetscReal a,epsilon;
401: PetscScalar *outptr;
402: const PetscScalar *inptr;
403: PetscInt i,j,len;
404: PetscErrorCode ierr;
405: IS from,to;
406: PetscInt *idx;
407: VecScatter scatter;
408: Vec tmp_in,tmp_out;
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: /* Get the length of parallel vector */
426: VecGetSize(globalin,&len);
428: /* Set the index sets */
429: PetscMalloc1(len,&idx);
430: for (i=0; i<len; i++) idx[i]=i;
432: /* Create local sequential vectors */
433: VecCreateSeq(PETSC_COMM_SELF,len,&tmp_in);
434: VecDuplicate(tmp_in,&tmp_out);
436: /* Create scatter context */
437: ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&from);
438: ISCreateGeneral(PETSC_COMM_SELF,len,idx,PETSC_COPY_VALUES,&to);
439: VecScatterCreate(globalin,from,tmp_in,to,&scatter);
440: VecScatterBegin(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
441: VecScatterEnd(scatter,globalin,tmp_in,INSERT_VALUES,SCATTER_FORWARD);
442: VecScatterDestroy(&scatter);
444: /*Extract income array - include ghost points */
445: VecGetArrayRead(tmp_in,&inptr);
447: /* Extract outcome array*/
448: VecGetArray(tmp_out,&outptr);
450: outptr[0] = xc*inptr[0]+xr*inptr[1]+yr*inptr[m];
451: outptr[m-1] = 2.0*xl*inptr[m-2]+xc*inptr[m-1]+yr*inptr[m-1+m];
452: for (i=1; i<m-1; i++) {
453: outptr[i] = xc*inptr[i]+xl*inptr[i-1]+xr*inptr[i+1]+yr*inptr[i+m];
454: }
456: for (j=1; j<n-1; j++) {
457: outptr[j*m] = xc*inptr[j*m]+xr*inptr[j*m+1]+ yl*inptr[j*m-m]+yr*inptr[j*m+m];
458: outptr[j*m+m-1] = xc*inptr[j*m+m-1]+2.0*xl*inptr[j*m+m-1-1]+ yl*inptr[j*m+m-1-m]+yr*inptr[j*m+m-1+m];
459: for (i=1; i<m-1; i++) {
460: outptr[j*m+i] = xc*inptr[j*m+i]+xl*inptr[j*m+i-1]+xr*inptr[j*m+i+1]+yl*inptr[j*m+i-m]+yr*inptr[j*m+i+m];
461: }
462: }
464: outptr[mn-m] = xc*inptr[mn-m]+xr*inptr[mn-m+1]+2.0*yl*inptr[mn-m-m];
465: outptr[mn-1] = 2.0*xl*inptr[mn-2]+xc*inptr[mn-1]+2.0*yl*inptr[mn-1-m];
466: for (i=1; i<m-1; i++) {
467: outptr[mn-m+i] = xc*inptr[mn-m+i]+xl*inptr[mn-m+i-1]+xr*inptr[mn-m+i+1]+2*yl*inptr[mn-m+i-m];
468: }
470: VecRestoreArrayRead(tmp_in,&inptr);
471: VecRestoreArray(tmp_out,&outptr);
473: VecScatterCreate(tmp_out,from,globalout,to,&scatter);
474: VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
475: VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
477: /* Destroy idx aand scatter */
478: VecDestroy(&tmp_in);
479: VecDestroy(&tmp_out);
480: ISDestroy(&from);
481: ISDestroy(&to);
482: VecScatterDestroy(&scatter);
484: PetscFree(idx);
485: return(0);
486: }
490: PetscErrorCode PostStep(TS ts)
491: {
493: PetscReal t;
496: TSGetTime(ts,&t);
497: PetscPrintf(PETSC_COMM_SELF," PostStep, t: %g\n",(double)t);
498: return(0);
499: }