xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex4.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Chemo-taxis Problems from Mathematical Biology.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown      Page 18, Chemo-taxis Problems from Mathematical Biology
6c4762a1bSJed Brown 
7c4762a1bSJed Brown         rho_t =
8c4762a1bSJed Brown         c_t   =
9c4762a1bSJed Brown 
10c4762a1bSJed Brown      Further discusson on Page 134 and in 2d on Page 409
11c4762a1bSJed Brown */
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown 
15c4762a1bSJed Brown    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
16c4762a1bSJed Brown    Include "petscts.h" so that we can use SNES solvers.  Note that this
17c4762a1bSJed Brown    file automatically includes:
18c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
19c4762a1bSJed Brown      petscmat.h - matrices
20c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
21c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
22c4762a1bSJed Brown      petscksp.h   - linear solvers
23c4762a1bSJed Brown */
24c4762a1bSJed Brown 
25c4762a1bSJed Brown #include <petscdm.h>
26c4762a1bSJed Brown #include <petscdmda.h>
27c4762a1bSJed Brown #include <petscts.h>
28c4762a1bSJed Brown 
29c4762a1bSJed Brown typedef struct {
30c4762a1bSJed Brown   PetscScalar rho, c;
31c4762a1bSJed Brown } Field;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown typedef struct {
34c4762a1bSJed Brown   PetscScalar epsilon, delta, alpha, beta, gamma, kappa, lambda, mu, cstar;
35c4762a1bSJed Brown   PetscBool   upwind;
36c4762a1bSJed Brown } AppCtx;
37c4762a1bSJed Brown 
38c4762a1bSJed Brown /*
39c4762a1bSJed Brown    User-defined routines
40c4762a1bSJed Brown */
41c4762a1bSJed Brown extern PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *), InitialConditions(DM, Vec);
42c4762a1bSJed Brown 
43*9371c9d4SSatish Balay int main(int argc, char **argv) {
44c4762a1bSJed Brown   TS     ts; /* nonlinear solver */
45c4762a1bSJed Brown   Vec    U;  /* solution, residual vectors */
46c4762a1bSJed Brown   DM     da;
47c4762a1bSJed Brown   AppCtx appctx;
48c4762a1bSJed Brown 
49c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50c4762a1bSJed Brown      Initialize program
51c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52327415f7SBarry Smith   PetscFunctionBeginUser;
539566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   appctx.epsilon = 1.0e-3;
56c4762a1bSJed Brown   appctx.delta   = 1.0;
57c4762a1bSJed Brown   appctx.alpha   = 10.0;
58c4762a1bSJed Brown   appctx.beta    = 4.0;
59c4762a1bSJed Brown   appctx.gamma   = 1.0;
60c4762a1bSJed Brown   appctx.kappa   = .75;
61c4762a1bSJed Brown   appctx.lambda  = 1.0;
62c4762a1bSJed Brown   appctx.mu      = 100.;
63c4762a1bSJed Brown   appctx.cstar   = .2;
64c4762a1bSJed Brown   appctx.upwind  = PETSC_TRUE;
65c4762a1bSJed Brown 
669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-delta", &appctx.delta, NULL));
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-upwind", &appctx.upwind, NULL));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
71c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
729566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, 8, 2, 1, NULL, &da));
739566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
749566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
759566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 0, "rho"));
769566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da, 1, "c"));
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
80c4762a1bSJed Brown      vectors that are the same types
81c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
829566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da, &U));
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
85c4762a1bSJed Brown      Create timestepping solver context
86c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
879566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
889566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
899566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, da));
909566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
919566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &appctx));
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94c4762a1bSJed Brown      Set initial conditions
95c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
969566063dSJacob Faibussowitsch   PetscCall(InitialConditions(da, U));
979566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown      Set solver options
101c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1029566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .0001));
1039566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.0));
1049566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1059566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
106c4762a1bSJed Brown 
107c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108c4762a1bSJed Brown      Solve nonlinear system
109c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1109566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
111c4762a1bSJed Brown 
112c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
114c4762a1bSJed Brown      are no longer needed.
115c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1169566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1179566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1189566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
119c4762a1bSJed Brown 
1209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
121b122ec5aSJacob Faibussowitsch   return 0;
122c4762a1bSJed Brown }
123c4762a1bSJed Brown /* ------------------------------------------------------------------- */
124c4762a1bSJed Brown /*
125c4762a1bSJed Brown    IFunction - Evaluates nonlinear function, F(U).
126c4762a1bSJed Brown 
127c4762a1bSJed Brown    Input Parameters:
128c4762a1bSJed Brown .  ts - the TS context
129c4762a1bSJed Brown .  U - input vector
130c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
131c4762a1bSJed Brown 
132c4762a1bSJed Brown    Output Parameter:
133c4762a1bSJed Brown .  F - function vector
134c4762a1bSJed Brown  */
135*9371c9d4SSatish Balay PetscErrorCode IFunction(TS ts, PetscReal ftime, Vec U, Vec Udot, Vec F, void *ptr) {
136c4762a1bSJed Brown   AppCtx     *appctx = (AppCtx *)ptr;
137c4762a1bSJed Brown   DM          da;
138c4762a1bSJed Brown   PetscInt    i, Mx, xs, xm;
139c4762a1bSJed Brown   PetscReal   hx, sx;
140c4762a1bSJed Brown   PetscScalar rho, c, rhoxx, cxx, cx, rhox, kcxrhox;
141c4762a1bSJed Brown   Field      *u, *f, *udot;
142c4762a1bSJed Brown   Vec         localU;
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   PetscFunctionBegin;
1459566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &da));
1469566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da, &localU));
1479566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
148c4762a1bSJed Brown 
149*9371c9d4SSatish Balay   hx = 1.0 / (PetscReal)(Mx - 1);
150*9371c9d4SSatish Balay   sx = 1.0 / (hx * hx);
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /*
153c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
154c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
155c4762a1bSJed Brown      By placing code between these two statements, computations can be
156c4762a1bSJed Brown      done while messages are in transition.
157c4762a1bSJed Brown   */
1589566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da, U, INSERT_VALUES, localU));
1599566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da, U, INSERT_VALUES, localU));
160c4762a1bSJed Brown 
161c4762a1bSJed Brown   /*
162c4762a1bSJed Brown      Get pointers to vector data
163c4762a1bSJed Brown   */
1649566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, localU, &u));
1659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da, Udot, &udot));
1669566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayWrite(da, F, &f));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /*
169c4762a1bSJed Brown      Get local grid boundaries
170c4762a1bSJed Brown   */
1719566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   if (!xs) {
174c4762a1bSJed Brown     f[0].rho = udot[0].rho; /* u[0].rho - 0.0; */
175c4762a1bSJed Brown     f[0].c   = udot[0].c;   /* u[0].c   - 1.0; */
176c4762a1bSJed Brown     xs++;
177c4762a1bSJed Brown     xm--;
178c4762a1bSJed Brown   }
179c4762a1bSJed Brown   if (xs + xm == Mx) {
180c4762a1bSJed Brown     f[Mx - 1].rho = udot[Mx - 1].rho; /* u[Mx-1].rho - 1.0; */
181c4762a1bSJed Brown     f[Mx - 1].c   = udot[Mx - 1].c;   /* u[Mx-1].c   - 0.0;  */
182c4762a1bSJed Brown     xm--;
183c4762a1bSJed Brown   }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   /*
186c4762a1bSJed Brown      Compute function over the locally owned part of the grid
187c4762a1bSJed Brown   */
188c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
189c4762a1bSJed Brown     rho   = u[i].rho;
190c4762a1bSJed Brown     rhoxx = (-2.0 * rho + u[i - 1].rho + u[i + 1].rho) * sx;
191c4762a1bSJed Brown     c     = u[i].c;
192c4762a1bSJed Brown     cxx   = (-2.0 * c + u[i - 1].c + u[i + 1].c) * sx;
193c4762a1bSJed Brown 
194c4762a1bSJed Brown     if (!appctx->upwind) {
195c4762a1bSJed Brown       rhox    = .5 * (u[i + 1].rho - u[i - 1].rho) / hx;
196c4762a1bSJed Brown       cx      = .5 * (u[i + 1].c - u[i - 1].c) / hx;
197c4762a1bSJed Brown       kcxrhox = appctx->kappa * (cxx * rho + cx * rhox);
198c4762a1bSJed Brown     } else {
199c4762a1bSJed Brown       kcxrhox = appctx->kappa * ((u[i + 1].c - u[i].c) * u[i + 1].rho - (u[i].c - u[i - 1].c) * u[i].rho) * sx;
200c4762a1bSJed Brown     }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown     f[i].rho = udot[i].rho - appctx->epsilon * rhoxx + kcxrhox - appctx->mu * PetscAbsScalar(rho) * (1.0 - rho) * PetscMax(0, PetscRealPart(c - appctx->cstar)) + appctx->beta * rho;
203c4762a1bSJed Brown     f[i].c   = udot[i].c - appctx->delta * cxx + appctx->lambda * c + appctx->alpha * rho * c / (appctx->gamma + c);
204c4762a1bSJed Brown   }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /*
207c4762a1bSJed Brown      Restore vectors
208c4762a1bSJed Brown   */
2099566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, localU, &u));
2109566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da, Udot, &udot));
2119566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayWrite(da, F, &f));
2129566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da, &localU));
213c4762a1bSJed Brown   PetscFunctionReturn(0);
214c4762a1bSJed Brown }
215c4762a1bSJed Brown 
216c4762a1bSJed Brown /* ------------------------------------------------------------------- */
217*9371c9d4SSatish Balay PetscErrorCode InitialConditions(DM da, Vec U) {
218c4762a1bSJed Brown   PetscInt  i, xs, xm, Mx;
219c4762a1bSJed Brown   Field    *u;
220c4762a1bSJed Brown   PetscReal hx, x;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   PetscFunctionBegin;
2239566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(da, PETSC_IGNORE, &Mx, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE, PETSC_IGNORE));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   hx = 1.0 / (PetscReal)(Mx - 1);
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /*
228c4762a1bSJed Brown      Get pointers to vector data
229c4762a1bSJed Brown   */
2309566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayWrite(da, U, &u));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /*
233c4762a1bSJed Brown      Get local grid boundaries
234c4762a1bSJed Brown   */
2359566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da, &xs, NULL, NULL, &xm, NULL, NULL));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /*
238c4762a1bSJed Brown      Compute function over the locally owned part of the grid
239c4762a1bSJed Brown   */
240c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
241c4762a1bSJed Brown     x = i * hx;
242af7a0b9fSSatish Balay     if (i < Mx - 1) u[i].rho = 0.0;
243c4762a1bSJed Brown     else u[i].rho = 1.0;
244c4762a1bSJed Brown     u[i].c = PetscCosReal(.5 * PETSC_PI * x);
245c4762a1bSJed Brown   }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /*
248c4762a1bSJed Brown      Restore vectors
249c4762a1bSJed Brown   */
2509566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayWrite(da, U, &u));
251c4762a1bSJed Brown   PetscFunctionReturn(0);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown 
254c4762a1bSJed Brown /*TEST
255c4762a1bSJed Brown 
256c4762a1bSJed Brown    test:
257c4762a1bSJed Brown       args: -pc_type lu -da_refine 2  -ts_view  -ts_monitor -ts_max_time 1
258c4762a1bSJed Brown       requires: double
259c4762a1bSJed Brown 
260c4762a1bSJed Brown    test:
261c4762a1bSJed Brown      suffix: 2
262c4762a1bSJed Brown      args:  -pc_type lu -da_refine 2  -ts_view  -ts_monitor_draw_solution -ts_monitor -ts_max_time 1
263c4762a1bSJed Brown      requires: x double
264c4762a1bSJed Brown      output_file: output/ex4_1.out
265c4762a1bSJed Brown 
266c4762a1bSJed Brown TEST*/
267