Actual source code: ex22.c

petsc-master 2016-09-25
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);if (ierr) return ierr;

 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);
 61:   DMSetFromOptions(da);
 62:   DMSetUp(da);

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

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

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

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

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

114:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115:      Set runtime options
116:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
117:   TSSetFromOptions(ts);

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

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

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

152:   TSGetDM(ts,&da);
153:   DMDAGetLocalInfo(da,&info);

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

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

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

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

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

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

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

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

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

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

244:   TSGetDM(ts,&da);
245:   DMDAGetLocalInfo(da,&info);

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

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

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

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

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

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

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

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

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