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