xref: /petsc/src/ts/tutorials/ex12.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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   PetscErrorCode       ierr;
38c4762a1bSJed Brown   DM                   da;
39c4762a1bSJed Brown   PetscReal            ftime;
40c4762a1bSJed Brown   SNES                 ts_snes;
41c4762a1bSJed Brown   PetscBool            usemonitor = PETSC_TRUE;
42c4762a1bSJed Brown   PetscViewerAndFormat *vf;
43c4762a1bSJed Brown 
44c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45c4762a1bSJed Brown      Initialize program
46c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
47c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-usemonitor",&usemonitor,NULL));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
52c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,NULL,&da));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,0,"u"));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDASetFieldName(da,1,"v"));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
61c4762a1bSJed Brown      vectors that are the same types
62c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
64*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4762a1bSJed Brown      Create timestepping solver context
68c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
70*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
71*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,FormFunction,da));
73c4762a1bSJed Brown 
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,1.0));
75c4762a1bSJed Brown   if (usemonitor) {
76*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,MyTSMonitor,0,0));
77c4762a1bSJed Brown   }
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80c4762a1bSJed Brown      Customize nonlinear solver
81c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
83*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSNES(ts,&ts_snes));
84c4762a1bSJed Brown   if (usemonitor) {
85*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
86*5f80ce2aSJacob Faibussowitsch     CHKERRQ(SNESMonitorSet(ts_snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void *))MySNESMonitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
87c4762a1bSJed Brown   }
88c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
89c4762a1bSJed Brown      Set initial conditions
90c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
91*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialSolution(da,x));
92*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.0001));
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Set runtime options
98c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102c4762a1bSJed Brown      Solve nonlinear system
103c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(x,NULL,"-final_sol"));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
111c4762a1bSJed Brown      are no longer needed.
112c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
113*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
114*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
115*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
116*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   ierr = PetscFinalize();
119c4762a1bSJed Brown   return ierr;
120c4762a1bSJed Brown }
121c4762a1bSJed Brown /* ------------------------------------------------------------------- */
122c4762a1bSJed Brown /*
123c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
124c4762a1bSJed Brown 
125c4762a1bSJed Brown    Input Parameters:
126c4762a1bSJed Brown .  ts - the TS context
127c4762a1bSJed Brown .  X - input vector
128c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
129c4762a1bSJed Brown 
130c4762a1bSJed Brown    Output Parameter:
131c4762a1bSJed Brown .  F - function vector
132c4762a1bSJed Brown  */
133c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
134c4762a1bSJed Brown {
135c4762a1bSJed Brown   DM             da = (DM)ptr;
136c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
137c4762a1bSJed Brown   PetscReal      hx,hy,/*hxdhy,hydhx,*/ sx,sy;
138c4762a1bSJed Brown   PetscScalar    u,uxx,uyy,v,***x,***f;
139c4762a1bSJed Brown   Vec            localX;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   PetscFunctionBeginUser;
142*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
146c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
147c4762a1bSJed Brown   /*hxdhy  = hx/hy;*/
148c4762a1bSJed Brown   /*hydhx  = hy/hx;*/
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /*
151c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
152c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
153c4762a1bSJed Brown      By placing code between these two statements, computations can be
154c4762a1bSJed Brown      done while messages are in transition.
155c4762a1bSJed Brown   */
156*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /*
160c4762a1bSJed Brown      Get pointers to vector data
161c4762a1bSJed Brown   */
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayDOF(da,localX,&x));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayDOF(da,F,&f));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /*
166c4762a1bSJed Brown      Get local grid boundaries
167c4762a1bSJed Brown   */
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /*
171c4762a1bSJed Brown      Compute function over the locally owned part of the grid
172c4762a1bSJed Brown   */
173c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
174c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
175c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
176c4762a1bSJed Brown         f[j][i][0] = x[j][i][0];
177c4762a1bSJed Brown         f[j][i][1] = x[j][i][1];
178c4762a1bSJed Brown         continue;
179c4762a1bSJed Brown       }
180c4762a1bSJed Brown       u          = x[j][i][0];
181c4762a1bSJed Brown       v          = x[j][i][1];
182c4762a1bSJed Brown       uxx        = (-2.0*u + x[j][i-1][0] + x[j][i+1][0])*sx;
183c4762a1bSJed Brown       uyy        = (-2.0*u + x[j-1][i][0] + x[j+1][i][0])*sy;
184c4762a1bSJed Brown       f[j][i][0] = v;
185c4762a1bSJed Brown       f[j][i][1] = uxx + uyy;
186c4762a1bSJed Brown     }
187c4762a1bSJed Brown   }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /*
190c4762a1bSJed Brown      Restore vectors
191c4762a1bSJed Brown   */
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayDOF(da,localX,&x));
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayDOF(da,F,&f));
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(11.0*ym*xm));
196c4762a1bSJed Brown   PetscFunctionReturn(0);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown /* ------------------------------------------------------------------- */
200c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U)
201c4762a1bSJed Brown {
202c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
203c4762a1bSJed Brown   PetscScalar    ***u;
204c4762a1bSJed Brown   PetscReal      hx,hy,x,y,r;
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   PetscFunctionBeginUser;
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
210c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1);
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /*
213c4762a1bSJed Brown      Get pointers to vector data
214c4762a1bSJed Brown   */
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayDOF(da,U,&u));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /*
218c4762a1bSJed Brown      Get local grid boundaries
219c4762a1bSJed Brown   */
220*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   /*
223c4762a1bSJed Brown      Compute function over the locally owned part of the grid
224c4762a1bSJed Brown   */
225c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
226c4762a1bSJed Brown     y = j*hy;
227c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
228c4762a1bSJed Brown       x = i*hx;
229c4762a1bSJed Brown       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
230c4762a1bSJed Brown       if (r < .125) {
231c4762a1bSJed Brown         u[j][i][0] = PetscExpReal(-30.0*r*r*r);
232c4762a1bSJed Brown         u[j][i][1] = 0.0;
233c4762a1bSJed Brown       } else {
234c4762a1bSJed Brown         u[j][i][0] = 0.0;
235c4762a1bSJed Brown         u[j][i][1] = 0.0;
236c4762a1bSJed Brown       }
237c4762a1bSJed Brown     }
238c4762a1bSJed Brown   }
239c4762a1bSJed Brown 
240c4762a1bSJed Brown   /*
241c4762a1bSJed Brown      Restore vectors
242c4762a1bSJed Brown   */
243*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayDOF(da,U,&u));
244c4762a1bSJed Brown   PetscFunctionReturn(0);
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown PetscErrorCode MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
248c4762a1bSJed Brown {
249c4762a1bSJed Brown   PetscReal      norm;
250c4762a1bSJed Brown   MPI_Comm       comm;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   PetscFunctionBeginUser;
253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(v,NORM_2,&norm));
254*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)ts,&comm));
255c4762a1bSJed Brown   if (step > -1) { /* -1 is used to indicate an interpolated value */
256*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscPrintf(comm,"timestep %D time %g norm %g\n",step,(double)ptime,(double)norm));
257c4762a1bSJed Brown   }
258c4762a1bSJed Brown   PetscFunctionReturn(0);
259c4762a1bSJed Brown }
260c4762a1bSJed Brown 
261c4762a1bSJed Brown /*
262c4762a1bSJed Brown    MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
263c4762a1bSJed Brown    Input Parameters:
264c4762a1bSJed Brown      snes - the SNES context
265c4762a1bSJed Brown      its - iteration number
266c4762a1bSJed Brown      fnorm - 2-norm function value (may be estimated)
267c4762a1bSJed Brown      ctx - optional user-defined context for private data for the
268c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
269c4762a1bSJed Brown  */
270c4762a1bSJed Brown PetscErrorCode MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *vf)
271c4762a1bSJed Brown {
272c4762a1bSJed Brown   PetscFunctionBeginUser;
273*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESMonitorDefaultShort(snes,its,fnorm,vf));
274c4762a1bSJed Brown   PetscFunctionReturn(0);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown /*TEST
277c4762a1bSJed Brown 
278c4762a1bSJed Brown     test:
279c4762a1bSJed 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
280c4762a1bSJed Brown       requires: !single
281c4762a1bSJed Brown 
282c4762a1bSJed Brown     test:
283c4762a1bSJed Brown       suffix: 2
284c4762a1bSJed Brown       args: -da_grid_x 20 -ts_max_time 0.11 -ts_dt 1e-1 -ts_type glle -ts_monitor -ksp_monitor_short
285c4762a1bSJed Brown       requires: !single
286c4762a1bSJed Brown 
287c4762a1bSJed Brown     test:
288c4762a1bSJed Brown       suffix: glvis_da_2d_vect
289c4762a1bSJed 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
290c4762a1bSJed Brown       requires: !single
291c4762a1bSJed Brown 
292c4762a1bSJed Brown     test:
293c4762a1bSJed Brown       suffix: glvis_da_2d_vect_ll
294c4762a1bSJed 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
295c4762a1bSJed Brown       requires: !single
296c4762a1bSJed Brown 
297c4762a1bSJed Brown TEST*/
298