Actual source code: ex17.c
petsc-3.4.5 2014-06-29
1: static const char help[] = "Time-dependent PDE in 1d. Simplified from ex15.c for illustrating how to solve DAEs. \n";
2: /*
3: u_t = uxx
4: 0 < x < 1;
5: At t=0: u(x) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5)) < .125
6: u(x) = 0.0 if r >= .125
9: Boundary conditions:
10: Dirichlet BC:
11: At x=0, x=1, u = 0.0
13: Neumann BC:
14: At x=0, x=1: du(x,t)/dx = 0
16: mpiexec -n 2 ./ex17 -da_grid_x 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
17: ./ex17 -da_grid_x 40 -monitor_solution
18: ./ex17 -da_grid_x 100 -ts_type theta -ts_theta_theta 0.5 # Midpoint is not L-stable
19: ./ex17 -jac_type fd_coloring -da_grid_x 500 -boundary 1
20: ./ex17 -da_grid_x 100 -ts_type gl -ts_adapt_type none -ts_max_steps 2
21: */
23: #include <petscdmda.h>
24: #include <petscts.h>
26: typedef enum {JACOBIAN_ANALYTIC,JACOBIAN_FD_COLORING,JACOBIAN_FD_FULL} JacobianType;
27: static const char *const JacobianTypes[] = {"analytic","fd_coloring","fd_full","JacobianType","fd_",0};
29: /*
30: User-defined data structures and routines
31: */
32: typedef struct {
33: PetscReal c;
34: PetscInt boundary; /* Type of boundary condition */
35: PetscBool viewJacobian;
36: } AppCtx;
38: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
39: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
40: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
44: int main(int argc,char **argv)
45: {
46: TS ts; /* nonlinear solver */
47: Vec u; /* solution, residual vectors */
48: Mat J; /* Jacobian matrix */
49: PetscInt maxsteps = 1000; /* iterations for convergence */
50: PetscInt nsteps;
51: PetscReal vmin,vmax,norm;
53: DM da;
54: PetscReal ftime,dt;
55: AppCtx user; /* user-defined work context */
56: JacobianType jacType;
58: PetscInitialize(&argc,&argv,(char*)0,help);
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create distributed array (DMDA) to manage parallel grid and vectors
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-11,1,1,NULL,&da);
65: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66: Extract global vectors from DMDA;
67: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68: DMCreateGlobalVector(da,&u);
70: /* Initialize user application context */
71: user.c = -30.0;
72: user.boundary = 0; /* 0: Dirichlet BC; 1: Neumann BC */
73: user.viewJacobian = PETSC_FALSE;
75: PetscOptionsGetInt(NULL,"-boundary",&user.boundary,NULL);
76: PetscOptionsHasName(NULL,"-viewJacobian",&user.viewJacobian);
78: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79: Create timestepping solver context
80: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
81: TSCreate(PETSC_COMM_WORLD,&ts);
82: TSSetProblemType(ts,TS_NONLINEAR);
83: TSSetType(ts,TSTHETA);
84: TSThetaSetTheta(ts,1.0); /* Make the Theta method behave like backward Euler */
85: TSSetIFunction(ts,NULL,FormIFunction,&user);
87: DMCreateMatrix(da,MATAIJ,&J);
88: jacType = JACOBIAN_ANALYTIC; /* use user-provide Jacobian */
90: TSSetDM(ts,da); /* Use TSGetDM() to access. Setting here allows easy use of geometric multigrid. */
92: ftime = 1.0;
93: TSSetDuration(ts,maxsteps,ftime);
95: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96: Set initial conditions
97: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98: FormInitialSolution(ts,u,&user);
99: TSSetSolution(ts,u);
100: dt = .01;
101: TSSetInitialTimeStep(ts,0.0,dt);
104: /* Use slow fd Jacobian or fast fd Jacobian with colorings.
105: Note: this requirs snes which is not created until TSSetUp()/TSSetFromOptions() is called */
106: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for Jacobian evaluation",NULL);
107: PetscOptionsEnum("-jac_type","Type of Jacobian","",JacobianTypes,(PetscEnum)jacType,(PetscEnum*)&jacType,0);
108: PetscOptionsEnd();
109: if (jacType == JACOBIAN_ANALYTIC) {
110: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
111: } else if (jacType == JACOBIAN_FD_COLORING) {
112: SNES snes;
113: TSGetSNES(ts,&snes);
114: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0);
115: } else if (jacType == JACOBIAN_FD_FULL) {
116: SNES snes;
117: TSGetSNES(ts,&snes);
118: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,&user);
119: }
121: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: Set runtime options
123: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
124: TSSetFromOptions(ts);
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Integrate ODE system
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129: TSSolve(ts,u);
131: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132: Compute diagnostics of the solution
133: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134: VecNorm(u,NORM_1,&norm);
135: VecMax(u,NULL,&vmax);
136: VecMin(u,NULL,&vmin);
137: TSGetTimeStepNumber(ts,&nsteps);
138: TSGetTime(ts,&ftime);
139: PetscPrintf(PETSC_COMM_WORLD,"timestep %D: time %G, solution norm %G, max %G, min %G\n",nsteps,ftime,norm,vmax,vmin);
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Free work space.
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: MatDestroy(&J);
145: VecDestroy(&u);
146: TSDestroy(&ts);
147: DMDestroy(&da);
148: PetscFinalize();
149: return(0);
150: }
151: /* ------------------------------------------------------------------- */
154: static PetscErrorCode FormIFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
155: {
156: AppCtx *user=(AppCtx*)ptr;
157: DM da;
159: PetscInt i,Mx,xs,xm;
160: PetscReal hx,sx;
161: PetscScalar *u,*udot,*f;
162: Vec localU;
165: TSGetDM(ts,&da);
166: DMGetLocalVector(da,&localU);
167: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
168: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
170: hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
172: /*
173: Scatter ghost points to local vector,using the 2-step process
174: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
175: By placing code between these two statements, computations can be
176: done while messages are in transition.
177: */
178: DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
179: DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
181: /* Get pointers to vector data */
182: DMDAVecGetArray(da,localU,&u);
183: DMDAVecGetArray(da,Udot,&udot);
184: DMDAVecGetArray(da,F,&f);
186: /* Get local grid boundaries */
187: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
189: /* Compute function over the locally owned part of the grid */
190: for (i=xs; i<xs+xm; i++) {
191: if (user->boundary == 0) { /* Dirichlet BC */
192: if (i == 0 || i == Mx-1) f[i] = u[i]; /* F = U */
193: else f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
194: } else { /* Neumann BC */
195: if (i == 0) f[i] = u[0] - u[1];
196: else if (i == Mx-1) f[i] = u[i] - u[i-1];
197: else f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
198: }
199: }
201: /* Restore vectors */
202: DMDAVecRestoreArray(da,localU,&u);
203: DMDAVecRestoreArray(da,Udot,&udot);
204: DMDAVecRestoreArray(da,F,&f);
205: DMRestoreLocalVector(da,&localU);
206: return(0);
207: }
209: /* --------------------------------------------------------------------- */
210: /*
211: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
212: */
215: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ctx)
216: {
218: PetscInt i,rstart,rend,Mx;
219: PetscReal hx,sx;
220: AppCtx *user = (AppCtx*)ctx;
221: DM da;
222: MatStencil col[3],row;
223: PetscInt nc;
224: PetscScalar vals[3];
227: TSGetDM(ts,&da);
228: MatGetOwnershipRange(*Jpre,&rstart,&rend);
229: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
230: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
231: hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
232: for (i=rstart; i<rend; i++) {
233: nc = 0;
234: row.i = i;
235: if (user->boundary == 0 && (i == 0 || i == Mx-1)) {
236: col[nc].i = i; vals[nc++] = 1.0;
237: } else if (user->boundary > 0 && i == 0) { /* Left Neumann */
238: col[nc].i = i; vals[nc++] = 1.0;
239: col[nc].i = i+1; vals[nc++] = -1.0;
240: } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
241: col[nc].i = i-1; vals[nc++] = -1.0;
242: col[nc].i = i; vals[nc++] = 1.0;
243: } else { /* Interior */
244: col[nc].i = i-1; vals[nc++] = -1.0*sx;
245: col[nc].i = i; vals[nc++] = 2.0*sx + a;
246: col[nc].i = i+1; vals[nc++] = -1.0*sx;
247: }
248: MatSetValuesStencil(*Jpre,1,&row,nc,col,vals,INSERT_VALUES);
249: }
251: MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);
252: MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);
253: if (*J != *Jpre) {
254: MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
255: MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
256: }
257: if (user->viewJacobian) {
258: PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n");
259: MatView(*Jpre,PETSC_VIEWER_STDOUT_WORLD);
260: }
261: return(0);
262: }
264: /* ------------------------------------------------------------------- */
267: PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr)
268: {
269: AppCtx *user=(AppCtx*)ptr;
270: PetscReal c =user->c;
271: DM da;
273: PetscInt i,xs,xm,Mx;
274: PetscScalar *u;
275: PetscReal hx,x,r;
278: TSGetDM(ts,&da);
279: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
280: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
282: hx = 1.0/(PetscReal)(Mx-1);
284: /* Get pointers to vector data */
285: DMDAVecGetArray(da,U,&u);
287: /* Get local grid boundaries */
288: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
290: /* Compute function over the locally owned part of the grid */
291: for (i=xs; i<xs+xm; i++) {
292: x = i*hx;
293: r = PetscSqrtScalar((x-.5)*(x-.5));
294: if (r < .125) u[i] = PetscExpScalar(c*r*r*r);
295: else u[i] = 0.0;
296: }
298: /* Restore vectors */
299: DMDAVecRestoreArray(da,U,&u);
300: return(0);
301: }