xref: /petsc/src/ts/tutorials/power_grid/ex7.c (revision 0e02354e4efb6ae7920bb0e7d191735ccbbd1e1e)
1 
2 static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n";
3 /*
4    p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy
5    xmin < x < xmax, ymin < y < ymax;
6 
7    Boundary conditions Neumman using mirror values
8 
9    Note that x_t and y_t in the above are given functions of x and y; they are not derivatives of x and y.
10    x_t = (y - ws)  y_t = (ws/2H)*(Pm - Pmax*sin(x))
11 
12 */
13 
14 #include <petscdm.h>
15 #include <petscdmda.h>
16 #include <petscts.h>
17 
18 /*
19    User-defined data structures and routines
20 */
21 typedef struct {
22   PetscScalar ws;   /* Synchronous speed */
23   PetscScalar H;    /* Inertia constant */
24   PetscScalar D;    /* Damping constant */
25   PetscScalar Pmax; /* Maximum power output of generator */
26   PetscScalar PM_min; /* Mean mechanical power input */
27   PetscScalar lambda; /* correlation time */
28   PetscScalar q;      /* noise strength */
29   PetscScalar mux;    /* Initial average angle */
30   PetscScalar sigmax; /* Standard deviation of initial angle */
31   PetscScalar muy;    /* Average speed */
32   PetscScalar sigmay; /* standard deviation of initial speed */
33   PetscScalar rho;    /* Cross-correlation coefficient */
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   PetscScalar disper_coe; /* Dispersion coefficient */
41   DM          da;
42   PetscInt    st_width; /* Stencil width */
43   DMBoundaryType bx; /* x boundary type */
44   DMBoundaryType by; /* y boundary type */
45   PetscBool        nonoiseinitial;
46 } AppCtx;
47 
48 PetscErrorCode Parameter_settings(AppCtx*);
49 PetscErrorCode ini_bou(Vec,AppCtx*);
50 PetscErrorCode IFunction(TS,PetscReal,Vec,Vec,Vec,void*);
51 PetscErrorCode IJacobian(TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
52 PetscErrorCode PostStep(TS);
53 
54 int main(int argc, char **argv)
55 {
56   Vec            x;  /* Solution vector */
57   TS             ts;   /* Time-stepping context */
58   AppCtx         user; /* Application context */
59   PetscViewer    viewer;
60 
61   PetscCall(PetscInitialize(&argc,&argv,"petscopt_ex7", help));
62 
63   /* Get physics and time parameters */
64   PetscCall(Parameter_settings(&user));
65   /* Create a 2D DA with dof = 1 */
66   PetscCall(DMDACreate2d(PETSC_COMM_WORLD,user.bx,user.by,DMDA_STENCIL_STAR,4,4,PETSC_DECIDE,PETSC_DECIDE,1,user.st_width,NULL,NULL,&user.da));
67   PetscCall(DMSetFromOptions(user.da));
68   PetscCall(DMSetUp(user.da));
69   /* Set x and y coordinates */
70   PetscCall(DMDASetUniformCoordinates(user.da,user.xmin,user.xmax,user.ymin,user.ymax,0,0));
71   PetscCall(DMDASetCoordinateName(user.da,0,"X - the angle"));
72   PetscCall(DMDASetCoordinateName(user.da,1,"Y - the speed"));
73 
74   /* Get global vector x from DM  */
75   PetscCall(DMCreateGlobalVector(user.da,&x));
76 
77   PetscCall(ini_bou(x,&user));
78   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"ini_x",FILE_MODE_WRITE,&viewer));
79   PetscCall(VecView(x,viewer));
80   PetscCall(PetscViewerDestroy(&viewer));
81 
82   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
83   PetscCall(TSSetDM(ts,user.da));
84   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
85   PetscCall(TSSetType(ts,TSARKIMEX));
86   PetscCall(TSSetIFunction(ts,NULL,IFunction,&user));
87   /*  PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobian,&user));  */
88   PetscCall(TSSetApplicationContext(ts,&user));
89   PetscCall(TSSetTimeStep(ts,.005));
90   PetscCall(TSSetFromOptions(ts));
91   PetscCall(TSSetPostStep(ts,PostStep));
92   PetscCall(TSSolve(ts,x));
93 
94   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"fin_x",FILE_MODE_WRITE,&viewer));
95   PetscCall(VecView(x,viewer));
96   PetscCall(PetscViewerDestroy(&viewer));
97 
98   PetscCall(VecDestroy(&x));
99   PetscCall(DMDestroy(&user.da));
100   PetscCall(TSDestroy(&ts));
101   PetscCall(PetscFinalize());
102   return 0;
103 }
104 
105 PetscErrorCode PostStep(TS ts)
106 {
107   Vec            X,gc;
108   AppCtx         *user;
109   PetscScalar    sum = 0,asum;
110   PetscReal      t,**p;
111   DMDACoor2d     **coors;
112   DM             cda;
113   PetscInt       i,j,xs,ys,xm,ym;
114 
115   PetscFunctionBegin;
116   PetscCall(TSGetApplicationContext(ts,&user));
117   PetscCall(TSGetTime(ts,&t));
118   PetscCall(TSGetSolution(ts,&X));
119 
120   PetscCall(DMGetCoordinateDM(user->da,&cda));
121   PetscCall(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0));
122   PetscCall(DMGetCoordinates(user->da,&gc));
123   PetscCall(DMDAVecGetArrayRead(cda,gc,&coors));
124   PetscCall(DMDAVecGetArrayRead(user->da,X,&p));
125   for (i=xs; i < xs+xm; i++) {
126     for (j=ys; j < ys+ym; j++) {
127       if (coors[j][i].y < 5) sum += p[j][i];
128     }
129   }
130   PetscCall(DMDAVecRestoreArrayRead(cda,gc,&coors));
131   PetscCall(DMDAVecRestoreArrayRead(user->da,X,&p));
132   PetscCallMPI(MPI_Allreduce(&sum,&asum,1,MPIU_SCALAR,MPIU_SUM,PetscObjectComm((PetscObject)ts)));
133   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"sum(p)*dw*dtheta at t = %f = %f\n",(double)t,(double)(asum)));
134   if (sum  < 1.0e-2) {
135     PetscCall(TSSetConvergedReason(ts,TS_CONVERGED_USER));
136     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"Exiting TS as the integral of PDF is almost zero\n"));
137   }
138   PetscFunctionReturn(0);
139 }
140 
141 PetscErrorCode ini_bou(Vec X,AppCtx* user)
142 {
143   DM             cda;
144   DMDACoor2d     **coors;
145   PetscScalar    **p;
146   Vec            gc;
147   PetscInt       i,j;
148   PetscInt       xs,ys,xm,ym,M,N;
149   PetscScalar    xi,yi;
150   PetscScalar    sigmax=user->sigmax,sigmay=user->sigmay;
151   PetscScalar    rho   =user->rho;
152   PetscScalar    muy=user->muy,mux;
153   PetscMPIInt    rank;
154   PetscScalar    sum;
155 
156   PetscFunctionBeginUser;
157   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
158   PetscCall(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
159   user->dx = (user->xmax - user->xmin)/(M-1); user->dy = (user->ymax - user->ymin)/(N-1);
160   PetscCall(DMGetCoordinateDM(user->da,&cda));
161   PetscCall(DMGetCoordinates(user->da,&gc));
162   PetscCall(DMDAVecGetArray(cda,gc,&coors));
163   PetscCall(DMDAVecGetArray(user->da,X,&p));
164   PetscCall(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0));
165 
166   /* mux and muy need to be grid points in the x and y-direction otherwise the solution goes unstable
167      muy is set by choosing the y domain, no. of grid points along y-direction so that muy is a grid point
168      in the y-direction. We only modify mux here
169   */
170   mux = user->mux = coors[0][M/2+10].x; /* For -pi < x < pi, this should be some angle between 0 and pi/2 */
171   if (user->nonoiseinitial) {
172     for (i=xs; i < xs+xm; i++) {
173       for (j=ys; j < ys+ym; j++) {
174         xi = coors[j][i].x; yi = coors[j][i].y;
175         if ((xi == mux) && (yi == muy)) {
176           p[j][i] = 1.0;
177         }
178       }
179     }
180   } else {
181     /* Change PM_min accordingly */
182     user->PM_min = user->Pmax*PetscSinScalar(mux);
183     for (i=xs; i < xs+xm; i++) {
184       for (j=ys; j < ys+ym; j++) {
185         xi = coors[j][i].x; yi = coors[j][i].y;
186         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)));
187       }
188     }
189   }
190   PetscCall(DMDAVecRestoreArray(cda,gc,&coors));
191   PetscCall(DMDAVecRestoreArray(user->da,X,&p));
192   PetscCall(VecSum(X,&sum));
193   PetscCall(VecScale(X,1.0/sum));
194   PetscFunctionReturn(0);
195 }
196 
197 /* First advection term */
198 PetscErrorCode adv1(PetscScalar **p,PetscScalar y,PetscInt i,PetscInt j,PetscInt M,PetscScalar *p1,AppCtx *user)
199 {
200   PetscScalar f,fpos,fneg;
201   PetscFunctionBegin;
202   f   =  (y - user->ws);
203   fpos = PetscMax(f,0);
204   fneg = PetscMin(f,0);
205   if (user->st_width == 1) {
206     *p1 = fpos*(p[j][i] - p[j][i-1])/user->dx + fneg*(p[j][i+1] - p[j][i])/user->dx;
207   } else if (user->st_width == 2) {
208     *p1 = fpos*(3*p[j][i] - 4*p[j][i-1] + p[j][i-2])/(2*user->dx) + fneg*(-p[j][i+2] + 4*p[j][i+1] - 3*p[j][i])/(2*user->dx);
209   } else if (user->st_width == 3) {
210     *p1 = fpos*(2*p[j][i+1] + 3*p[j][i] - 6*p[j][i-1] + p[j][i-2])/(6*user->dx) + fneg*(-p[j][i+2] + 6*p[j][i+1] - 3*p[j][i] - 2*p[j][i-1])/(6*user->dx);
211   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No support for wider stencils");
212   PetscFunctionReturn(0);
213 }
214 
215 /* Second advection term */
216 PetscErrorCode adv2(PetscScalar **p,PetscScalar x,PetscInt i,PetscInt j,PetscInt N,PetscScalar *p2,AppCtx *user)
217 {
218   PetscScalar f,fpos,fneg;
219   PetscFunctionBegin;
220   f   = (user->ws/(2*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(x));
221   fpos = PetscMax(f,0);
222   fneg = PetscMin(f,0);
223   if (user->st_width == 1) {
224     *p2 = fpos*(p[j][i] - p[j-1][i])/user->dy + fneg*(p[j+1][i] - p[j][i])/user->dy;
225   } else if (user->st_width ==2) {
226     *p2 = fpos*(3*p[j][i] - 4*p[j-1][i] + p[j-2][i])/(2*user->dy) + fneg*(-p[j+2][i] + 4*p[j+1][i] - 3*p[j][i])/(2*user->dy);
227   } else if (user->st_width == 3) {
228     *p2 = fpos*(2*p[j+1][i] + 3*p[j][i] - 6*p[j-1][i] + p[j-2][i])/(6*user->dy) + fneg*(-p[j+2][i] + 6*p[j+1][i] - 3*p[j][i] - 2*p[j-1][i])/(6*user->dy);
229   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No support for wider stencils");
230   PetscFunctionReturn(0);
231 }
232 
233 /* Diffusion term */
234 PetscErrorCode diffuse(PetscScalar **p,PetscInt i,PetscInt j,PetscReal t,PetscScalar *p_diff,AppCtx * user)
235 {
236   PetscFunctionBeginUser;
237   if (user->st_width == 1) {
238     *p_diff = user->disper_coe*((p[j-1][i] - 2*p[j][i] + p[j+1][i])/(user->dy*user->dy));
239   } else if (user->st_width == 2) {
240     *p_diff = user->disper_coe*((-p[j-2][i] + 16*p[j-1][i] - 30*p[j][i] + 16*p[j+1][i] - p[j+2][i])/(12.0*user->dy*user->dy));
241   } else if (user->st_width == 3) {
242     *p_diff = user->disper_coe*((2*p[j-3][i] - 27*p[j-2][i] + 270*p[j-1][i] - 490*p[j][i] + 270*p[j+1][i] - 27*p[j+2][i] + 2*p[j+3][i])/(180.0*user->dy*user->dy));
243   } else SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"No support for wider stencils");
244   PetscFunctionReturn(0);
245 }
246 
247 PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
248 {
249   AppCtx         *user=(AppCtx*)ctx;
250   DM             cda;
251   DMDACoor2d     **coors;
252   PetscScalar    **p,**f,**pdot;
253   PetscInt       i,j;
254   PetscInt       xs,ys,xm,ym,M,N;
255   Vec            localX,gc,localXdot;
256   PetscScalar    p_adv1 = 0.0,p_adv2 = 0.0,p_diff; /* initialize to prevent incorrect compiler warnings */
257 
258   PetscFunctionBeginUser;
259   PetscCall(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
260   PetscCall(DMGetCoordinateDM(user->da,&cda));
261   PetscCall(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0));
262 
263   PetscCall(DMGetLocalVector(user->da,&localX));
264   PetscCall(DMGetLocalVector(user->da,&localXdot));
265 
266   PetscCall(DMGlobalToLocalBegin(user->da,X,INSERT_VALUES,localX));
267   PetscCall(DMGlobalToLocalEnd(user->da,X,INSERT_VALUES,localX));
268   PetscCall(DMGlobalToLocalBegin(user->da,Xdot,INSERT_VALUES,localXdot));
269   PetscCall(DMGlobalToLocalEnd(user->da,Xdot,INSERT_VALUES,localXdot));
270 
271   PetscCall(DMGetCoordinatesLocal(user->da,&gc));
272 
273   PetscCall(DMDAVecGetArrayRead(cda,gc,&coors));
274   PetscCall(DMDAVecGetArrayRead(user->da,localX,&p));
275   PetscCall(DMDAVecGetArrayRead(user->da,localXdot,&pdot));
276   PetscCall(DMDAVecGetArray(user->da,F,&f));
277 
278   user->disper_coe = PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda));
279   for (i=xs; i < xs+xm; i++) {
280     for (j=ys; j < ys+ym; j++) {
281       PetscCall(adv1(p,coors[j][i].y,i,j,M,&p_adv1,user));
282       PetscCall(adv2(p,coors[j][i].x,i,j,N,&p_adv2,user));
283       PetscCall(diffuse(p,i,j,t,&p_diff,user));
284       f[j][i] = -p_adv1 - p_adv2  + p_diff - pdot[j][i];
285     }
286   }
287   PetscCall(DMDAVecRestoreArrayRead(user->da,localX,&p));
288   PetscCall(DMDAVecRestoreArrayRead(user->da,localX,&pdot));
289   PetscCall(DMRestoreLocalVector(user->da,&localX));
290   PetscCall(DMRestoreLocalVector(user->da,&localXdot));
291   PetscCall(DMDAVecRestoreArray(user->da,F,&f));
292   PetscCall(DMDAVecRestoreArrayRead(cda,gc,&coors));
293 
294   PetscFunctionReturn(0);
295 }
296 
297 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat J,Mat Jpre,void *ctx)
298 {
299   AppCtx         *user=(AppCtx*)ctx;
300   DM             cda;
301   DMDACoor2d     **coors;
302   PetscInt       i,j;
303   PetscInt       xs,ys,xm,ym,M,N;
304   Vec            gc;
305   PetscScalar    val[5],xi,yi;
306   MatStencil     row,col[5];
307   PetscScalar    c1,c3,c5,c1pos,c1neg,c3pos,c3neg;
308 
309   PetscFunctionBeginUser;
310   PetscCall(DMDAGetInfo(user->da,NULL,&M,&N,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL,NULL));
311   PetscCall(DMGetCoordinateDM(user->da,&cda));
312   PetscCall(DMDAGetCorners(cda,&xs,&ys,0,&xm,&ym,0));
313 
314   PetscCall(DMGetCoordinatesLocal(user->da,&gc));
315   PetscCall(DMDAVecGetArrayRead(cda,gc,&coors));
316   for (i=xs; i < xs+xm; i++) {
317     for (j=ys; j < ys+ym; j++) {
318       PetscInt nc = 0;
319       xi = coors[j][i].x; yi = coors[j][i].y;
320       row.i = i; row.j = j;
321       c1        = (yi-user->ws)/user->dx;
322       c1pos    = PetscMax(c1,0);
323       c1neg    = PetscMin(c1,0);
324       c3        = (user->ws/(2.0*user->H))*(user->PM_min - user->Pmax*PetscSinScalar(xi))/user->dy;
325       c3pos    = PetscMax(c3,0);
326       c3neg    = PetscMin(c3,0);
327       c5        = (PetscPowScalar((user->lambda*user->ws)/(2*user->H),2)*user->q*(1.0-PetscExpScalar(-t/user->lambda)))/(user->dy*user->dy);
328       col[nc].i = i-1; col[nc].j = j;   val[nc++] = c1pos;
329       col[nc].i = i+1; col[nc].j = j;   val[nc++] = -c1neg;
330       col[nc].i = i;   col[nc].j = j-1; val[nc++] = c3pos + c5;
331       col[nc].i = i;   col[nc].j = j+1; val[nc++] = -c3neg + c5;
332       col[nc].i = i;   col[nc].j = j;   val[nc++] = -c1pos + c1neg -c3pos + c3neg -2*c5 -a;
333       PetscCall(MatSetValuesStencil(Jpre,1,&row,nc,col,val,INSERT_VALUES));
334     }
335   }
336   PetscCall(DMDAVecRestoreArrayRead(cda,gc,&coors));
337 
338   PetscCall(MatAssemblyBegin(Jpre,MAT_FINAL_ASSEMBLY));
339   PetscCall(MatAssemblyEnd(Jpre,MAT_FINAL_ASSEMBLY));
340   if (J != Jpre) {
341     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
342     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
343   }
344   PetscFunctionReturn(0);
345 }
346 
347 PetscErrorCode Parameter_settings(AppCtx *user)
348 {
349   PetscBool      flg;
350 
351   PetscFunctionBeginUser;
352 
353   /* Set default parameters */
354   user->ws     = 1.0;
355   user->H      = 5.0;
356   user->Pmax   = 2.1;
357   user->PM_min = 1.0;
358   user->lambda = 0.1;
359   user->q      = 1.0;
360   user->mux    = PetscAsinScalar(user->PM_min/user->Pmax);
361   user->sigmax = 0.1;
362   user->sigmay = 0.1;
363   user->rho    = 0.0;
364   user->xmin   = -PETSC_PI;
365   user->xmax   =  PETSC_PI;
366   user->bx     = DM_BOUNDARY_PERIODIC;
367   user->by     = DM_BOUNDARY_MIRROR;
368   user->nonoiseinitial = PETSC_FALSE;
369 
370   /*
371      ymin of -3 seems to let the unstable solution move up and leave a zero in its wake
372      with an ymin of -1 the wake is never exactly zero
373   */
374   user->ymin   = -3.0;
375   user->ymax   = 10.0;
376   user->st_width = 1;
377 
378   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-ws",&user->ws,&flg));
379   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-Inertia",&user->H,&flg));
380   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-Pmax",&user->Pmax,&flg));
381   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-PM_min",&user->PM_min,&flg));
382   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-lambda",&user->lambda,&flg));
383   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-q",&user->q,&flg));
384   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-mux",&user->mux,&flg));
385   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-sigmax",&user->sigmax,&flg));
386   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-muy",&user->muy,&flg));
387   if (flg == 0) {
388     user->muy = user->ws;
389   }
390   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-sigmay",&user->sigmay,&flg));
391   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-rho",&user->rho,&flg));
392   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-xmin",&user->xmin,&flg));
393   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-xmax",&user->xmax,&flg));
394   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-ymin",&user->ymin,&flg));
395   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-ymax",&user->ymax,&flg));
396   PetscCall(PetscOptionsGetInt(NULL,NULL,"-stencil_width",&user->st_width,&flg));
397   PetscCall(PetscOptionsGetEnum(NULL,NULL,"-bx",DMBoundaryTypes,(PetscEnum*)&user->bx,&flg));
398   PetscCall(PetscOptionsGetEnum(NULL,NULL,"-by",DMBoundaryTypes,(PetscEnum*)&user->by,&flg));
399   PetscCall(PetscOptionsGetBool(NULL,NULL,"-nonoiseinitial",&user->nonoiseinitial,&flg));
400 
401   PetscFunctionReturn(0);
402 }
403 
404 /*TEST
405 
406    build:
407       requires: !complex !single
408 
409    test:
410       args: -ts_max_steps 2
411       localrunfiles: petscopt_ex7
412 
413    test:
414       suffix: 2
415       args: -ts_max_steps 2 -snes_mf_operator
416       output_file: output/ex7_1.out
417       localrunfiles: petscopt_ex7
418       timeoutfactor: 2
419 
420    test:
421       suffix: 3
422       args: -ts_max_steps 2 -snes_mf -pc_type none
423       output_file: output/ex7_1.out
424       localrunfiles: petscopt_ex7
425       timeoutfactor: 2
426 
427 TEST*/
428