Actual source code: ex6.c

petsc-3.7.4 2016-10-02
Report Typos and Errors
  1: static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n";
  2: /*
  3:    p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy
  4:    xmin < x < xmax, ymin < y < ymax;
  5:    x_t = (y - ws)  y_t = (ws/2H)*(Pm - Pmax*sin(x))

  7:    Boundary conditions: -bc_type 0 => Zero dirichlet boundary
  8:                         -bc_type 1 => Steady state boundary condition
  9:    Steady state boundary condition found by setting p_t = 0
 10: */

 12: #include <petscdm.h>
 13: #include <petscdmda.h>
 14: #include <petscts.h>

 16: /*
 17:    User-defined data structures and routines
 18: */
 19: typedef struct {
 20:   PetscScalar ws;   /* Synchronous speed */
 21:   PetscScalar H;    /* Inertia constant */
 22:   PetscScalar D;    /* Damping constant */
 23:   PetscScalar Pmax; /* Maximum power output of generator */
 24:   PetscScalar PM_min; /* Mean mechanical power input */
 25:   PetscScalar lambda; /* correlation time */
 26:   PetscScalar q;      /* noise strength */
 27:   PetscScalar mux;    /* Initial average angle */
 28:   PetscScalar sigmax; /* Standard deviation of initial angle */
 29:   PetscScalar muy;    /* Average speed */
 30:   PetscScalar sigmay; /* standard deviation of initial speed */
 31:   PetscScalar rho;    /* Cross-correlation coefficient */
 32:   PetscScalar t0;     /* Initial time */
 33:   PetscScalar tmax;   /* Final time */
 34:   PetscScalar xmin;   /* left boundary of angle */
 35:   PetscScalar xmax;   /* right boundary of angle */
 36:   PetscScalar ymin;   /* bottom boundary of speed */
 37:   PetscScalar ymax;   /* top boundary of speed */
 38:   PetscScalar dx;     /* x step size */
 39:   PetscScalar dy;     /* y step size */
 40:   PetscInt    bc; /* Boundary conditions */
 41:   PetscScalar disper_coe; /* Dispersion coefficient */
 42:   DM          da;
 43: } AppCtx;

 45: PetscErrorCode Parameter_settings(AppCtx*);
 46: PetscErrorCode ini_bou(Vec,AppCtx*);
 47: PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
 48: PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
 49: PetscErrorCode PostStep(TS);

 53: int main(int argc, char **argv)
 54: {
 56:   Vec            x;  /* Solution vector */
 57:   TS             ts;   /* Time-stepping context */
 58:   AppCtx         user; /* Application context */
 59:   Mat            J;
 60:   PetscViewer    viewer;

 62:   PetscInitialize(&argc,&argv,"petscopt_ex6", help);

 64:   /* Get physics and time parameters */
 65:   Parameter_settings(&user);
 66:   /* Create a 2D DA with dof = 1 */
 67:   DMDACreate2d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,DM_BOUNDARY_NONE,DMDA_STENCIL_STAR,-4,-4,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,NULL,&user.da);
 68:   /* Set x and y coordinates */
 69:   DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0.0,1.0);

 71:   /* Get global vector x from DM  */
 72:   DMCreateGlobalVector(user.da,&x);

 74:   ini_bou(x,&user);
 75:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer);
 76:   VecView(x,viewer);
 77:   PetscViewerDestroy(&viewer);

 79:   /* Get Jacobian matrix structure from the da */
 80:   DMSetMatType(user.da,MATAIJ);
 81:   DMCreateMatrix(user.da,&J);

 83:   TSCreate(PETSC_COMM_WORLD,&ts);
 84:   TSSetProblemType(ts,TS_NONLINEAR);
 85:   TSSetIFunction(ts,NULL,IFunction,&user);
 86:   TSSetIJacobian(ts,J,J,IJacobian,&user);
 87:   TSSetApplicationContext(ts,&user);
 88:   TSSetDuration(ts,PETSC_DEFAULT,user.tmax);
 89:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
 90:   TSSetInitialTimeStep(ts,user.t0,.005);
 91:   TSSetFromOptions(ts);
 92:   TSSetPostStep(ts,PostStep);
 93:   TSSolve(ts,x);

 95:   PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer);
 96:   VecView(x,viewer);
 97:   PetscViewerDestroy(&viewer);

 99:   VecDestroy(&x);
