xref: /petsc/src/ts/tutorials/power_grid/ex7.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Time-dependent PDE in 2d for calculating joint PDF. \n";
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    p_t = -x_t*p_x -y_t*p_y + f(t)*p_yy
5c4762a1bSJed Brown    xmin < x < xmax, ymin < y < ymax;
6c4762a1bSJed Brown 
7c4762a1bSJed Brown    Boundary conditions Neumman using mirror values
8c4762a1bSJed Brown 
9c4762a1bSJed 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.
10c4762a1bSJed Brown    x_t = (y - ws)  y_t = (ws/2H)*(Pm - Pmax*sin(x))
11c4762a1bSJed Brown 
12c4762a1bSJed Brown */
13c4762a1bSJed Brown 
14c4762a1bSJed Brown #include <petscdm.h>
15c4762a1bSJed Brown #include <petscdmda.h>
16c4762a1bSJed Brown #include <petscts.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown /*
19c4762a1bSJed Brown    User-defined data structures and routines
20c4762a1bSJed Brown */
21c4762a1bSJed Brown typedef struct {
22c4762a1bSJed Brown   PetscScalar    ws;         /* Synchronous speed */
23c4762a1bSJed Brown   PetscScalar    H;          /* Inertia constant */
24c4762a1bSJed Brown   PetscScalar    D;          /* Damping constant */
25c4762a1bSJed Brown   PetscScalar    Pmax;       /* Maximum power output of generator */
26c4762a1bSJed Brown   PetscScalar    PM_min;     /* Mean mechanical power input */
27c4762a1bSJed Brown   PetscScalar    lambda;     /* correlation time */
28c4762a1bSJed Brown   PetscScalar    q;          /* noise strength */
29c4762a1bSJed Brown   PetscScalar    mux;        /* Initial average angle */
30c4762a1bSJed Brown   PetscScalar    sigmax;     /* Standard deviation of initial angle */
31c4762a1bSJed Brown   PetscScalar    muy;        /* Average speed */
32c4762a1bSJed Brown   PetscScalar    sigmay;     /* standard deviation of initial speed */
33c4762a1bSJed Brown   PetscScalar    rho;        /* Cross-correlation coefficient */
34c4762a1bSJed Brown   PetscScalar    xmin;       /* left boundary of angle */
35c4762a1bSJed Brown   PetscScalar    xmax;       /* right boundary of angle */
36c4762a1bSJed Brown   PetscScalar    ymin;       /* bottom boundary of speed */
37c4762a1bSJed Brown   PetscScalar    ymax;       /* top boundary of speed */
38c4762a1bSJed Brown   PetscScalar    dx;         /* x step size */
39c4762a1bSJed Brown   PetscScalar    dy;         /* y step size */
40c4762a1bSJed Brown   PetscScalar    disper_coe; /* Dispersion coefficient */
41c4762a1bSJed Brown   DM             da;
42c4762a1bSJed Brown   PetscInt       st_width; /* Stencil width */
43c4762a1bSJed Brown   DMBoundaryType bx;       /* x boundary type */
44c4762a1bSJed Brown   DMBoundaryType by;       /* y boundary type */
45c4762a1bSJed Brown   PetscBool      nonoiseinitial;
46c4762a1bSJed Brown } AppCtx;
47c4762a1bSJed Brown 
48c4762a1bSJed Brown PetscErrorCode Parameter_settings(AppCtx *);
49c4762a1bSJed Brown PetscErrorCode ini_bou(Vec, AppCtx *);
50c4762a1bSJed Brown PetscErrorCode IFunction(TS, PetscReal, Vec, Vec, Vec, void *);
51c4762a1bSJed Brown PetscErrorCode IJacobian(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
52c4762a1bSJed Brown PetscErrorCode PostStep(TS);
53c4762a1bSJed Brown 
54d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
55d71ae5a4SJacob Faibussowitsch {
56c4762a1bSJed Brown   Vec         x;    /* Solution vector */
57c4762a1bSJed Brown   TS          ts;   /* Time-stepping context */
58c4762a1bSJed Brown   AppCtx      user; /* Application context */
59c4762a1bSJed Brown   PetscViewer viewer;
60c4762a1bSJed Brown 
61327415f7SBarry Smith   PetscFunctionBeginUser;
629566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, "petscopt_ex7", help));
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   /* Get physics and time parameters */
659566063dSJacob Faibussowitsch   PetscCall(Parameter_settings(&user));
66c4762a1bSJed Brown   /* Create a 2D DA with dof = 1 */
679566063dSJacob 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));
689566063dSJacob Faibussowitsch   PetscCall(DMSetFromOptions(user.da));
699566063dSJacob Faibussowitsch   PetscCall(DMSetUp(user.da));
70c4762a1bSJed Brown   /* Set x and y coordinates */
719566063dSJacob Faibussowitsch   PetscCall(DMDASetUniformCoordinates(user.da, user.xmin, user.xmax, user.ymin, user.ymax, 0, 0));
729566063dSJacob Faibussowitsch   PetscCall(DMDASetCoordinateName(user.da, 0, "X - the angle"));
739566063dSJacob Faibussowitsch   PetscCall(DMDASetCoordinateName(user.da, 1, "Y - the speed"));
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   /* Get global vector x from DM  */
769566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(user.da, &x));
77c4762a1bSJed Brown 
789566063dSJacob Faibussowitsch   PetscCall(ini_bou(x, &user));
799566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "ini_x", FILE_MODE_WRITE, &viewer));
809566063dSJacob Faibussowitsch   PetscCall(VecView(x, viewer));
819566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
849566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, user.da));
859566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
869566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSARKIMEX));
879566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
889566063dSJacob Faibussowitsch   /*  PetscCall(TSSetIJacobian(ts,NULL,NULL,IJacobian,&user));  */
899566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(ts, &user));
909566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .005));
919566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
929566063dSJacob Faibussowitsch   PetscCall(TSSetPostStep(ts, PostStep));
939566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "fin_x", FILE_MODE_WRITE, &viewer));
969566063dSJacob Faibussowitsch   PetscCall(VecView(x, viewer));
979566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
98c4762a1bSJed Brown 
999566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1009566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&user.da));
1019566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
1029566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
103b122ec5aSJacob Faibussowitsch   return 0;
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStep(TS ts)
107d71ae5a4SJacob Faibussowitsch {
108c4762a1bSJed Brown   Vec          X, gc;
109c4762a1bSJed Brown   AppCtx      *user;
110c4762a1bSJed Brown   PetscScalar  sum = 0, asum;
111c4762a1bSJed Brown   PetscReal    t, **p;
112c4762a1bSJed Brown   DMDACoor2d **coors;
113c4762a1bSJed Brown   DM           cda;
114c4762a1bSJed Brown   PetscInt     i, j, xs, ys, xm, ym;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBegin;
1179566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
1189566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &t));
1199566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &X));
120c4762a1bSJed Brown 
1219566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(user->da, &cda));
1229566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
1239566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(user->da, &gc));
1249566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
1259566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(user->da, X, &p));
126c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
127c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
128c4762a1bSJed Brown       if (coors[j][i].y < 5) sum += p[j][i];
129c4762a1bSJed Brown     }
130c4762a1bSJed Brown   }
1319566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
1329566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(user->da, X, &p));
1339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Allreduce(&sum, &asum, 1, MPIU_SCALAR, MPIU_SUM, PetscObjectComm((PetscObject)ts)));
1349566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "sum(p)*dw*dtheta at t = %f = %f\n", (double)t, (double)(asum)));
135c4762a1bSJed Brown   if (sum < 1.0e-2) {
1369566063dSJacob Faibussowitsch     PetscCall(TSSetConvergedReason(ts, TS_CONVERGED_USER));
1379566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Exiting TS as the integral of PDF is almost zero\n"));
138c4762a1bSJed Brown   }
139*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
140c4762a1bSJed Brown }
141c4762a1bSJed Brown 
142d71ae5a4SJacob Faibussowitsch PetscErrorCode ini_bou(Vec X, AppCtx *user)
143d71ae5a4SJacob Faibussowitsch {
144c4762a1bSJed Brown   DM            cda;
145c4762a1bSJed Brown   DMDACoor2d  **coors;
146c4762a1bSJed Brown   PetscScalar **p;
147c4762a1bSJed Brown   Vec           gc;
148c4762a1bSJed Brown   PetscInt      i, j;
149c4762a1bSJed Brown   PetscInt      xs, ys, xm, ym, M, N;
150c4762a1bSJed Brown   PetscScalar   xi, yi;
151c4762a1bSJed Brown   PetscScalar   sigmax = user->sigmax, sigmay = user->sigmay;
152c4762a1bSJed Brown   PetscScalar   rho = user->rho;
153c4762a1bSJed Brown   PetscScalar   muy = user->muy, mux;
154c4762a1bSJed Brown   PetscMPIInt   rank;
155c4762a1bSJed Brown   PetscScalar   sum;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBeginUser;
1589566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
1599566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
1609371c9d4SSatish Balay   user->dx = (user->xmax - user->xmin) / (M - 1);
1619371c9d4SSatish Balay   user->dy = (user->ymax - user->ymin) / (N - 1);
1629566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(user->da, &cda));
1639566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinates(user->da, &gc));
1649566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(cda, gc, &coors));
1659566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->da, X, &p));
1669566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* mux and muy need to be grid points in the x and y-direction otherwise the solution goes unstable
169c4762a1bSJed Brown      muy is set by choosing the y domain, no. of grid points along y-direction so that muy is a grid point
170c4762a1bSJed Brown      in the y-direction. We only modify mux here
171c4762a1bSJed Brown   */
172c4762a1bSJed Brown   mux = user->mux = coors[0][M / 2 + 10].x; /* For -pi < x < pi, this should be some angle between 0 and pi/2 */
173c4762a1bSJed Brown   if (user->nonoiseinitial) {
174c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
175c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
1769371c9d4SSatish Balay         xi = coors[j][i].x;
1779371c9d4SSatish Balay         yi = coors[j][i].y;
178ad540459SPierre Jolivet         if ((xi == mux) && (yi == muy)) p[j][i] = 1.0;
179c4762a1bSJed Brown       }
180c4762a1bSJed Brown     }
181c4762a1bSJed Brown   } else {
182c4762a1bSJed Brown     /* Change PM_min accordingly */
183c4762a1bSJed Brown     user->PM_min = user->Pmax * PetscSinScalar(mux);
184c4762a1bSJed Brown     for (i = xs; i < xs + xm; i++) {
185c4762a1bSJed Brown       for (j = ys; j < ys + ym; j++) {
1869371c9d4SSatish Balay         xi = coors[j][i].x;
1879371c9d4SSatish Balay         yi = coors[j][i].y;
188c4762a1bSJed Brown         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)));
189c4762a1bSJed Brown       }
190c4762a1bSJed Brown     }
191c4762a1bSJed Brown   }
1929566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(cda, gc, &coors));
1939566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->da, X, &p));
1949566063dSJacob Faibussowitsch   PetscCall(VecSum(X, &sum));
1959566063dSJacob Faibussowitsch   PetscCall(VecScale(X, 1.0 / sum));
196*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
197c4762a1bSJed Brown }
198c4762a1bSJed Brown 
199c4762a1bSJed Brown /* First advection term */
200d71ae5a4SJacob Faibussowitsch PetscErrorCode adv1(PetscScalar **p, PetscScalar y, PetscInt i, PetscInt j, PetscInt M, PetscScalar *p1, AppCtx *user)
201d71ae5a4SJacob Faibussowitsch {
202c4762a1bSJed Brown   PetscScalar f, fpos, fneg;
203c4762a1bSJed Brown   PetscFunctionBegin;
204c4762a1bSJed Brown   f    = (y - user->ws);
205c4762a1bSJed Brown   fpos = PetscMax(f, 0);
206c4762a1bSJed Brown   fneg = PetscMin(f, 0);
207c4762a1bSJed Brown   if (user->st_width == 1) {
208c4762a1bSJed Brown     *p1 = fpos * (p[j][i] - p[j][i - 1]) / user->dx + fneg * (p[j][i + 1] - p[j][i]) / user->dx;
209c4762a1bSJed Brown   } else if (user->st_width == 2) {
210c4762a1bSJed 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);
211c4762a1bSJed Brown   } else if (user->st_width == 3) {
212c4762a1bSJed 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);
213c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils");
214*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
215c4762a1bSJed Brown }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown /* Second advection term */
218d71ae5a4SJacob Faibussowitsch PetscErrorCode adv2(PetscScalar **p, PetscScalar x, PetscInt i, PetscInt j, PetscInt N, PetscScalar *p2, AppCtx *user)
219d71ae5a4SJacob Faibussowitsch {
220c4762a1bSJed Brown   PetscScalar f, fpos, fneg;
221c4762a1bSJed Brown   PetscFunctionBegin;
222c4762a1bSJed Brown   f    = (user->ws / (2 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(x));
223c4762a1bSJed Brown   fpos = PetscMax(f, 0);
224c4762a1bSJed Brown   fneg = PetscMin(f, 0);
225c4762a1bSJed Brown   if (user->st_width == 1) {
226c4762a1bSJed Brown     *p2 = fpos * (p[j][i] - p[j - 1][i]) / user->dy + fneg * (p[j + 1][i] - p[j][i]) / user->dy;
227c4762a1bSJed Brown   } else if (user->st_width == 2) {
228c4762a1bSJed 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);
229c4762a1bSJed Brown   } else if (user->st_width == 3) {
230c4762a1bSJed 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);
231c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils");
232*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
233c4762a1bSJed Brown }
234c4762a1bSJed Brown 
235c4762a1bSJed Brown /* Diffusion term */
236d71ae5a4SJacob Faibussowitsch PetscErrorCode diffuse(PetscScalar **p, PetscInt i, PetscInt j, PetscReal t, PetscScalar *p_diff, AppCtx *user)
237d71ae5a4SJacob Faibussowitsch {
238c4762a1bSJed Brown   PetscFunctionBeginUser;
239c4762a1bSJed Brown   if (user->st_width == 1) {
240c4762a1bSJed Brown     *p_diff = user->disper_coe * ((p[j - 1][i] - 2 * p[j][i] + p[j + 1][i]) / (user->dy * user->dy));
241c4762a1bSJed Brown   } else if (user->st_width == 2) {
242c4762a1bSJed 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));
243c4762a1bSJed Brown   } else if (user->st_width == 3) {
244c4762a1bSJed 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));
245c4762a1bSJed Brown   } else SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "No support for wider stencils");
246*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
247c4762a1bSJed Brown }
248c4762a1bSJed Brown 
249d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
250d71ae5a4SJacob Faibussowitsch {
251c4762a1bSJed Brown   AppCtx       *user = (AppCtx *)ctx;
252c4762a1bSJed Brown   DM            cda;
253c4762a1bSJed Brown   DMDACoor2d  **coors;
254c4762a1bSJed Brown   PetscScalar **p, **f, **pdot;
255c4762a1bSJed Brown   PetscInt      i, j;
256c4762a1bSJed Brown   PetscInt      xs, ys, xm, ym, M, N;
257c4762a1bSJed Brown   Vec           localX, gc, localXdot;
258bbadc9eeSBarry Smith   PetscScalar   p_adv1 = 0.0, p_adv2 = 0.0, p_diff; /* initialize to prevent incorrect compiler warnings */
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   PetscFunctionBeginUser;
2619566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
2629566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(user->da, &cda));
2639566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
264c4762a1bSJed Brown 
2659566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->da, &localX));
2669566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(user->da, &localXdot));
267c4762a1bSJed Brown 
2689566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->da, X, INSERT_VALUES, localX));
2699566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->da, X, INSERT_VALUES, localX));
2709566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(user->da, Xdot, INSERT_VALUES, localXdot));
2719566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(user->da, Xdot, INSERT_VALUES, localXdot));
272c4762a1bSJed Brown 
2739566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(user->da, &gc));
274c4762a1bSJed Brown 
2759566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
2769566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(user->da, localX, &p));
2779566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(user->da, localXdot, &pdot));
2789566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArray(user->da, F, &f));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   user->disper_coe = PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda));
281c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
282c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
2839566063dSJacob Faibussowitsch       PetscCall(adv1(p, coors[j][i].y, i, j, M, &p_adv1, user));
2849566063dSJacob Faibussowitsch       PetscCall(adv2(p, coors[j][i].x, i, j, N, &p_adv2, user));
2859566063dSJacob Faibussowitsch       PetscCall(diffuse(p, i, j, t, &p_diff, user));
286c4762a1bSJed Brown       f[j][i] = -p_adv1 - p_adv2 + p_diff - pdot[j][i];
287c4762a1bSJed Brown     }
288c4762a1bSJed Brown   }
2899566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &p));
2909566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(user->da, localX, &pdot));
2919566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->da, &localX));
2929566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(user->da, &localXdot));
2939566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArray(user->da, F, &f));
2949566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
295c4762a1bSJed Brown 
296*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
297c4762a1bSJed Brown }
298c4762a1bSJed Brown 
299d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat J, Mat Jpre, void *ctx)
300d71ae5a4SJacob Faibussowitsch {
301c4762a1bSJed Brown   AppCtx      *user = (AppCtx *)ctx;
302c4762a1bSJed Brown   DM           cda;
303c4762a1bSJed Brown   DMDACoor2d **coors;
304c4762a1bSJed Brown   PetscInt     i, j;
305c4762a1bSJed Brown   PetscInt     xs, ys, xm, ym, M, N;
306c4762a1bSJed Brown   Vec          gc;
307c4762a1bSJed Brown   PetscScalar  val[5], xi, yi;
308c4762a1bSJed Brown   MatStencil   row, col[5];
309c4762a1bSJed Brown   PetscScalar  c1, c3, c5, c1pos, c1neg, c3pos, c3neg;
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   PetscFunctionBeginUser;
3129566063dSJacob Faibussowitsch   PetscCall(DMDAGetInfo(user->da, NULL, &M, &N, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL, NULL));
3139566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinateDM(user->da, &cda));
3149566063dSJacob Faibussowitsch   PetscCall(DMDAGetCorners(cda, &xs, &ys, 0, &xm, &ym, 0));
315c4762a1bSJed Brown 
3169566063dSJacob Faibussowitsch   PetscCall(DMGetCoordinatesLocal(user->da, &gc));
3179566063dSJacob Faibussowitsch   PetscCall(DMDAVecGetArrayRead(cda, gc, &coors));
318c4762a1bSJed Brown   for (i = xs; i < xs + xm; i++) {
319c4762a1bSJed Brown     for (j = ys; j < ys + ym; j++) {
320c4762a1bSJed Brown       PetscInt nc = 0;
3219371c9d4SSatish Balay       xi          = coors[j][i].x;
3229371c9d4SSatish Balay       yi          = coors[j][i].y;
3239371c9d4SSatish Balay       row.i       = i;
3249371c9d4SSatish Balay       row.j       = j;
325c4762a1bSJed Brown       c1          = (yi - user->ws) / user->dx;
326c4762a1bSJed Brown       c1pos       = PetscMax(c1, 0);
327c4762a1bSJed Brown       c1neg       = PetscMin(c1, 0);
328c4762a1bSJed Brown       c3          = (user->ws / (2.0 * user->H)) * (user->PM_min - user->Pmax * PetscSinScalar(xi)) / user->dy;
329c4762a1bSJed Brown       c3pos       = PetscMax(c3, 0);
330c4762a1bSJed Brown       c3neg       = PetscMin(c3, 0);
331c4762a1bSJed Brown       c5          = (PetscPowScalar((user->lambda * user->ws) / (2 * user->H), 2) * user->q * (1.0 - PetscExpScalar(-t / user->lambda))) / (user->dy * user->dy);
3329371c9d4SSatish Balay       col[nc].i   = i - 1;
3339371c9d4SSatish Balay       col[nc].j   = j;
3349371c9d4SSatish Balay       val[nc++]   = c1pos;
3359371c9d4SSatish Balay       col[nc].i   = i + 1;
3369371c9d4SSatish Balay       col[nc].j   = j;
3379371c9d4SSatish Balay       val[nc++]   = -c1neg;
3389371c9d4SSatish Balay       col[nc].i   = i;
3399371c9d4SSatish Balay       col[nc].j   = j - 1;
3409371c9d4SSatish Balay       val[nc++]   = c3pos + c5;
3419371c9d4SSatish Balay       col[nc].i   = i;
3429371c9d4SSatish Balay       col[nc].j   = j + 1;
3439371c9d4SSatish Balay       val[nc++]   = -c3neg + c5;
3449371c9d4SSatish Balay       col[nc].i   = i;
3459371c9d4SSatish Balay       col[nc].j   = j;
3469371c9d4SSatish Balay       val[nc++]   = -c1pos + c1neg - c3pos + c3neg - 2 * c5 - a;
3479566063dSJacob Faibussowitsch       PetscCall(MatSetValuesStencil(Jpre, 1, &row, nc, col, val, INSERT_VALUES));
348c4762a1bSJed Brown     }
349c4762a1bSJed Brown   }
3509566063dSJacob Faibussowitsch   PetscCall(DMDAVecRestoreArrayRead(cda, gc, &coors));
351c4762a1bSJed Brown 
3529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(Jpre, MAT_FINAL_ASSEMBLY));
3539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(Jpre, MAT_FINAL_ASSEMBLY));
354c4762a1bSJed Brown   if (J != Jpre) {
3559566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
3569566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
357c4762a1bSJed Brown   }
358*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
359c4762a1bSJed Brown }
360c4762a1bSJed Brown 
361d71ae5a4SJacob Faibussowitsch PetscErrorCode Parameter_settings(AppCtx *user)
362d71ae5a4SJacob Faibussowitsch {
363c4762a1bSJed Brown   PetscBool flg;
364c4762a1bSJed Brown 
365c4762a1bSJed Brown   PetscFunctionBeginUser;
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   /* Set default parameters */
368c4762a1bSJed Brown   user->ws             = 1.0;
369c4762a1bSJed Brown   user->H              = 5.0;
370c4762a1bSJed Brown   user->Pmax           = 2.1;
371c4762a1bSJed Brown   user->PM_min         = 1.0;
372c4762a1bSJed Brown   user->lambda         = 0.1;
373c4762a1bSJed Brown   user->q              = 1.0;
374c4762a1bSJed Brown   user->mux            = PetscAsinScalar(user->PM_min / user->Pmax);
375c4762a1bSJed Brown   user->sigmax         = 0.1;
376c4762a1bSJed Brown   user->sigmay         = 0.1;
377c4762a1bSJed Brown   user->rho            = 0.0;
378c4762a1bSJed Brown   user->xmin           = -PETSC_PI;
379c4762a1bSJed Brown   user->xmax           = PETSC_PI;
380c4762a1bSJed Brown   user->bx             = DM_BOUNDARY_PERIODIC;
381c4762a1bSJed Brown   user->by             = DM_BOUNDARY_MIRROR;
382c4762a1bSJed Brown   user->nonoiseinitial = PETSC_FALSE;
383c4762a1bSJed Brown 
384c4762a1bSJed Brown   /*
385c4762a1bSJed Brown      ymin of -3 seems to let the unstable solution move up and leave a zero in its wake
386c4762a1bSJed Brown      with an ymin of -1 the wake is never exactly zero
387c4762a1bSJed Brown   */
388c4762a1bSJed Brown   user->ymin     = -3.0;
389c4762a1bSJed Brown   user->ymax     = 10.0;
390c4762a1bSJed Brown   user->st_width = 1;
391c4762a1bSJed Brown 
3929566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ws", &user->ws, &flg));
3939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Inertia", &user->H, &flg));
3949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-Pmax", &user->Pmax, &flg));
3959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-PM_min", &user->PM_min, &flg));
3969566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-lambda", &user->lambda, &flg));
3979566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-q", &user->q, &flg));
3989566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-mux", &user->mux, &flg));
3999566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmax", &user->sigmax, &flg));
4009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-muy", &user->muy, &flg));
401ad540459SPierre Jolivet   if (flg == 0) user->muy = user->ws;
4029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-sigmay", &user->sigmay, &flg));
4039566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-rho", &user->rho, &flg));
4049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmin", &user->xmin, &flg));
4059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-xmax", &user->xmax, &flg));
4069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymin", &user->ymin, &flg));
4079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-ymax", &user->ymax, &flg));
4089566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-stencil_width", &user->st_width, &flg));
4099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-bx", DMBoundaryTypes, (PetscEnum *)&user->bx, &flg));
4109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetEnum(NULL, NULL, "-by", DMBoundaryTypes, (PetscEnum *)&user->by, &flg));
4119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-nonoiseinitial", &user->nonoiseinitial, &flg));
412c4762a1bSJed Brown 
413*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
414c4762a1bSJed Brown }
415c4762a1bSJed Brown 
416c4762a1bSJed Brown /*TEST
417c4762a1bSJed Brown 
418c4762a1bSJed Brown    build:
419c4762a1bSJed Brown       requires: !complex !single
420c4762a1bSJed Brown 
421c4762a1bSJed Brown    test:
422c4762a1bSJed Brown       args: -ts_max_steps 2
423c4762a1bSJed Brown       localrunfiles: petscopt_ex7
424c4762a1bSJed Brown 
425c4762a1bSJed Brown    test:
426c4762a1bSJed Brown       suffix: 2
427c4762a1bSJed Brown       args: -ts_max_steps 2 -snes_mf_operator
428c4762a1bSJed Brown       output_file: output/ex7_1.out
429c4762a1bSJed Brown       localrunfiles: petscopt_ex7
430c4762a1bSJed Brown       timeoutfactor: 2
431c4762a1bSJed Brown 
432c4762a1bSJed Brown    test:
433c4762a1bSJed Brown       suffix: 3
434c4762a1bSJed Brown       args: -ts_max_steps 2 -snes_mf -pc_type none
435c4762a1bSJed Brown       output_file: output/ex7_1.out
436c4762a1bSJed Brown       localrunfiles: petscopt_ex7
437c4762a1bSJed Brown       timeoutfactor: 2
438c4762a1bSJed Brown 
439c4762a1bSJed Brown TEST*/
440