xref: /petsc/src/ts/tutorials/ex17.c (revision c3e4dd79e17d3f0a51a524340fae4661630fd09e)
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
7 
8    Boundary conditions:
9    Dirichlet BC:
10    At x=0, x=1, u = 0.0
11 
12    Neumann BC:
13    At x=0, x=1: du(x,t)/dx = 0
14 
15    mpiexec -n 2 ./ex17 -da_grid_x 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
16          ./ex17 -da_grid_x 40 -monitor_solution
17          ./ex17 -da_grid_x 100  -ts_type theta -ts_theta_theta 0.5 # Midpoint is not L-stable
18          ./ex17 -jac_type fd_coloring  -da_grid_x 500 -boundary 1
19          ./ex17 -da_grid_x 100  -ts_type gl -ts_adapt_type none -ts_max_steps 2
20 */
21 
22 #include <petscdm.h>
23 #include <petscdmda.h>
24 #include <petscts.h>
25 
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};
28 
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;
37 
38 static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
39 static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
40 static PetscErrorCode FormInitialSolution(TS,Vec,void*);
41 
42 int main(int argc,char **argv)
43 {
44   TS             ts;                   /* nonlinear solver */
45   Vec            u;                    /* solution, residual vectors */
46   Mat            J;                    /* Jacobian matrix */
47   PetscInt       nsteps;
48   PetscReal      vmin,vmax,norm;
49   DM             da;
50   PetscReal      ftime,dt;
51   AppCtx         user;              /* user-defined work context */
52   JacobianType   jacType;
53 
54   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
55 
56   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57      Create distributed array (DMDA) to manage parallel grid and vectors
58   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59   PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,11,1,1,NULL,&da));
60   PetscCall(DMSetFromOptions(da));
61   PetscCall(DMSetUp(da));
62 
63   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64      Extract global vectors from DMDA;
65    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
66   PetscCall(DMCreateGlobalVector(da,&u));
67 
68   /* Initialize user application context */
69   user.c            = -30.0;
70   user.boundary     = 0;  /* 0: Dirichlet BC; 1: Neumann BC */
71   user.viewJacobian = PETSC_FALSE;
72 
73   PetscCall(PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL));
74   PetscCall(PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian));
75 
76   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77      Create timestepping solver context
78      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
80   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
81   PetscCall(TSSetType(ts,TSTHETA));
82   PetscCall(TSThetaSetTheta(ts,1.0)); /* Make the Theta method behave like backward Euler */
83   PetscCall(TSSetIFunction(ts,NULL,FormIFunction,&user));
84 
85   PetscCall(DMSetMatType(da,MATAIJ));
86   PetscCall(DMCreateMatrix(da,&J));
87   jacType = JACOBIAN_ANALYTIC; /* use user-provide Jacobian */
88 
89   PetscCall(TSSetDM(ts,da)); /* Use TSGetDM() to access. Setting here allows easy use of geometric multigrid. */
90 
91   ftime = 1.0;
92   PetscCall(TSSetMaxTime(ts,ftime));
93   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
94 
95   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96      Set initial conditions
97    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98   PetscCall(FormInitialSolution(ts,u,&user));
99   PetscCall(TSSetSolution(ts,u));
100   dt   = .01;
101   PetscCall(TSSetTimeStep(ts,dt));
102 
103   /* Use slow fd Jacobian or fast fd Jacobian with colorings.
104      Note: this requirs snes which is not created until TSSetUp()/TSSetFromOptions() is called */
105   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Options for Jacobian evaluation",NULL);
106   PetscCall(PetscOptionsEnum("-jac_type","Type of Jacobian","",JacobianTypes,(PetscEnum)jacType,(PetscEnum*)&jacType,0));
107   PetscOptionsEnd();
108   if (jacType == JACOBIAN_ANALYTIC) {
109     PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,&user));
110   } else if (jacType == JACOBIAN_FD_COLORING) {
111     SNES snes;
112     PetscCall(TSGetSNES(ts,&snes));
113     PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefaultColor,0));
114   } else if (jacType == JACOBIAN_FD_FULL) {
115     SNES snes;
116     PetscCall(TSGetSNES(ts,&snes));
117     PetscCall(SNESSetJacobian(snes,J,J,SNESComputeJacobianDefault,&user));
118   }
119 
120   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121      Set runtime options
122    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123   PetscCall(TSSetFromOptions(ts));
124 
125   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
126      Integrate ODE system
127      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
128   PetscCall(TSSolve(ts,u));
129 
130   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131    Compute diagnostics of the solution
132    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133   PetscCall(VecNorm(u,NORM_1,&norm));
134   PetscCall(VecMax(u,NULL,&vmax));
135   PetscCall(VecMin(u,NULL,&vmin));
136   PetscCall(TSGetStepNumber(ts,&nsteps));
137   PetscCall(TSGetTime(ts,&ftime));
138   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"timestep %" PetscInt_FMT ": time %g, solution norm %g, max %g, min %g\n",nsteps,(double)ftime,(double)norm,(double)vmax,(double)vmin));
139 
140   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141      Free work space.
142    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143   PetscCall(MatDestroy(&J));
144   PetscCall(VecDestroy(&u));
145   PetscCall(TSDestroy(&ts));
146   PetscCall(DMDestroy(&da));
147   PetscCall(PetscFinalize());
148   return 0;
149 }
150 /* ------------------------------------------------------------------- */
151 static PetscErrorCode FormIFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
152 {
153   AppCtx         *user=(AppCtx*)ptr;
154   DM             da;
155   PetscInt       i,Mx,xs,xm;
156   PetscReal      hx,sx;
157   PetscScalar    *u,*udot,*f;
158   Vec            localU;
159 
160   PetscFunctionBeginUser;
161   PetscCall(TSGetDM(ts,&da));
162   PetscCall(DMGetLocalVector(da,&localU));
163   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
164 
165   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
166 
167   /*
168      Scatter ghost points to local vector,using the 2-step process
169         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
170      By placing code between these two statements, computations can be
171      done while messages are in transition.
172   */
173   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
174   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
175 
176   /* Get pointers to vector data */
177   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
178   PetscCall(DMDAVecGetArrayRead(da,Udot,&udot));
179   PetscCall(DMDAVecGetArray(da,F,&f));
180 
181   /* Get local grid boundaries */
182   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
183 
184   /* Compute function over the locally owned part of the grid */
185   for (i=xs; i<xs+xm; i++) {
186     if (user->boundary == 0) { /* Dirichlet BC */
187       if (i == 0 || i == Mx-1) f[i] = u[i]; /* F = U */
188       else                     f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
189     } else { /* Neumann BC */
190       if (i == 0)         f[i] = u[0] - u[1];
191       else if (i == Mx-1) f[i] = u[i] - u[i-1];
192       else                f[i] = udot[i] + (2.*u[i] - u[i-1] - u[i+1])*sx;
193     }
194   }
195 
196   /* Restore vectors */
197   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
198   PetscCall(DMDAVecRestoreArrayRead(da,Udot,&udot));
199   PetscCall(DMDAVecRestoreArray(da,F,&f));
200   PetscCall(DMRestoreLocalVector(da,&localU));
201   PetscFunctionReturn(0);
202 }
203 
204 /* --------------------------------------------------------------------- */
205 /*
206   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
207 */
208 PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
209 {
210   PetscInt       i,rstart,rend,Mx;
211   PetscReal      hx,sx;
212   AppCtx         *user = (AppCtx*)ctx;
213   DM             da;
214   MatStencil     col[3],row;
215   PetscInt       nc;
216   PetscScalar    vals[3];
217 
218   PetscFunctionBeginUser;
219   PetscCall(TSGetDM(ts,&da));
220   PetscCall(MatGetOwnershipRange(Jpre,&rstart,&rend));
221   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
222   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
223   for (i=rstart; i<rend; i++) {
224     nc    = 0;
225     row.i = i;
226     if (user->boundary == 0 && (i == 0 || i == Mx-1)) {
227       col[nc].i = i; vals[nc++] = 1.0;
228     } else if (user->boundary > 0 && i == 0) { /* Left Neumann */
229       col[nc].i = i;   vals[nc++] = 1.0;
230       col[nc].i = i+1; vals[nc++] = -1.0;
231     } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
232       col[nc].i = i-1; vals[nc++] = -1.0;
233       col[nc].i = i;   vals[nc++] = 1.0;
234     } else {                    /* Interior */
235       col[nc].i = i-1; vals[nc++] = -1.0*sx;
236       col[nc].i = i;   vals[nc++] = 2.0*sx + a;
237       col[nc].i = i+1; vals[nc++] = -1.0*sx;
238     }
239     PetscCall(MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES));
240   }
241 
242   PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
243   PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
244   if (J != Jpre) {
245     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
246     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
247   }
248   if (user->viewJacobian) {
249     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Jpre:\n"));
250     PetscCall(MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD));
251   }
252   PetscFunctionReturn(0);
253 }
254 
255 /* ------------------------------------------------------------------- */
256 PetscErrorCode FormInitialSolution(TS ts,Vec U,void *ptr)
257 {
258   AppCtx         *user=(AppCtx*)ptr;
259   PetscReal      c    =user->c;
260   DM             da;
261   PetscInt       i,xs,xm,Mx;
262   PetscScalar    *u;
263   PetscReal      hx,x,r;
264 
265   PetscFunctionBeginUser;
266   PetscCall(TSGetDM(ts,&da));
267   PetscCall(DMDAGetInfo(da,PETSC_IGNORE,&Mx,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE,PETSC_IGNORE));
268 
269   hx = 1.0/(PetscReal)(Mx-1);
270 
271   /* Get pointers to vector data */
272   PetscCall(DMDAVecGetArray(da,U,&u));
273 
274   /* Get local grid boundaries */
275   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
276 
277   /* Compute function over the locally owned part of the grid */
278   for (i=xs; i<xs+xm; i++) {
279     x = i*hx;
280     r = PetscSqrtReal((x-.5)*(x-.5));
281     if (r < .125) u[i] = PetscExpReal(c*r*r*r);
282     else          u[i] = 0.0;
283   }
284 
285   /* Restore vectors */
286   PetscCall(DMDAVecRestoreArray(da,U,&u));
287   PetscFunctionReturn(0);
288 }
289 
290 /*TEST
291 
292     test:
293       requires: !single
294       args: -da_grid_x 40 -ts_max_steps 2 -snes_monitor_short -ksp_monitor_short -ts_monitor
295 
296     test:
297       suffix: 2
298       requires: !single
299       args: -da_grid_x 100 -ts_type theta -ts_theta_theta 0.5
300 
301 TEST*/
302