xref: /petsc/src/ts/tutorials/advection-diffusion-reaction/ex4.c (revision af7a0b9f4186d892160c7d2f0bec9dd32b55d72c)
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   PetscErrorCode ierr;
48c4762a1bSJed Brown   DM             da;
49c4762a1bSJed Brown   AppCtx         appctx;
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52c4762a1bSJed Brown      Initialize program
53c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
54c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
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 
67c4762a1bSJed Brown   ierr = PetscOptionsGetScalar(NULL,NULL,"-delta",&appctx.delta,NULL);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = PetscOptionsGetBool(NULL,NULL,"-upwind",&appctx.upwind,NULL);CHKERRQ(ierr);
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71c4762a1bSJed Brown      Create distributed array (DMDA) to manage parallel grid and vectors
72c4762a1bSJed Brown   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73c4762a1bSJed Brown   ierr = DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE,8,2,1,NULL,&da);CHKERRQ(ierr);
74c4762a1bSJed Brown   ierr = DMSetFromOptions(da);CHKERRQ(ierr);
75c4762a1bSJed Brown   ierr = DMSetUp(da);CHKERRQ(ierr);
76c4762a1bSJed Brown   ierr = DMDASetFieldName(da,0,"rho");CHKERRQ(ierr);
77c4762a1bSJed Brown   ierr = DMDASetFieldName(da,1,"c");CHKERRQ(ierr);
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    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
83c4762a1bSJed Brown   ierr = DMCreateGlobalVector(da,&U);CHKERRQ(ierr);
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
86c4762a1bSJed Brown      Create timestepping solver context
87c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
88c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
89c4762a1bSJed Brown   ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr);
90c4762a1bSJed Brown   ierr = TSSetDM(ts,da);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
92c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,IFunction,&appctx);CHKERRQ(ierr);
93c4762a1bSJed Brown 
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
96c4762a1bSJed Brown      Set initial conditions
97c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
98c4762a1bSJed Brown   ierr = InitialConditions(da,U);CHKERRQ(ierr);
99c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
102c4762a1bSJed Brown      Set solver options
103c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
104c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.0001);CHKERRQ(ierr);
105c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr);
106c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
107c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110c4762a1bSJed Brown      Solve nonlinear system
111c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
112c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
113c4762a1bSJed Brown 
114c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
116c4762a1bSJed Brown      are no longer needed.
117c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
118c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
120c4762a1bSJed Brown   ierr = DMDestroy(&da);CHKERRQ(ierr);
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   ierr = PetscFinalize();
123c4762a1bSJed Brown   return ierr;
124c4762a1bSJed Brown }
125c4762a1bSJed Brown /* ------------------------------------------------------------------- */
126c4762a1bSJed Brown /*
127c4762a1bSJed Brown    IFunction - Evaluates nonlinear function, F(U).
128c4762a1bSJed Brown 
129c4762a1bSJed Brown    Input Parameters:
130c4762a1bSJed Brown .  ts - the TS context
131c4762a1bSJed Brown .  U - input vector
132c4762a1bSJed Brown .  ptr - optional user-defined context, as set by SNESSetFunction()
133c4762a1bSJed Brown 
134c4762a1bSJed Brown    Output Parameter:
135c4762a1bSJed Brown .  F - function vector
136c4762a1bSJed Brown  */
137c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
138c4762a1bSJed Brown {
139c4762a1bSJed Brown   AppCtx         *appctx = (AppCtx*)ptr;
140c4762a1bSJed Brown   DM             da;
141c4762a1bSJed Brown   PetscErrorCode ierr;
142c4762a1bSJed Brown   PetscInt       i,Mx,xs,xm;
143c4762a1bSJed Brown   PetscReal      hx,sx;
144c4762a1bSJed Brown   PetscScalar    rho,c,rhoxx,cxx,cx,rhox,kcxrhox;
145c4762a1bSJed Brown   Field          *u,*f,*udot;
146c4762a1bSJed Brown   Vec            localU;
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   PetscFunctionBegin;
149c4762a1bSJed Brown   ierr = TSGetDM(ts,&da);CHKERRQ(ierr);
150c4762a1bSJed Brown   ierr = DMGetLocalVector(da,&localU);CHKERRQ(ierr);
151c4762a1bSJed 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);
152c4762a1bSJed Brown 
153c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /*
156c4762a1bSJed Brown      Scatter ghost points to local vector,using the 2-step process
157c4762a1bSJed Brown         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
158c4762a1bSJed Brown      By placing code between these two statements, computations can be
159c4762a1bSJed Brown      done while messages are in transition.
160c4762a1bSJed Brown   */
161c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
162c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);CHKERRQ(ierr);
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   /*
165c4762a1bSJed Brown      Get pointers to vector data
166c4762a1bSJed Brown   */
167c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,localU,&u);CHKERRQ(ierr);
168c4762a1bSJed Brown   ierr = DMDAVecGetArrayRead(da,Udot,&udot);CHKERRQ(ierr);
169303a5415SBarry Smith   ierr = DMDAVecGetArrayWrite(da,F,&f);CHKERRQ(ierr);
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /*
172c4762a1bSJed Brown      Get local grid boundaries
173c4762a1bSJed Brown   */
174c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   if (!xs) {
177c4762a1bSJed Brown     f[0].rho = udot[0].rho; /* u[0].rho - 0.0; */
178c4762a1bSJed Brown     f[0].c   = udot[0].c; /* u[0].c   - 1.0; */
179c4762a1bSJed Brown     xs++;
180c4762a1bSJed Brown     xm--;
181c4762a1bSJed Brown   }
182c4762a1bSJed Brown   if (xs+xm == Mx) {
183c4762a1bSJed Brown     f[Mx-1].rho = udot[Mx-1].rho; /* u[Mx-1].rho - 1.0; */
184c4762a1bSJed Brown     f[Mx-1].c   = udot[Mx-1].c;  /* u[Mx-1].c   - 0.0;  */
185c4762a1bSJed Brown     xm--;
186c4762a1bSJed Brown   }
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /*
189c4762a1bSJed Brown      Compute function over the locally owned part of the grid
190c4762a1bSJed Brown   */
191c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
192c4762a1bSJed Brown     rho   = u[i].rho;
193c4762a1bSJed Brown     rhoxx = (-2.0*rho + u[i-1].rho + u[i+1].rho)*sx;
194c4762a1bSJed Brown     c     = u[i].c;
195c4762a1bSJed Brown     cxx   = (-2.0*c + u[i-1].c + u[i+1].c)*sx;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown     if (!appctx->upwind) {
198c4762a1bSJed Brown       rhox    = .5*(u[i+1].rho - u[i-1].rho)/hx;
199c4762a1bSJed Brown       cx      = .5*(u[i+1].c - u[i-1].c)/hx;
200c4762a1bSJed Brown       kcxrhox = appctx->kappa*(cxx*rho + cx*rhox);
201c4762a1bSJed Brown     } else {
202c4762a1bSJed 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;
203c4762a1bSJed Brown     }
204c4762a1bSJed Brown 
205c4762a1bSJed 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;
206c4762a1bSJed Brown     f[i].c   = udot[i].c - appctx->delta*cxx + appctx->lambda*c + appctx->alpha*rho*c/(appctx->gamma + c);
207c4762a1bSJed Brown   }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /*
210c4762a1bSJed Brown      Restore vectors
211c4762a1bSJed Brown   */
212c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,localU,&u);CHKERRQ(ierr);
213c4762a1bSJed Brown   ierr = DMDAVecRestoreArrayRead(da,Udot,&udot);CHKERRQ(ierr);
214303a5415SBarry Smith   ierr = DMDAVecRestoreArrayWrite(da,F,&f);CHKERRQ(ierr);
215c4762a1bSJed Brown   ierr = DMRestoreLocalVector(da,&localU);CHKERRQ(ierr);
216c4762a1bSJed Brown   PetscFunctionReturn(0);
217c4762a1bSJed Brown }
218c4762a1bSJed Brown 
219c4762a1bSJed Brown /* ------------------------------------------------------------------- */
220c4762a1bSJed Brown PetscErrorCode InitialConditions(DM da,Vec U)
221c4762a1bSJed Brown {
222c4762a1bSJed Brown   PetscErrorCode ierr;
223c4762a1bSJed Brown   PetscInt       i,xs,xm,Mx;
224c4762a1bSJed Brown   Field          *u;
225c4762a1bSJed Brown   PetscReal      hx,x;
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   PetscFunctionBegin;
228c4762a1bSJed 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);
229c4762a1bSJed Brown 
230c4762a1bSJed Brown   hx = 1.0/(PetscReal)(Mx-1);
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /*
233c4762a1bSJed Brown      Get pointers to vector data
234c4762a1bSJed Brown   */
235303a5415SBarry Smith   ierr = DMDAVecGetArrayWrite(da,U,&u);CHKERRQ(ierr);
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /*
238c4762a1bSJed Brown      Get local grid boundaries
239c4762a1bSJed Brown   */
240c4762a1bSJed Brown   ierr = DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);CHKERRQ(ierr);
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /*
243c4762a1bSJed Brown      Compute function over the locally owned part of the grid
244c4762a1bSJed Brown   */
245c4762a1bSJed Brown   for (i=xs; i<xs+xm; i++) {
246c4762a1bSJed Brown     x = i*hx;
247*af7a0b9fSSatish Balay     if (i < Mx-1) u[i].rho = 0.0;
248c4762a1bSJed Brown     else          u[i].rho = 1.0;
249c4762a1bSJed Brown     u[i].c = PetscCosReal(.5*PETSC_PI*x);
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /*
253c4762a1bSJed Brown      Restore vectors
254c4762a1bSJed Brown   */
255303a5415SBarry Smith   ierr = DMDAVecRestoreArrayWrite(da,U,&u);CHKERRQ(ierr);
256c4762a1bSJed Brown   PetscFunctionReturn(0);
257c4762a1bSJed Brown }
258c4762a1bSJed Brown 
259c4762a1bSJed Brown 
260c4762a1bSJed Brown 
261c4762a1bSJed Brown 
262c4762a1bSJed Brown /*TEST
263c4762a1bSJed Brown 
264c4762a1bSJed Brown    test:
265c4762a1bSJed Brown       args: -pc_type lu -da_refine 2  -ts_view  -ts_monitor -ts_max_time 1
266c4762a1bSJed Brown       requires: double
267c4762a1bSJed Brown 
268c4762a1bSJed Brown    test:
269c4762a1bSJed Brown      suffix: 2
270c4762a1bSJed Brown      args:  -pc_type lu -da_refine 2  -ts_view  -ts_monitor_draw_solution -ts_monitor -ts_max_time 1
271c4762a1bSJed Brown      requires: x double
272c4762a1bSJed Brown      output_file: output/ex4_1.out
273c4762a1bSJed Brown 
274c4762a1bSJed Brown TEST*/
275