Actual source code: ex4.c
petsc-3.5.4 2015-05-23
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,"-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: TSSetSolution(ts,global);
100: /* set user provided RHSFunction and RHSJacobian */
101: TSSetRHSFunction(ts,NULL,RHSFunction,&data);
102: MatCreate(PETSC_COMM_WORLD,&J);
103: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,mn,mn);
104: MatSetFromOptions(J);
105: MatSeqAIJSetPreallocation(J,5,NULL);
106: MatMPIAIJSetPreallocation(J,5,NULL,5,NULL);
108: PetscOptionsHasName(NULL,"-ts_fd",&flg);
109: if (!flg) {
110: TSSetRHSJacobian(ts,J,J,RHSJacobian,&data);
111: } else {
112: TSGetSNES(ts,&snes);
113: PetscOptionsHasName(NULL,"-fd_color",&fd_jacobian_coloring);
114: if (fd_jacobian_coloring) { /* Use finite differences with coloring */
115: /* Get data structure of J */
116: PetscBool pc_diagonal;
117: PetscOptionsHasName(NULL,"-pc_diagonal",&pc_diagonal);
118: if (pc_diagonal) { /* the preconditioner of J is a diagonal matrix */
119: PetscInt rstart,rend,i;
120: PetscScalar zero=0.0;
121: MatGetOwnershipRange(J,&rstart,&rend);
122: for (i=rstart; i<rend; i++) {
123: MatSetValues(J,1,&i,1,&i,&zero,INSERT_VALUES);
124: }
125: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
126: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
127: } else {
128: /* Fill the structure using the expensive SNESComputeJacobianDefault. Temporarily set up the TS so we can call this function */
129: TSSetType(ts,TSBEULER);
130: TSSetUp(ts);
131: SNESComputeJacobianDefault(snes,x,J,J,ts);
132: }
134: /* create coloring context */
135: MatColoringCreate(J,&mc);
136: MatColoringSetType(mc,MATCOLORINGSL);
137: MatColoringSetFromOptions(mc);
138: MatColoringApply(mc,&iscoloring);
139: MatColoringDestroy(&mc);
140: MatFDColoringCreate(J,iscoloring,&matfdcoloring);
141: MatFDColoringSetFunction(matfdcoloring,(PetscErrorCode (*)(void))SNESTSFormFunction,ts);
142: MatFDColoringSetFromOptions(matfdcoloring);
143: MatFDColoringSetUp(J,iscoloring,matfdcoloring);
144: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,matfdcoloring);
145: ISColoringDestroy(&iscoloring);
146: } else { /* Use finite differences (slow) */
147: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,NULL);
148: }
149: }
151: /* Pick up a Petsc preconditioner */
152: /* one can always set method or preconditioner during the run time */
153: TSGetSNES(ts,&snes);
154: SNESGetKSP(snes,&ksp);
155: KSPGetPC(ksp,&pc);
156: PCSetType(pc,PCJACOBI);
158: TSSetFromOptions(ts);
159: TSSetUp(ts);
161: /* Test TSSetPostStep() */
162: PetscOptionsHasName(NULL,"-test_PostStep",&flg);
163: if (flg) {
164: TSSetPostStep(ts,PostStep);
165: }
167: PetscOptionsGetInt(NULL,"-NOUT",&NOUT,NULL);
168: for (iout=1; iout<=NOUT; iout++) {
169: TSSetDuration(ts,time_steps,iout*ftime_original/NOUT);
170: TSSolve(ts,global);
171: TSGetSolveTime(ts,&ftime);
172: TSSetInitialTimeStep(ts,ftime,dt);
173: }
174: /* Interpolate solution at tfinal */
175: TSGetSolution(ts,&global);
176: TSInterpolate(ts,ftime_original,global);
178: PetscOptionsHasName(NULL,"-matlab_view",&flg);
179: if (flg) { /* print solution into a MATLAB file */
180: PetscViewerASCIIOpen(PETSC_COMM_WORLD,"out.m",&viewfile);
181: PetscViewerSetFormat(viewfile,PETSC_VIEWER_ASCII_MATLAB);
182: VecView(global,viewfile);
183: PetscViewerDestroy(&viewfile);
184: }
186: /* display solver info for Sundials */
187: TSGetType(ts,&tstype);
188: PetscObjectTypeCompare((PetscObject)ts,TSSUNDIALS,&sundials);
189: if (sundials) {
190: PetscViewerStringOpen(PETSC_COMM_WORLD,tsinfo,120,&viewer);
191: TSView(ts,viewer);
192: PetscViewerDestroy(&viewer);
193: PetscViewerStringOpen(PETSC_COMM_WORLD,pcinfo,120,&viewer);
194: PCView(pc,viewer);
195: PetscPrintf(PETSC_COMM_WORLD,"%d Procs,%s TSType, %s Preconditioner\n",size,tsinfo,pcinfo);
196: PetscViewerDestroy(&viewer);
197: }
199: /* free the memories */
200: TSDestroy(&ts);
201: VecDestroy(&global);
202: VecDestroy(&x);
203: MatDestroy(&J);
204: if (fd_jacobian_coloring) {MatFDColoringDestroy(&matfdcoloring);}
205: PetscFinalize();
206: return 0;
207: }
209: /* -------------------------------------------------------------------*/
210: /* the initial function */
211: PetscReal f_ini(PetscReal x,PetscReal y)
212: {
213: PetscReal f;
215: f=PetscExpReal(-20.0*(PetscPowRealInt(x-0.5,2)+PetscPowRealInt(y-0.5,2)));
216: return f;
217: }
221: PetscErrorCode Initial(Vec global,void *ctx)
222: {
223: Data *data = (Data*)ctx;
224: PetscInt m,row,col;
225: PetscReal x,y,dx,dy;
226: PetscScalar *localptr;
227: PetscInt i,mybase,myend,locsize;
231: /* make the local copies of parameters */
232: m = data->m;
233: dx = data->dx;
234: dy = data->dy;
236: /* determine starting point of each processor */
237: VecGetOwnershipRange(global,&mybase,&myend);
238: VecGetLocalSize(global,&locsize);
240: /* Initialize the array */
241: VecGetArray(global,&localptr);
243: for (i=0; i<locsize; i++) {
244: row = 1+(mybase+i)-((mybase+i)/m)*m;
245: col = (mybase+i)/m+1;
246: x = dx*row;
247: y = dy*col;
248: localptr[i] = f_ini(x,y);
249: }
251: VecRestoreArray(global,&localptr);
252: return(0);
253: }
257: PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal time,Vec global,void *ctx)
258: {
259: VecScatter scatter;
260: IS from,to;
261: PetscInt i,n,*idx,nsteps,maxsteps;
262: Vec tmp_vec;
264: PetscScalar *tmp;
265: PetscReal maxtime;
266: Data *data = (Data*)ctx;
267: PetscReal tfinal = data->tfinal;
270: if (time > tfinal) return(0);
272: TSGetTimeStepNumber(ts,&nsteps);
273: /* display output at selected time steps */
274: TSGetDuration(ts, &maxsteps, &maxtime);
275: if (nsteps % 10 != 0 && time < maxtime) return(0);
277: /* Get the size of the vector */
278: VecGetSize(global,&n);
280: /* Set the index sets */
281: PetscMalloc1(n,&idx);
282: for (i=0; i<n; i++) idx[i]=i;
284: /* Create local sequential vectors */
285: VecCreateSeq(PETSC_COMM_SELF,n,&tmp_vec);
287: /* Create scatter context */
288: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&from);
289: ISCreateGeneral(PETSC_COMM_SELF,n,idx,PETSC_COPY_VALUES,&to);
290: VecScatterCreate(global,from,tmp_vec,to,&scatter);
291: VecScatterBegin(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);
292: VecScatterEnd(scatter,global,tmp_vec,INSERT_VALUES,SCATTER_FORWARD);
294: VecGetArray(tmp_vec,&tmp);
295: PetscPrintf(PETSC_COMM_WORLD,"At t[%D] =%14.2e u= %14.2e at the center \n",nsteps,(double)time,(double)PetscRealPart(tmp[n/2]));
296: VecRestoreArray(tmp_vec,&tmp);
298: PetscFree(idx);
299: ISDestroy(&from);
300: ISDestroy(&to);
301: VecScatterDestroy(&scatter);
302: VecDestroy(&tmp_vec);
303: return(0);
304: }
308: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec x,Mat A,Mat BB,void *ptr)
309: {
310: Data *data = (Data*)ptr;
311: PetscScalar v[5];
312: PetscInt idx[5],i,j,row;
314: PetscInt m,n,mn;
315: PetscReal dx,dy,a,epsilon,xc,xl,xr,yl,yr;
318: m = data->m;
319: n = data->n;
320: mn = m*n;
321: dx = data->dx;
322: dy = data->dy;
323: a = data->a;
324: epsilon = data->epsilon;
326: xc = -2.0*epsilon*(1.0/(dx*dx)+1.0/(dy*dy));
327: xl = 0.5*a/dx+epsilon/(dx*dx);
328: xr = -0.5*a/dx+epsilon/(dx*dx);
329: yl = 0.5*a/dy+epsilon/(dy*dy);
330: yr = -0.5*a/dy+epsilon/(dy*dy);
332: row = 0;
333: v[0] = xc; v[1] = xr; v[2] = yr;
334: idx[0] = 0; idx[1] = 2; idx[2] = m;
335: MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
337: row = m-1;
338: v[0] = 2.0*xl; v[1] = xc; v[2] = yr;
339: idx[0] = m-2; idx[1] = m-1; idx[2] = m-1+m;
340: MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
342: for (i=1; i<m-1; i++) {
343: row = i;
344: v[0] = xl; v[1] = xc; v[2] = xr; v[3] = yr;
345: idx[0] = i-1; idx[1] = i; idx[2] = i+1; idx[3] = i+m;
346: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
347: }
349: for (j=1; j<n-1; j++) {
350: row = j*m;
351: v[0] = xc; v[1] = xr; v[2] = yl; v[3] = yr;
352: idx[0] = j*m; idx[1] = j*m; idx[2] = j*m-m; idx[3] = j*m+m;
353: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
355: row = j*m+m-1;
356: v[0] = xc; v[1] = 2.0*xl; v[2] = yl; v[3] = yr;
357: 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;
358: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
360: for (i=1; i<m-1; i++) {
361: row = j*m+i;
362: v[0] = xc; v[1] = xl; v[2] = xr; v[3] = yl; v[4]=yr;
363: idx[0] = j*m+i; idx[1] = j*m+i-1; idx[2] = j*m+i+1; idx[3] = j*m+i-m;
364: idx[4] = j*m+i+m;
365: MatSetValues(A,1,&row,5,idx,v,INSERT_VALUES);
366: }
367: }
369: row = mn-m;
370: v[0] = xc; v[1] = xr; v[2] = 2.0*yl;
371: idx[0] = mn-m; idx[1] = mn-m+1; idx[2] = mn-m-m;
372: MatSetValues(A,1,&row,3,idx,v,INSERT_VALUES);
374: row = mn-1;
375: v[0] = xc; v[1] = 2.0*xl; v[2] = 2.0*yl;
376: idx[0] = mn-1; idx[1] = mn-2; idx[2] = mn-1-m;
377: MatSetValues(A,1,&i,3,idx,v,INSERT_VALUES);
379: for (i=1; i<m-1; i++) {
380: row = mn-m+i;
381: v[0] = xl; v[1] = xc; v[2] = xr; v[3] = 2.0*yl;
382: idx[0] = mn-m+i-1; idx[1] = mn-m+i; idx[2] = mn-m+i+1; idx[3] = mn-m+i-m;
383: MatSetValues(A,1,&row,4,idx,v,INSERT_VALUES);
384: }
386: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
387: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
389: return(0);
390: }
392: /* globalout = -a*(u_x+u_y) + epsilon*(u_xx+u_yy) */
395: PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec globalin,Vec globalout,void *ctx)
396: {
397: Data *data = (Data*)ctx;
398: PetscInt m,n,mn;
399: PetscReal dx,dy;
400: PetscReal xc,xl,xr,yl,yr;
401: PetscReal a,epsilon;
402: PetscScalar *inptr,*outptr;
403: PetscInt i,j,len;
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: VecGetArray(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]
454: +yr*inptr[i+m];
455: }
457: for (j=1; j<n-1; j++) {
458: outptr[j*m] = xc*inptr[j*m]+xr*inptr[j*m+1]+
459: yl*inptr[j*m-m]+yr*inptr[j*m+m];
460: outptr[j*m+m-1] = xc*inptr[j*m+m-1]+2.0*xl*inptr[j*m+m-1-1]+
461: yl*inptr[j*m+m-1-m]+yr*inptr[j*m+m-1+m];
462: for (i=1; i<m-1; i++) {
463: outptr[j*m+i] = xc*inptr[j*m+i]+xl*inptr[j*m+i-1]+xr*inptr[j*m+i+1]
464: +yl*inptr[j*m+i-m]+yr*inptr[j*m+i+m];
465: }
466: }
468: outptr[mn-m] = xc*inptr[mn-m]+xr*inptr[mn-m+1]+2.0*yl*inptr[mn-m-m];
469: outptr[mn-1] = 2.0*xl*inptr[mn-2]+xc*inptr[mn-1]+2.0*yl*inptr[mn-1-m];
470: for (i=1; i<m-1; i++) {
471: outptr[mn-m+i] = xc*inptr[mn-m+i]+xl*inptr[mn-m+i-1]+xr*inptr[mn-m+i+1]
472: +2*yl*inptr[mn-m+i-m];
473: }
475: VecRestoreArray(tmp_in,&inptr);
476: VecRestoreArray(tmp_out,&outptr);
478: VecScatterCreate(tmp_out,from,globalout,to,&scatter);
479: VecScatterBegin(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
480: VecScatterEnd(scatter,tmp_out,globalout,INSERT_VALUES,SCATTER_FORWARD);
482: /* Destroy idx aand scatter */
483: VecDestroy(&tmp_in);
484: VecDestroy(&tmp_out);
485: ISDestroy(&from);
486: ISDestroy(&to);
487: VecScatterDestroy(&scatter);
489: PetscFree(idx);
490: return(0);
491: }
495: PetscErrorCode PostStep(TS ts)
496: {
498: PetscReal t;
501: TSGetTime(ts,&t);
502: PetscPrintf(PETSC_COMM_SELF," PostStep, t: %g\n",(double)t);
503: return(0);
504: }