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