xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex4.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
43c4762a1bSJed Brown int main(int argc,char **argv)
44c4762a1bSJed Brown {
45c4762a1bSJed Brown   TS             ts;                  /* nonlinear solver */
46c4762a1bSJed Brown   Vec            U;                   /* solution, residual vectors */
47c4762a1bSJed Brown   DM             da;
48c4762a1bSJed Brown   AppCtx         appctx;
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
51c4762a1bSJed Brown      Initialize program
52c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
53*327415f7SBarry Smith   PetscFunctionBeginUser;
549566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
55c4762a1bSJed Brown 
56c4762a1bSJed Brown   appctx.epsilon = 1.0e-3;
57c4762a1bSJed Brown   appctx.delta   = 1.0;
58c4762a1bSJed Brown   appctx.alpha   = 10.0;
59c4762a1bSJed Brown   appctx.beta    = 4.0;
60c4762a1bSJed Brown   appctx.gamma   = 1.0;
61c4762a1bSJed Brown   appctx.kappa   = .75;
62c4762a1bSJed Brown   appctx.lambda  = 1.0;
63c4762a1bSJed Brown   appctx.mu      = 100.;
64c4762a1bSJed Brown   appctx.cstar   = .2;
65c4762a1bSJed Brown   appctx.upwind  = PETSC_TRUE;
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-delta",&appctx.delta,NULL));
689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-upwind",&appctx.upwind,NULL));
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
72c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
739566063dSJacob Faibussowitsch   PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE,8,2,1,NULL,&da));
749566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(da));
759566063dSJacob Faibussowitsch   PetscCall(DMSetUp(da));
769566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,0,"rho"));
779566063dSJacob Faibussowitsch   PetscCall(DMDASetFieldName(da,1,"c"));
78c4762a1bSJed Brown 
79c4762a1bSJed Brown   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
80c4762a1bSJed Brown      Extract global vectors from DMDA; then duplicate for remaining
81c4762a1bSJed Brown      vectors that are the same types
82c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
839566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(da,&U));
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86c4762a1bSJed Brown      Create timestepping solver context
87c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
889566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
899566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSROSW));
909566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,da));
919566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
929566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,IFunction,&appctx));
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
95c4762a1bSJed Brown      Set initial conditions
96c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
979566063dSJacob Faibussowitsch   PetscCall(InitialConditions(da,U));
989566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,U));
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101c4762a1bSJed Brown      Set solver options
102c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1039566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.0001));
1049566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,1.0));
1059566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
1069566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109c4762a1bSJed Brown      Solve nonlinear system
110c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1119566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
115c4762a1bSJed Brown      are no longer needed.
116c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
1189566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1199566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&da));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
122b122ec5aSJacob Faibussowitsch   return 0;
123c4762a1bSJed Brown }
124c4762a1bSJed Brown /* ------------------------------------------------------------------- */
125c4762a1bSJed Brown /*
126c4762a1bSJed Brown    IFunction - Evaluates nonlinear function, F(U).
127c4762a1bSJed Brown 
128c4762a1bSJed Brown    Input Parameters:
129c4762a1bSJed Brown .  ts - the TS context
130c4762a1bSJed Brown .  U - input vector
131c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
132c4762a1bSJed Brown 
133c4762a1bSJed Brown    Output Parameter:
134c4762a1bSJed Brown .  F - function vector
135c4762a1bSJed Brown  */
136c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
137c4762a1bSJed Brown {
138c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
139c4762a1bSJed Brown   DM             da;
140c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
141c4762a1bSJed Brown   PetscReal      hx,sx;
142c4762a1bSJed Brown   PetscScalar    rho,c,rhoxx,cxx,cx,rhox,kcxrhox;
143c4762a1bSJed Brown   Field          *u,*f,*udot;
144c4762a1bSJed Brown   Vec            localU;
145c4762a1bSJed Brown 
146c4762a1bSJed Brown   PetscFunctionBegin;
1479566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts,&da));
1489566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(da,&localU));
1499566063dSJacob 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));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   /*
154c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
155c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
156c4762a1bSJed Brown      By placing code between these two statements, computations can be
157c4762a1bSJed Brown      done while messages are in transition.
158c4762a1bSJed Brown   */
1599566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU));
1609566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /*
163c4762a1bSJed Brown      Get pointers to vector data
164c4762a1bSJed Brown   */
1659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,localU,&u));
1669566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(da,Udot,&udot));
1679566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayWrite(da,F,&f));
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /*
170c4762a1bSJed Brown      Get local grid boundaries
171c4762a1bSJed Brown   */
1729566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   if (!xs) {
175c4762a1bSJed Brown     f[0].rho = udot[0].rho; /* u[0].rho - 0.0; */
176c4762a1bSJed Brown     f[0].c   = udot[0].c; /* u[0].c   - 1.0; */
177c4762a1bSJed Brown     xs++;
178c4762a1bSJed Brown     xm--;
179c4762a1bSJed Brown   }
180c4762a1bSJed Brown   if (xs+xm == Mx) {
181c4762a1bSJed Brown     f[Mx-1].rho = udot[Mx-1].rho; /* u[Mx-1].rho - 1.0; */
182c4762a1bSJed Brown     f[Mx-1].c   = udot[Mx-1].c;  /* u[Mx-1].c   - 0.0;  */
183c4762a1bSJed Brown     xm--;
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /*
187c4762a1bSJed Brown      Compute function over the locally owned part of the grid
188c4762a1bSJed Brown   */
189c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
190c4762a1bSJed Brown     rho   = u[i].rho;
191c4762a1bSJed Brown     rhoxx = (-2.0*rho + u[i-1].rho + u[i+1].rho)*sx;
192c4762a1bSJed Brown     c     = u[i].c;
193c4762a1bSJed Brown     cxx   = (-2.0*c + u[i-1].c + u[i+1].c)*sx;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown     if (!appctx->upwind) {
196c4762a1bSJed Brown       rhox    = .5*(u[i+1].rho - u[i-1].rho)/hx;
197c4762a1bSJed Brown       cx      = .5*(u[i+1].c - u[i-1].c)/hx;
198c4762a1bSJed Brown       kcxrhox = appctx->kappa*(cxx*rho + cx*rhox);
199c4762a1bSJed Brown     } else {
200c4762a1bSJed 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;
201c4762a1bSJed Brown     }
202c4762a1bSJed Brown 
203c4762a1bSJed 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;
204c4762a1bSJed Brown     f[i].c   = udot[i].c - appctx->delta*cxx + appctx->lambda*c + appctx->alpha*rho*c/(appctx->gamma + c);
205c4762a1bSJed Brown   }
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   /*
208c4762a1bSJed Brown      Restore vectors
209c4762a1bSJed Brown   */
2109566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,localU,&u));
2119566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(da,Udot,&udot));
2129566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayWrite(da,F,&f));
2139566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(da,&localU));
214c4762a1bSJed Brown   PetscFunctionReturn(0);
215c4762a1bSJed Brown }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown /* ------------------------------------------------------------------- */
218c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U)
219c4762a1bSJed Brown {
220c4762a1bSJed Brown   PetscInt       i,xs,xm,Mx;
221c4762a1bSJed Brown   Field          *u;
222c4762a1bSJed Brown   PetscReal      hx,x;
223c4762a1bSJed Brown 
224c4762a1bSJed Brown   PetscFunctionBegin;
2259566063dSJacob 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));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /*
230c4762a1bSJed Brown      Get pointers to vector data
231c4762a1bSJed Brown   */
2329566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayWrite(da,U,&u));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   /*
235c4762a1bSJed Brown      Get local grid boundaries
236c4762a1bSJed Brown   */
2379566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL));
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /*
240c4762a1bSJed Brown      Compute function over the locally owned part of the grid
241c4762a1bSJed Brown   */
242c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
243c4762a1bSJed Brown     x = i*hx;
244af7a0b9fSSatish Balay     if (i < Mx-1) u[i].rho = 0.0;
245c4762a1bSJed Brown     else          u[i].rho = 1.0;
246c4762a1bSJed Brown     u[i].c = PetscCosReal(.5*PETSC_PI*x);
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /*
250c4762a1bSJed Brown      Restore vectors
251c4762a1bSJed Brown   */
2529566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayWrite(da,U,&u));
253c4762a1bSJed Brown   PetscFunctionReturn(0);
254c4762a1bSJed Brown }
255c4762a1bSJed Brown 
256c4762a1bSJed Brown /*TEST
257c4762a1bSJed Brown 
258c4762a1bSJed Brown    test:
259c4762a1bSJed Brown       args: -pc_type lu -da_refine 2  -ts_view  -ts_monitor -ts_max_time 1
260c4762a1bSJed Brown       requires: double
261c4762a1bSJed Brown 
262c4762a1bSJed Brown    test:
263c4762a1bSJed Brown      suffix: 2
264c4762a1bSJed Brown      args:  -pc_type lu -da_refine 2  -ts_view  -ts_monitor_draw_solution -ts_monitor -ts_max_time 1
265c4762a1bSJed Brown      requires: x double
266c4762a1bSJed Brown      output_file: output/ex4_1.out
267c4762a1bSJed Brown 
268c4762a1bSJed Brown TEST*/
269