xref: /petsc/src/ts/tutorials/ex25.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1a80e93e8SEmil static const char help[] = "Time-dependent Brusselator reaction-diffusion PDE in 1d formulated as a PDAE. Demonstrates solving PDEs with algebraic constraints (PDAE).\n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    u_t - alpha u_xx = A + u^2 v - (B+1) u
4c4762a1bSJed Brown    v_t - alpha v_xx = B u - u^2 v
5c4762a1bSJed Brown    0 < x < 1;
6c4762a1bSJed Brown    A = 1, B = 3, alpha = 1/50
7c4762a1bSJed Brown 
8c4762a1bSJed Brown    Initial conditions:
9c4762a1bSJed Brown    u(x,0) = 1 + sin(2 pi x)
10c4762a1bSJed Brown    v(x,0) = 3
11c4762a1bSJed Brown 
12c4762a1bSJed Brown    Boundary conditions:
13c4762a1bSJed Brown    u(0,t) = u(1,t) = 1
14c4762a1bSJed Brown    v(0,t) = v(1,t) = 3
15c4762a1bSJed Brown */
16c4762a1bSJed Brown 
17c4762a1bSJed Brown #include <petscdm.h>
18c4762a1bSJed Brown #include <petscdmda.h>
19c4762a1bSJed Brown #include <petscts.h>
20c4762a1bSJed Brown 
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscScalar u, v;
23c4762a1bSJed Brown } Field;
24c4762a1bSJed Brown 
25c4762a1bSJed Brown typedef struct _User *User;
26c4762a1bSJed Brown struct _User {
27c4762a1bSJed Brown   PetscReal A, B;          /* Reaction coefficients */
28c4762a1bSJed Brown   PetscReal alpha;         /* Diffusion coefficient */
29c4762a1bSJed Brown   PetscReal uleft, uright; /* Dirichlet boundary conditions */
30c4762a1bSJed Brown   PetscReal vleft, vright; /* Dirichlet boundary conditions */
31c4762a1bSJed Brown };
32c4762a1bSJed Brown 
33c4762a1bSJed Brown static PetscErrorCode FormRHSFunction(TS, PetscReal, Vec, Vec, void *);
34c4762a1bSJed Brown static PetscErrorCode FormIFunction(TS, PetscReal, Vec, Vec, Vec, void *);
35c4762a1bSJed Brown static PetscErrorCode FormIJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
36c4762a1bSJed Brown static PetscErrorCode FormInitialSolution(TS, Vec, void *);
37c4762a1bSJed Brown 
38*9371c9d4SSatish Balay int main(int argc, char **argv) {
39c4762a1bSJed Brown   TS                ts; /* nonlinear solver */
40c4762a1bSJed Brown   Vec               X;  /* solution, residual vectors */
41c4762a1bSJed Brown   Mat               J;  /* Jacobian matrix */
42c4762a1bSJed Brown   PetscInt          steps, mx;
43c4762a1bSJed Brown   DM                da;
44c4762a1bSJed Brown   PetscReal         ftime, hx, dt;
45c4762a1bSJed Brown   struct _User      user; /* user-defined work context */
46c4762a1bSJed Brown   TSConvergedReason reason;
47c4762a1bSJed Brown 
48327415f7SBarry Smith   PetscFunctionBeginUser;
499566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
50c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
52c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
539566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 11, 2, 2, NULL, &da));
549566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
559566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Extract global vectors from DMDA;
59c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
609566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &X));
61c4762a1bSJed Brown 
62c4762a1bSJed Brown   /* Initialize user application context */
63d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Advection-reaction options", "");
64c4762a1bSJed Brown   {
65c4762a1bSJed Brown     user.A      = 1;
66c4762a1bSJed Brown     user.B      = 3;
67c4762a1bSJed Brown     user.alpha  = 0.02;
68c4762a1bSJed Brown     user.uleft  = 1;
69c4762a1bSJed Brown     user.uright = 1;
70c4762a1bSJed Brown     user.vleft  = 3;
71c4762a1bSJed Brown     user.vright = 3;
729566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-A", "Reaction rate", "", user.A, &user.A, NULL));
739566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-B", "Reaction rate", "", user.B, &user.B, NULL));
749566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-alpha", "Diffusion coefficient", "", user.alpha, &user.alpha, NULL));
759566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-uleft", "Dirichlet boundary condition", "", user.uleft, &user.uleft, NULL));
769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-uright", "Dirichlet boundary condition", "", user.uright, &user.uright, NULL));
779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-vleft", "Dirichlet boundary condition", "", user.vleft, &user.vleft, NULL));
789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-vright", "Dirichlet boundary condition", "", user.vright, &user.vright, NULL));
79c4762a1bSJed Brown   }
80d0609cedSBarry Smith   PetscOptionsEnd();
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83c4762a1bSJed Brown      Create timestepping solver context
84c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
859566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
869566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
879566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSARKIMEX));
889566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts, TS_EQ_DAE_IMPLICIT_INDEX1));
899566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, FormRHSFunction, &user));
909566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, FormIFunction, &user));
919566063dSJacob Faibussowitsch   PetscCall(DMSetMatType(da, MATAIJ));
929566063dSJacob Faibussowitsch   PetscCall(DMCreateMatrix(da, &J));
939566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, J, J, FormIJacobian, &user));
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   ftime = 10.0;
969566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
979566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown      Set initial conditions
101c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1029566063dSJacob Faibussowitsch   PetscCall(FormInitialSolution(ts, X, &user));
1039566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
1049566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &mx));
105c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(mx / 2 - 1);
106c4762a1bSJed Brown   dt = 0.4 * PetscSqr(hx) / user.alpha; /* Diffusive stability limit */
1079566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown      Set runtime options
111c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1129566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown      Solve nonlinear system
116c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1179566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
1189566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
1199566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
1209566063dSJacob Faibussowitsch   PetscCall(TSGetConvergedReason(ts, &reason));
12163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%s at time %g after %" PetscInt_FMT " steps\n", TSConvergedReasons[reason], (double)ftime, steps));
122c4762a1bSJed Brown 
123c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
124c4762a1bSJed Brown      Free work space.
125c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1269566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
1279566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
1289566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1299566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
1309566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
131b122ec5aSJacob Faibussowitsch   return 0;
132c4762a1bSJed Brown }
133c4762a1bSJed Brown 
134*9371c9d4SSatish Balay static PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ptr) {
135c4762a1bSJed Brown   User          user = (User)ptr;
136c4762a1bSJed Brown   DM            da;
137c4762a1bSJed Brown   DMDALocalInfo info;
138c4762a1bSJed Brown   PetscInt      i;
139c4762a1bSJed Brown   Field        *x, *xdot, *f;
140c4762a1bSJed Brown   PetscReal     hx;
141c4762a1bSJed Brown   Vec           Xloc;
142c4762a1bSJed Brown 
143c4762a1bSJed Brown   PetscFunctionBeginUser;
1449566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1459566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
146c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info.mx - 1);
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /*
149c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
150c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
151c4762a1bSJed Brown      By placing code between these two statements, computations can be
152c4762a1bSJed Brown      done while messages are in transition.
153c4762a1bSJed Brown   */
1549566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &Xloc));
1559566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, X, INSERT_VALUES, Xloc));
1569566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, X, INSERT_VALUES, Xloc));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* Get pointers to vector data */
1599566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Xloc, &x));
1609566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Xdot, &xdot));
1619566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
164c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
165c4762a1bSJed Brown     if (i == 0) {
166c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uleft);
167c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vleft);
168c4762a1bSJed Brown     } else if (i == info.mx - 1) {
169c4762a1bSJed Brown       f[i].u = hx * (x[i].u - user->uright);
170c4762a1bSJed Brown       f[i].v = hx * (x[i].v - user->vright);
171c4762a1bSJed Brown     } else {
172c4762a1bSJed Brown       f[i].u = hx * xdot[i].u - user->alpha * (x[i - 1].u - 2. * x[i].u + x[i + 1].u) / hx;
173c4762a1bSJed Brown       f[i].v = hx * xdot[i].v - user->alpha * (x[i - 1].v - 2. * x[i].v + x[i + 1].v) / hx;
174c4762a1bSJed Brown     }
175c4762a1bSJed Brown   }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* Restore vectors */
1789566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Xloc, &x));
1799566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Xdot, &xdot));
1809566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
1819566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &Xloc));
182c4762a1bSJed Brown   PetscFunctionReturn(0);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185*9371c9d4SSatish Balay static PetscErrorCode FormRHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ptr) {
186c4762a1bSJed Brown   User          user = (User)ptr;
187c4762a1bSJed Brown   DM            da;
188c4762a1bSJed Brown   DMDALocalInfo info;
189c4762a1bSJed Brown   PetscInt      i;
190c4762a1bSJed Brown   PetscReal     hx;
191c4762a1bSJed Brown   Field        *x, *f;
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   PetscFunctionBeginUser;
1949566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1959566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
196c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info.mx - 1);
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* Get pointers to vector data */
1999566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, X, &x));
2009566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, F, &f));
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
203c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
204c4762a1bSJed Brown     PetscScalar u = x[i].u, v = x[i].v;
205c4762a1bSJed Brown     f[i].u = hx * (user->A + u * u * v - (user->B + 1) * u);
206c4762a1bSJed Brown     f[i].v = hx * (user->B * u - u * u * v);
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* Restore vectors */
2109566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, X, &x));
2119566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, F, &f));
212c4762a1bSJed Brown   PetscFunctionReturn(0);
213c4762a1bSJed Brown }
214c4762a1bSJed Brown 
215c4762a1bSJed Brown /* --------------------------------------------------------------------- */
216c4762a1bSJed Brown /*
217c4762a1bSJed Brown   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
218c4762a1bSJed Brown */
219*9371c9d4SSatish Balay PetscErrorCode FormIJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ptr) {
220c4762a1bSJed Brown   User          user = (User)ptr;
221c4762a1bSJed Brown   DMDALocalInfo info;
222c4762a1bSJed Brown   PetscInt      i;
223c4762a1bSJed Brown   PetscReal     hx;
224c4762a1bSJed Brown   DM            da;
225c4762a1bSJed Brown   Field        *x, *xdot;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   PetscFunctionBeginUser;
2289566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2299566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
230c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info.mx - 1);
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* Get pointers to vector data */
2339566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, X, &x));
2349566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Xdot, &xdot));
235c4762a1bSJed Brown 
236c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
237c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
238c4762a1bSJed Brown     if (i == 0 || i == info.mx - 1) {
239c4762a1bSJed Brown       const PetscInt    row = i, col = i;
240*9371c9d4SSatish Balay       const PetscScalar vals[2][2] = {
241*9371c9d4SSatish Balay         {hx, 0 },
242*9371c9d4SSatish Balay         {0,  hx}
243*9371c9d4SSatish Balay       };
2449566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 1, &col, &vals[0][0], INSERT_VALUES));
245c4762a1bSJed Brown     } else {
246c4762a1bSJed Brown       const PetscInt    row = i, col[] = {i - 1, i, i + 1};
247c4762a1bSJed Brown       const PetscScalar dxxL = -user->alpha / hx, dxx0 = 2. * user->alpha / hx, dxxR = -user->alpha / hx;
248*9371c9d4SSatish Balay       const PetscScalar vals[2][3][2] = {
249*9371c9d4SSatish Balay         {{dxxL, 0}, {a * hx + dxx0, 0}, {dxxR, 0}},
250*9371c9d4SSatish Balay         {{0, dxxL}, {0, a * hx + dxx0}, {0, dxxR}}
251*9371c9d4SSatish Balay       };
2529566063dSJacob Faibussowitsch       PetscCall(MatSetValuesBlocked(Jpre, 1, &row, 3, col, &vals[0][0][0], INSERT_VALUES));
253c4762a1bSJed Brown     }
254c4762a1bSJed Brown   }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* Restore vectors */
2579566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, X, &x));
2589566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Xdot, &xdot));
259c4762a1bSJed Brown 
2609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
2619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
262c4762a1bSJed Brown   if (J != Jpre) {
2639566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2649566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
265c4762a1bSJed Brown   }
266c4762a1bSJed Brown   PetscFunctionReturn(0);
267c4762a1bSJed Brown }
268c4762a1bSJed Brown 
269*9371c9d4SSatish Balay PetscErrorCode FormInitialSolution(TS ts, Vec X, void *ctx) {
270c4762a1bSJed Brown   User          user = (User)ctx;
271c4762a1bSJed Brown   DM            da;
272c4762a1bSJed Brown   PetscInt      i;
273c4762a1bSJed Brown   DMDALocalInfo info;
274c4762a1bSJed Brown   Field        *x;
275c4762a1bSJed Brown   PetscReal     hx;
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   PetscFunctionBeginUser;
2789566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
2799566063dSJacob Faibussowitsch   PetscCall(DMDAGetLocalInfo(da, &info));
280c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(info.mx - 1);
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   /* Get pointers to vector data */
2839566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(da, X, &x));
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   /* Compute function over the locally owned part of the grid */
286c4762a1bSJed Brown   for (i = info.xs; i < info.xs + info.xm; i++) {
287c4762a1bSJed Brown     PetscReal xi = i * hx;
288c4762a1bSJed Brown     x[i].u       = user->uleft * (1. - xi) + user->uright * xi + PetscSinReal(2. * PETSC_PI * xi);
289c4762a1bSJed Brown     x[i].v       = user->vleft * (1. - xi) + user->vright * xi;
290c4762a1bSJed Brown   }
2919566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(da, X, &x));
292c4762a1bSJed Brown   PetscFunctionReturn(0);
293c4762a1bSJed Brown }
294c4762a1bSJed Brown 
295c4762a1bSJed Brown /*TEST
296c4762a1bSJed Brown 
297c4762a1bSJed Brown     test:
298c4762a1bSJed Brown       args: -nox -da_grid_x 20 -ts_monitor_draw_solution -ts_type rosw -ts_rosw_type 2p -ts_dt 5e-2 -ts_adapt_type none
299c4762a1bSJed Brown 
300c4762a1bSJed Brown TEST*/
301