xref: /petsc/src/ts/tutorials/power_grid/ex7.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
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   Vec         x;    /* Solution vector */
56   TS          ts;   /* Time-stepping context */
57   AppCtx      user; /* Application context */
58   PetscViewer viewer;
59 
60   PetscFunctionBeginUser;
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   Vec          X, gc;
107   AppCtx      *user;
108   PetscScalar  sum = 0, asum;
109   PetscReal    t, **p;
110   DMDACoor2d **coors;
111   DM           cda;
112   PetscInt     i, j, xs, ys, xm, ym;
113 
114   PetscFunctionBegin;
115   PetscCall(TSGetApplicationContext(ts, &user));
116   PetscCall(TSGetTime(ts, &t));
117   PetscCall(TSGetSolution(ts, &X));
118 
119   PetscCall(DMGetCoordinateDM(user->da, &cda));
120   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
121   PetscCall(DMGetCoordinates(user->da, &gc));
122   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
123   PetscCall(DMDAVecGetArrayRead(user->da, X, &p));
124   for (i = xs; i < xs + xm; i++) {
125     for (j = ys; j < ys + ym; j++) {
126       if (coors[j][i].y < 5) sum += p[j][i];
127     }
128   }
129   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
130   PetscCall(DMDAVecRestoreArrayRead(user->da, X, &p));
131   PetscCallMPI(MPI_Allreduce(&sum, &asum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)ts)));
132   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p)*dw*dtheta at t = %f = %f\n", (double)t, (double)(asum)));
133   if (sum < 1.0e-2) {
134     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
135     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Exiting TS as the integral of PDF is almost zero\n"));
136   }
137   PetscFunctionReturn(0);
138 }
139 
140 PetscErrorCode ini_bou(Vec X, AppCtx *user) {
141   DM            cda;
142   DMDACoor2d  **coors;
143   PetscScalar **p;
144   Vec           gc;
145   PetscInt      i, j;
146   PetscInt      xs, ys, xm, ym, M, N;
147   PetscScalar   xi, yi;
148   PetscScalar   sigmax = user->sigmax, sigmay = user->sigmay;
149   PetscScalar   rho = user->rho;
150   PetscScalar   muy = user->muy, mux;
151   PetscMPIInt   rank;
152   PetscScalar   sum;
153 
154   PetscFunctionBeginUser;
155   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
156   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
157   user->dx = (user->xmax - user->xmin) / (M - 1);
158   user->dy = (user->ymax - user->ymin) / (N - 1);
159   PetscCall(DMGetCoordinateDM(user->da, &cda));
160   PetscCall(DMGetCoordinates(user->da, &gc));
161   PetscCall(DMDAVecGetArray(cda, gc, &coors));
162   PetscCall(DMDAVecGetArray(user->da, X, &p));
163   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
164 
165   /* mux and muy need to be grid points in the x and y-direction otherwise the solution goes unstable
166      muy is set by choosing the y domain, no. of grid points along y-direction so that muy is a grid point
167      in the y-direction. We only modify mux here
168   */
169   mux = user->mux = coors[0][M / 2 + 10].x; /* For -pi < x < pi, this should be some angle between 0 and pi/2 */
170   if (user->nonoiseinitial) {
171     for (i = xs; i < xs + xm; i++) {
172       for (j = ys; j < ys + ym; j++) {
173         xi = coors[j][i].x;
174         yi = coors[j][i].y;
175         if ((xi == mux) && (yi == muy)) { p[j][i] = 1.0; }
176       }
177     }
178   } else {
179     /* Change PM_min accordingly */
180     user->PM_min = user->Pmax * PetscSinScalar(mux);
181     for (i = xs; i < xs + xm; i++) {
182       for (j = ys; j < ys + ym; j++) {
183         xi      = coors[j][i].x;
184         yi      = coors[j][i].y;
185         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)));
186       }
187     }
188   }
189   PetscCall(DMDAVecRestoreArray(cda, gc, &coors));
190   PetscCall(DMDAVecRestoreArray(user->da, X, &p));
191   PetscCall(VecSum(X, &sum));
192   PetscCall(VecScale(X, 1.0 / sum));
193   PetscFunctionReturn(0);
194 }
195 
196 /* First advection term */
197 PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user) {
198   PetscScalar f, fpos, fneg;
199   PetscFunctionBegin;
200   f    = (y - user->ws);
201   fpos = PetscMax(f, 0);
202   fneg = PetscMin(f, 0);
203   if (user->st_width == 1) {
204     *p1 = fpos * (p[j][i] - p[j][i - 1]) / user->dx + fneg * (p[j][i + 1] - p[j][i]) / user->dx;
205   } else if (user->st_width == 2) {
206     *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);
207   } else if (user->st_width == 3) {
208     *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);
209   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils");
210   PetscFunctionReturn(0);
211 }
212 
213 /* Second advection term */
214 PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user) {
215   PetscScalar f, fpos, fneg;
216   PetscFunctionBegin;
217   f    = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x));
218   fpos = PetscMax(f, 0);
219   fneg = PetscMin(f, 0);
220   if (user->st_width == 1) {
221     *p2 = fpos * (p[j][i] - p[j - 1][i]) / user->dy + fneg * (p[j + 1][i] - p[j][i]) / user->dy;
222   } else if (user->st_width == 2) {
223     *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);
224   } else if (user->st_width == 3) {
225     *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);
226   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils");
227   PetscFunctionReturn(0);
228 }
229 
230 /* Diffusion term */
231 PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user) {
232   PetscFunctionBeginUser;
233   if (user->st_width == 1) {
234     *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy));
235   } else if (user->st_width == 2) {
236     *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));
237   } else if (user->st_width == 3) {
238     *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));
239   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils");
240   PetscFunctionReturn(0);
241 }
242 
243 PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
244   AppCtx       *user = (AppCtx *)ctx;
245   DM            cda;
246   DMDACoor2d  **coors;
247   PetscScalar **p, **f, **pdot;
248   PetscInt      i, j;
249   PetscInt      xs, ys, xm, ym, M, N;
250   Vec           localX, gc, localXdot;
251   PetscScalar   p_adv1 = 0.0, p_adv2 = 0.0, p_diff; /* initialize to prevent incorrect compiler warnings */
252 
253   PetscFunctionBeginUser;
254   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
255   PetscCall(DMGetCoordinateDM(user->da, &cda));
256   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
257 
258   PetscCall(DMGetLocalVector(user->da, &localX));
259   PetscCall(DMGetLocalVector(user->da, &localXdot));
260 
261   PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX));
262   PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX));
263   PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot));
264   PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot));
265 
266   PetscCall(DMGetCoordinatesLocal(user->da, &gc));
267 
268   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
269   PetscCall(DMDAVecGetArrayRead(user->da, localX, &p));
270   PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot));
271   PetscCall(DMDAVecGetArray(user->da, F, &f));
272 
273   user->disper_coe = PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda));
274   for (i = xs; i < xs + xm; i++) {
275     for (j = ys; j < ys + ym; j++) {
276       PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user));
277       PetscCall(adv2(p, coors[j][i].x, i, j, N, &p_adv2, user));
278       PetscCall(diffuse(p, i, j, t, &p_diff, user));
279       f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i];
280     }
281   }
282   PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p));
283   PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot));
284   PetscCall(DMRestoreLocalVector(user->da, &localX));
285   PetscCall(DMRestoreLocalVector(user->da, &localXdot));
286   PetscCall(DMDAVecRestoreArray(user->da, F, &f));
287   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
288 
289   PetscFunctionReturn(0);
290 }
291 
292 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx) {
293   AppCtx      *user = (AppCtx *)ctx;
294   DM           cda;
295   DMDACoor2d **coors;
296   PetscInt     i, j;
297   PetscInt     xs, ys, xm, ym, M, N;
298   Vec          gc;
299   PetscScalar  val[5], xi, yi;
300   MatStencil   row, col[5];
301   PetscScalar  c1, c3, c5, c1pos, c1neg, c3pos, c3neg;
302 
303   PetscFunctionBeginUser;
304   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
305   PetscCall(DMGetCoordinateDM(user->da, &cda));
306   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
307 
308   PetscCall(DMGetCoordinatesLocal(user->da, &gc));
309   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
310   for (i = xs; i < xs + xm; i++) {
311     for (j = ys; j < ys + ym; j++) {
312       PetscInt nc = 0;
313       xi          = coors[j][i].x;
314       yi          = coors[j][i].y;
315       row.i       = i;
316       row.j       = j;
317       c1          = (yi - user->ws) / user->dx;
318       c1pos       = PetscMax(c1, 0);
319       c1neg       = PetscMin(c1, 0);
320       c3          = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi)) / user->dy;
321       c3pos       = PetscMax(c3, 0);
322       c3neg       = PetscMin(c3, 0);
323       c5          = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy);
324       col[nc].i   = i - 1;
325       col[nc].j   = j;
326       val[nc++]   = c1pos;
327       col[nc].i   = i + 1;
328       col[nc].j   = j;
329       val[nc++]   = -c1neg;
330       col[nc].i   = i;
331       col[nc].j   = j - 1;
332       val[nc++]   = c3pos + c5;
333       col[nc].i   = i;
334       col[nc].j   = j + 1;
335       val[nc++]   = -c3neg + c5;
336       col[nc].i   = i;
337       col[nc].j   = j;
338       val[nc++]   = -c1pos + c1neg - c3pos + c3neg - 2 * c5 - a;
339       PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES));
340     }
341   }
342   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
343 
344   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
345   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
346   if (J != Jpre) {
347     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
348     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
349   }
350   PetscFunctionReturn(0);
351 }
352 
353 PetscErrorCode Parameter_settings(AppCtx *user) {
354   PetscBool flg;
355 
356   PetscFunctionBeginUser;
357 
358   /* Set default parameters */
359   user->ws             = 1.0;
360   user->H              = 5.0;
361   user->Pmax           = 2.1;
362   user->PM_min         = 1.0;
363   user->lambda         = 0.1;
364   user->q              = 1.0;
365   user->mux            = PetscAsinScalar(user->PM_min / user->Pmax);
366   user->sigmax         = 0.1;
367   user->sigmay         = 0.1;
368   user->rho            = 0.0;
369   user->xmin           = -PETSC_PI;
370   user->xmax           = PETSC_PI;
371   user->bx             = DM_BOUNDARY_PERIODIC;
372   user->by             = DM_BOUNDARY_MIRROR;
373   user->nonoiseinitial = PETSC_FALSE;
374 
375   /*
376      ymin of -3 seems to let the unstable solution move up and leave a zero in its wake
377      with an ymin of -1 the wake is never exactly zero
378   */
379   user->ymin     = -3.0;
380   user->ymax     = 10.0;
381   user->st_width = 1;
382 
383   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg));
384   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg));
385   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg));
386   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg));
387   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg));
388   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg));
389   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg));
390   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmax", &user->sigmax, &flg));
391   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg));
392   if (flg == 0) { user->muy = user->ws; }
393   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmay", &user->sigmay, &flg));
394   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-rho", &user->rho, &flg));
395   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg));
396   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg));
397   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg));
398   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg));
399   PetscCall(PetscOptionsGetInt(NULL, NULL, "-stencil_width", &user->st_width, &flg));
400   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-bx", DMBoundaryTypes, (PetscEnum *)&user->bx, &flg));
401   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-by", DMBoundaryTypes, (PetscEnum *)&user->by, &flg));
402   PetscCall(PetscOptionsGetBool(NULL, NULL, "-nonoiseinitial", &user->nonoiseinitial, &flg));
403 
404   PetscFunctionReturn(0);
405 }
406 
407 /*TEST
408 
409    build:
410       requires: !complex !single
411 
412    test:
413       args: -ts_max_steps 2
414       localrunfiles: petscopt_ex7
415 
416    test:
417       suffix: 2
418       args: -ts_max_steps 2 -snes_mf_operator
419       output_file: output/ex7_1.out
420       localrunfiles: petscopt_ex7
421       timeoutfactor: 2
422 
423    test:
424       suffix: 3
425       args: -ts_max_steps 2 -snes_mf -pc_type none
426       output_file: output/ex7_1.out
427       localrunfiles: petscopt_ex7
428       timeoutfactor: 2
429 
430 TEST*/
431