xref: /petsc/src/ts/tutorials/ex15.c (revision 76d901e46dda72c1afe96306c7cb4731c47d4e87)
1 
2 static char help[] = "Time-dependent PDE in 2d. Modified from ex13.c for illustrating how to solve DAEs. \n";
3 /*
4    u_t = uxx + uyy
5    0 < x < 1, 0 < y < 1;
6    At t=0: u(x,y) = exp(c*r*r*r), if r=PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5)) < .125
7            u(x,y) = 0.0           if r >= .125
8 
9    Boundary conditions:
10    Drichlet BC:
11    At x=0, x=1, y=0, y=1: u = 0.0
12 
13    Neumann BC:
14    At x=0, x=1: du(x,y,t)/dx = 0
15    At y=0, y=1: du(x,y,t)/dy = 0
16 
17    mpiexec -n 2 ./ex15 -da_grid_x 40 -da_grid_y 40 -ts_max_steps 2 -snes_monitor -ksp_monitor
18          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -ts_monitor_draw_solution
19          ./ex15 -da_grid_x 40 -da_grid_y 40  -draw_pause .1 -boundary 1 -Jtype 2 -nstencilpts 9
20 
21 */
22 
23 #include <petscdm.h>
24 #include <petscdmda.h>
25 #include <petscts.h>
26 
27 /*
28    User-defined data structures and routines
29 */
30 
31 /* AppCtx: used by FormIFunction() and FormIJacobian() */
32 typedef struct {
33   DM        da;
34   PetscInt  nstencilpts;         /* number of stencil points: 5 or 9 */
35   PetscReal c;
36   PetscInt  boundary;            /* Type of boundary condition */
37   PetscBool viewJacobian;
38 } AppCtx;
39 
40 extern PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
41 extern PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
42 extern PetscErrorCode FormInitialSolution(Vec,void*);
43 
44 int main(int argc,char **argv)
45 {
46   TS             ts;                   /* nonlinear solver */
47   Vec            u,r;                  /* solution, residual vectors */
48   Mat            J,Jmf = NULL;   /* Jacobian matrices */
49   DM             da;
50   PetscReal      dt;
51   AppCtx         user;              /* user-defined work context */
52   SNES           snes;
53   PetscInt       Jtype; /* Jacobian type
54                             0: user provide Jacobian;
55                             1: slow finite difference;
56                             2: fd with coloring; */
57 
58   PetscFunctionBeginUser;
59   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
60   /* Initialize user application context */
61   user.da           = NULL;
62   user.nstencilpts  = 5;
63   user.c            = -30.0;
64   user.boundary     = 0;  /* 0: Drichlet BC; 1: Neumann BC */
65   user.viewJacobian = PETSC_FALSE;
66 
67   PetscCall(PetscOptionsGetInt(NULL,NULL,"-nstencilpts",&user.nstencilpts,NULL));
68   PetscCall(PetscOptionsGetInt(NULL,NULL,"-boundary",&user.boundary,NULL));
69   PetscCall(PetscOptionsHasName(NULL,NULL,"-viewJacobian",&user.viewJacobian));
70 
71   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
72      Create distributed array (DMDA) to manage parallel grid and vectors
73   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74   if (user.nstencilpts == 5) {
75     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
76   } else if (user.nstencilpts == 9) {
77     PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_BOX,11,11,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
78   } else SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"nstencilpts %" PetscInt_FMT " is not supported",user.nstencilpts);
79   PetscCall(DMSetFromOptions(da));
80   PetscCall(DMSetUp(da));
81   user.da = da;
82 
83   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
84      Extract global vectors from DMDA;
85    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
86   PetscCall(DMCreateGlobalVector(da,&u));
87   PetscCall(VecDuplicate(u,&r));
88 
89   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
90      Create timestepping solver context
91      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
92   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
93   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
94   PetscCall(TSSetType(ts,TSBEULER));
95   PetscCall(TSSetDM(ts,da));
96   PetscCall(TSSetIFunction(ts,r,FormIFunction,&user));
97   PetscCall(TSSetMaxTime(ts,1.0));
98   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
99 
100   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101      Set initial conditions
102    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103   PetscCall(FormInitialSolution(u,&user));
104   PetscCall(TSSetSolution(ts,u));
105   dt   = .01;
106   PetscCall(TSSetTimeStep(ts,dt));
107 
108   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109    Set Jacobian evaluation routine
110   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111   PetscCall(DMSetMatType(da,MATAIJ));
112   PetscCall(DMCreateMatrix(da,&J));
113   Jtype = 0;
114   PetscCall(PetscOptionsGetInt(NULL,NULL, "-Jtype",&Jtype,NULL));
115   if (Jtype == 0) { /* use user provided Jacobian evaluation routine */
116     PetscCheck(user.nstencilpts == 5,PETSC_COMM_WORLD,PETSC_ERR_SUP,"user Jacobian routine FormIJacobian() does not support nstencilpts=%" PetscInt_FMT,user.nstencilpts);
117     PetscCall(TSSetIJacobian(ts,J,J,FormIJacobian,&user));
118   } else { /* use finite difference Jacobian J as preconditioner and '-snes_mf_operator' for Mat*vec */
119     PetscCall(TSGetSNES(ts,&snes));
120     PetscCall(MatCreateSNESMF(snes,&Jmf));
121     if (Jtype == 1) { /* slow finite difference J; */
122       PetscCall(SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefault,NULL));
123     } else if (Jtype == 2) { /* Use coloring to compute  finite difference J efficiently */
124       PetscCall(SNESSetJacobian(snes,Jmf,J,SNESComputeJacobianDefaultColor,0));
125     } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP,"Jtype is not supported");
126   }
127 
128   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129    Sets various TS parameters from user options
130    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131   PetscCall(TSSetFromOptions(ts));
132 
133   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
134      Solve nonlinear system
135      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
136   PetscCall(TSSolve(ts,u));
137 
138   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139      Free work space.
140    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141   PetscCall(MatDestroy(&J));
142   PetscCall(MatDestroy(&Jmf));
143   PetscCall(VecDestroy(&u));
144   PetscCall(VecDestroy(&r));
145   PetscCall(TSDestroy(&ts));
146   PetscCall(DMDestroy(&da));
147 
148   PetscCall(PetscFinalize());
149   return 0;
150 }
151 
152 /* --------------------------------------------------------------------- */
153 /*
154   FormIFunction = Udot - RHSFunction
155 */
156 PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx)
157 {
158   AppCtx         *user=(AppCtx*)ctx;
159   DM             da   = (DM)user->da;
160   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
161   PetscReal      hx,hy,sx,sy;
162   PetscScalar    u,uxx,uyy,**uarray,**f,**udot;
163   Vec            localU;
164 
165   PetscFunctionBeginUser;
166   PetscCall(DMGetLocalVector(da,&localU));
167   PetscCall(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));
168 
169   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
170   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
171   PetscCheck(user->nstencilpts != 9 || hx == hy,PETSC_COMM_WORLD,PETSC_ERR_SUP,"hx must equal hy when nstencilpts = 9 for this example");
172 
173   /*
174      Scatter ghost points to local vector,using the 2-step process
175         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
176      By placing code between these two statements, computations can be
177      done while messages are in transition.
178   */
179   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
180   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
181 
182   /* Get pointers to vector data */
183   PetscCall(DMDAVecGetArrayRead(da,localU,&uarray));
184   PetscCall(DMDAVecGetArray(da,F,&f));
185   PetscCall(DMDAVecGetArray(da,Udot,&udot));
186 
187   /* Get local grid boundaries */
188   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
189 
190   /* Compute function over the locally owned part of the grid */
191   for (j=ys; j<ys+ym; j++) {
192     for (i=xs; i<xs+xm; i++) {
193       /* Boundary conditions */
194       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
195         if (user->boundary == 0) { /* Drichlet BC */
196           f[j][i] = uarray[j][i]; /* F = U */
197         } else {                  /* Neumann BC */
198           if (i == 0 && j == 0) {              /* SW corner */
199             f[j][i] = uarray[j][i] - uarray[j+1][i+1];
200           } else if (i == Mx-1 && j == 0) {    /* SE corner */
201             f[j][i] = uarray[j][i] - uarray[j+1][i-1];
202           } else if (i == 0 && j == My-1) {    /* NW corner */
203             f[j][i] = uarray[j][i] - uarray[j-1][i+1];
204           } else if (i == Mx-1 && j == My-1) { /* NE corner */
205             f[j][i] = uarray[j][i] - uarray[j-1][i-1];
206           } else if (i == 0) {                  /* Left */
207             f[j][i] = uarray[j][i] - uarray[j][i+1];
208           } else if (i == Mx-1) {               /* Right */
209             f[j][i] = uarray[j][i] - uarray[j][i-1];
210           } else if (j == 0) {                 /* Bottom */
211             f[j][i] = uarray[j][i] - uarray[j+1][i];
212           } else if (j == My-1) {               /* Top */
213             f[j][i] = uarray[j][i] - uarray[j-1][i];
214           }
215         }
216       } else { /* Interior */
217         u = uarray[j][i];
218         /* 5-point stencil */
219         uxx = (-2.0*u + uarray[j][i-1] + uarray[j][i+1]);
220         uyy = (-2.0*u + uarray[j-1][i] + uarray[j+1][i]);
221         if (user->nstencilpts == 9) {
222           /* 9-point stencil: assume hx=hy */
223           uxx = 2.0*uxx/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0;
224           uyy = 2.0*uyy/3.0 + (0.5*(uarray[j-1][i-1]+uarray[j-1][i+1]+uarray[j+1][i-1]+uarray[j+1][i+1]) - 2.0*u)/6.0;
225         }
226         f[j][i] = udot[j][i] - (uxx*sx + uyy*sy);
227       }
228     }
229   }
230 
231   /* Restore vectors */
232   PetscCall(DMDAVecRestoreArrayRead(da,localU,&uarray));
233   PetscCall(DMDAVecRestoreArray(da,F,&f));
234   PetscCall(DMDAVecRestoreArray(da,Udot,&udot));
235   PetscCall(DMRestoreLocalVector(da,&localU));
236   PetscCall(PetscLogFlops(11.0*ym*xm));
237   PetscFunctionReturn(0);
238 }
239 
240 /* --------------------------------------------------------------------- */
241 /*
242   FormIJacobian() - Compute IJacobian = dF/dU + a dF/dUdot
243   This routine is not used with option '-use_coloring'
244 */
245 PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat J,Mat Jpre,void *ctx)
246 {
247   PetscInt       i,j,Mx,My,xs,ys,xm,ym,nc;
248   AppCtx         *user = (AppCtx*)ctx;
249   DM             da    = (DM)user->da;
250   MatStencil     col[5],row;
251   PetscScalar    vals[5],hx,hy,sx,sy;
252 
253   PetscFunctionBeginUser;
254   PetscCall(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));
255   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
256 
257   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
258   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
259 
260   for (j=ys; j<ys+ym; j++) {
261     for (i=xs; i<xs+xm; i++) {
262       nc    = 0;
263       row.j = j; row.i = i;
264       if (user->boundary == 0 && (i == 0 || i == Mx-1 || j == 0 || j == My-1)) {
265         col[nc].j = j; col[nc].i = i; vals[nc++] = 1.0;
266 
267       } else if (user->boundary > 0 && i == 0) {  /* Left Neumann */
268         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
269         col[nc].j = j; col[nc].i = i+1; vals[nc++] = -1.0;
270       } else if (user->boundary > 0 && i == Mx-1) { /* Right Neumann */
271         col[nc].j = j; col[nc].i = i;   vals[nc++] = 1.0;
272         col[nc].j = j; col[nc].i = i-1; vals[nc++] = -1.0;
273       } else if (user->boundary > 0 && j == 0) {  /* Bottom Neumann */
274         col[nc].j = j;   col[nc].i = i; vals[nc++] = 1.0;
275         col[nc].j = j+1; col[nc].i = i; vals[nc++] = -1.0;
276       } else if (user->boundary > 0 && j == My-1) { /* Top Neumann */
277         col[nc].j = j;   col[nc].i = i;  vals[nc++] = 1.0;
278         col[nc].j = j-1; col[nc].i = i;  vals[nc++] = -1.0;
279       } else {   /* Interior */
280         col[nc].j = j-1; col[nc].i = i;   vals[nc++] = -sy;
281         col[nc].j = j;   col[nc].i = i-1; vals[nc++] = -sx;
282         col[nc].j = j;   col[nc].i = i;   vals[nc++] = 2.0*(sx + sy) + a;
283         col[nc].j = j;   col[nc].i = i+1; vals[nc++] = -sx;
284         col[nc].j = j+1; col[nc].i = i;   vals[nc++] = -sy;
285       }
286       PetscCall(MatSetValuesStencil(Jpre,1,&row,nc,col,vals,INSERT_VALUES));
287     }
288   }
289   PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
290   PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
291   if (J != Jpre) {
292     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
293     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
294   }
295 
296   if (user->viewJacobian) {
297     PetscCall(PetscPrintf(PetscObjectComm((PetscObject)Jpre),"Jpre:\n"));
298     PetscCall(MatView(Jpre,PETSC_VIEWER_STDOUT_WORLD));
299   }
300   PetscFunctionReturn(0);
301 }
302 
303 /* ------------------------------------------------------------------- */
304 PetscErrorCode FormInitialSolution(Vec U,void *ptr)
305 {
306   AppCtx         *user=(AppCtx*)ptr;
307   DM             da   =user->da;
308   PetscReal      c    =user->c;
309   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
310   PetscScalar    **u;
311   PetscReal      hx,hy,x,y,r;
312 
313   PetscFunctionBeginUser;
314   PetscCall(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));
315 
316   hx = 1.0/(PetscReal)(Mx-1);
317   hy = 1.0/(PetscReal)(My-1);
318 
319   /* Get pointers to vector data */
320   PetscCall(DMDAVecGetArray(da,U,&u));
321 
322   /* Get local grid boundaries */
323   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
324 
325   /* Compute function over the locally owned part of the grid */
326   for (j=ys; j<ys+ym; j++) {
327     y = j*hy;
328     for (i=xs; i<xs+xm; i++) {
329       x = i*hx;
330       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
331       if (r < .125) u[j][i] = PetscExpReal(c*r*r*r);
332       else u[j][i] = 0.0;
333     }
334   }
335 
336   /* Restore vectors */
337   PetscCall(DMDAVecRestoreArray(da,U,&u));
338   PetscFunctionReturn(0);
339 }
340 
341 /*TEST
342 
343     test:
344       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -ts_monitor
345 
346     test:
347       suffix: 2
348       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 2 -ts_monitor
349 
350     test:
351       suffix: 3
352       requires: !single
353       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
354 
355     test:
356       suffix: 4
357       requires: !single
358       nsize: 2
359       args: -da_grid_x 20 -da_grid_y 20 -boundary 1 -ts_max_steps 10 -ts_monitor
360 
361     test:
362       suffix: 5
363       nsize: 1
364       args: -da_grid_x 20 -da_grid_y 20 -boundary 0 -ts_max_steps 10 -Jtype 1 -ts_monitor
365 
366 TEST*/
367