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