Actual source code: ex13.c
petsc-3.8.4 2018-03-24
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 -snes_monitor -ksp_monitor
11: mpiexec -n 1 ./ex13 -snes_fd_color -ts_monitor_draw_solution
12: mpiexec -n 2 ./ex13 -ts_type sundials -ts_monitor
13: */
15: #include <petscdm.h>
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,void*);
28: extern PetscErrorCode FormInitialSolution(DM,Vec,void*);
30: int main(int argc,char **argv)
31: {
32: TS ts; /* nonlinear solver */
33: Vec u,r; /* solution, residual vector */
34: Mat J; /* Jacobian matrix */
35: PetscInt steps; /* iterations for convergence */
37: DM da;
38: PetscReal ftime,dt;
39: AppCtx user; /* user-defined work context */
41: PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
42: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43: Create distributed array (DMDA) to manage parallel grid and vectors
44: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
45: DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da);
46: DMSetFromOptions(da);
47: DMSetUp(da);
49: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50: Extract global vectors from DMDA;
51: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: DMCreateGlobalVector(da,&u);
53: VecDuplicate(u,&r);
55: /* Initialize user application context */
56: user.c = -30.0;
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Create timestepping solver context
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: TSCreate(PETSC_COMM_WORLD,&ts);
62: TSSetDM(ts,da);
63: TSSetType(ts,TSBEULER);
64: TSSetRHSFunction(ts,r,RHSFunction,&user);
66: /* Set Jacobian */
67: DMSetMatType(da,MATAIJ);
68: DMCreateMatrix(da,&J);
69: TSSetRHSJacobian(ts,J,J,RHSJacobian,NULL);
71: ftime = 1.0;
72: TSSetMaxTime(ts,ftime);
73: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Set initial conditions
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
78: FormInitialSolution(da,u,&user);
79: dt = .01;
80: TSSetTimeStep(ts,dt);
82: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83: Set runtime options
84: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
85: TSSetFromOptions(ts);
87: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
88: Solve nonlinear system
89: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
90: TSSolve(ts,u);
91: TSGetSolveTime(ts,&ftime);
92: TSGetStepNumber(ts,&steps);
94: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95: Free work space.
96: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: MatDestroy(&J);
98: VecDestroy(&u);
99: VecDestroy(&r);
100: TSDestroy(&ts);
101: DMDestroy(&da);
103: PetscFinalize();
104: return ierr;
105: }
106: /* ------------------------------------------------------------------- */
107: /*
108: RHSFunction - Evaluates nonlinear function, F(u).
110: Input Parameters:
111: . ts - the TS context
112: . U - input vector
113: . ptr - optional user-defined context, as set by TSSetFunction()
115: Output Parameter:
116: . F - function vector
117: */
118: PetscErrorCode RHSFunction(TS ts,PetscReal ftime,Vec U,Vec F,void *ptr)
119: {
120: /* PETSC_UNUSED AppCtx *user=(AppCtx*)ptr; */
121: DM da;
123: PetscInt i,j,Mx,My,xs,ys,xm,ym;
124: PetscReal two = 2.0,hx,hy,sx,sy;
125: PetscScalar u,uxx,uyy,**uarray,**f;
126: Vec localU;
129: TSGetDM(ts,&da);
130: DMGetLocalVector(da,&localU);
131: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
133: hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
134: hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
136: /*
137: Scatter ghost points to local vector,using the 2-step process
138: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
139: By placing code between these two statements, computations can be
140: done while messages are in transition.
141: */
142: DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
143: DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
145: /* Get pointers to vector data */
146: DMDAVecGetArrayRead(da,localU,&uarray);
147: DMDAVecGetArray(da,F,&f);
149: /* Get local grid boundaries */
150: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
152: /* Compute function over the locally owned part of the grid */
153: for (j=ys; j<ys+ym; j++) {
154: for (i=xs; i<xs+xm; i++) {
155: if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
156: f[j][i] = uarray[j][i];
157: continue;
158: }
159: u = uarray[j][i];
160: uxx = (-two*u + uarray[j][i-1] + uarray[j][i+1])*sx;
161: uyy = (-two*u + uarray[j-1][i] + uarray[j+1][i])*sy;
162: f[j][i] = uxx + uyy;
163: }
164: }
166: /* Restore vectors */
167: DMDAVecRestoreArrayRead(da,localU,&uarray);
168: DMDAVecRestoreArray(da,F,&f);
169: DMRestoreLocalVector(da,&localU);
170: PetscLogFlops(11.0*ym*xm);
171: return(0);
172: }
174: /* --------------------------------------------------------------------- */
175: /*
176: RHSJacobian - User-provided routine to compute the Jacobian of
177: the nonlinear right-hand-side function of the ODE.
179: Input Parameters:
180: ts - the TS context
181: t - current time
182: U - global input vector
183: dummy - optional user-defined context, as set by TSetRHSJacobian()
185: Output Parameters:
186: J - Jacobian matrix
187: Jpre - optionally different preconditioning matrix
188: str - flag indicating matrix structure
189: */
190: PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat J,Mat Jpre,void *ctx)
191: {
193: DM da;
194: DMDALocalInfo info;
195: PetscInt i,j;
196: PetscReal hx,hy,sx,sy;
199: TSGetDM(ts,&da);
200: DMDAGetLocalInfo(da,&info);
201: hx = 1.0/(PetscReal)(info.mx-1); sx = 1.0/(hx*hx);
202: hy = 1.0/(PetscReal)(info.my-1); sy = 1.0/(hy*hy);
203: for (j=info.ys; j<info.ys+info.ym; j++) {
204: for (i=info.xs; i<info.xs+info.xm; i++) {
205: PetscInt nc = 0;
206: MatStencil row,col[5];
207: PetscScalar val[5];
208: row.i = i; row.j = j;
209: if (i == 0 || j == 0 || i == info.mx-1 || j == info.my-1) {
210: col[nc].i = i; col[nc].j = j; val[nc++] = 1.0;
211: } else {
212: col[nc].i = i-1; col[nc].j = j; val[nc++] = sx;
213: col[nc].i = i+1; col[nc].j = j; val[nc++] = sx;
214: col[nc].i = i; col[nc].j = j-1; val[nc++] = sy;
215: col[nc].i = i; col[nc].j = j+1; val[nc++] = sy;
216: col[nc].i = i; col[nc].j = j; val[nc++] = -2*sx - 2*sy;
217: }
218: MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);
219: }
220: }
221: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
222: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
223: if (J != Jpre) {
224: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
225: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
226: }
227: return(0);
228: }
230: /* ------------------------------------------------------------------- */
231: PetscErrorCode FormInitialSolution(DM da,Vec U,void* ptr)
232: {
233: AppCtx *user=(AppCtx*)ptr;
234: PetscReal c=user->c;
236: PetscInt i,j,xs,ys,xm,ym,Mx,My;
237: PetscScalar **u;
238: PetscReal hx,hy,x,y,r;
241: DMDAGetInfo(da,PETSC_IGNORE,&Mx,&My,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
243: hx = 1.0/(PetscReal)(Mx-1);
244: hy = 1.0/(PetscReal)(My-1);
246: /* Get pointers to vector data */
247: DMDAVecGetArray(da,U,&u);
249: /* Get local grid boundaries */
250: DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL);
252: /* Compute function over the locally owned part of the grid */
253: for (j=ys; j<ys+ym; j++) {
254: y = j*hy;
255: for (i=xs; i<xs+xm; i++) {
256: x = i*hx;
257: r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
258: if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
259: else u[j][i] = 0.0;
260: }
261: }
263: /* Restore vectors */
264: DMDAVecRestoreArray(da,U,&u);
265: return(0);
266: }