100:   MatDestroy(&J);
101:   DMDestroy(&user.da);
102:   TSDestroy(&ts);
103:   PetscFinalize();
104:   return 0;
105: }

109: PetscErrorCode PostStep(TS ts)
110: {
112:   Vec            X;
113:   AppCtx         *user;
114:   PetscScalar    sum;
115:   PetscReal      t;
117:   TSGetApplicationContext(ts,&user);
118:   TSGetTime(ts,&t);
119:   TSGetSolution(ts,&X);
120:   VecSum(X,&sum);
121:   PetscPrintf(PETSC_COMM_WORLD,"sum(p)*dw*dtheta at t = %3.2f = %3.6f\n",t,sum*user->dx*user->dy);
122:   return(0);
123: }

127: PetscErrorCode ini_bou(Vec X,AppCtx* user)
128: {
130:   DM             cda;
131:   DMDACoor2d     **coors;
132:   PetscScalar    **p;
133:   Vec            gc;
134:   PetscInt       i,j;
135:   PetscInt       xs,ys,xm,ym,M,N;
136:   PetscScalar    xi,yi;
137:   PetscScalar    sigmax=user->sigmax,sigmay=user->sigmay;
138:   PetscScalar    rho   =user->rho;
139:   PetscScalar    mux   =user->mux,muy=user->muy;
140:   PetscMPIInt    rank;

143:   MPI_Comm_rank(PETSC_COMM_WORLD,&rank);
144:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
145:   user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
146:   DMGetCoordinateDM(user->da,&cda);
147:   DMGetCoordinates(user->da,&gc);
148:   DMDAVecGetArray(cda,gc,&coors);
149:   DMDAVecGetArray(user->da,X,&p);
150:   DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);
151:   for (i=xs; i < xs+xm; i++) {
152:     for (j=ys; j < ys+ym; j++) {
153:       xi = coors[j][i].x; yi = coors[j][i].y;
154:       if (i == 0 || j == 0 || i == M-1 || j == N-1) p[j][i] = 0.0;
155:       else p[j][i] = (0.5/(PETSC_PI*sigmax*sigmay*PetscSqrtScalar(1.0-rho*rho)))*PetscExpScalar(-0.5/(1-rho*rho)*(PetscPowScalar((xi-mux)/sigmax,2) + PetscPowScalar((yi-muy)/sigmay,2) - 2*rho*(xi-mux)*(yi-muy)/(sigmax*sigmay)));
156:     }
157:   }
158:   /*  p[N/2+N%2][M/2+M%2] = 1/(user->dx*user->dy); */

160:   DMDAVecRestoreArray(cda,gc,&coors);
161:   DMDAVecRestoreArray(user->da,X,&p);
162:   return(0);
163: }

165: /* First advection term */
168: PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
169: {
170:   PetscScalar f;
171:   /*  PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */
173:   /*  if (i > 2 && i < M-2) {
174:     v1 = (y-user->ws)*(p[j][i-2] - p[j][i-3])/user->dx;
175:     v2 = (y-user->ws)*(p[j][i-1] - p[j][i-2])/user->dx;
176:     v3 = (y-user->ws)*(p[j][i] - p[j][i-1])/user->dx;
177:     v4 = (y-user->ws)*(p[j][i+1] - p[j][i])/user->dx;
178:     v5 = (y-user->ws)*(p[j][i+1] - p[j][i+2])/user->dx;

180:     s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3;
181:     s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4;
182:     s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5;

184:     *p1 = 0.1*s1 + 0.6*s2 + 0.3*s3;
185:     } else *p1 = 0.0; */
186:   f   =  (y - user->ws);
187:   *p1 = f*(p[j][i+1] - p[j][i-1])/(2*user->dx);
188:   return(0);
189: }

