Actual source code: ex4.c

petsc-3.8.4 2018-03-24
Report Typos and Errors

  2: static char help[] = "Chemo-taxis Problems from Mathematical Biology.\n";

  4: /*
  5:      Page 18, Chemo-taxis Problems from Mathematical Biology

  7:         rho_t =
  8:         c_t   =

 10:      Further discusson on Page 134 and in 2d on Page 409
 11: */

 13: /*

 15:    Include "petscdmda.h" so that we can use distributed arrays (DMDAs).
 16:    Include "petscts.h" so that we can use SNES solvers.  Note that this
 17:    file automatically includes:
 18:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 19:      petscmat.h - matrices
 20:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 21:      petscviewer.h - viewers               petscpc.h  - preconditioners
 22:      petscksp.h   - linear solvers
 23: */
 24: #include <petscdm.h>
 25: #include <petscdmda.h>
 26: #include <petscts.h>

 28: typedef struct {
 29:   PetscScalar rho,c;
 30: } Field;

 32: typedef struct {
 33:   PetscScalar epsilon,delta,alpha,beta,gamma,kappa,lambda,mu,cstar;
 34:   PetscBool   upwind;
 35: } AppCtx;

 37: /*
 38:    User-defined routines
 39: */
 40: extern PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*),InitialConditions(DM,Vec);

 42: int main(int argc,char **argv)
 43: {
 44:   TS             ts;                  /* nonlinear solver */
 45:   Vec            U;                   /* solution, residual vectors */
 47:   DM             da;
 48:   AppCtx         appctx;

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:      Initialize program
 52:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 53:   PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;

 55:   appctx.epsilon = 1.0e-3;
 56:   appctx.delta   = 1.0;
 57:   appctx.alpha   = 10.0;
 58:   appctx.beta    = 4.0;
 59:   appctx.gamma   = 1.0;
 60:   appctx.kappa   = .75;
 61:   appctx.lambda  = 1.0;
 62:   appctx.mu      = 100.;
 63:   appctx.cstar   = .2;
 64:   appctx.upwind  = PETSC_TRUE;

 66:   PetscOptionsGetScalar(NULL,NULL,"-delta",&appctx.delta,NULL);
 67:   PetscOptionsGetBool(NULL,NULL,"-upwind",&appctx.upwind,NULL);

 69:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 70:      Create distributed array (DMDA) to manage parallel grid and vectors
 71:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 72:   DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE,-8,2,1,NULL,&da);
 73:   DMSetFromOptions(da);
 74:   DMSetUp(da);
 75:   DMDASetFieldName(da,0,"rho");
 76:   DMDASetFieldName(da,1,"c");

 78:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 79:      Extract global vectors from DMDA; then duplicate for remaining
 80:      vectors that are the same types
 81:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 82:   DMCreateGlobalVector(da,&U);

 84:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 85:      Create timestepping solver context
 86:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 87:   TSCreate(PETSC_COMM_WORLD,&ts);
 88:   TSSetType(ts,TSROSW);
 89:   TSSetDM(ts,da);
 90:   TSSetProblemType(ts,TS_NONLINEAR);
 91:   TSSetIFunction(ts,NULL,IFunction,&appctx);


 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:      Set initial conditions
 96:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   InitialConditions(da,U);
 98:   TSSetSolution(ts,U);

100:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101:      Set solver options
102:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103:   TSSetTimeStep(ts,.0001);
104:   TSSetMaxTime(ts,1.0);
105:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
106:   TSSetFromOptions(ts);

108:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
109:      Solve nonlinear system
110:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
111:   TSSolve(ts,U);

113:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
114:      Free work space.  All PETSc objects should be destroyed when they
115:      are no longer needed.
116:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117:   VecDestroy(&U);
118:   TSDestroy(&ts);
119:   DMDestroy(&da);

