xref: /petsc/src/ts/tutorials/ex12.c (revision 1baa6e3354bfe224b33a0c158f482508792a8a6e)
1c4762a1bSJed Brown static char help[] = "Nonlinear, time-dependent PDE in 2d.\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown   Solves the equation
4c4762a1bSJed Brown 
5c4762a1bSJed Brown     u_tt - \Delta u = 0
6c4762a1bSJed Brown 
7c4762a1bSJed Brown   which we split into two first-order systems
8c4762a1bSJed Brown 
9c4762a1bSJed Brown     u_t -     v    = 0    so that     F(u,v).u = v
10c4762a1bSJed Brown     v_t - \Delta u = 0    so that     F(u,v).v = Delta u
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
13c4762a1bSJed Brown    Include "petscts.h" so that we can use SNES solvers.  Note that this
14c4762a1bSJed Brown    file automatically includes:
15c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
16c4762a1bSJed Brown      petscmat.h - matrices
17c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
18c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
19c4762a1bSJed Brown      petscksp.h   - linear solvers
20c4762a1bSJed Brown */
21c4762a1bSJed Brown #include <petscdm.h>
22c4762a1bSJed Brown #include <petscdmda.h>
23c4762a1bSJed Brown #include <petscts.h>
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    User-defined routines
27c4762a1bSJed Brown */
28c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec);
29c4762a1bSJed Brown extern PetscErrorCode MyTSMonitor(TS,PetscInt,PetscReal,Vec,void*);
30c4762a1bSJed Brown extern PetscErrorCode MySNESMonitor(SNES,PetscInt,PetscReal,PetscViewerAndFormat*);
31c4762a1bSJed Brown 
32c4762a1bSJed Brown int main(int argc,char **argv)
33c4762a1bSJed Brown {
34c4762a1bSJed Brown   TS                   ts;                         /* nonlinear solver */
35c4762a1bSJed Brown   Vec                  x,r;                        /* solution, residual vectors */
36c4762a1bSJed Brown   PetscInt             steps;                      /* iterations for convergence */
37c4762a1bSJed Brown   DM                   da;
38c4762a1bSJed Brown   PetscReal            ftime;
39c4762a1bSJed Brown   SNES                 ts_snes;
40c4762a1bSJed Brown   PetscBool            usemonitor = PETSC_TRUE;
41c4762a1bSJed Brown   PetscViewerAndFormat *vf;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
44c4762a1bSJed Brown      Initialize program
45c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
469566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
479566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-usemonitor",&usemonitor,NULL));
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
51c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
529566063dSJacob Faibussowitsch   PetscCall(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da));
539566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
549566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
559566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,0,"u"));
569566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,1,"v"));
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
60c4762a1bSJed Brown      vectors that are the same types
61c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
629566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&x));
639566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4762a1bSJed Brown      Create timestepping solver context
67c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
689566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
699566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,da));
709566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
719566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,FormFunction,da));
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,1.0));
74*1baa6e33SBarry Smith   if (usemonitor) PetscCall(TSMonitorSet(ts,MyTSMonitor,0,0));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77c4762a1bSJed Brown      Customize nonlinear solver
78c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
799566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSBEULER));
809566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&ts_snes));
81c4762a1bSJed Brown   if (usemonitor) {
829566063dSJacob Faibussowitsch     PetscCall(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
839566063dSJacob Faibussowitsch     PetscCall(SNESMonitorSet(ts_snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void *))MySNESMonitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86c4762a1bSJed Brown      Set initial conditions
87c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
889566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(da,x));
899566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.0001));
909566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
919566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,x));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown      Set runtime options
95c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
969566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99c4762a1bSJed Brown      Solve nonlinear system
100c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1019566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
1029566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
1039566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
1049566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(x,NULL,"-final_sol"));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
107c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
108c4762a1bSJed Brown      are no longer needed.
109c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1119566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
1129566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1139566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
114c4762a1bSJed Brown 
1159566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
116b122ec5aSJacob Faibussowitsch   return 0;
117c4762a1bSJed Brown }
118c4762a1bSJed Brown /* ------------------------------------------------------------------- */
119c4762a1bSJed Brown /*
120c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
121c4762a1bSJed Brown 
122c4762a1bSJed Brown    Input Parameters:
123c4762a1bSJed Brown .  ts - the TS context
124c4762a1bSJed Brown .  X - input vector
125c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
126c4762a1bSJed Brown 
127c4762a1bSJed Brown    Output Parameter:
128c4762a1bSJed Brown .  F - function vector
129c4762a1bSJed Brown  */
130c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
131c4762a1bSJed Brown {
132c4762a1bSJed Brown   DM             da = (DM)ptr;
133c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
134c4762a1bSJed Brown   PetscReal      hx,hy,/*hxdhy,hydhx,*/ sx,sy;
135c4762a1bSJed Brown   PetscScalar    u,uxx,uyy,v,***x,***f;
136c4762a1bSJed Brown   Vec            localX;
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   PetscFunctionBeginUser;
1399566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localX));
1409566063dSJacob Faibussowitsch   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));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
143c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
144c4762a1bSJed Brown   /*hxdhy  = hx/hy;*/
145c4762a1bSJed Brown   /*hydhx  = hy/hx;*/
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /*
148c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
149c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
150c4762a1bSJed Brown      By placing code between these two statements, computations can be
151c4762a1bSJed Brown      done while messages are in transition.
152c4762a1bSJed Brown   */
1539566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
1549566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /*
157c4762a1bSJed Brown      Get pointers to vector data
158c4762a1bSJed Brown   */
1599566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(da,localX,&x));
1609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(da,F,&f));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /*
163c4762a1bSJed Brown      Get local grid boundaries
164c4762a1bSJed Brown   */
1659566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   /*
168c4762a1bSJed Brown      Compute function over the locally owned part of the grid
169c4762a1bSJed Brown   */
170c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
171c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
172c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
173c4762a1bSJed Brown         f[j][i][0] = x[j][i][0];
174c4762a1bSJed Brown         f[j][i][1] = x[j][i][1];
175c4762a1bSJed Brown         continue;
176c4762a1bSJed Brown       }
177c4762a1bSJed Brown       u          = x[j][i][0];
178c4762a1bSJed Brown       v          = x[j][i][1];
179c4762a1bSJed Brown       uxx        = (-2.0*u + x[j][i-1][0] + x[j][i+1][0])*sx;
180c4762a1bSJed Brown       uyy        = (-2.0*u + x[j-1][i][0] + x[j+1][i][0])*sy;
181c4762a1bSJed Brown       f[j][i][0] = v;
182c4762a1bSJed Brown       f[j][i][1] = uxx + uyy;
183c4762a1bSJed Brown     }
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /*
187c4762a1bSJed Brown      Restore vectors
188c4762a1bSJed Brown   */
1899566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da,localX,&x));
1909566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da,F,&f));
1919566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localX));
1929566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(11.0*ym*xm));
193c4762a1bSJed Brown   PetscFunctionReturn(0);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown /* ------------------------------------------------------------------- */
197c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U)
198c4762a1bSJed Brown {
199c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
200c4762a1bSJed Brown   PetscScalar    ***u;
201c4762a1bSJed Brown   PetscReal      hx,hy,x,y,r;
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   PetscFunctionBeginUser;
2049566063dSJacob Faibussowitsch   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));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
207c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1);
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /*
210c4762a1bSJed Brown      Get pointers to vector data
211c4762a1bSJed Brown   */
2129566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayDOF(da,U,&u));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /*
215c4762a1bSJed Brown      Get local grid boundaries
216c4762a1bSJed Brown   */
2179566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /*
220c4762a1bSJed Brown      Compute function over the locally owned part of the grid
221c4762a1bSJed Brown   */
222c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
223c4762a1bSJed Brown     y = j*hy;
224c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
225c4762a1bSJed Brown       x = i*hx;
226c4762a1bSJed Brown       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
227c4762a1bSJed Brown       if (r < .125) {
228c4762a1bSJed Brown         u[j][i][0] = PetscExpReal(-30.0*r*r*r);
229c4762a1bSJed Brown         u[j][i][1] = 0.0;
230c4762a1bSJed Brown       } else {
231c4762a1bSJed Brown         u[j][i][0] = 0.0;
232c4762a1bSJed Brown         u[j][i][1] = 0.0;
233c4762a1bSJed Brown       }
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /*
238c4762a1bSJed Brown      Restore vectors
239c4762a1bSJed Brown   */
2409566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayDOF(da,U,&u));
241c4762a1bSJed Brown   PetscFunctionReturn(0);
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown PetscErrorCode MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
245c4762a1bSJed Brown {
246c4762a1bSJed Brown   PetscReal      norm;
247c4762a1bSJed Brown   MPI_Comm       comm;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   PetscFunctionBeginUser;
2509566063dSJacob Faibussowitsch   PetscCall(VecNorm(v,NORM_2,&norm));
2519566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)ts,&comm));
252c4762a1bSJed Brown   if (step > -1) { /* -1 is used to indicate an interpolated value */
25363a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(comm,"timestep %" PetscInt_FMT " time %g norm %g\n",step,(double)ptime,(double)norm));
254c4762a1bSJed Brown   }
255c4762a1bSJed Brown   PetscFunctionReturn(0);
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown /*
259c4762a1bSJed Brown    MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
260c4762a1bSJed Brown    Input Parameters:
261c4762a1bSJed Brown      snes - the SNES context
262c4762a1bSJed Brown      its - iteration number
263c4762a1bSJed Brown      fnorm - 2-norm function value (may be estimated)
264c4762a1bSJed Brown      ctx - optional user-defined context for private data for the
265c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
266c4762a1bSJed Brown  */
267c4762a1bSJed Brown PetscErrorCode MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *vf)
268c4762a1bSJed Brown {
269c4762a1bSJed Brown   PetscFunctionBeginUser;
2709566063dSJacob Faibussowitsch   PetscCall(SNESMonitorDefaultShort(snes,its,fnorm,vf));
271c4762a1bSJed Brown   PetscFunctionReturn(0);
272c4762a1bSJed Brown }
273c4762a1bSJed Brown /*TEST
274c4762a1bSJed Brown 
275c4762a1bSJed Brown     test:
276c4762a1bSJed Brown       args: -da_grid_x 20 -ts_max_time 3 -ts_dt 1e-1 -ts_theta_initial_guess_extrapolate 0 -ts_monitor -ksp_monitor_short
277c4762a1bSJed Brown       requires: !single
278c4762a1bSJed Brown 
279c4762a1bSJed Brown     test:
280c4762a1bSJed Brown       suffix: 2
281c4762a1bSJed Brown       args: -da_grid_x 20 -ts_max_time 0.11 -ts_dt 1e-1 -ts_type glle -ts_monitor -ksp_monitor_short
282c4762a1bSJed Brown       requires: !single
283c4762a1bSJed Brown 
284c4762a1bSJed Brown     test:
285c4762a1bSJed Brown       suffix: glvis_da_2d_vect
286c4762a1bSJed Brown       args: -usemonitor 0 -da_grid_x 20 -ts_max_time 0.3 -ts_dt 1e-1 -ts_type glle -final_sol glvis: -viewer_glvis_dm_da_bs 2,0
287c4762a1bSJed Brown       requires: !single
288c4762a1bSJed Brown 
289c4762a1bSJed Brown     test:
290c4762a1bSJed Brown       suffix: glvis_da_2d_vect_ll
291c4762a1bSJed Brown       args: -usemonitor 0 -da_grid_x 20 -ts_max_time 0.3 -ts_dt 1e-1 -ts_type glle -final_sol glvis: -viewer_glvis_dm_da_bs 2,0 -viewer_glvis_dm_da_ll
292c4762a1bSJed Brown       requires: !single
293c4762a1bSJed Brown 
294c4762a1bSJed Brown TEST*/
295