Actual source code: ex22.c

petsc-3.7.5 2017-01-01
Report Typos and Errors
  1: static const char help[] = "Time-dependent advection-reaction PDE in 1d, demonstrates IMEX methods.\n";
  2: /*
  3:    u_t + a1*u_x = -k1*u + k2*v + s1
  4:    v_t + a2*v_x = k1*u - k2*v + s2
  5:    0 < x < 1;
  6:    a1 = 1, k1 = 10^6, s1 = 0,
  7:    a2 = 0, k2 = 2*k1, s2 = 1

  9:    Initial conditions:
 10:    u(x,0) = 1 + s2*x
 11:    v(x,0) = k0/k1*u(x,0) + s1/k1

 13:    Upstream boundary conditions:
 14:    u(0,t) = 1-sin(12*t)^4

 16:    Useful command-line parameters:
 17:    -ts_arkimex_fully_implicit - treats advection implicitly, only relevant with -ts_type arkimex (default)
 18:    -ts_type rosw - use Rosenbrock-W method (linearly implicit IMEX, amortizes assembly and preconditioner setup)
 19: */

 21: #include <petscdm.h>
 22: #include <petscdmda.h>
 23: #include <petscts.h>

 25: typedef PetscScalar Field[2];

 27: typedef struct _User *User;
 28: struct _User {
 29:   PetscReal a[2];              /* Advection speeds */
 30:   PetscReal k[2];              /* Reaction coefficients */
 31:   PetscReal s[2];              /* Source terms */
 32: };

 34: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
 35: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 36: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
 37: static PetscErrorCode FormInitialSolution(TS,Vec,void*);

 41: int main(int argc,char **argv)
 42: {
 43:   TS                ts;         /* time integrator */
 44:   SNES              snes;       /* nonlinear solver */
 45:   SNESLineSearch    linesearch; /* line search */
 46:   Vec               X;          /* solution, residual vectors */
 47:   Mat               J;          /* Jacobian matrix */
 48:   PetscInt          steps,maxsteps,mx;
 49:   PetscErrorCode    ierr;
 50:   DM                da;
 51:   PetscReal         ftime,dt;
 52:   struct _User      user;       /* user-defined work context */
 53:   TSConvergedReason reason;

 55:   PetscInitialize(&argc,&argv,(char*)0,help);

 57:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 58:      Create distributed array (DMDA) to manage parallel grid and vectors
 59:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 60:   DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,-11,2,2,NULL,&da);

 62:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 63:      Extract global vectors from DMDA;
 64:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 65:   DMCreateGlobalVector(da,&X);

 67:   /* Initialize user application context */
 68:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Advection-reaction options","");
 69:   {
 70:     user.a[0] = 1;           PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],NULL);
 71:     user.a[1] = 0;           PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],NULL);
 72:     user.k[0] = 1e6;         PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],NULL);
 73:     user.k[1] = 2*user.k[0]; PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],NULL);
 74:     user.s[0] = 0;           PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],NULL);
 75:     user.s[1] = 1;           PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],NULL);
 76:   }
 77:   PetscOptionsEnd();

 79:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 80:      Create timestepping solver context
 81:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 82:   TSCreate(PETSC_COMM_WORLD,&ts);
 83:   TSSetDM(ts,da);
 84:   TSSetType(ts,TSARKIMEX);
 85:   TSSetRHSFunction(ts,NULL,FormRHSFunction,&user);
 86:   TSSetIFunction(ts,NULL,FormIFunction,&user);
 87:   DMSetMatType(da,MATAIJ);
 88:   DMCreateMatrix(da,&J);
 89:   TSSetIJacobian(ts,J,J,FormIJacobian,&user);

 91:   /* A line search in the nonlinear solve can fail due to ill-conditioning unless an absolute tolerance is set. Since
 92:    * this problem is linear, we deactivate the line search. For a linear problem, it is usually recommended to also use
 93:    * SNESSetType(snes,SNESKSPONLY). */
 94:   TSGetSNES(ts,&snes);
 95:   SNESGetLineSearch(snes,&linesearch);
 96:   SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);

 98:   ftime    = .1;
 99:   maxsteps = 10000;
100:   TSSetDuration(ts,maxsteps,ftime);
101:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
102: 
103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:      Set initial conditions
105:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
106:   FormInitialSolution(ts,X,&user);
107:   TSSetSolution(ts,X);
108:   VecGetSize(X,&mx);
109:   dt   = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
110:   TSSetInitialTimeStep(ts,0.0,dt);

112:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113:      Set runtime options
114:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115:   TSSetFromOptions(ts);

117:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118:      Solve nonlinear system
119:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120:   TSSolve(ts,X);
121:   TSGetSolveTime(ts,&ftime);
122:   TSGetTimeStepNumber(ts,&steps);
123:   TSGetConvergedReason(ts,&reason);
124:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %g after %D steps\n",TSConvergedReasons[reason],(double)ftime,steps);

126:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127:      Free work space.
128:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129:   MatDestroy(&J);
130:   VecDestroy(&X);
131:   TSDestroy(&ts);
132:   DMDestroy(&da);
133:   PetscFinalize();
134:   return 0;
135: }

