Actual source code: ex17.c
petsc-3.6.1 2015-08-06
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 <petscdm.h>
24: #include <petscdmda.h>
25: #include <petscts.h>
27: typedef enum {JACOBIAN_ANALYTIC,JACOBIAN_FD_COLORING,JACOBIAN_FD_FULL} JacobianType;
28: static const char *const JacobianTypes[] = {"analytic","fd_coloring","fd_full","JacobianType","fd_",0};
30: /*
31: User-defined data structures and routines
32: */
33: typedef struct {
34: PetscReal c;
35: PetscInt boundary; /* Type of boundary condition */
36: PetscBool viewJacobian;
37: } AppCtx;
39: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
40: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
41: static PetscErrorCode FormInitialSolution(TS,Vec,void*);
45: int main(int argc,char **argv)
46: {
47: TS ts; /* nonlinear solver */
48: Vec u; /* solution, residual vectors */
49: Mat J; /* Jacobian matrix */
50: PetscInt maxsteps = 1000; /* iterations for convergence */
51: PetscInt nsteps;
52: PetscReal vmin,vmax,norm;
54: DM da;
55: PetscReal ftime,dt;
56: AppCtx user; /* user-defined work context */
57: JacobianType jacType;
59: PetscInitialize(&argc,&argv,(char*)0,help);
61: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
62: Create distributed array (DMDA) to manage parallel grid and vectors
63: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,-11,1,1,NULL,&da);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Extract global vectors from DMDA;
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69: DMCreateGlobalVector(da,&u);
71: /* Initialize user application context */
72: user.c = -30.0;
73: user.boundary = 0; /* 0: Dirichlet BC; 1: Neumann BC */
74: user.viewJacobian = PETSC_FALSE;
76: PetscOptionsGetInt(NULL,"-boundary",&user.boundary,NULL);
77: PetscOptionsHasName(NULL,"-viewJacobian",&user.viewJacobian);
79: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80: Create timestepping solver context
81: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82: TSCreate(PETSC_COMM_WORLD,&ts);
83: TSSetProblemType(ts,TS_NONLINEAR);
84: TSSetType(ts,TSTHETA);
85: TSThetaSetTheta(ts,1.0); /* Make the Theta method behave like backward Euler */
86: TSSetIFunction(ts,NULL,FormIFunction,&user);
88: DMSetMatType(da,MATAIJ);
89: DMCreateMatrix(da,&J);
90: jacType = JACOBIAN_ANALYTIC; /* use user-provide Jacobian */
92: TSSetDM(ts,da); /* Use TSGetDM() to access. Setting here allows easy use of geometric multigrid. */
94: ftime = 1.0;
95: TSSetDuration(ts,maxsteps,ftime);
97: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
98: Set initial conditions
99: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: FormInitialSolution(ts,u,&user);
101: TSSetSolution(ts,u);
102: dt = .01;
103: TSSetInitialTimeStep(ts,0.0,dt);
106: /* Use slow fd Jacobian or fast fd Jacobian with colorings.
107: Note: this requirs snes which is not created until TSSetUp()/TSSetFromOptions() is called */
108: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for Jacobian evaluation",NULL);
109: PetscOptionsEnum("-jac_type","Type of Jacobian","",JacobianTypes,(PetscEnum)jacType,(PetscEnum*)&jacType,0);
110: PetscOptionsEnd();
111: if (jacType == JACOBIAN_ANALYTIC) {
112: TSSetIJacobian(ts,J,J,FormIJacobian,&user);
113: } else if (jacType == JACOBIAN_FD_COLORING) {
114: SNES snes;
115: TSGetSNES(ts,&snes);
116: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0);
117: } else if (jacType == JACOBIAN_FD_FULL) {
118: SNES snes;
119: TSGetSNES(ts,&snes);
120: SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,&user);
121: }
123: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124: Set runtime options
125: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
126: TSSetFromOptions(ts);
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Integrate ODE system
130: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131: TSSolve(ts,u);
133: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134: Compute diagnostics of the solution
135: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136: VecNorm(u,NORM_1,&norm);
137: VecMax(u,NULL,&vmax);
138: VecMin(u,NULL,&vmin);
139: TSGetTimeStepNumber(ts,&nsteps);
140: TSGetTime(ts,&ftime);
141: PetscPrintf(PETSC_COMM_WORLD,"timestep %D: time %g, solution norm %g, max %g, min %g\n",nsteps,(double)ftime,(double)norm,(double)vmax,(double)vmin);
143: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
144: Free work space.
145: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146: MatDestroy(&J);
147: VecDestroy(&u);
148: TSDestroy(&ts);
149: DMDestroy(&da);
150: PetscFinalize();
151: return(0);
152: }
153: /* ------------------------------------------------------------------- */
156: static PetscErrorCode FormIFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
157: {
158: AppCtx *user=(AppCtx*)ptr;
159: DM da;
161: PetscInt i,Mx,xs,xm;
162: PetscReal hx,sx;
163: PetscScalar *u,*udot,*f;
164: Vec localU;
167: TSGetDM(ts,&da);
168: DMGetLocalVector(da,&localU);
169: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
170: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
172: hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
174: /*
175: Scatter ghost points to local vector,using the 2-step process
176: DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
177: By placing code between these two statements, computations can be
178: done while messages are in transition.
179: */
180: DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
181: DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);
183: /* Get pointers to vector data */
184: DMDAVecGetArrayRead(da,localU,&u);
185: DMDAVecGetArrayRead(da,Udot,&udot);
186: DMDAVecGetArray(da,F,&f);
188: /* Get local grid boundaries */
189: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
191: /* Compute function over the locally owned part of the grid */
192: for (i=xs; i<xs+xm; i++) {
193: if (user->boundary == 0) { /* Dirichlet BC */
194: if (i == 0 || i == Mx-1) f[i] = u[i]; /* F = U */
195: else f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
196: } else { /* Neumann BC */
197: if (i == 0) f[i] = u[0] - u[1];
198: else if (i == Mx-1) f[i] = u[i] - u[i-1];
199: else f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
200: }
201: }
203: /* Restore vectors */
204: DMDAVecRestoreArrayRead(da,localU,&u);
205: DMDAVecRestoreArrayRead(da,Udot,&udot);
206: DMDAVecRestoreArray(da,F,&f);
207: DMRestoreLocalVector(da,&localU);
208: return(0);
209: }
211: /* --------------------------------------------------------------------- */
212: /*
213: IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
214: */
217: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
218: {
220: PetscInt i,rstart,rend,Mx;
221: PetscReal hx,sx;
222: AppCtx *user = (AppCtx*)ctx;
223: DM da;
224: MatStencil col[3],row;
225: PetscInt nc;
226: PetscScalar vals[3];
229: TSGetDM(ts,&da);
230: MatGetOwnershipRange(Jpre,&rstart,&rend);
231: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
232: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
233: hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
234: for (i=rstart; i<rend; i++) {
235: nc = 0;
236: row.i = i;
237: if (user->boundary == 0 && (i == 0 || i == Mx-1)) {
238: col[nc].i = i; vals[nc++] = 1.0;
239: } else if (user->boundary > 0 && i == 0) { /* Left Neumann */
240: col[nc].i = i; vals[nc++] = 1.0;
241: col[nc].i = i+1; vals[nc++] = -1.0;
242: } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
243: col[nc].i = i-1; vals[nc++] = -1.0;
244: col[nc].i = i; vals[nc++] = 1.0;
245: } else { /* Interior */
246: col[nc].i = i-1; vals[nc++] = -1.0*sx;
247: col[nc].i = i; vals[nc++] = 2.0*sx + a;
248: col[nc].i = i+1; vals[nc++] = -1.0*sx;
249: }
250: MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES);
251: }
253: MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
254: MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
255: if (J != Jpre) {
256: MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
257: MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
258: }
259: if (user->viewJacobian) {
260: PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n");
261: MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD);
262: }
263: return(0);
264: }
266: /* ------------------------------------------------------------------- */
269: PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr)
270: {
271: AppCtx *user=(AppCtx*)ptr;
272: PetscReal c =user->c;
273: DM da;
275: PetscInt i,xs,xm,Mx;
276: PetscScalar *u;
277: PetscReal hx,x,r;
280: TSGetDM(ts,&da);
281: DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,
282: PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE);
284: hx = 1.0/(PetscReal)(Mx-1);
286: /* Get pointers to vector data */
287: DMDAVecGetArray(da,U,&u);
289: /* Get local grid boundaries */
290: DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);
292: /* Compute function over the locally owned part of the grid */
293: for (i=xs; i<xs+xm; i++) {
294: x = i*hx;
295: r = PetscSqrtReal((x-.5)*(x-.5));
296: if (r < .125) u[i] = PetscExpReal(c*r*r*r);
297: else u[i] = 0.0;
298: }
300: /* Restore vectors */
301: DMDAVecRestoreArray(da,U,&u);
302: return(0);
303: }