xref: /petsc/src/ts/tutorials/power_grid/ex9opt.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*F
5c4762a1bSJed Brown 
6c4762a1bSJed Brown \begin{eqnarray}
7c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8c4762a1bSJed Brown                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9c4762a1bSJed Brown \end{eqnarray}
10c4762a1bSJed Brown 
11c4762a1bSJed Brown   Ensemble of initial conditions
12c4762a1bSJed Brown    ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3      -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   Fault at .1 seconds
15c4762a1bSJed Brown    ./ex2           -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   Initial conditions same as when fault is ended
18c4762a1bSJed Brown    ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05  -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
19c4762a1bSJed Brown 
20c4762a1bSJed Brown F*/
21c4762a1bSJed Brown 
22c4762a1bSJed Brown /*
23c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
24c4762a1bSJed Brown    file automatically includes:
25c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
26c4762a1bSJed Brown      petscmat.h - matrices
27c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
28c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
29c4762a1bSJed Brown      petscksp.h   - linear solvers
30c4762a1bSJed Brown */
31c4762a1bSJed Brown 
32c4762a1bSJed Brown #include <petsctao.h>
33c4762a1bSJed Brown #include <petscts.h>
34c4762a1bSJed Brown 
35c4762a1bSJed Brown typedef struct {
36c4762a1bSJed Brown   TS          ts;
37c4762a1bSJed Brown   PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c;
38c4762a1bSJed Brown   PetscInt    beta;
39c4762a1bSJed Brown   PetscReal   tf, tcl, dt;
40c4762a1bSJed Brown } AppCtx;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
43c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *);
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /*
46c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
47c4762a1bSJed Brown */
48*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) {
49c4762a1bSJed Brown   PetscScalar       *f, Pmax;
50c4762a1bSJed Brown   const PetscScalar *u;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   PetscFunctionBegin;
53c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
56c4762a1bSJed Brown   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
57c4762a1bSJed Brown   else Pmax = ctx->Pmax;
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   f[0] = ctx->omega_b * (u[1] - ctx->omega_s);
60c4762a1bSJed Brown   f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H);
61c4762a1bSJed Brown 
629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
64c4762a1bSJed Brown   PetscFunctionReturn(0);
65c4762a1bSJed Brown }
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /*
68c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
69c4762a1bSJed Brown */
70*9371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) {
71c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
72c4762a1bSJed Brown   PetscScalar        J[2][2], Pmax;
73c4762a1bSJed Brown   const PetscScalar *u;
74c4762a1bSJed Brown 
75c4762a1bSJed Brown   PetscFunctionBegin;
769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
77c4762a1bSJed Brown   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
78c4762a1bSJed Brown   else Pmax = ctx->Pmax;
79c4762a1bSJed Brown 
80*9371c9d4SSatish Balay   J[0][0] = 0;
81*9371c9d4SSatish Balay   J[0][1] = ctx->omega_b;
82*9371c9d4SSatish Balay   J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H);
83*9371c9d4SSatish Balay   J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H);
84c4762a1bSJed Brown 
859566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
90c4762a1bSJed Brown   if (A != B) {
919566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown   }
94c4762a1bSJed Brown   PetscFunctionReturn(0);
95c4762a1bSJed Brown }
96c4762a1bSJed Brown 
97*9371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0) {
98c4762a1bSJed Brown   PetscInt    row[] = {0, 1}, col[] = {0};
99c4762a1bSJed Brown   PetscScalar J[2][1];
100c4762a1bSJed Brown   AppCtx     *ctx = (AppCtx *)ctx0;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBeginUser;
103c4762a1bSJed Brown   J[0][0] = 0;
104c4762a1bSJed Brown   J[1][0] = ctx->omega_s / (2.0 * ctx->H);
1059566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
1069566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
108c4762a1bSJed Brown   PetscFunctionReturn(0);
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111*9371c9d4SSatish Balay static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) {
112c4762a1bSJed Brown   PetscScalar       *r;
113c4762a1bSJed Brown   const PetscScalar *u;
114c4762a1bSJed Brown 
115c4762a1bSJed Brown   PetscFunctionBegin;
1169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R, &r));
1182f613bf5SBarry Smith   r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta);
1199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R, &r));
1209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
121c4762a1bSJed Brown   PetscFunctionReturn(0);
122c4762a1bSJed Brown }
123c4762a1bSJed Brown 
124*9371c9d4SSatish Balay static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) {
125c4762a1bSJed Brown   PetscScalar        ru[1];
126c4762a1bSJed Brown   const PetscScalar *u;
127c4762a1bSJed Brown   PetscInt           row[] = {0}, col[] = {0};
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   PetscFunctionBegin;
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1312f613bf5SBarry Smith   ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1);
1329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
1339566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES));
1349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY));
1359566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY));
136c4762a1bSJed Brown   PetscFunctionReturn(0);
137c4762a1bSJed Brown }
138c4762a1bSJed Brown 
139*9371c9d4SSatish Balay static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx) {
140c4762a1bSJed Brown   PetscFunctionBegin;
1419566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
1429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY));
1439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY));
144c4762a1bSJed Brown   PetscFunctionReturn(0);
145c4762a1bSJed Brown }
146c4762a1bSJed Brown 
147*9371c9d4SSatish Balay PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) {
148c4762a1bSJed Brown   PetscScalar       *y, sensip;
149c4762a1bSJed Brown   const PetscScalar *x;
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   PetscFunctionBegin;
1529566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(lambda, &x));
1539566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu, &y));
154c4762a1bSJed Brown   sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0];
155c4762a1bSJed Brown   y[0]   = sensip;
1569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu, &y));
1579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(lambda, &x));
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161*9371c9d4SSatish Balay int main(int argc, char **argv) {
162c4762a1bSJed Brown   Vec          p;
163c4762a1bSJed Brown   PetscScalar *x_ptr;
164c4762a1bSJed Brown   PetscMPIInt  size;
165c4762a1bSJed Brown   AppCtx       ctx;
166c4762a1bSJed Brown   Vec          lowerb, upperb;
167c4762a1bSJed Brown   Tao          tao;
168c4762a1bSJed Brown   KSP          ksp;
169c4762a1bSJed Brown   PC           pc;
170c4762a1bSJed Brown   Vec          U, lambda[1], mu[1];
171c4762a1bSJed Brown   Mat          A;    /* Jacobian matrix */
172c4762a1bSJed Brown   Mat          Jacp; /* Jacobian matrix */
173c4762a1bSJed Brown   Mat          DRDU, DRDP;
174c4762a1bSJed Brown   PetscInt     n = 2;
175c4762a1bSJed Brown   TS           quadts;
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown      Initialize program
179c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180327415f7SBarry Smith   PetscFunctionBeginUser;
1819566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
182c4762a1bSJed Brown   PetscFunctionBeginUser;
1839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1843c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown     Set runtime options
188c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
190c4762a1bSJed Brown   {
191c4762a1bSJed Brown     ctx.beta    = 2;
192c4762a1bSJed Brown     ctx.c       = PetscRealConstant(10000.0);
193c4762a1bSJed Brown     ctx.u_s     = PetscRealConstant(1.0);
194c4762a1bSJed Brown     ctx.omega_s = PetscRealConstant(1.0);
195c4762a1bSJed Brown     ctx.omega_b = PetscRealConstant(120.0) * PETSC_PI;
196c4762a1bSJed Brown     ctx.H       = PetscRealConstant(5.0);
1979566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
198c4762a1bSJed Brown     ctx.D = PetscRealConstant(5.0);
1999566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
200c4762a1bSJed Brown     ctx.E    = PetscRealConstant(1.1378);
201c4762a1bSJed Brown     ctx.V    = PetscRealConstant(1.0);
202c4762a1bSJed Brown     ctx.X    = PetscRealConstant(0.545);
203c4762a1bSJed Brown     ctx.Pmax = ctx.E * ctx.V / ctx.X;
2049566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
205c4762a1bSJed Brown     ctx.Pm = PetscRealConstant(1.0194);
2069566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
207c4762a1bSJed Brown     ctx.tf  = PetscRealConstant(0.1);
208c4762a1bSJed Brown     ctx.tcl = PetscRealConstant(0.2);
2099566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
2109566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
211c4762a1bSJed Brown   }
212d0609cedSBarry Smith   PetscOptionsEnd();
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215c4762a1bSJed Brown     Create necessary matrix and vectors
216c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2179566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2189566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
2199566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
2209566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
222c4762a1bSJed Brown 
2239566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
224c4762a1bSJed Brown 
2259566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
2269566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1));
2279566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jacp));
2289566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jacp));
2299566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP));
2309566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDP));
2319566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU));
2329566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDU));
233c4762a1bSJed Brown 
234c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235c4762a1bSJed Brown      Create timestepping solver context
236c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2379566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ctx.ts));
2389566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ctx.ts, TS_NONLINEAR));
2399566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ctx.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
2409566063dSJacob Faibussowitsch   PetscCall(TSSetType(ctx.ts, TSRK));
2419566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ctx.ts, NULL, (TSRHSFunction)RHSFunction, &ctx));
2429566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ctx.ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx));
2439566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ctx.ts, TS_EXACTFINALTIME_MATCHSTEP));
244c4762a1bSJed Brown 
2459566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &lambda[0], NULL));
2469566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp, &mu[0], NULL));
2479566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ctx.ts, 1, lambda, mu));
2489566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ctx.ts, Jacp, RHSJacobianP, &ctx));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
251c4762a1bSJed Brown      Set solver options
252c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2539566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ctx.ts, PetscRealConstant(1.0)));
2549566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ctx.ts, PetscRealConstant(0.01)));
2559566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ctx.ts));
256c4762a1bSJed Brown 
2579566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ctx.ts, &ctx.dt)); /* save the stepsize */
258c4762a1bSJed Brown 
2599566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(ctx.ts, PETSC_TRUE, &quadts));
2609566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, &ctx));
2619566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobian)DRDUJacobianTranspose, &ctx));
2629566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianP)DRDPJacobianTranspose, &ctx));
2639566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ctx.ts, U));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
2669566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
2679566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBLMVM));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /*
270c4762a1bSJed Brown      Optimization starts
271c4762a1bSJed Brown   */
272c4762a1bSJed Brown   /* Set initial solution guess */
2739566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p));
2749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(p, &x_ptr));
275c4762a1bSJed Brown   x_ptr[0] = ctx.Pm;
2769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(p, &x_ptr));
277c4762a1bSJed Brown 
2789566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, p));
279c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
2809566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx));
2819566063dSJacob Faibussowitsch   PetscCall(TaoSetGradient(tao, NULL, FormGradient, (void *)&ctx));
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /* Set bounds for the optimization */
2849566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p, &lowerb));
2859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p, &upperb));
2869566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lowerb, &x_ptr));
287c4762a1bSJed Brown   x_ptr[0] = 0.;
2889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lowerb, &x_ptr));
2899566063dSJacob Faibussowitsch   PetscCall(VecGetArray(upperb, &x_ptr));
290c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(1.1);
2919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(upperb, &x_ptr));
2929566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   /* Check for any TAO command line options */
2959566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
2969566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
297c4762a1bSJed Brown   if (ksp) {
2989566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
2999566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCNONE));
300c4762a1bSJed Brown   }
301c4762a1bSJed Brown 
302c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
3039566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
304c4762a1bSJed Brown 
3059566063dSJacob Faibussowitsch   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
306c4762a1bSJed Brown   /* Free TAO data structures */
3079566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
3089566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&p));
3099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lowerb));
3109566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&upperb));
311c4762a1bSJed Brown 
3129566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ctx.ts));
3139566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3149566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3159566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jacp));
3169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDU));
3179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDP));
3189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3209566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
321b122ec5aSJacob Faibussowitsch   return 0;
322c4762a1bSJed Brown }
323c4762a1bSJed Brown 
324c4762a1bSJed Brown /* ------------------------------------------------------------------ */
325c4762a1bSJed Brown /*
326c4762a1bSJed Brown    FormFunction - Evaluates the function
327c4762a1bSJed Brown 
328c4762a1bSJed Brown    Input Parameters:
329c4762a1bSJed Brown    tao - the Tao context
330c4762a1bSJed Brown    X   - the input vector
331a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
332c4762a1bSJed Brown 
333c4762a1bSJed Brown    Output Parameters:
334c4762a1bSJed Brown    f   - the newly evaluated function
335c4762a1bSJed Brown */
336*9371c9d4SSatish Balay PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0) {
337c4762a1bSJed Brown   AppCtx      *ctx = (AppCtx *)ctx0;
338c4762a1bSJed Brown   TS           ts  = ctx->ts;
339c4762a1bSJed Brown   Vec          U; /* solution will be stored here */
340c4762a1bSJed Brown   PetscScalar *u;
341c4762a1bSJed Brown   PetscScalar *x_ptr;
342c4762a1bSJed Brown   Vec          q;
343c4762a1bSJed Brown 
3449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
345c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
3469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   /* reset time */
3499566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
350c4762a1bSJed Brown   /* reset step counter, this is critical for adjoint solver */
3519566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
352c4762a1bSJed Brown   /* reset step size, the step size becomes negative after TSAdjointSolve */
3539566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, ctx->dt));
354c4762a1bSJed Brown   /* reinitialize the integral value */
3559566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &q));
3569566063dSJacob Faibussowitsch   PetscCall(VecSet(q, 0.0));
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
359c4762a1bSJed Brown      Set initial conditions
360c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3619566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
3629566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
363c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
364c4762a1bSJed Brown   u[1] = PetscRealConstant(1.0);
3659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
366c4762a1bSJed Brown 
367c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
368c4762a1bSJed Brown      Solve nonlinear system
369c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3709566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
3719566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &q));
3729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(q, &x_ptr));
373c4762a1bSJed Brown   *f = -ctx->Pm + x_ptr[0];
3749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(q, &x_ptr));
375c4762a1bSJed Brown   return 0;
376c4762a1bSJed Brown }
377c4762a1bSJed Brown 
378*9371c9d4SSatish Balay PetscErrorCode FormGradient(Tao tao, Vec P, Vec G, void *ctx0) {
379c4762a1bSJed Brown   AppCtx      *ctx = (AppCtx *)ctx0;
380c4762a1bSJed Brown   TS           ts  = ctx->ts;
381c4762a1bSJed Brown   Vec          U; /* solution will be stored here */
382c4762a1bSJed Brown   PetscReal    ftime;
383c4762a1bSJed Brown   PetscInt     steps;
384c4762a1bSJed Brown   PetscScalar *u;
385c4762a1bSJed Brown   PetscScalar *x_ptr, *y_ptr;
386c4762a1bSJed Brown   Vec         *lambda, q, *mu;
387c4762a1bSJed Brown 
3889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr));
389c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
3909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr));
391c4762a1bSJed Brown 
392c4762a1bSJed Brown   /* reset time */
3939566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, 0.0));
394c4762a1bSJed Brown   /* reset step counter, this is critical for adjoint solver */
3959566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ts, 0));
396c4762a1bSJed Brown   /* reset step size, the step size becomes negative after TSAdjointSolve */
3979566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, ctx->dt));
398c4762a1bSJed Brown   /* reinitialize the integral value */
3999566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &q));
4009566063dSJacob Faibussowitsch   PetscCall(VecSet(q, 0.0));
401c4762a1bSJed Brown 
402c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
403c4762a1bSJed Brown      Set initial conditions
404c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4059566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts, &U));
4069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
407c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
408c4762a1bSJed Brown   u[1] = PetscRealConstant(1.0);
4099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
410c4762a1bSJed Brown 
411f32d6360SSatish Balay   /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */
4129566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
4139566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
416c4762a1bSJed Brown      Solve nonlinear system
417c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4189566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
419c4762a1bSJed Brown 
4209566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
4219566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
422c4762a1bSJed Brown 
423c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
424c4762a1bSJed Brown      Adjoint model starts here
425c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4269566063dSJacob Faibussowitsch   PetscCall(TSGetCostGradients(ts, NULL, &lambda, &mu));
427c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
4289566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0], &y_ptr));
429*9371c9d4SSatish Balay   y_ptr[0] = 0.0;
430*9371c9d4SSatish Balay   y_ptr[1] = 0.0;
4319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0], &y_ptr));
4329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0], &x_ptr));
433c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(-1.0);
4349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0], &x_ptr));
435c4762a1bSJed Brown 
4369566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
4379566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts, &q));
4389566063dSJacob Faibussowitsch   PetscCall(ComputeSensiP(lambda[0], mu[0], ctx));
4399566063dSJacob Faibussowitsch   PetscCall(VecCopy(mu[0], G));
440c4762a1bSJed Brown   return 0;
441c4762a1bSJed Brown }
442c4762a1bSJed Brown 
443c4762a1bSJed Brown /*TEST
444c4762a1bSJed Brown 
445c4762a1bSJed Brown    build:
446c4762a1bSJed Brown       requires: !complex
447c4762a1bSJed Brown 
448c4762a1bSJed Brown    test:
449c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason
450c4762a1bSJed Brown 
451c4762a1bSJed Brown    test:
452c4762a1bSJed Brown       suffix: 2
453c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason -tao_test_gradient
454c4762a1bSJed Brown 
455c4762a1bSJed Brown TEST*/
456