121:   PetscFinalize();
122:   return ierr;
123: }
124: /* ------------------------------------------------------------------- */
125: /*
126:    IFunction - Evaluates nonlinear function, F(U).

128:    Input Parameters:
129: .  ts - the TS context
130: .  U - input vector
131: .  ptr - optional user-defined context, as set by SNESSetFunction()

133:    Output Parameter:
134: .  F - function vector
135:  */
136: PetscErrorCode IFunction(TS ts,PetscReal ftime,Vec U,Vec Udot,Vec F,void *ptr)
137: {
138:   AppCtx         *appctx = (AppCtx*)ptr;
139:   DM             da;
141:   PetscInt       i,Mx,xs,xm;
142:   PetscReal      hx,sx;
143:   PetscScalar    rho,c,rhoxx,cxx,cx,rhox,kcxrhox;
144:   Field          *u,*f,*udot;
145:   Vec            localU;

148:   TSGetDM(ts,&da);
149:   DMGetLocalVector(da,&localU);
150:   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);

152:   hx = 1.0/(PetscReal)(Mx-1); sx = 1.0/(hx*hx);

154:   /*
155:      Scatter ghost points to local vector,using the 2-step process
156:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
157:      By placing code between these two statements, computations can be
158:      done while messages are in transition.
159:   */
160:   DMGlobalToLocalBegin(da,U,INSERT_VALUES,localU);
161:   DMGlobalToLocalEnd(da,U,INSERT_VALUES,localU);

163:   /*
164:      Get pointers to vector data
165:   */
166:   DMDAVecGetArrayRead(da,localU,&u);
167:   DMDAVecGetArrayRead(da,Udot,&udot);
168:   DMDAVecGetArray(da,F,&f);

170:   /*
171:      Get local grid boundaries
172:   */
173:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

175:   if (!xs) {
176:     f[0].rho = udot[0].rho; /* u[0].rho - 0.0; */
177:     f[0].c   = udot[0].c; /* u[0].c   - 1.0; */
178:     xs++;
179:     xm--;
180:   }
181:   if (xs+xm == Mx) {
182:     f[Mx-1].rho = udot[Mx-1].rho; /* u[Mx-1].rho - 1.0; */
183:     f[Mx-1].c   = udot[Mx-1].c;  /* u[Mx-1].c   - 0.0;  */
184:     xm--;
185:   }

187:   /*
188:      Compute function over the locally owned part of the grid
189:   */
190:   for (i=xs; i<xs+xm; i++) {
191:     rho   = u[i].rho;
192:     rhoxx = (-2.0*rho + u[i-1].rho + u[i+1].rho)*sx;
193:     c     = u[i].c;
194:     cxx   = (-2.0*c + u[i-1].c + u[i+1].c)*sx;

196:     if (!appctx->upwind) {
197:       rhox    = .5*(u[i+1].rho - u[i-1].rho)/hx;
198:       cx      = .5*(u[i+1].c - u[i-1].c)/hx;
199:       kcxrhox = appctx->kappa*(cxx*rho + cx*rhox);
200:     } else {
201:       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;
202:     }

204:     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;
205:     f[i].c   = udot[i].c - appctx->delta*cxx + appctx->lambda*c + appctx->alpha*rho*c/(appctx->gamma + c);
206:   }

208:   /*
209:      Restore vectors
210:   */
211:   DMDAVecRestoreArrayRead(da,localU,&u);
212:   DMDAVecRestoreArrayRead(da,Udot,&udot);
213:   DMDAVecRestoreArray(da,F,&f);
214:   DMRestoreLocalVector(da,&localU);
215:   return(0);
216: }

218: /* ------------------------------------------------------------------- */
219: PetscErrorCode InitialConditions(DM da,Vec U)
220: {
222:   PetscInt       i,xs,xm,Mx;
223:   Field          *u;
224:   PetscReal      hx,x;

227:   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);

229:   hx = 1.0/(PetscReal)(Mx-1);

231:   /*
232:      Get pointers to vector data
233:   */
234:   DMDAVecGetArray(da,U,&u);

236:   /*
237:      Get local grid boundaries
238:   */
239:   DMDAGetCorners(da,&xs,NULL,NULL,&xm,NULL,NULL);

241:   /*
242:      Compute function over the locally owned part of the grid
243:   */
244:   for (i=xs; i<xs+xm; i++) {
245:     x = i*hx;
246:     if (x < 1.0) u[i].rho = 0.0;
247:     else         u[i].rho = 1.0;
248:     u[i].c = PetscCosReal(.5*PETSC_PI*x);
249:   }

251:   /*
252:      Restore vectors
253:   */
254:   DMDAVecRestoreArray(da,U,&u);
255:   return(0);
256: }