139: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
140: {
141:   User           user = (User)ptr;
142:   DM             da;
143:   DMDALocalInfo  info;
144:   PetscInt       i;
145:   Field          *f;
146:   const Field    *x,*xdot;

150:   TSGetDM(ts,&da);
151:   DMDAGetLocalInfo(da,&info);

153:   /* Get pointers to vector data */
154:   DMDAVecGetArrayRead(da,X,(void*)&x);
155:   DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);
156:   DMDAVecGetArray(da,F,&f);

158:   /* Compute function over the locally owned part of the grid */
159:   for (i=info.xs; i<info.xs+info.xm; i++) {
160:     f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0];
161:     f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1];
162:   }

164:   /* Restore vectors */
165:   DMDAVecRestoreArrayRead(da,X,(void*)&x);
166:   DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);
167:   DMDAVecRestoreArray(da,F,&f);
168:   return(0);
169: }

173: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
174: {
175:   User           user = (User)ptr;
176:   DM             da;
177:   Vec            Xloc;
178:   DMDALocalInfo  info;
179:   PetscInt       i,j;
180:   PetscReal      hx;
181:   Field          *f;
182:   const Field    *x;

186:   TSGetDM(ts,&da);
187:   DMDAGetLocalInfo(da,&info);
188:   hx   = 1.0/(PetscReal)info.mx;

190:   /*
191:      Scatter ghost points to local vector,using the 2-step process
192:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
193:      By placing code between these two statements, computations can be
194:      done while messages are in transition.
195:   */
196:   DMGetLocalVector(da,&Xloc);
197:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
198:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);

200:   /* Get pointers to vector data */
201:   DMDAVecGetArrayRead(da,Xloc,(void*)&x);
202:   DMDAVecGetArray(da,F,&f);

204:   /* Compute function over the locally owned part of the grid */
205:   for (i=info.xs; i<info.xs+info.xm; i++) {
206:     const PetscReal *a  = user->a;
207:     PetscReal       u0t[2];
208:     u0t[0] = 1.0 - PetscPowRealInt(PetscSinReal(12*t),4);
209:     u0t[1] = 0.0;
210:     for (j=0; j<2; j++) {
211:       if (i == 0)              f[i][j] = a[j]/hx*(1./3*u0t[j] + 0.5*x[i][j] - x[i+1][j] + 1./6*x[i+2][j]);
212:       else if (i == 1)         f[i][j] = a[j]/hx*(-1./12*u0t[j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
213:       else if (i == info.mx-2) f[i][j] = a[j]/hx*(-1./6*x[i-2][j] + x[i-1][j] - 0.5*x[i][j] - 1./3*x[i+1][j]);
214:       else if (i == info.mx-1) f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]);
215:       else                     f[i][j] = a[j]/hx*(-1./12*x[i-2][j] + 2./3*x[i-1][j] - 2./3*x[i+1][j] + 1./12*x[i+2][j]);
216:     }
217:   }

219:   /* Restore vectors */
220:   DMDAVecRestoreArrayRead(da,Xloc,(void*)&x);
221:   DMDAVecRestoreArray(da,F,&f);
222:   DMRestoreLocalVector(da,&Xloc);
223:   return(0);
224: }

226: /* --------------------------------------------------------------------- */
227: /*
228:   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
229: */
232: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ptr)
233: {
234:   User           user = (User)ptr;
236:   DMDALocalInfo  info;
237:   PetscInt       i;
238:   DM             da;
239:   const Field    *x,*xdot;

242:   TSGetDM(ts,&da);
243:   DMDAGetLocalInfo(da,&info);

245:   /* Get pointers to vector data */
246:   DMDAVecGetArrayRead(da,X,(void*)&x);
247:   DMDAVecGetArrayRead(da,Xdot,(void*)&xdot);

249:   /* Compute function over the locally owned part of the grid */
250:   for (i=info.xs; i<info.xs+info.xm; i++) {
251:     const PetscReal *k = user->k;
252:     PetscScalar     v[2][2];

254:     v[0][0] = a + k[0]; v[0][1] =  -k[1];
255:     v[1][0] =    -k[0]; v[1][1] = a+k[1];
256:     MatSetValuesBlocked(Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES);
257:   }

259:   /* Restore vectors */
260:   DMDAVecRestoreArrayRead(da,X,(void*)&x);
261:   DMDAVecRestoreArrayRead(da,Xdot,(void*)&xdot);

263:   MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
264:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
265:   if (J != Jpre) {
266:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
267:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
268:   }
269:   return(0);
270: }

274: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
275: {
276:   User           user = (User)ctx;
277:   DM             da;
278:   PetscInt       i;
279:   DMDALocalInfo  info;
280:   Field          *x;
281:   PetscReal      hx;

285:   TSGetDM(ts,&da);
286:   DMDAGetLocalInfo(da,&info);
287:   hx   = 1.0/(PetscReal)info.mx;

289:   /* Get pointers to vector data */
290:   DMDAVecGetArray(da,X,&x);

292:   /* Compute function over the locally owned part of the grid */
293:   for (i=info.xs; i<info.xs+info.xm; i++) {
294:     PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0;
295:     x[i][0] = 1 + user->s[1]*r;
296:     x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik;
297:   }
298:   DMDAVecRestoreArray(da,X,&x);
299:   return(0);
300: }