191: /* Second advection term */
194: PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
195: {
196:   PetscScalar f;
197:   /*  PetscScalar v1,v2,v3,v4,v5,s1,s2,s3; */
199:   /*  if (j > 2 && j < N-2) {
200:     v1 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-2][i] - p[j-3][i])/user->dy;
201:     v2 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j-1][i] - p[j-2][i])/user->dy;
202:     v3 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j][i] - p[j-1][i])/user->dy;
203:     v4 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+1][i] - p[j][i])/user->dy;
204:     v5 = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*sin(x))*(p[j+2][i] - p[j+1][i])/user->dy;

206:     s1 = v1/3.0 - (7.0/6.0)*v2 + (11.0/6.0)*v3;
207:     s2 =-v2/6.0 + (5.0/6.0)*v3 + (1.0/3.0)*v4;
208:     s3 = v3/3.0 + (5.0/6.0)*v4 - (1.0/6.0)*v5;

210:     *p2 = 0.1*s1 + 0.6*s2 + 0.3*s3;
211:     } else *p2 = 0.0; */
212:   f   = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x));
213:   *p2 = f*(p[j+1][i] - p[j-1][i])/(2*user->dy);
214:   return(0);
215: }

217: /* Diffusion term */
220: PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
221: {

224:   *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
225:   return(0);
226: }

230: PetscErrorCode BoundaryConditions(PetscScalar **p,DMDACoor2d **coors,PetscInt i,PetscInt j,PetscInt M, PetscInt N,PetscScalar **f,AppCtx *user)
231: {
232:   PetscScalar fwc,fthetac;
233:   PetscScalar w=coors[j][i].y,theta=coors[j][i].x;

236:   if (user->bc == 0) { /* Natural boundary condition */
237:     f[j][i] = p[j][i];
238:   } else { /* Steady state boundary condition */
239:     fthetac = user->ws/(2*user->H)*(user->PM_min - user->Pmax*PetscSinScalar(theta));
240:     fwc = (w*w/2.0 - user->ws*w);
241:     if (i == 0 && j == 0) { /* left bottom corner */
242:       f[j][i] = fwc*(p[j][i+1] - p[j][i])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j][i])/user->dy;
243:     } else if (i == 0 && j == N-1) { /* right bottom corner */
244:       f[j][i] = fwc*(p[j][i+1] - p[j][i])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j][i] - p[j-1][i])/user->dy;
245:     } else if (i == M-1 && j == 0) { /* left top corner */
246:       f[j][i] = fwc*(p[j][i] - p[j][i-1])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j][i])/user->dy;
247:     } else if (i == M-1 && j == N-1) { /* right top corner */
248:       f[j][i] = fwc*(p[j][i] - p[j][i-1])/user->dx + fthetac*p[j][i] - user->disper_coe*(p[j][i] - p[j-1][i])/user->dy;
249:     } else if (i == 0) { /* Bottom edge */
250:       f[j][i] = fwc*(p[j][i+1] - p[j][i])/(user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j-1][i])/(2*user->dy);
251:     } else if (i == M-1) { /* Top edge */
252:       f[j][i] = fwc*(p[j][i] - p[j][i-1])/(user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j-1][i])/(2*user->dy);
253:     } else if (j == 0) { /* Left edge */
254:       f[j][i] = fwc*(p[j][i+1] - p[j][i-1])/(2*user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j+1][i] - p[j][i])/(user->dy);
255:     } else if (j == N-1) { /* Right edge */
256:       f[j][i] = fwc*(p[j][i+1] - p[j][i-1])/(2*user->dx) + fthetac*p[j][i] - user->disper_coe*(p[j][i] - p[j-1][i])/(user->dy);
257:     }
258:   }
259:   return(0);
260: }

264: PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
265: {
267:   AppCtx         *user=(AppCtx*)ctx;
268:   DM             cda;
269:   DMDACoor2d     **coors;
270:   PetscScalar    **p,**f,**pdot;
271:   PetscInt       i,j;
272:   PetscInt       xs,ys,xm,ym,M,N;
273:   Vec            localX,gc,localXdot;
274:   PetscScalar    p_adv1,p_adv2,p_diff;

277:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
278:   DMGetCoordinateDM(user->da,&cda);
279:   DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);

281:   DMGetLocalVector(user->da,&localX);
282:   DMGetLocalVector(user->da,&localXdot);

