xref: /petsc/src/ts/tutorials/power_grid/ex8.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n";
2c4762a1bSJed Brown /*
3c4762a1bSJed Brown    p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy
4c4762a1bSJed Brown    xmin < x < xmax, ymin < y < ymax;
5c4762a1bSJed Brown 
6c4762a1bSJed Brown    Boundary conditions:
7c4762a1bSJed Brown    Zero dirichlet in y using ghosted values
8c4762a1bSJed Brown    Periodic in x
9c4762a1bSJed Brown 
10c4762a1bSJed Brown    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.
11c4762a1bSJed Brown    x_t = (y - ws)
12c4762a1bSJed Brown    y_t = (ws/2H)*(Pm - Pmax*sin(x) - D*(w - ws))
13c4762a1bSJed Brown 
14c4762a1bSJed Brown    In this example, we can see the effect of a fault, that zeroes the electrical power output
15c4762a1bSJed Brown    Pmax*sin(x), on the PDF. The fault on/off times can be controlled by options -tf and -tcl respectively.
16c4762a1bSJed Brown 
17c4762a1bSJed Brown */
18c4762a1bSJed Brown 
19c4762a1bSJed Brown #include <petscdm.h>
20c4762a1bSJed Brown #include <petscdmda.h>
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown /*
24c4762a1bSJed Brown    User-defined data structures and routines
25c4762a1bSJed Brown */
26c4762a1bSJed Brown typedef struct {
27c4762a1bSJed Brown   PetscScalar    ws;           /* Synchronous speed */
28c4762a1bSJed Brown   PetscScalar    H;            /* Inertia constant */
29c4762a1bSJed Brown   PetscScalar    D;            /* Damping constant */
30c4762a1bSJed Brown   PetscScalar    Pmax, Pmax_s; /* Maximum power output of generator */
31c4762a1bSJed Brown   PetscScalar    PM_min;       /* Mean mechanical power input */
32c4762a1bSJed Brown   PetscScalar    lambda;       /* correlation time */
33c4762a1bSJed Brown   PetscScalar    q;            /* noise strength */
34c4762a1bSJed Brown   PetscScalar    mux;          /* Initial average angle */
35c4762a1bSJed Brown   PetscScalar    sigmax;       /* Standard deviation of initial angle */
36c4762a1bSJed Brown   PetscScalar    muy;          /* Average speed */
37c4762a1bSJed Brown   PetscScalar    sigmay;       /* standard deviation of initial speed */
38c4762a1bSJed Brown   PetscScalar    rho;          /* Cross-correlation coefficient */
39c4762a1bSJed Brown   PetscScalar    xmin;         /* left boundary of angle */
40c4762a1bSJed Brown   PetscScalar    xmax;         /* right boundary of angle */
41c4762a1bSJed Brown   PetscScalar    ymin;         /* bottom boundary of speed */
42c4762a1bSJed Brown   PetscScalar    ymax;         /* top boundary of speed */
43c4762a1bSJed Brown   PetscScalar    dx;           /* x step size */
44c4762a1bSJed Brown   PetscScalar    dy;           /* y step size */
45c4762a1bSJed Brown   PetscScalar    disper_coe;   /* Dispersion coefficient */
46c4762a1bSJed Brown   DM             da;
47c4762a1bSJed Brown   PetscInt       st_width; /* Stencil width */
48c4762a1bSJed Brown   DMBoundaryType bx;       /* x boundary type */
49c4762a1bSJed Brown   DMBoundaryType by;       /* y boundary type */
50c4762a1bSJed Brown   PetscReal      tf, tcl;  /* Fault incidence and clearing times */
51c4762a1bSJed Brown } AppCtx;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *);
54c4762a1bSJed Brown PetscErrorCode ini_bou(Vec, AppCtx *);
55c4762a1bSJed Brown PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *);
56c4762a1bSJed Brown PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
57c4762a1bSJed Brown PetscErrorCode PostStep(TS);
58c4762a1bSJed Brown 
59d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
60d71ae5a4SJacob Faibussowitsch {
61c4762a1bSJed Brown   Vec         x;    /* Solution vector */
62c4762a1bSJed Brown   TS          ts;   /* Time-stepping context */
63c4762a1bSJed Brown   AppCtx      user; /* Application context */
64c4762a1bSJed Brown   PetscViewer viewer;
65c4762a1bSJed Brown 
66327415f7SBarry Smith   PetscFunctionBeginUser;
679566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, "petscopt_ex8", help));
68c4762a1bSJed Brown 
69c4762a1bSJed Brown   /* Get physics and time parameters */
709566063dSJacob Faibussowitsch   PetscCall(Parameter_settings(&user));
71c4762a1bSJed Brown   /* Create a 2D DA with dof = 1 */
729566063dSJacob Faibussowitsch   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));
739566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.da));
749566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.da));
75c4762a1bSJed Brown   /* Set x and y coordinates */
769566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(user.da, user.xmin, user.xmax, user.ymin, user.ymax, 0, 0));
779566063dSJacob Faibussowitsch   PetscCall(DMDASetCoordinateName(user.da, 0, "X - the angle"));
789566063dSJacob Faibussowitsch   PetscCall(DMDASetCoordinateName(user.da, 1, "Y - the speed"));
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   /* Get global vector x from DM  */
819566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.da, &x));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(ini_bou(x, &user));
849566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "ini_x", FILE_MODE_WRITE, &viewer));
859566063dSJacob Faibussowitsch   PetscCall(VecView(x, viewer));
869566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
899566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, user.da));
909566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
919566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSARKIMEX));
929566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
939566063dSJacob Faibussowitsch   /*  PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobian,&user)); */
949566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(ts, &user));
959566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .005));
969566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
979566063dSJacob Faibussowitsch   PetscCall(TSSetPostStep(ts, PostStep));
989566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
99c4762a1bSJed Brown 
1009566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "fin_x", FILE_MODE_WRITE, &viewer));
1019566063dSJacob Faibussowitsch   PetscCall(VecView(x, viewer));
1029566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
103c4762a1bSJed Brown 
1049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1059566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.da));
1069566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1079566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
108b122ec5aSJacob Faibussowitsch   return 0;
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStep(TS ts)
112d71ae5a4SJacob Faibussowitsch {
113c4762a1bSJed Brown   Vec         X;
114c4762a1bSJed Brown   AppCtx     *user;
115c4762a1bSJed Brown   PetscReal   t;
116c4762a1bSJed Brown   PetscScalar asum;
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   PetscFunctionBegin;
1199566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
1209566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
1219566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &X));
122c4762a1bSJed Brown   /*
123c4762a1bSJed Brown   if (t >= .2) {
1249566063dSJacob Faibussowitsch     PetscCall(TSGetSolution(ts,&X));
1259566063dSJacob Faibussowitsch     PetscCall(VecView(X,PETSC_VIEWER_BINARY_WORLD));
126c4762a1bSJed Brown     exit(0);
127c4762a1bSJed Brown      results in initial conditions after fault in binaryoutput
128c4762a1bSJed Brown   }*/
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   if ((t > user->tf) && (t < user->tcl)) user->Pmax = 0.0; /* A short-circuit that drives the electrical power output (Pmax*sin(delta)) to zero */
131c4762a1bSJed Brown   else user->Pmax = user->Pmax_s;
132c4762a1bSJed Brown 
1339566063dSJacob Faibussowitsch   PetscCall(VecSum(X, &asum));
1349566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p) at t = %f = %f\n", (double)t, (double)(asum)));
135*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
138d71ae5a4SJacob Faibussowitsch PetscErrorCode ini_bou(Vec X, AppCtx *user)
139d71ae5a4SJacob Faibussowitsch {
140c4762a1bSJed Brown   DM            cda;
141c4762a1bSJed Brown   DMDACoor2d  **coors;
142c4762a1bSJed Brown   PetscScalar **p;
143c4762a1bSJed Brown   Vec           gc;
144c4762a1bSJed Brown   PetscInt      M, N, Ir, J;
145c4762a1bSJed Brown   PetscMPIInt   rank;
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   PetscFunctionBeginUser;
1489566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1499566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
1509371c9d4SSatish Balay   user->dx = (user->xmax - user->xmin) / (M - 1);
1519371c9d4SSatish Balay   user->dy = (user->ymax - user->ymin) / (N - 1);
1529566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(user->da, &cda));
1539566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(user->da, &gc));
1549566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
1559566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->da, X, &p));
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   /* Point mass at (mux,muy) */
15863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Original user->mux = %f, user->muy = %f\n", (double)PetscRealPart(user->mux), (double)PetscRealPart(user->muy)));
1599566063dSJacob Faibussowitsch   PetscCall(DMDAGetLogicalCoordinate(user->da, user->mux, user->muy, 0.0, &Ir, &J, NULL, &user->mux, &user->muy, NULL));
160c4762a1bSJed Brown   user->PM_min = user->Pmax * PetscSinScalar(user->mux);
1619371c9d4SSatish Balay   PetscCall(
1629371c9d4SSatish Balay     PetscPrintf(PETSC_COMM_WORLD, "Corrected user->mux = %f, user->muy = %f user->PM_min = %f,user->dx = %f\n", (double)PetscRealPart(user->mux), (double)PetscRealPart(user->muy), (double)PetscRealPart(user->PM_min), (double)PetscRealPart(user->dx)));
163ad540459SPierre Jolivet   if (Ir > -1 && J > -1) p[J][Ir] = 1.0;
164c4762a1bSJed Brown 
1659566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
1669566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->da, X, &p));
167*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170c4762a1bSJed Brown /* First advection term */
171d71ae5a4SJacob Faibussowitsch PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user)
172d71ae5a4SJacob Faibussowitsch {
173c4762a1bSJed Brown   PetscScalar f, fpos, fneg;
174c4762a1bSJed Brown   PetscFunctionBegin;
175c4762a1bSJed Brown   f    = (y - user->ws);
176c4762a1bSJed Brown   fpos = PetscMax(f, 0);
177c4762a1bSJed Brown   fneg = PetscMin(f, 0);
178c4762a1bSJed Brown   if (user->st_width == 1) {
179c4762a1bSJed Brown     *p1 = fpos * (p[j][i] - p[j][i - 1]) / user->dx + fneg * (p[j][i + 1] - p[j][i]) / user->dx;
180c4762a1bSJed Brown   } else if (user->st_width == 2) {
181c4762a1bSJed Brown     *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);
182c4762a1bSJed Brown   } else if (user->st_width == 3) {
183c4762a1bSJed Brown     *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);
184c4762a1bSJed Brown   }
185c4762a1bSJed Brown   /* *p1 = f*(p[j][i+1] - p[j][i-1])/user->dx;*/
186*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
189c4762a1bSJed Brown /* Second advection term */
190d71ae5a4SJacob Faibussowitsch PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscScalar y, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user)
191d71ae5a4SJacob Faibussowitsch {
192c4762a1bSJed Brown   PetscScalar f, fpos, fneg;
193c4762a1bSJed Brown   PetscFunctionBegin;
194c4762a1bSJed Brown   f    = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x) - user->D * (y - user->ws));
195c4762a1bSJed Brown   fpos = PetscMax(f, 0);
196c4762a1bSJed Brown   fneg = PetscMin(f, 0);
197c4762a1bSJed Brown   if (user->st_width == 1) {
198c4762a1bSJed Brown     *p2 = fpos * (p[j][i] - p[j - 1][i]) / user->dy + fneg * (p[j + 1][i] - p[j][i]) / user->dy;
199c4762a1bSJed Brown   } else if (user->st_width == 2) {
200c4762a1bSJed Brown     *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);
201c4762a1bSJed Brown   } else if (user->st_width == 3) {
202c4762a1bSJed Brown     *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);
203c4762a1bSJed Brown   }
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /* *p2 = f*(p[j+1][i] - p[j-1][i])/user->dy;*/
206*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207c4762a1bSJed Brown }
208c4762a1bSJed Brown 
209c4762a1bSJed Brown /* Diffusion term */
210d71ae5a4SJacob Faibussowitsch PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user)
211d71ae5a4SJacob Faibussowitsch {
212c4762a1bSJed Brown   PetscFunctionBeginUser;
213c4762a1bSJed Brown   if (user->st_width == 1) {
214c4762a1bSJed Brown     *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy));
215c4762a1bSJed Brown   } else if (user->st_width == 2) {
216c4762a1bSJed Brown     *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));
217c4762a1bSJed Brown   } else if (user->st_width == 3) {
218c4762a1bSJed Brown     *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));
219c4762a1bSJed Brown   }
220*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
221c4762a1bSJed Brown }
222c4762a1bSJed Brown 
223d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
224d71ae5a4SJacob Faibussowitsch {
225c4762a1bSJed Brown   AppCtx       *user = (AppCtx *)ctx;
226c4762a1bSJed Brown   DM            cda;
227c4762a1bSJed Brown   DMDACoor2d  **coors;
228c4762a1bSJed Brown   PetscScalar **p, **f, **pdot;
229c4762a1bSJed Brown   PetscInt      i, j;
230c4762a1bSJed Brown   PetscInt      xs, ys, xm, ym, M, N;
231c4762a1bSJed Brown   Vec           localX, gc, localXdot;
232c4762a1bSJed Brown   PetscScalar   p_adv1 = 0.0, p_adv2 = 0.0, p_diff = 0;
233c4762a1bSJed Brown   PetscScalar   diffuse1, gamma;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   PetscFunctionBeginUser;
2369566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
2379566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(user->da, &cda));
2389566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
239c4762a1bSJed Brown 
2409566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->da, &localX));
2419566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->da, &localXdot));
242c4762a1bSJed Brown 
2439566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX));
2449566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX));
2459566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot));
2469566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot));
247c4762a1bSJed Brown 
2489566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(user->da, &gc));
249c4762a1bSJed Brown 
2509566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
2519566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(user->da, localX, &p));
2529566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot));
2539566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->da, F, &f));
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   gamma            = user->D * user->ws / (2 * user->H);
256c4762a1bSJed Brown   diffuse1         = user->lambda * user->lambda * user->q / (user->lambda * gamma + 1) * (1.0 - PetscExpScalar(-t * (gamma + 1.0) / user->lambda));
257c4762a1bSJed Brown   user->disper_coe = user->ws * user->ws / (4 * user->H * user->H) * diffuse1;
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
260c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
2619566063dSJacob Faibussowitsch       PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user));
2629566063dSJacob Faibussowitsch       PetscCall(adv2(p, coors[j][i].x, coors[j][i].y, i, j, N, &p_adv2, user));
2639566063dSJacob Faibussowitsch       PetscCall(diffuse(p, i, j, t, &p_diff, user));
264c4762a1bSJed Brown       f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i];
265c4762a1bSJed Brown     }
266c4762a1bSJed Brown   }
2679566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p));
2689566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot));
2699566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->da, &localX));
2709566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->da, &localXdot));
2719566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->da, F, &f));
2729566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
273c4762a1bSJed Brown 
274*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
275c4762a1bSJed Brown }
276c4762a1bSJed Brown 
277d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx)
278d71ae5a4SJacob Faibussowitsch {
279c4762a1bSJed Brown   AppCtx      *user = (AppCtx *)ctx;
280c4762a1bSJed Brown   DM           cda;
281c4762a1bSJed Brown   DMDACoor2d **coors;
282c4762a1bSJed Brown   PetscInt     i, j;
283c4762a1bSJed Brown   PetscInt     xs, ys, xm, ym, M, N;
284c4762a1bSJed Brown   Vec          gc;
285c4762a1bSJed Brown   PetscScalar  val[5], xi, yi;
286c4762a1bSJed Brown   MatStencil   row, col[5];
287c4762a1bSJed Brown   PetscScalar  c1, c3, c5, c1pos, c1neg, c3pos, c3neg;
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   PetscFunctionBeginUser;
2909566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
2919566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(user->da, &cda));
2929566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(user->da, &gc));
2959566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
296c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
297c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
298c4762a1bSJed Brown       PetscInt nc = 0;
2999371c9d4SSatish Balay       xi          = coors[j][i].x;
3009371c9d4SSatish Balay       yi          = coors[j][i].y;
3019371c9d4SSatish Balay       row.i       = i;
3029371c9d4SSatish Balay       row.j       = j;
303c4762a1bSJed Brown       c1          = (yi - user->ws) / user->dx;
304c4762a1bSJed Brown       c1pos       = PetscMax(c1, 0);
305c4762a1bSJed Brown       c1neg       = PetscMin(c1, 0);
306c4762a1bSJed Brown       c3          = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi) - user->D * (yi - user->ws)) / user->dy;
307c4762a1bSJed Brown       c3pos       = PetscMax(c3, 0);
308c4762a1bSJed Brown       c3neg       = PetscMin(c3, 0);
309c4762a1bSJed Brown       c5          = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy);
3109371c9d4SSatish Balay       col[nc].i   = i - 1;
3119371c9d4SSatish Balay       col[nc].j   = j;
3129371c9d4SSatish Balay       val[nc++]   = c1pos;
3139371c9d4SSatish Balay       col[nc].i   = i + 1;
3149371c9d4SSatish Balay       col[nc].j   = j;
3159371c9d4SSatish Balay       val[nc++]   = -c1neg;
3169371c9d4SSatish Balay       col[nc].i   = i;
3179371c9d4SSatish Balay       col[nc].j   = j - 1;
3189371c9d4SSatish Balay       val[nc++]   = c3pos + c5;
3199371c9d4SSatish Balay       col[nc].i   = i;
3209371c9d4SSatish Balay       col[nc].j   = j + 1;
3219371c9d4SSatish Balay       val[nc++]   = -c3neg + c5;
3229371c9d4SSatish Balay       col[nc].i   = i;
3239371c9d4SSatish Balay       col[nc].j   = j;
3249371c9d4SSatish Balay       val[nc++]   = -c1pos + c1neg - c3pos + c3neg - 2 * c5 - a;
3259566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES));
326c4762a1bSJed Brown     }
327c4762a1bSJed Brown   }
3289566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
329c4762a1bSJed Brown 
3309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
3319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
332c4762a1bSJed Brown   if (J != Jpre) {
3339566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
3349566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
335c4762a1bSJed Brown   }
336*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
337c4762a1bSJed Brown }
338c4762a1bSJed Brown 
339d71ae5a4SJacob Faibussowitsch PetscErrorCode Parameter_settings(AppCtx *user)
340d71ae5a4SJacob Faibussowitsch {
341c4762a1bSJed Brown   PetscBool flg;
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   PetscFunctionBeginUser;
344c4762a1bSJed Brown   /* Set default parameters */
345c4762a1bSJed Brown   user->ws   = 1.0;
346c4762a1bSJed Brown   user->H    = 5.0;
347c4762a1bSJed Brown   user->D    = 0.0;
348c4762a1bSJed Brown   user->Pmax = user->Pmax_s = 2.1;
349c4762a1bSJed Brown   user->PM_min              = 1.0;
350c4762a1bSJed Brown   user->lambda              = 0.1;
351c4762a1bSJed Brown   user->q                   = 1.0;
352c4762a1bSJed Brown   user->mux                 = PetscAsinScalar(user->PM_min / user->Pmax);
353c4762a1bSJed Brown   user->sigmax              = 0.1;
354c4762a1bSJed Brown   user->sigmay              = 0.1;
355c4762a1bSJed Brown   user->rho                 = 0.0;
356c4762a1bSJed Brown   user->xmin                = -PETSC_PI;
357c4762a1bSJed Brown   user->xmax                = PETSC_PI;
358c4762a1bSJed Brown   user->bx                  = DM_BOUNDARY_PERIODIC;
359c4762a1bSJed Brown   user->by                  = DM_BOUNDARY_GHOSTED;
360c4762a1bSJed Brown   user->tf = user->tcl = -1;
361c4762a1bSJed Brown   user->ymin           = -2.0;
362c4762a1bSJed Brown   user->ymax           = 2.0;
363c4762a1bSJed Brown   user->st_width       = 1;
364c4762a1bSJed Brown 
3659566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg));
3669566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg));
3679566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-D", &user->D, &flg));
3689566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg));
3699566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg));
3709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg));
3719566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg));
3729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg));
3739566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg));
3745f80ce2aSJacob Faibussowitsch   if (flg == 0) user->muy = user->ws;
3759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg));
3769566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg));
3779566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg));
3789566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg));
3799566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-stencil_width", &user->st_width, &flg));
3809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-bx", DMBoundaryTypes, (PetscEnum *)&user->bx, &flg));
3819566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-by", DMBoundaryTypes, (PetscEnum *)&user->by, &flg));
3829566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tf", &user->tf, &flg));
3839566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-tcl", &user->tcl, &flg));
384*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
385c4762a1bSJed Brown }
386c4762a1bSJed Brown 
387c4762a1bSJed Brown /*TEST
388c4762a1bSJed Brown 
389c4762a1bSJed Brown    build:
390c4762a1bSJed Brown       requires: !complex x
391c4762a1bSJed Brown 
392c4762a1bSJed Brown    test:
393c4762a1bSJed Brown       args: -ts_max_steps 1
394c4762a1bSJed Brown       localrunfiles: petscopt_ex8
395c4762a1bSJed Brown       timeoutfactor: 3
396c4762a1bSJed Brown 
397c4762a1bSJed Brown TEST*/
398