Actual source code: ex13.c
petsc-3.4.5 2014-06-29
3: static char help[] = "Time-dependent PDE in 2d. Simplified from ex7.c for illustrating how to use TS on a structured domain. \n";
4: /*
5: u_t = uxx + uyy
6: 0 < x < 1, 0 < y < 1;
7: At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125
8: u(x,y) = 0.0 if r >= .125
10: mpiexec -n 2 ./ex13 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -use_coloring -snes_monitor -ksp_monitor
11: ./ex13 -use_coloring
12: ./ex13 -use_coloring -draw_pause -1
13: mpiexec -n 2 ./ex13 -ts_type sundials -ts_sundials_monitor_steps
14: */
16: #include <petscdmda.h>
17: #include <petscts.h>
19: /*
20: User-defined data structures and routines
21: */
22: typedef struct {
23: PetscReal c;
24: } AppCtx;
26: extern PetscErrorCode RHSFunction(TS,PetscReal,Vec,Vec,void*);
27: extern PetscErrorCode RHSJacobian(TS,PetscReal,Vec,Mat*,Mat*,MatStructure*,void*);
28: extern PetscErrorCode FormInitialSolution(DM,Vec,void*);
32: int main(int argc,char **argv)
33: {
34: TS ts; /* nonlinear solver */
35: Vec u,r; /* solution, residual vector */
36: Mat J; /* Jacobian matrix */
37: PetscInt steps,maxsteps = 1000; /* iterations for convergence */
39: DM da;
40: PetscReal ftime,dt;
41: AppCtx user; /* user-defined work context */
43: PetscInitialize(&argc,&argv,(char*)0,help);
44: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45: Create distributed array (DMDA) to manage parallel grid and vectors
46: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47: DMDACreate2d(PETSC_COMM_WORLD, DMDA_BOUNDARY_NONE, DMDA_BOUNDARY_NONE,DMDA_STENCIL_STAR,-8,-8,PETSC_DECIDE,PETSC_DECIDE,
48: 1,1,NULL,NULL,&da);
50: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51: Extract global vectors from DMDA;
52: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53: DMCreateGlobalVector(da,&u);
54: VecDuplicate(u,&r);
56: /* Initialize user application context */
57: user.c = -30.0;
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create timestepping solver context
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: TSCreate(PETSC_COMM_WORLD,&ts);
63: TSSetDM(ts,da);
64: TSSetType(ts,TSBEULER);
65: TSSetRHSFunction(ts,r,RHSFunction,&user);
67: /* Set Jacobian */
68: DMCreateMatrix(da,MATAIJ,&J);
69: TSSetRHSJacobian(ts,J,J,RHSJacobian,NULL);
71: ftime = 1.0;
72: TSSetDuration(ts,maxsteps,ftime);
74: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75: Set initial conditions
76: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
77: FormInitialSolution(da,u,&user);
78: dt = .01;
79: TSSetInitialTimeStep(ts,0.0,dt);
81: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
82: Set runtime options
83: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
84: TSSetFromOptions(ts);
86: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
87: Solve nonlinear system
88: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
89: TSSolve(ts,u);
90: TSGetSolveTime(ts,&ftime);
91: TSGetTimeStepNumber(ts,&steps);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Free work space.
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: MatDestroy(&J);
97: VecDestroy(&u);
98: VecDestroy(&r);
99: TSDestroy(&ts);
100: DMDestroy(&da);
102: PetscFinalize();
103: return(0);
104: }
105: /* ------------------------------------------------------------------- */
108: /*
109: RHSFunction - Evaluates nonlinear function, F(u).
111: Input Parameters:
112: . ts - the TS context
113: . U - input vector
114: . ptr - optional user-defined context, as set by TSSetFunction()
116: Output Parameter:
117: . F - function vector
118: */
119: PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
120: {
121: /* PETSC_UNUSED AppCtx *user=(AppCtx*)ptr; */
122: DM da;
124: PetscInt i,j,Mx,My,xs,ys,xm,ym;
125: PetscReal two = 2.0,hx,hy,sx,sy;
126: PetscScalar u,uxx,uyy,**uarray,**f;
127: Vec localU;
130: TSGetDM(ts,&da);
131: DMGetLocalVector(da,&localU);
132: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
133: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
135: hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
136: hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
138: /*
139: Scatter ghost points to local vector,using the 2-step process
140: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
141: By placing code between these two statements, computations can be
142: done while messages are in transition.
143: */
144: DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
145: DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
147: /* Get pointers to vector data */
148: DMDAVecGetArray(da,localU,&uarray);
149: DMDAVecGetArray(da,F,&f);
151: /* Get local grid boundaries */
152: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
154: /* Compute function over the locally owned part of the grid */
155: for (j=ys; j<ys+ym; j++) {
156: for (i=xs; i<xs+xm; i++) {
157: if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
158: f[j][i] = uarray[j][i];
159: continue;
160: }
161: u = uarray[j][i];
162: uxx = (-two*u + uarray[j][i-1] + uarray[j][i+1])*sx;
163: uyy = (-two*u + uarray[j-1][i] + uarray[j+1][i])*sy;
164: f[j][i] = uxx + uyy;
165: }
166: }
168: /* Restore vectors */
169: DMDAVecRestoreArray(da,localU,&uarray);
170: DMDAVecRestoreArray(da,F,&f);
171: DMRestoreLocalVector(da,&localU);
172: PetscLogFlops(11.0*ym*xm);
173: return(0);
174: }
176: /* --------------------------------------------------------------------- */
179: /*
180: RHSJacobian - User-provided routine to compute the Jacobian of
181: the nonlinear right-hand-side function of the ODE.
183: Input Parameters:
184: ts - the TS context
185: t - current time
186: U - global input vector
187: dummy - optional user-defined context, as set by TSetRHSJacobian()
189: Output Parameters:
190: J - Jacobian matrix
191: Jpre - optionally different preconditioning matrix
192: str - flag indicating matrix structure
193: */
194: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat *J,Mat *Jpre,MatStructure *str,void *ctx)
195: {
197: DM da;
198: DMDALocalInfo info;
199: PetscInt i,j;
200: PetscReal hx,hy,sx,sy;
203: TSGetDM(ts,&da);
204: DMDAGetLocalInfo(da,&info);
205: hx = 1.0/(PetscReal)(info.mx-1); sx = 1.0/(hx*hx);
206: hy = 1.0/(PetscReal)(info.my-1); sy = 1.0/(hy*hy);
207: for (j=info.ys; j<info.ys+info.ym; j++) {
208: for (i=info.xs; i<info.xs+info.xm; i++) {
209: PetscInt nc = 0;
210: MatStencil row,col[5];
211: PetscScalar val[5];
212: row.i = i; row.j = j;
213: if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
214: col[nc].i = i; col[nc].j = j; val[nc++] = 1.0;
215: } else {
216: col[nc].i = i-1; col[nc].j = j; val[nc++] = sx;
217: col[nc].i = i+1; col[nc].j = j; val[nc++] = sx;
218: col[nc].i = i; col[nc].j = j-1; val[nc++] = sy;
219: col[nc].i = i; col[nc].j = j+1; val[nc++] = sy;
220: col[nc].i = i; col[nc].j = j; val[nc++] = -2*sx - 2*sy;
221: }
222: MatSetValuesStencil(*Jpre,1,&row,nc,col,val,INSERT_VALUES);
223: }
224: }
225: MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);
226: MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);
227: if (*J != *Jpre) {
228: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
229: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
230: }
231: return(0);
232: }
234: /* ------------------------------------------------------------------- */
237: PetscErrorCode FormInitialSolution(DM da,Vec U,void* ptr)
238: {
239: AppCtx *user=(AppCtx*)ptr;
240: PetscReal c=user->c;
242: PetscInt i,j,xs,ys,xm,ym,Mx,My;
243: PetscScalar **u;
244: PetscReal hx,hy,x,y,r;
247: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
248: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
250: hx = 1.0/(PetscReal)(Mx-1);
251: hy = 1.0/(PetscReal)(My-1);
253: /* Get pointers to vector data */
254: DMDAVecGetArray(da,U,&u);
256: /* Get local grid boundaries */
257: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
259: /* Compute function over the locally owned part of the grid */
260: for (j=ys; j<ys+ym; j++) {
261: y = j*hy;
262: for (i=xs; i<xs+xm; i++) {
263: x = i*hx;
264: r = PetscSqrtScalar((x-.5)*(x-.5) + (y-.5)*(y-.5));
265: if (r < .125) u[j][i] = PetscExpScalar(c*r*r*r);
266: else u[j][i] = 0.0;
267: }
268: }
270: /* Restore vectors */
271: DMDAVecRestoreArray(da,U,&u);
272: return(0);
273: }