284:   DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX);
285:   DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX);
286:   DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot);
287:   DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot);

289:   DMGetCoordinatesLocal(user->da,&gc);

291:   DMDAVecGetArrayRead(cda,gc,&coors);
292:   DMDAVecGetArrayRead(user->da,localX,&p);
293:   DMDAVecGetArrayRead(user->da,localXdot,&pdot);
294:   DMDAVecGetArray(user->da,F,&f);

296:   user->disper_coe = PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda));
297:   for (i=xs; i < xs+xm; i++) {
298:     for (j=ys; j < ys+ym; j++) {
299:       if (i == 0 || j == 0 || i == M-1 || j == N-1) {
300:         BoundaryConditions(p,coors,i,j,M,N,f,user);
301:       } else {
302:         adv1(p,coors[j][i].y,i,j,M,&p_adv1,user);
303:         adv2(p,coors[j][i].x,i,j,N,&p_adv2,user);
304:         diffuse(p,i,j,t,&p_diff,user);
305:         f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i];
306:       }
307:     }
308:   }
309:   DMDAVecRestoreArrayRead(user->da,localX,&p);
310:   DMDAVecRestoreArrayRead(user->da,localX,&pdot);
311:   DMRestoreLocalVector(user->da,&localX);
312:   DMRestoreLocalVector(user->da,&localXdot);
313:   DMDAVecRestoreArray(user->da,F,&f);
314:   DMDAVecRestoreArrayRead(cda,gc,&coors);

316:   return(0);
317: }

321: PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
322: {
324:   AppCtx         *user=(AppCtx*)ctx;
325:   DM             cda;
326:   DMDACoor2d     **coors;
327:   PetscInt       i,j;
328:   PetscInt       xs,ys,xm,ym,M,N;
329:   Vec            gc;
330:   PetscScalar    val[5],xi,yi;
331:   MatStencil     row,col[5];
332:   PetscScalar    c1,c3,c5;

335:   DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL);
336:   DMGetCoordinateDM(user->da,&cda);
337:   DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0);

