xref: /petsc/src/ts/tutorials/ex7.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Nonlinear, time-dependent PDE in 2d.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
6c4762a1bSJed Brown    Include "petscts.h" so that we can use SNES solvers.  Note that this
7c4762a1bSJed Brown    file automatically includes:
8c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
9c4762a1bSJed Brown      petscmat.h - matrices
10c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
11c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
12c4762a1bSJed Brown      petscksp.h   - linear solvers
13c4762a1bSJed Brown */
14c4762a1bSJed Brown #include <petscdm.h>
15c4762a1bSJed Brown #include <petscdmda.h>
16c4762a1bSJed Brown #include <petscts.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown    User-defined routines
20c4762a1bSJed Brown */
21c4762a1bSJed Brown extern PetscErrorCode FormFunction(TS,PetscReal,Vec,Vec,void*),FormInitialSolution(DM,Vec);
22c4762a1bSJed Brown extern PetscErrorCode MyTSMonitor(TS,PetscInt,PetscReal,Vec,void*);
23c4762a1bSJed Brown extern PetscErrorCode MySNESMonitor(SNES,PetscInt,PetscReal,PetscViewerAndFormat*);
24c4762a1bSJed Brown 
25c4762a1bSJed Brown int main(int argc,char **argv)
26c4762a1bSJed Brown {
27c4762a1bSJed Brown   TS                   ts;                         /* time integrator */
28c4762a1bSJed Brown   SNES                 snes;
29c4762a1bSJed Brown   Vec                  x,r;                        /* solution, residual vectors */
30c4762a1bSJed Brown   PetscErrorCode       ierr;
31c4762a1bSJed Brown   DM                   da;
32c4762a1bSJed Brown   PetscViewerAndFormat *vf;
33c4762a1bSJed Brown 
34c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35c4762a1bSJed Brown      Initialize program
36c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
38c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
40c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDACreate2d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,8,8,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&da));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetFromOptions(da));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(da));
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
46c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
47c4762a1bSJed Brown      vectors that are the same types
48c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMCreateGlobalVector(da,&x));
50*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
53c4762a1bSJed Brown      Create timestepping solver context
54c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
56*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
57*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,FormFunction,da));
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60c4762a1bSJed Brown      Create matrix data structure; set Jacobian evaluation routine
61c4762a1bSJed Brown 
62c4762a1bSJed Brown      Set Jacobian matrix data structure and default Jacobian evaluation
63c4762a1bSJed Brown      routine. User can override with:
64c4762a1bSJed Brown      -snes_mf : matrix-free Newton-Krylov method with no preconditioning
65c4762a1bSJed Brown                 (unless user explicitly sets preconditioner)
66c4762a1bSJed Brown      -snes_mf_operator : form preconditioning matrix as set by the user,
67c4762a1bSJed Brown                          but use matrix-free approx for Jacobian-vector
68c4762a1bSJed Brown                          products within Newton-Krylov method
69c4762a1bSJed Brown 
70c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
71c4762a1bSJed Brown 
72*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,1.0));
73*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
74*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSMonitorSet(ts,MyTSMonitor,PETSC_VIEWER_STDOUT_WORLD,NULL));
75*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,da));
76c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77c4762a1bSJed Brown      Customize nonlinear solver
78c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
80*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSNES(ts,&snes));
81*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscViewerAndFormatCreate(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_DEFAULT,&vf));
82*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESMonitorSet(snes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))MySNESMonitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85c4762a1bSJed Brown      Set initial conditions
86c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
87*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialSolution(da,x));
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.0001));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
90c4762a1bSJed Brown 
91c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
92c4762a1bSJed Brown      Set runtime options
93c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97c4762a1bSJed Brown      Solve nonlinear system
98c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
103c4762a1bSJed Brown      are no longer needed.
104c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
105*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
106*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
107*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
108*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDestroy(&da));
109c4762a1bSJed Brown 
110c4762a1bSJed Brown   ierr = PetscFinalize();
111c4762a1bSJed Brown   return ierr;
112c4762a1bSJed Brown }
113c4762a1bSJed Brown /* ------------------------------------------------------------------- */
114c4762a1bSJed Brown /*
115c4762a1bSJed Brown    FormFunction - Evaluates nonlinear function, F(x).
116c4762a1bSJed Brown 
117c4762a1bSJed Brown    Input Parameters:
118c4762a1bSJed Brown .  ts - the TS context
119c4762a1bSJed Brown .  X - input vector
120c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
121c4762a1bSJed Brown 
122c4762a1bSJed Brown    Output Parameter:
123c4762a1bSJed Brown .  F - function vector
124c4762a1bSJed Brown  */
125c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal ftime,Vec X,Vec F,void *ptr)
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   DM             da;
128c4762a1bSJed Brown   PetscInt       i,j,Mx,My,xs,ys,xm,ym;
129c4762a1bSJed Brown   PetscReal      two = 2.0,hx,hy,sx,sy;
130c4762a1bSJed Brown   PetscScalar    u,uxx,uyy,**x,**f;
131c4762a1bSJed Brown   Vec            localX;
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   PetscFunctionBeginUser;
134*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetDM(ts,&da));
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGetLocalVector(da,&localX));
136*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));
137c4762a1bSJed Brown 
138c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
139c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1); sy = 1.0/(hy*hy);
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /*
142c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
143c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
144c4762a1bSJed Brown      By placing code between these two statements, computations can be
145c4762a1bSJed Brown      done while messages are in transition.
146c4762a1bSJed Brown   */
147*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalBegin(da,X,INSERT_VALUES,localX));
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMGlobalToLocalEnd(da,X,INSERT_VALUES,localX));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /*
151c4762a1bSJed Brown      Get pointers to vector data
152c4762a1bSJed Brown   */
153*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArrayRead(da,localX,&x));
154*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,F,&f));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   /*
157c4762a1bSJed Brown      Get local grid boundaries
158c4762a1bSJed Brown   */
159*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown      Compute function over the locally owned part of the grid
163c4762a1bSJed Brown   */
164c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
165c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
166c4762a1bSJed Brown       if (i == 0 || j == 0 || i == Mx-1 || j == My-1) {
167c4762a1bSJed Brown         f[j][i] = x[j][i];
168c4762a1bSJed Brown         continue;
169c4762a1bSJed Brown       }
170c4762a1bSJed Brown       u   = x[j][i];
171c4762a1bSJed Brown       uxx = (two*u - x[j][i-1] - x[j][i+1])*sx;
172c4762a1bSJed Brown       uyy = (two*u - x[j-1][i] - x[j+1][i])*sy;
173c4762a1bSJed Brown       /*      f[j][i] = -(uxx + uyy); */
174c4762a1bSJed Brown       f[j][i] = -u*(uxx + uyy) - (4.0 - 1.0)*((x[j][i+1] - x[j][i-1])*(x[j][i+1] - x[j][i-1])*.25*sx +
175c4762a1bSJed Brown                                               (x[j+1][i] - x[j-1][i])*(x[j+1][i] - x[j-1][i])*.25*sy);
176c4762a1bSJed Brown     }
177c4762a1bSJed Brown   }
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /*
180c4762a1bSJed Brown      Restore vectors
181c4762a1bSJed Brown   */
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArrayRead(da,localX,&x));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,F,&f));
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMRestoreLocalVector(da,&localX));
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscLogFlops(11.0*ym*xm));
186c4762a1bSJed Brown   PetscFunctionReturn(0);
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown /* ------------------------------------------------------------------- */
190c4762a1bSJed Brown PetscErrorCode FormInitialSolution(DM da,Vec U)
191c4762a1bSJed Brown {
192c4762a1bSJed Brown   PetscInt       i,j,xs,ys,xm,ym,Mx,My;
193c4762a1bSJed Brown   PetscScalar    **u;
194c4762a1bSJed Brown   PetscReal      hx,hy,x,y,r;
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   PetscFunctionBeginUser;
197*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));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
200c4762a1bSJed Brown   hy = 1.0/(PetscReal)(My-1);
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /*
203c4762a1bSJed Brown      Get pointers to vector data
204c4762a1bSJed Brown   */
205*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecGetArray(da,U,&u));
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /*
208c4762a1bSJed Brown      Get local grid boundaries
209c4762a1bSJed Brown   */
210*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAGetCorners(da,&xs,&ys,NULL,&xm,&ym,NULL));
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /*
213c4762a1bSJed Brown      Compute function over the locally owned part of the grid
214c4762a1bSJed Brown   */
215c4762a1bSJed Brown   for (j=ys; j<ys+ym; j++) {
216c4762a1bSJed Brown     y = j*hy;
217c4762a1bSJed Brown     for (i=xs; i<xs+xm; i++) {
218c4762a1bSJed Brown       x = i*hx;
219c4762a1bSJed Brown       r = PetscSqrtReal((x-.5)*(x-.5) + (y-.5)*(y-.5));
220c4762a1bSJed Brown       if (r < .125) u[j][i] = PetscExpReal(-30.0*r*r*r);
221c4762a1bSJed Brown       else          u[j][i] = 0.0;
222c4762a1bSJed Brown     }
223c4762a1bSJed Brown   }
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /*
226c4762a1bSJed Brown      Restore vectors
227c4762a1bSJed Brown   */
228*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMDAVecRestoreArray(da,U,&u));
229c4762a1bSJed Brown   PetscFunctionReturn(0);
230c4762a1bSJed Brown }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown PetscErrorCode MyTSMonitor(TS ts,PetscInt step,PetscReal ptime,Vec v,void *ctx)
233c4762a1bSJed Brown {
234c4762a1bSJed Brown   PetscReal      norm;
235c4762a1bSJed Brown   MPI_Comm       comm;
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   PetscFunctionBeginUser;
238c4762a1bSJed Brown   if (step < 0) PetscFunctionReturn(0); /* step of -1 indicates an interpolated solution */
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecNorm(v,NORM_2,&norm));
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscObjectGetComm((PetscObject)ts,&comm));
241*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(comm,"timestep %D time %g norm %g\n",step,(double)ptime,(double)norm));
242c4762a1bSJed Brown   PetscFunctionReturn(0);
243c4762a1bSJed Brown }
244c4762a1bSJed Brown 
245c4762a1bSJed Brown /*
246c4762a1bSJed Brown    MySNESMonitor - illustrate how to set user-defined monitoring routine for SNES.
247c4762a1bSJed Brown    Input Parameters:
248c4762a1bSJed Brown      snes - the SNES context
249c4762a1bSJed Brown      its - iteration number
250c4762a1bSJed Brown      fnorm - 2-norm function value (may be estimated)
251c4762a1bSJed Brown      ctx - optional user-defined context for private data for the
252c4762a1bSJed Brown          monitor routine, as set by SNESMonitorSet()
253c4762a1bSJed Brown  */
254c4762a1bSJed Brown PetscErrorCode MySNESMonitor(SNES snes,PetscInt its,PetscReal fnorm,PetscViewerAndFormat *vf)
255c4762a1bSJed Brown {
256c4762a1bSJed Brown   PetscFunctionBeginUser;
257*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESMonitorDefaultShort(snes,its,fnorm,vf));
258c4762a1bSJed Brown   PetscFunctionReturn(0);
259c4762a1bSJed Brown }
260c4762a1bSJed Brown 
261c4762a1bSJed Brown /*TEST
262c4762a1bSJed Brown 
263c4762a1bSJed Brown     test:
264c4762a1bSJed Brown       args: -ts_max_steps 5
265c4762a1bSJed Brown 
266c4762a1bSJed Brown     test:
267c4762a1bSJed Brown       suffix: 2
268c4762a1bSJed Brown       args: -ts_max_steps 5  -snes_mf_operator
269c4762a1bSJed Brown 
270c4762a1bSJed Brown     test:
271c4762a1bSJed Brown       suffix: 3
272c4762a1bSJed Brown       args: -ts_max_steps 5  -snes_mf -pc_type none
273c4762a1bSJed Brown 
274c4762a1bSJed Brown TEST*/
275