xref: /petsc/src/ts/tutorials/ex1.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266) !
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] ="Solves the time independent Bratu problem using pseudo-timestepping.";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown    Concepts: TS^pseudo-timestepping
6ee300463SSatish Balay    Concepts: TS^pseudo-timestepping
7c4762a1bSJed Brown    Concepts: TS^nonlinear problems
8c4762a1bSJed Brown    Processors: 1
9c4762a1bSJed Brown 
10c4762a1bSJed Brown */
11c4762a1bSJed Brown 
12c4762a1bSJed Brown /* ------------------------------------------------------------------------
13c4762a1bSJed Brown 
14c4762a1bSJed Brown     This code demonstrates how one may solve a nonlinear problem
15c4762a1bSJed Brown     with pseudo-timestepping. In this simple example, the pseudo-timestep
16c4762a1bSJed Brown     is the same for all grid points, i.e., this is equivalent to using
17c4762a1bSJed Brown     the backward Euler method with a variable timestep.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown     Note: This example does not require pseudo-timestepping since it
20c4762a1bSJed Brown     is an easy nonlinear problem, but it is included to demonstrate how
21c4762a1bSJed Brown     the pseudo-timestepping may be done.
22c4762a1bSJed Brown 
23c4762a1bSJed Brown     See snes/tutorials/ex4.c[ex4f.F] and
24c4762a1bSJed Brown     snes/tutorials/ex5.c[ex5f.F] where the problem is described
25c4762a1bSJed Brown     and solved using Newton's method alone.
26c4762a1bSJed Brown 
27c4762a1bSJed Brown   ----------------------------------------------------------------------------- */
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown     Include "petscts.h" to use the PETSc timestepping routines. Note that
30c4762a1bSJed Brown     this file automatically includes "petscsys.h" and other lower-level
31c4762a1bSJed Brown     PETSc include files.
32c4762a1bSJed Brown */
33c4762a1bSJed Brown #include <petscts.h>
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown   Create an application context to contain data needed by the
37c4762a1bSJed Brown   application-provided call-back routines, FormJacobian() and
38c4762a1bSJed Brown   FormFunction().
39c4762a1bSJed Brown */
40c4762a1bSJed Brown typedef struct {
41c4762a1bSJed Brown   PetscReal param;          /* test problem parameter */
42c4762a1bSJed Brown   PetscInt  mx;             /* Discretization in x-direction */
43c4762a1bSJed Brown   PetscInt  my;             /* Discretization in y-direction */
44c4762a1bSJed Brown } AppCtx;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown /*
47c4762a1bSJed Brown    User-defined routines
48c4762a1bSJed Brown */
49c4762a1bSJed Brown extern PetscErrorCode  FormJacobian(TS,PetscReal,Vec,Mat,Mat,void*), FormFunction(TS,PetscReal,Vec,Vec,void*), FormInitialGuess(Vec,AppCtx*);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown int main(int argc,char **argv)
52c4762a1bSJed Brown {
53c4762a1bSJed Brown   TS             ts;                 /* timestepping context */
54c4762a1bSJed Brown   Vec            x,r;               /* solution, residual vectors */
55c4762a1bSJed Brown   Mat            J;                  /* Jacobian matrix */
56c4762a1bSJed Brown   AppCtx         user;               /* user-defined work context */
57c4762a1bSJed Brown   PetscInt       its,N;                /* iterations for convergence */
58c4762a1bSJed Brown   PetscErrorCode ierr;
59c4762a1bSJed Brown   PetscReal      param_max = 6.81,param_min = 0.,dt;
60c4762a1bSJed Brown   PetscReal      ftime;
61c4762a1bSJed Brown   PetscMPIInt    size;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
64*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
653c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only");
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   user.mx    = 4;
68c4762a1bSJed Brown   user.my    = 4;
69c4762a1bSJed Brown   user.param = 6.0;
70c4762a1bSJed Brown 
71c4762a1bSJed Brown   /*
72c4762a1bSJed Brown      Allow user to set the grid dimensions and nonlinearity parameter at run-time
73c4762a1bSJed Brown   */
74c4762a1bSJed Brown   PetscOptionsGetInt(NULL,NULL,"-mx",&user.mx,NULL);
75c4762a1bSJed Brown   PetscOptionsGetInt(NULL,NULL,"-my",&user.my,NULL);
76c4762a1bSJed Brown   N  = user.mx*user.my;
77c4762a1bSJed Brown   dt = .5/PetscMax(user.mx,user.my);
78c4762a1bSJed Brown   PetscOptionsGetReal(NULL,NULL,"-param",&user.param,NULL);
793c633725SBarry Smith   PetscCheck(user.param < param_max && user.param >= param_min,PETSC_COMM_WORLD,PETSC_ERR_ARG_OUTOFRANGE,"Parameter is out of range");
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   /*
82c4762a1bSJed Brown       Create vectors to hold the solution and function value
83c4762a1bSJed Brown   */
84*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_SELF,N,&x));
85*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(x,&r));
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   /*
88c4762a1bSJed Brown     Create matrix to hold Jacobian. Preallocate 5 nonzeros per row
89c4762a1bSJed Brown     in the sparse matrix. Note that this is not the optimal strategy; see
90c4762a1bSJed Brown     the Performance chapter of the users manual for information on
91c4762a1bSJed Brown     preallocating memory in sparse matrices.
92c4762a1bSJed Brown   */
93*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,5,0,&J));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /*
96c4762a1bSJed Brown      Create timestepper context
97c4762a1bSJed Brown   */
98*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
99*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /*
102c4762a1bSJed Brown      Tell the timestepper context where to compute solutions
103c4762a1bSJed Brown   */
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,x));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   /*
107c4762a1bSJed Brown      Provide the call-back for the nonlinear function we are
108c4762a1bSJed Brown      evaluating. Thus whenever the timestepping routines need the
109c4762a1bSJed Brown      function they will call this routine. Note the final argument
110c4762a1bSJed Brown      is the application context used by the call-back functions.
111c4762a1bSJed Brown   */
112*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,FormFunction,&user));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /*
115c4762a1bSJed Brown      Set the Jacobian matrix and the function used to compute
116c4762a1bSJed Brown      Jacobians.
117c4762a1bSJed Brown   */
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,J,J,FormJacobian,&user));
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /*
121c4762a1bSJed Brown        Form the initial guess for the problem
122c4762a1bSJed Brown   */
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(FormInitialGuess(x,&user));
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /*
126c4762a1bSJed Brown        This indicates that we are using pseudo timestepping to
127c4762a1bSJed Brown      find a steady state solution to the nonlinear problem.
128c4762a1bSJed Brown   */
129*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSPSEUDO));
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /*
132c4762a1bSJed Brown        Set the initial time to start at (this is arbitrary for
133c4762a1bSJed Brown      steady state problems); and the initial timestep given above
134c4762a1bSJed Brown   */
135*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,dt));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /*
138c4762a1bSJed Brown       Set a large number of timesteps and final duration time
139c4762a1bSJed Brown      to insure convergence to steady state.
140c4762a1bSJed Brown   */
141*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxSteps(ts,10000));
142*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,1e12));
143*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /*
146c4762a1bSJed Brown       Use the default strategy for increasing the timestep
147c4762a1bSJed Brown   */
148*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSPseudoSetTimeStep(ts,TSPseudoTimeStepDefault,0));
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   /*
151c4762a1bSJed Brown       Set any additional options from the options database. This
152c4762a1bSJed Brown      includes all options for the nonlinear and linear solvers used
153c4762a1bSJed Brown      internally the timestepping routines.
154c4762a1bSJed Brown   */
155*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
156c4762a1bSJed Brown 
157*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetUp(ts));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /*
160c4762a1bSJed Brown       Perform the solve. This is where the timestepping takes place.
161c4762a1bSJed Brown   */
162*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
163*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
164c4762a1bSJed Brown 
165c4762a1bSJed Brown   /*
166c4762a1bSJed Brown       Get the number of steps
167c4762a1bSJed Brown   */
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&its));
169c4762a1bSJed Brown 
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"Number of pseudo timesteps = %D final time %4.2e\n",its,(double)ftime));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /*
173c4762a1bSJed Brown      Free the data structures constructed above
174c4762a1bSJed Brown   */
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
176*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
177*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&J));
178*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
179c4762a1bSJed Brown   ierr = PetscFinalize();
180c4762a1bSJed Brown   return ierr;
181c4762a1bSJed Brown }
182c4762a1bSJed Brown /* ------------------------------------------------------------------ */
183c4762a1bSJed Brown /*           Bratu (Solid Fuel Ignition) Test Problem                 */
184c4762a1bSJed Brown /* ------------------------------------------------------------------ */
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /* --------------------  Form initial approximation ----------------- */
187c4762a1bSJed Brown 
188c4762a1bSJed Brown PetscErrorCode FormInitialGuess(Vec X,AppCtx *user)
189c4762a1bSJed Brown {
190c4762a1bSJed Brown   PetscInt       i,j,row,mx,my;
191c4762a1bSJed Brown   PetscReal      one = 1.0,lambda;
192c4762a1bSJed Brown   PetscReal      temp1,temp,hx,hy;
193c4762a1bSJed Brown   PetscScalar    *x;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   mx     = user->mx;
196c4762a1bSJed Brown   my     = user->my;
197c4762a1bSJed Brown   lambda = user->param;
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   hx = one / (PetscReal)(mx-1);
200c4762a1bSJed Brown   hy = one / (PetscReal)(my-1);
201c4762a1bSJed Brown 
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(X,&x));
203c4762a1bSJed Brown   temp1 = lambda/(lambda + one);
204c4762a1bSJed Brown   for (j=0; j<my; j++) {
205c4762a1bSJed Brown     temp = (PetscReal)(PetscMin(j,my-j-1))*hy;
206c4762a1bSJed Brown     for (i=0; i<mx; i++) {
207c4762a1bSJed Brown       row = i + j*mx;
208c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
209c4762a1bSJed Brown         x[row] = 0.0;
210c4762a1bSJed Brown         continue;
211c4762a1bSJed Brown       }
212c4762a1bSJed Brown       x[row] = temp1*PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i,mx-i-1))*hx,temp));
213c4762a1bSJed Brown     }
214c4762a1bSJed Brown   }
215*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(X,&x));
216c4762a1bSJed Brown   return 0;
217c4762a1bSJed Brown }
218c4762a1bSJed Brown /* --------------------  Evaluate Function F(x) --------------------- */
219c4762a1bSJed Brown 
220c4762a1bSJed Brown PetscErrorCode FormFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
221c4762a1bSJed Brown {
222c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
223c4762a1bSJed Brown   PetscInt          i,j,row,mx,my;
224c4762a1bSJed Brown   PetscReal         two = 2.0,one = 1.0,lambda;
225c4762a1bSJed Brown   PetscReal         hx,hy,hxdhy,hydhx;
226c4762a1bSJed Brown   PetscScalar       ut,ub,ul,ur,u,uxx,uyy,sc,*f;
227c4762a1bSJed Brown   const PetscScalar *x;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   mx     = user->mx;
230c4762a1bSJed Brown   my     = user->my;
231c4762a1bSJed Brown   lambda = user->param;
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   hx    = one / (PetscReal)(mx-1);
234c4762a1bSJed Brown   hy    = one / (PetscReal)(my-1);
235c4762a1bSJed Brown   sc    = hx*hy;
236c4762a1bSJed Brown   hxdhy = hx/hy;
237c4762a1bSJed Brown   hydhx = hy/hx;
238c4762a1bSJed Brown 
239*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
241c4762a1bSJed Brown   for (j=0; j<my; j++) {
242c4762a1bSJed Brown     for (i=0; i<mx; i++) {
243c4762a1bSJed Brown       row = i + j*mx;
244c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
245c4762a1bSJed Brown         f[row] = x[row];
246c4762a1bSJed Brown         continue;
247c4762a1bSJed Brown       }
248c4762a1bSJed Brown       u      = x[row];
249c4762a1bSJed Brown       ub     = x[row - mx];
250c4762a1bSJed Brown       ul     = x[row - 1];
251c4762a1bSJed Brown       ut     = x[row + mx];
252c4762a1bSJed Brown       ur     = x[row + 1];
253c4762a1bSJed Brown       uxx    = (-ur + two*u - ul)*hydhx;
254c4762a1bSJed Brown       uyy    = (-ut + two*u - ub)*hxdhy;
255c4762a1bSJed Brown       f[row] = -uxx + -uyy + sc*lambda*PetscExpScalar(u);
256c4762a1bSJed Brown     }
257c4762a1bSJed Brown   }
258*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
259*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
260c4762a1bSJed Brown   return 0;
261c4762a1bSJed Brown }
262c4762a1bSJed Brown /* --------------------  Evaluate Jacobian F'(x) -------------------- */
263c4762a1bSJed Brown 
264c4762a1bSJed Brown /*
265c4762a1bSJed Brown    Calculate the Jacobian matrix J(X,t).
266c4762a1bSJed Brown 
267c4762a1bSJed Brown    Note: We put the Jacobian in the preconditioner storage B instead of J. This
268c4762a1bSJed Brown    way we can give the -snes_mf_operator flag to check our work. This replaces
269c4762a1bSJed Brown    J with a finite difference approximation, using our analytic Jacobian B for
270c4762a1bSJed Brown    the preconditioner.
271c4762a1bSJed Brown */
272c4762a1bSJed Brown PetscErrorCode FormJacobian(TS ts,PetscReal t,Vec X,Mat J,Mat B,void *ptr)
273c4762a1bSJed Brown {
274c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
275c4762a1bSJed Brown   PetscInt          i,j,row,mx,my,col[5];
276c4762a1bSJed Brown   PetscScalar       two = 2.0,one = 1.0,lambda,v[5],sc;
277c4762a1bSJed Brown   const PetscScalar *x;
278c4762a1bSJed Brown   PetscReal         hx,hy,hxdhy,hydhx;
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   mx     = user->mx;
281c4762a1bSJed Brown   my     = user->my;
282c4762a1bSJed Brown   lambda = user->param;
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   hx    = 1.0 / (PetscReal)(mx-1);
285c4762a1bSJed Brown   hy    = 1.0 / (PetscReal)(my-1);
286c4762a1bSJed Brown   sc    = hx*hy;
287c4762a1bSJed Brown   hxdhy = hx/hy;
288c4762a1bSJed Brown   hydhx = hy/hx;
289c4762a1bSJed Brown 
290*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
291c4762a1bSJed Brown   for (j=0; j<my; j++) {
292c4762a1bSJed Brown     for (i=0; i<mx; i++) {
293c4762a1bSJed Brown       row = i + j*mx;
294c4762a1bSJed Brown       if (i == 0 || j == 0 || i == mx-1 || j == my-1) {
295*5f80ce2aSJacob Faibussowitsch         CHKERRQ(MatSetValues(B,1,&row,1,&row,&one,INSERT_VALUES));
296c4762a1bSJed Brown         continue;
297c4762a1bSJed Brown       }
298c4762a1bSJed Brown       v[0] = hxdhy; col[0] = row - mx;
299c4762a1bSJed Brown       v[1] = hydhx; col[1] = row - 1;
300c4762a1bSJed Brown       v[2] = -two*(hydhx + hxdhy) + sc*lambda*PetscExpScalar(x[row]); col[2] = row;
301c4762a1bSJed Brown       v[3] = hydhx; col[3] = row + 1;
302c4762a1bSJed Brown       v[4] = hxdhy; col[4] = row + mx;
303*5f80ce2aSJacob Faibussowitsch       CHKERRQ(MatSetValues(B,1,&row,5,col,v,INSERT_VALUES));
304c4762a1bSJed Brown     }
305c4762a1bSJed Brown   }
306*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
307*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
308*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
309c4762a1bSJed Brown   if (J != B) {
310*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
311*5f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
312c4762a1bSJed Brown   }
313c4762a1bSJed Brown   return 0;
314c4762a1bSJed Brown }
315c4762a1bSJed Brown 
316c4762a1bSJed Brown /*TEST
317c4762a1bSJed Brown 
318c4762a1bSJed Brown     test:
319c4762a1bSJed Brown       args: -ksp_gmres_cgs_refinement_type refine_always -snes_type newtonls -ts_monitor_pseudo -snes_atol 1.e-7 -ts_pseudo_frtol 1.e-5 -ts_view draw:tikz:fig.tex
320c4762a1bSJed Brown 
321c4762a1bSJed Brown     test:
322c4762a1bSJed Brown       suffix: 2
323c4762a1bSJed Brown       args: -ts_monitor_pseudo -ts_pseudo_frtol 1.e-5
324c4762a1bSJed Brown 
325c4762a1bSJed Brown TEST*/
326