339:   DMGetCoordinatesLocal(user->da,&gc);
340:   DMDAVecGetArrayRead(cda,gc,&coors);
341:   for (i=xs; i < xs+xm; i++) {
342:     for (j=ys; j < ys+ym; j++) {
343:       PetscInt nc = 0;
344:       xi = coors[j][i].x; yi = coors[j][i].y;
345:       row.i = i; row.j = j;
346:       if (i == 0 || j == 0 || i == M-1 || j == N-1) {
347:         if (user->bc == 0) {
348:           col[nc].i = i; col[nc].j = j; val[nc++] = 1.0;
349:         } else {
350:           PetscScalar fthetac,fwc;
351:           fthetac = user->ws/(2*user->H)*(user->PM_min - user->Pmax*PetscSinScalar(xi));
352:           fwc     = (yi*yi/2.0 - user->ws*yi);
353:           if (i==0 && j==0) {
354:             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/user->dx;
355:             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy;
356:             col[nc].i = i;   col[nc].j = j;   val[nc++] = -fwc/user->dx + fthetac + user->disper_coe/user->dy;
357:           } else if (i==0 && j == N-1) {
358:             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/user->dx;
359:             col[nc].i = i;   col[nc].j = j-1; val[nc++] = user->disper_coe/user->dy;
360:             col[nc].i = i;   col[nc].j = j;   val[nc++] = -fwc/user->dx + fthetac - user->disper_coe/user->dy;
361:           } else if (i== M-1 && j == 0) {
362:             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/user->dx;
363:             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy;
364:             col[nc].i = i;   col[nc].j = j;   val[nc++] =  fwc/user->dx + fthetac + user->disper_coe/user->dy;
365:           } else if (i == M-1 && j == N-1) {
366:             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/user->dx;
367:             col[nc].i = i;   col[nc].j = j-1; val[nc++] =  user->disper_coe/user->dy;
368:             col[nc].i = i;   col[nc].j = j;   val[nc++] =  fwc/user->dx + fthetac - user->disper_coe/user->dy;
369:           } else if (i==0) {
370:             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/user->dx;
371:             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/(2*user->dy);
372:             col[nc].i = i;   col[nc].j = j-1; val[nc++] =  user->disper_coe/(2*user->dy);
373:             col[nc].i = i;   col[nc].j = j;   val[nc++] = -fwc/user->dx + fthetac;
374:           } else if (i == M-1) {
375:             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/user->dx;
376:             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/(2*user->dy);
377:             col[nc].i = i;   col[nc].j = j-1; val[nc++] =  user->disper_coe/(2*user->dy);
378:             col[nc].i = i;   col[nc].j = j;   val[nc++] = fwc/user->dx + fthetac;
379:           } else if (j==0) {
380:             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/(2*user->dx);
381:             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/(2*user->dx);
382:             col[nc].i = i;   col[nc].j = j+1; val[nc++] = -user->disper_coe/user->dy;
383:             col[nc].i = i;   col[nc].j = j;   val[nc++] = user->disper_coe/user->dy + fthetac;
384:           } else if (j == N-1) {
385:             col[nc].i = i+1; col[nc].j = j;   val[nc++] = fwc/(2*user->dx);
386:             col[nc].i = i-1; col[nc].j = j;   val[nc++] = -fwc/(2*user->dx);
387:             col[nc].i = i;   col[nc].j = j-1; val[nc++] = user->disper_coe/user->dy;
388:             col[nc].i = i;   col[nc].j = j;   val[nc++] = -user->disper_coe/user->dy + fthetac;
389:           }
390:         }
391:       } else {
392:         c1        = (yi-user->ws)/(2*user->dx);
393:         c3        = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi))/(2*user->dy);
394:         c5        = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
395:         col[nc].i = i-1; col[nc].j = j;   val[nc++] = c1;
396:         col[nc].i = i+1; col[nc].j = j;   val[nc++] = -c1;
397:         col[nc].i = i;   col[nc].j = j-1; val[nc++] = c3 + c5;
398:         col[nc].i = i;   col[nc].j = j+1; val[nc++] = -c3 + c5;
399:         col[nc].i = i;   col[nc].j = j;   val[nc++] = -2*c5 -a;
400:       }
401:       MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES);
402:     }
403:   }
404:   DMDAVecRestoreArrayRead(cda,gc,&coors);

406:    MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY);
407:   MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY);
408:   if (J != Jpre) {
409:     MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);
410:     MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);
411:   }
412:   return(0);
413: }



419: PetscErrorCode Parameter_settings(AppCtx *user)
420: {
422:   PetscBool      flg;


426:   /* Set default parameters */
427:   user->ws     = 1.0;
428:   user->H      = 5.0;  user->Pmax   = 2.1;
429:   user->PM_min = 1.0;  user->lambda = 0.1;
430:   user->q      = 1.0;  user->mux    = PetscAsinScalar(user->PM_min/user->Pmax);
431:   user->sigmax = 0.1;
432:   user->sigmay = 0.1;  user->rho  = 0.0;
433:   user->t0     = 0.0;  user->tmax = 2.0;
434:   user->xmin   = -1.0; user->xmax = 10.0;
435:   user->ymin   = -1.0; user->ymax = 10.0;
436:   user->bc     = 0;
437: 
438:   PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg);
439:   PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg);
440:   PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg);
441:   PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg);
442:   PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg);
443:   PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg);
444:   PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg);
445:   PetscOptionsGetScalar(NULL,NULL,"-sigmax",&user->sigmax,&flg);
446:   PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg);
447:   PetscOptionsGetScalar(NULL,NULL,"-sigmay",&user->sigmay,&flg);
448:   PetscOptionsGetScalar(NULL,NULL,"-rho",&user->rho,&flg);
449:   PetscOptionsGetScalar(NULL,NULL,"-t0",&user->t0,&flg);
450:   PetscOptionsGetScalar(NULL,NULL,"-tmax",&user->tmax,&flg);
451:   PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg);
452:   PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg);
453:   PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg);
454:   PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg);
455:   PetscOptionsGetInt(NULL,NULL,"-bc_type",&user->bc,&flg);
456:   user->muy = user->ws;
457:   return(0);
458: }