Actual source code: ex22.c

petsc-3.3-p7 2013-05-11
  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
 15: */

 17: #include <petscdmda.h>
 18: #include <petscts.h>

 20: typedef PetscScalar Field[2];

 22: typedef struct _User *User;
 23: struct _User {
 24:   PetscReal      a[2];         /* Advection speeds */
 25:   PetscReal      k[2];         /* Reaction coefficients */
 26:   PetscReal      s[2];         /* Source terms */
 27: };

 29: static PetscErrorCode FormRHSFunction(TS,PetscReal,Vec,Vec,void*);
 30: static PetscErrorCode FormIFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 31: static PetscErrorCode FormIJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat*,Mat*,MatStructure*,void*);
 32: static PetscErrorCode FormInitialSolution(TS,Vec,void*);

 36: int main(int argc,char **argv)
 37: {
 38:   TS                ts;         /* nonlinear solver */
 39:   Vec               X;          /* solution, residual vectors */
 40:   Mat               J;          /* Jacobian matrix */
 41:   PetscInt          steps,maxsteps,mx;
 42:   PetscErrorCode    ierr;
 43:   DM                da;
 44:   PetscReal         ftime,dt;
 45:   struct _User      user;       /* user-defined work context */
 46:   TSConvergedReason reason;

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

 50:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 51:      Create distributed array (DMDA) to manage parallel grid and vectors
 52:   - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 53:   DMDACreate1d(PETSC_COMM_WORLD,DMDA_BOUNDARY_NONE,-11,2,2,PETSC_NULL,&da);

 55:   /*  - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 56:      Extract global vectors from DMDA;
 57:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 58:   DMCreateGlobalVector(da,&X);

 60:   /* Initialize user application context */
 61:   PetscOptionsBegin(PETSC_COMM_WORLD,PETSC_NULL,"Advection-reaction options",""); {
 62:     user.a[0] = 1;           PetscOptionsReal("-a0","Advection rate 0","",user.a[0],&user.a[0],PETSC_NULL);
 63:     user.a[1] = 0;           PetscOptionsReal("-a1","Advection rate 1","",user.a[1],&user.a[1],PETSC_NULL);
 64:     user.k[0] = 1e6;         PetscOptionsReal("-k0","Reaction rate 0","",user.k[0],&user.k[0],PETSC_NULL);
 65:     user.k[1] = 2*user.k[0]; PetscOptionsReal("-k1","Reaction rate 1","",user.k[1],&user.k[1],PETSC_NULL);
 66:     user.s[0] = 0;           PetscOptionsReal("-s0","Source 0","",user.s[0],&user.s[0],PETSC_NULL);
 67:     user.s[1] = 1;           PetscOptionsReal("-s1","Source 1","",user.s[1],&user.s[1],PETSC_NULL);
 68:   } PetscOptionsEnd();

 70:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 71:      Create timestepping solver context
 72:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 73:   TSCreate(PETSC_COMM_WORLD,&ts);
 74:   TSSetDM(ts,da);
 75:   TSSetType(ts,TSARKIMEX);
 76:   TSSetRHSFunction(ts,PETSC_NULL,FormRHSFunction,&user);
 77:   TSSetIFunction(ts,PETSC_NULL,FormIFunction,&user);
 78:   DMCreateMatrix(da,MATAIJ,&J);
 79:   TSSetIJacobian(ts,J,J,FormIJacobian,&user);

 81:   ftime = 1.0;
 82:   maxsteps = 10000;
 83:   TSSetDuration(ts,maxsteps,ftime);

 85:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 86:      Set initial conditions
 87:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 88:   FormInitialSolution(ts,X,&user);
 89:   TSSetSolution(ts,X);
 90:   VecGetSize(X,&mx);
 91:   dt   = .1 * PetscMax(user.a[0],user.a[1]) / mx; /* Advective CFL, I don't know why it needs so much safety factor. */
 92:   TSSetInitialTimeStep(ts,0.0,dt);

 94:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 95:      Set runtime options
 96:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 97:   TSSetFromOptions(ts);

 99:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100:      Solve nonlinear system
101:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102:   TSSolve(ts,X,&ftime);
103:   TSGetTimeStepNumber(ts,&steps);
104:   TSGetConvergedReason(ts,&reason);
105:   PetscPrintf(PETSC_COMM_WORLD,"%s at time %G after %D steps\n",TSConvergedReasons[reason],ftime,steps);

107:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
108:      Free work space.
109:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
110:   MatDestroy(&J);
111:   VecDestroy(&X);
112:   TSDestroy(&ts);
113:   DMDestroy(&da);
114:   PetscFinalize();
115:   return 0;
116: }

