1c4762a1bSJed Brown static char help[] = "Solves the time independent Bratu problem using pseudo-timestepping.";
2c4762a1bSJed Brown
3c4762a1bSJed Brown /* ------------------------------------------------------------------------
4c4762a1bSJed Brown
5c4762a1bSJed Brown This code demonstrates how one may solve a nonlinear problem
6c4762a1bSJed Brown with pseudo-timestepping. In this simple example, the pseudo-timestep
7c4762a1bSJed Brown is the same for all grid points, i.e., this is equivalent to using
8c4762a1bSJed Brown the backward Euler method with a variable timestep.
9c4762a1bSJed Brown
10c4762a1bSJed Brown Note: This example does not require pseudo-timestepping since it
11c4762a1bSJed Brown is an easy nonlinear problem, but it is included to demonstrate how
12c4762a1bSJed Brown the pseudo-timestepping may be done.
13c4762a1bSJed Brown
14c4762a1bSJed Brown See snes/tutorials/ex4.c[ex4f.F] and
15c4762a1bSJed Brown snes/tutorials/ex5.c[ex5f.F] where the problem is described
16c4762a1bSJed Brown and solved using Newton's method alone.
17c4762a1bSJed Brown
18c4762a1bSJed Brown ----------------------------------------------------------------------------- */
19c4762a1bSJed Brown /*
20c4762a1bSJed Brown Include "petscts.h" to use the PETSc timestepping routines. Note that
21c4762a1bSJed Brown this file automatically includes "petscsys.h" and other lower-level
22c4762a1bSJed Brown PETSc include files.
23c4762a1bSJed Brown */
24c4762a1bSJed Brown #include <petscts.h>
25c4762a1bSJed Brown
26c4762a1bSJed Brown /*
27c4762a1bSJed Brown Create an application context to contain data needed by the
28c4762a1bSJed Brown application-provided call-back routines, FormJacobian() and
29c4762a1bSJed Brown FormFunction().
30c4762a1bSJed Brown */
31c4762a1bSJed Brown typedef struct {
32c4762a1bSJed Brown PetscReal param; /* test problem parameter */
33c4762a1bSJed Brown PetscInt mx; /* Discretization in x-direction */
34c4762a1bSJed Brown PetscInt my; /* Discretization in y-direction */
35c4762a1bSJed Brown } AppCtx;
36c4762a1bSJed Brown
37c4762a1bSJed Brown /*
38c4762a1bSJed Brown User-defined routines
39c4762a1bSJed Brown */
40c4762a1bSJed Brown extern PetscErrorCode FormJacobian(TS, PetscReal, Vec, Mat, Mat, void *), FormFunction(TS, PetscReal, Vec, Vec, void *), FormInitialGuess(Vec, AppCtx *);
41c4762a1bSJed Brown
main(int argc,char ** argv)42d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
43d71ae5a4SJacob Faibussowitsch {
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 */
65*3ba16761SJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-mx", &user.mx, NULL));
66*3ba16761SJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL, NULL, "-my", &user.my, NULL));
67c4762a1bSJed Brown N = user.mx * user.my;
68c4762a1bSJed Brown dt = .5 / PetscMax(user.mx, user.my);
69*3ba16761SJacob Faibussowitsch PetscCall(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
FormInitialGuess(Vec X,AppCtx * user)179d71ae5a4SJacob Faibussowitsch PetscErrorCode FormInitialGuess(Vec X, AppCtx *user)
180d71ae5a4SJacob Faibussowitsch {
181c4762a1bSJed Brown PetscInt i, j, row, mx, my;
182c4762a1bSJed Brown PetscReal one = 1.0, lambda;
183c4762a1bSJed Brown PetscReal temp1, temp, hx, hy;
184c4762a1bSJed Brown PetscScalar *x;
185c4762a1bSJed Brown
186*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
187c4762a1bSJed Brown mx = user->mx;
188c4762a1bSJed Brown my = user->my;
189c4762a1bSJed Brown lambda = user->param;
190c4762a1bSJed Brown
191c4762a1bSJed Brown hx = one / (PetscReal)(mx - 1);
192c4762a1bSJed Brown hy = one / (PetscReal)(my - 1);
193c4762a1bSJed Brown
1949566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x));
195c4762a1bSJed Brown temp1 = lambda / (lambda + one);
196c4762a1bSJed Brown for (j = 0; j < my; j++) {
197c4762a1bSJed Brown temp = (PetscReal)(PetscMin(j, my - j - 1)) * hy;
198c4762a1bSJed Brown for (i = 0; i < mx; i++) {
199c4762a1bSJed Brown row = i + j * mx;
200c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
201c4762a1bSJed Brown x[row] = 0.0;
202c4762a1bSJed Brown continue;
203c4762a1bSJed Brown }
204c4762a1bSJed Brown x[row] = temp1 * PetscSqrtReal(PetscMin((PetscReal)(PetscMin(i, mx - i - 1)) * hx, temp));
205c4762a1bSJed Brown }
206c4762a1bSJed Brown }
2079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(X, &x));
208*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
209c4762a1bSJed Brown }
210c4762a1bSJed Brown /* -------------------- Evaluate Function F(x) --------------------- */
211c4762a1bSJed Brown
FormFunction(TS ts,PetscReal t,Vec X,Vec F,void * ptr)212d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr)
213d71ae5a4SJacob Faibussowitsch {
214c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
215c4762a1bSJed Brown PetscInt i, j, row, mx, my;
216c4762a1bSJed Brown PetscReal two = 2.0, one = 1.0, lambda;
217c4762a1bSJed Brown PetscReal hx, hy, hxdhy, hydhx;
218c4762a1bSJed Brown PetscScalar ut, ub, ul, ur, u, uxx, uyy, sc, *f;
219c4762a1bSJed Brown const PetscScalar *x;
220c4762a1bSJed Brown
221*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
222c4762a1bSJed Brown mx = user->mx;
223c4762a1bSJed Brown my = user->my;
224c4762a1bSJed Brown lambda = user->param;
225c4762a1bSJed Brown
226c4762a1bSJed Brown hx = one / (PetscReal)(mx - 1);
227c4762a1bSJed Brown hy = one / (PetscReal)(my - 1);
228c4762a1bSJed Brown sc = hx * hy;
229c4762a1bSJed Brown hxdhy = hx / hy;
230c4762a1bSJed Brown hydhx = hy / hx;
231c4762a1bSJed Brown
2329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
2339566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
234c4762a1bSJed Brown for (j = 0; j < my; j++) {
235c4762a1bSJed Brown for (i = 0; i < mx; i++) {
236c4762a1bSJed Brown row = i + j * mx;
237c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
238c4762a1bSJed Brown f[row] = x[row];
239c4762a1bSJed Brown continue;
240c4762a1bSJed Brown }
241c4762a1bSJed Brown u = x[row];
242c4762a1bSJed Brown ub = x[row - mx];
243c4762a1bSJed Brown ul = x[row - 1];
244c4762a1bSJed Brown ut = x[row + mx];
245c4762a1bSJed Brown ur = x[row + 1];
246c4762a1bSJed Brown uxx = (-ur + two * u - ul) * hydhx;
247c4762a1bSJed Brown uyy = (-ut + two * u - ub) * hxdhy;
248c4762a1bSJed Brown f[row] = -uxx + -uyy + sc * lambda * PetscExpScalar(u);
249c4762a1bSJed Brown }
250c4762a1bSJed Brown }
2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
2529566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
253*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown /* -------------------- Evaluate Jacobian F'(x) -------------------- */
256c4762a1bSJed Brown
257c4762a1bSJed Brown /*
258c4762a1bSJed Brown Calculate the Jacobian matrix J(X,t).
259c4762a1bSJed Brown
260c4762a1bSJed Brown Note: We put the Jacobian in the preconditioner storage B instead of J. This
261c4762a1bSJed Brown way we can give the -snes_mf_operator flag to check our work. This replaces
262c4762a1bSJed Brown J with a finite difference approximation, using our analytic Jacobian B for
263c4762a1bSJed Brown the preconditioner.
264c4762a1bSJed Brown */
FormJacobian(TS ts,PetscReal t,Vec X,Mat J,Mat B,void * ptr)265d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian(TS ts, PetscReal t, Vec X, Mat J, Mat B, void *ptr)
266d71ae5a4SJacob Faibussowitsch {
267c4762a1bSJed Brown AppCtx *user = (AppCtx *)ptr;
268c4762a1bSJed Brown PetscInt i, j, row, mx, my, col[5];
269c4762a1bSJed Brown PetscScalar two = 2.0, one = 1.0, lambda, v[5], sc;
270c4762a1bSJed Brown const PetscScalar *x;
271c4762a1bSJed Brown PetscReal hx, hy, hxdhy, hydhx;
272c4762a1bSJed Brown
273*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser;
274c4762a1bSJed Brown mx = user->mx;
275c4762a1bSJed Brown my = user->my;
276c4762a1bSJed Brown lambda = user->param;
277c4762a1bSJed Brown
278c4762a1bSJed Brown hx = 1.0 / (PetscReal)(mx - 1);
279c4762a1bSJed Brown hy = 1.0 / (PetscReal)(my - 1);
280c4762a1bSJed Brown sc = hx * hy;
281c4762a1bSJed Brown hxdhy = hx / hy;
282c4762a1bSJed Brown hydhx = hy / hx;
283c4762a1bSJed Brown
2849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
285c4762a1bSJed Brown for (j = 0; j < my; j++) {
286c4762a1bSJed Brown for (i = 0; i < mx; i++) {
287c4762a1bSJed Brown row = i + j * mx;
288c4762a1bSJed Brown if (i == 0 || j == 0 || i == mx - 1 || j == my - 1) {
2899566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &row, 1, &row, &one, INSERT_VALUES));
290c4762a1bSJed Brown continue;
291c4762a1bSJed Brown }
2929371c9d4SSatish Balay v[0] = hxdhy;
2939371c9d4SSatish Balay col[0] = row - mx;
2949371c9d4SSatish Balay v[1] = hydhx;
2959371c9d4SSatish Balay col[1] = row - 1;
2969371c9d4SSatish Balay v[2] = -two * (hydhx + hxdhy) + sc * lambda * PetscExpScalar(x[row]);
2979371c9d4SSatish Balay col[2] = row;
2989371c9d4SSatish Balay v[3] = hydhx;
2999371c9d4SSatish Balay col[3] = row + 1;
3009371c9d4SSatish Balay v[4] = hxdhy;
3019371c9d4SSatish Balay col[4] = row + mx;
3029566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 1, &row, 5, col, v, INSERT_VALUES));
303c4762a1bSJed Brown }
304c4762a1bSJed Brown }
3059566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
3069566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
3079566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
308c4762a1bSJed Brown if (J != B) {
3099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
3109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
311c4762a1bSJed Brown }
312*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
313c4762a1bSJed Brown }
314c4762a1bSJed Brown
315c4762a1bSJed Brown /*TEST
316c4762a1bSJed Brown
317c4762a1bSJed Brown test:
318c4762a1bSJed 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
319c4762a1bSJed Brown
320c4762a1bSJed Brown test:
321c4762a1bSJed Brown suffix: 2
322c4762a1bSJed Brown args: -ts_monitor_pseudo -ts_pseudo_frtol 1.e-5
323c4762a1bSJed Brown
324c4762a1bSJed Brown TEST*/
325