Actual source code: ex6.c

petsc-3.6.4 2016-04-12
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:   TSSetInitialTimeStep(ts,user.t0,.005);
 90:   TSSetFromOptions(ts);
 91:   TSSetPostStep(ts,PostStep);
 92:   TSSolve(ts,x);

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

315:   return(0);
316: }

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

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

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

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



418: PetscErrorCode Parameter_settings(AppCtx *user)
419: {
421:   PetscBool      flg;


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