120: static PetscErrorCode FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ptr)
121: {
122:   User           user = (User)ptr;
123:   DM             da;
124:   DMDALocalInfo  info;
125:   PetscInt       i;
126:   Field          *x,*xdot,*f;

130:   TSGetDM(ts,&da);
131:   DMDAGetLocalInfo(da,&info);

133:   /* Get pointers to vector data */
134:   DMDAVecGetArray(da,X,&x);
135:   DMDAVecGetArray(da,Xdot,&xdot);
136:   DMDAVecGetArray(da,F,&f);

138:   /* Compute function over the locally owned part of the grid */
139:   for (i=info.xs; i<info.xs+info.xm; i++) {
140:     f[i][0] = xdot[i][0] + user->k[0]*x[i][0] - user->k[1]*x[i][1] - user->s[0];
141:     f[i][1] = xdot[i][1] - user->k[0]*x[i][0] + user->k[1]*x[i][1] - user->s[1];
142:   }

144:   /* Restore vectors */
145:   DMDAVecRestoreArray(da,X,&x);
146:   DMDAVecRestoreArray(da,Xdot,&xdot);
147:   DMDAVecRestoreArray(da,F,&f);
148:   return(0);
149: }

153: static PetscErrorCode FormRHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ptr)
154: {
155:   User           user = (User)ptr;
156:   DM             da;
157:   Vec            Xloc;
158:   DMDALocalInfo  info;
159:   PetscInt       i,j;
160:   PetscReal      hx;
161:   Field          *x,*f;

165:   TSGetDM(ts,&da);
166:   DMDAGetLocalInfo(da,&info);
167:   hx = 1.0/(PetscReal)info.mx;

169:   /*
170:      Scatter ghost points to local vector,using the 2-step process
171:         DMGlobalToLocalBegin(),DMGlobalToLocalEnd().
172:      By placing code between these two statements, computations can be
173:      done while messages are in transition.
174:   */
175:   DMGetLocalVector(da,&Xloc);
176:   DMGlobalToLocalBegin(da,X,INSERT_VALUES,Xloc);
177:   DMGlobalToLocalEnd(da,X,INSERT_VALUES,Xloc);

179:   /* Get pointers to vector data */
180:   DMDAVecGetArray(da,Xloc,&x);
181:   DMDAVecGetArray(da,F,&f);

183:   /* Compute function over the locally owned part of the grid */
184:   for (i=info.xs; i<info.xs+info.xm; i++) {
185:     const PetscReal *a = user->a;
186:     PetscReal u0t[2] = {1. - PetscPowScalar(sin(12*t),4.),0};
187:     for (j=0; j<2; j++) {
188:       if (i == 0) {
189:         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]);
190:       } else if (i == 1) {
191:         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]);
192:       } else if (i == info.mx-2) {
193:         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]);
194:       } else if (i == info.mx-1) {
195:         f[i][j] = a[j]/hx*(-x[i][j] + x[i-1][j]);
196:       } else {
197:         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]);
198:       }
199:     }
200:   }

202:   /* Restore vectors */
203:   DMDAVecRestoreArray(da,Xloc,&x);
204:   DMDAVecRestoreArray(da,F,&f);
205:   DMRestoreLocalVector(da,&Xloc);
206:   return(0);
207: }

209: /* --------------------------------------------------------------------- */
210: /*
211:   IJacobian - Compute IJacobian = dF/dU + a dF/dUdot
212: */
215: PetscErrorCode FormIJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat *J,Mat *Jpre,MatStructure *str,void *ptr)
216: {
217:   User           user = (User)ptr;
219:   DMDALocalInfo  info;
220:   PetscInt       i;
221:   DM             da;
222:   Field          *x,*xdot;

225:   TSGetDM(ts,&da);
226:   DMDAGetLocalInfo(da,&info);

228:   /* Get pointers to vector data */
229:   DMDAVecGetArray(da,X,&x);
230:   DMDAVecGetArray(da,Xdot,&xdot);

232:   /* Compute function over the locally owned part of the grid */
233:   for (i=info.xs; i<info.xs+info.xm; i++) {
234:     const PetscReal *k = user->k;
235:     PetscScalar v[2][2] = {{a + k[0], -k[1]},
236:                            {-k[0], a+k[1]}};
237:     MatSetValuesBlocked(*Jpre,1,&i,1,&i,&v[0][0],INSERT_VALUES);
238:   }

240:   /* Restore vectors */
241:   DMDAVecRestoreArray(da,X,&x);
242:   DMDAVecRestoreArray(da,Xdot,&xdot);

244:   MatAssemblyBegin(*Jpre,MAT_FINAL_ASSEMBLY);
245:   MatAssemblyEnd(*Jpre,MAT_FINAL_ASSEMBLY);
246:   if (*J != *Jpre) {
247:     MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);
248:     MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);
249:   }
250:   return(0);
251: }

255: PetscErrorCode FormInitialSolution(TS ts,Vec X,void *ctx)
256: {
257:   User user = (User)ctx;
258:   DM             da;
259:   PetscInt       i;
260:   DMDALocalInfo  info;
261:   Field          *x;
262:   PetscReal      hx;

266:   TSGetDM(ts,&da);
267:   DMDAGetLocalInfo(da,&info);
268:   hx = 1.0/(PetscReal)info.mx;

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

273:   /* Compute function over the locally owned part of the grid */
274:   for (i=info.xs; i<info.xs+info.xm; i++) {
275:     PetscReal r = (i+1)*hx,ik = user->k[1] != 0.0 ? 1.0/user->k[1] : 1.0;
276:     x[i][0] = 1 + user->s[1]*r;
277:     x[i][1] = user->k[0]*ik*x[i][0] + user->s[1]*ik;
278:   }
279:   DMDAVecRestoreArray(da,X,&x);
280:   return(0);
281: }