xref: /petsc/src/ts/tutorials/power_grid/ex9adj.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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    ./ex9 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   Fault at .1 seconds
15c4762a1bSJed Brown    ./ex9 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
16c4762a1bSJed Brown 
17c4762a1bSJed Brown   Initial conditions same as when fault is ended
18c4762a1bSJed Brown    ./ex9 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -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 <petscts.h>
33c4762a1bSJed Brown 
34c4762a1bSJed Brown typedef struct {
35c4762a1bSJed Brown   PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c;
36c4762a1bSJed Brown   PetscInt    beta;
37c4762a1bSJed Brown   PetscReal   tf,tcl;
38c4762a1bSJed Brown } AppCtx;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown PetscErrorCode PostStepFunction(TS ts)
41c4762a1bSJed Brown {
42c4762a1bSJed Brown   Vec               U;
43c4762a1bSJed Brown   PetscReal         t;
44c4762a1bSJed Brown   const PetscScalar *u;
45c4762a1bSJed Brown 
46c4762a1bSJed Brown   PetscFunctionBegin;
479566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts,&t));
489566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ts,&U));
499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
509566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",(double)t,(double)u[0]));
519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
52c4762a1bSJed Brown   PetscFunctionReturn(0);
53c4762a1bSJed Brown }
54c4762a1bSJed Brown 
55c4762a1bSJed Brown /*
56c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
57c4762a1bSJed Brown */
58c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
59c4762a1bSJed Brown {
60c4762a1bSJed Brown   PetscScalar       *f,Pmax;
61c4762a1bSJed Brown   const PetscScalar *u;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   PetscFunctionBegin;
64c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
659566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
669566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
67c4762a1bSJed 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 */
68c4762a1bSJed Brown   else Pmax = ctx->Pmax;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
71c4762a1bSJed Brown   f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);
72c4762a1bSJed Brown 
739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
75c4762a1bSJed Brown   PetscFunctionReturn(0);
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
78c4762a1bSJed Brown /*
79c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
80c4762a1bSJed Brown */
81c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
82c4762a1bSJed Brown {
83c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
84c4762a1bSJed Brown   PetscScalar       J[2][2],Pmax;
85c4762a1bSJed Brown   const PetscScalar *u;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   PetscFunctionBegin;
889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
89c4762a1bSJed 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 */
90c4762a1bSJed Brown   else Pmax = ctx->Pmax;
91c4762a1bSJed Brown 
92c4762a1bSJed Brown   J[0][0] = 0;                                  J[0][1] = ctx->omega_b;
93c4762a1bSJed Brown   J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H);  J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H);
94c4762a1bSJed Brown 
959566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
97c4762a1bSJed Brown 
989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
999566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
100c4762a1bSJed Brown   if (A != B) {
1019566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1029566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
103c4762a1bSJed Brown   }
104c4762a1bSJed Brown   PetscFunctionReturn(0);
105c4762a1bSJed Brown }
106c4762a1bSJed Brown 
107c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
108c4762a1bSJed Brown {
109c4762a1bSJed Brown   PetscInt       row[] = {0,1},col[]={0};
110c4762a1bSJed Brown   PetscScalar    J[2][1];
111c4762a1bSJed Brown   AppCtx         *ctx=(AppCtx*)ctx0;
112c4762a1bSJed Brown 
113c4762a1bSJed Brown   PetscFunctionBeginUser;
114c4762a1bSJed Brown   J[0][0] = 0;
115c4762a1bSJed Brown   J[1][0] = ctx->omega_s/(2.0*ctx->H);
1169566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
1179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1189566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
119c4762a1bSJed Brown   PetscFunctionReturn(0);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
123c4762a1bSJed Brown {
124c4762a1bSJed Brown   PetscScalar       *r;
125c4762a1bSJed Brown   const PetscScalar *u;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBegin;
1289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(R,&r));
1302f613bf5SBarry Smith   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(R,&r));
1329566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
133c4762a1bSJed Brown   PetscFunctionReturn(0);
134c4762a1bSJed Brown }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx)
137c4762a1bSJed Brown {
138c4762a1bSJed Brown   PetscScalar       ru[1];
139c4762a1bSJed Brown   const PetscScalar *u;
140c4762a1bSJed Brown   PetscInt          row[] = {0},col[] = {0};
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBegin;
1439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1442f613bf5SBarry Smith   ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
1459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
1469566063dSJacob Faibussowitsch   PetscCall(MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES));
1479566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY));
1489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY));
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx)
153c4762a1bSJed Brown {
154c4762a1bSJed Brown   PetscFunctionBegin;
1559566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(DRDP));
1569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY));
1579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY));
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
162c4762a1bSJed Brown {
163c4762a1bSJed Brown   PetscScalar       sensip;
164c4762a1bSJed Brown   const PetscScalar *x,*y;
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   PetscFunctionBegin;
1679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(lambda,&x));
1689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(mu,&y));
169c4762a1bSJed Brown   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
1709566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip));
1719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(lambda,&x));
1729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(mu,&y));
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown int main(int argc,char **argv)
177c4762a1bSJed Brown {
178c4762a1bSJed Brown   TS             ts,quadts;     /* ODE integrator */
179c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
180c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
181c4762a1bSJed Brown   Mat            Jacp;          /* Jacobian matrix */
182c4762a1bSJed Brown   Mat            DRDU,DRDP;
183c4762a1bSJed Brown   PetscMPIInt    size;
184c4762a1bSJed Brown   PetscInt       n = 2;
185c4762a1bSJed Brown   AppCtx         ctx;
186c4762a1bSJed Brown   PetscScalar    *u;
187c4762a1bSJed Brown   PetscReal      du[2] = {0.0,0.0};
188c4762a1bSJed Brown   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
189c4762a1bSJed Brown   PetscReal      ftime;
190c4762a1bSJed Brown   PetscInt       steps;
191c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr;
192c4762a1bSJed Brown   Vec            lambda[1],q,mu[1];
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195c4762a1bSJed Brown      Initialize program
196c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197*327415f7SBarry Smith   PetscFunctionBeginUser;
1989566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
1999566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
2003c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203c4762a1bSJed Brown     Create necessary matrix and vectors
204c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2059566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
2069566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
2079566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
2089566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2099566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
210c4762a1bSJed Brown 
2119566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&U,NULL));
212c4762a1bSJed Brown 
2139566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&Jacp));
2149566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
2159566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(Jacp));
2169566063dSJacob Faibussowitsch   PetscCall(MatSetUp(Jacp));
217c4762a1bSJed Brown 
2189566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP));
2199566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDP));
2209566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU));
2219566063dSJacob Faibussowitsch   PetscCall(MatSetUp(DRDU));
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224c4762a1bSJed Brown     Set runtime options
225c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
227c4762a1bSJed Brown   {
228c4762a1bSJed Brown     ctx.beta    = 2;
229c4762a1bSJed Brown     ctx.c       = 10000.0;
230c4762a1bSJed Brown     ctx.u_s     = 1.0;
231c4762a1bSJed Brown     ctx.omega_s = 1.0;
232c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
233c4762a1bSJed Brown     ctx.H       = 5.0;
2349566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL));
235c4762a1bSJed Brown     ctx.D       = 5.0;
2369566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL));
237c4762a1bSJed Brown     ctx.E       = 1.1378;
238c4762a1bSJed Brown     ctx.V       = 1.0;
239c4762a1bSJed Brown     ctx.X       = 0.545;
240c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
2419566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL));
242c4762a1bSJed Brown     ctx.Pm      = 1.1;
2439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL));
244c4762a1bSJed Brown     ctx.tf      = 0.1;
245c4762a1bSJed Brown     ctx.tcl     = 0.2;
2469566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL));
2479566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL));
2489566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL));
249c4762a1bSJed Brown     if (ensemble) {
250c4762a1bSJed Brown       ctx.tf      = -1;
251c4762a1bSJed Brown       ctx.tcl     = -1;
252c4762a1bSJed Brown     }
253c4762a1bSJed Brown 
2549566063dSJacob Faibussowitsch     PetscCall(VecGetArray(U,&u));
255c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
256c4762a1bSJed Brown     u[1] = 1.0;
2579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1));
258c4762a1bSJed Brown     n    = 2;
2599566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2));
260c4762a1bSJed Brown     u[0] += du[0];
261c4762a1bSJed Brown     u[1] += du[1];
2629566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(U,&u));
263c4762a1bSJed Brown     if (flg1 || flg2) {
264c4762a1bSJed Brown       ctx.tf      = -1;
265c4762a1bSJed Brown       ctx.tcl     = -1;
266c4762a1bSJed Brown     }
267c4762a1bSJed Brown   }
268d0609cedSBarry Smith   PetscOptionsEnd();
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271c4762a1bSJed Brown      Create timestepping solver context
272c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2739566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
2749566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
2759566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
2769566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSRK));
2779566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx));
2789566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx));
2799566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts));
2809566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx));
2819566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx));
2829566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx));
283c4762a1bSJed Brown 
284c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
285c4762a1bSJed Brown      Set initial conditions
286c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2879566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,U));
288c4762a1bSJed Brown 
289c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
290c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
291c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2929566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
293c4762a1bSJed Brown 
2949566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&lambda[0],NULL));
295c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
2969566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0],&y_ptr));
297c4762a1bSJed Brown   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
2989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0],&y_ptr));
299c4762a1bSJed Brown 
3009566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(Jacp,&mu[0],NULL));
3019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0],&x_ptr));
302c4762a1bSJed Brown   x_ptr[0] = -1.0;
3039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0],&x_ptr));
3049566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,lambda,mu));
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307c4762a1bSJed Brown      Set solver options
308c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3099566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,10.0));
3109566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
3119566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.01));
3129566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
313c4762a1bSJed Brown 
314c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
315c4762a1bSJed Brown      Solve nonlinear system
316c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
317c4762a1bSJed Brown   if (ensemble) {
318c4762a1bSJed Brown     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
3199566063dSJacob Faibussowitsch       PetscCall(VecGetArray(U,&u));
320c4762a1bSJed Brown       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
321c4762a1bSJed Brown       u[1] = ctx.omega_s;
322c4762a1bSJed Brown       u[0] += du[0];
323c4762a1bSJed Brown       u[1] += du[1];
3249566063dSJacob Faibussowitsch       PetscCall(VecRestoreArray(U,&u));
3259566063dSJacob Faibussowitsch       PetscCall(TSSetTimeStep(ts,.01));
3269566063dSJacob Faibussowitsch       PetscCall(TSSolve(ts,U));
327c4762a1bSJed Brown     }
328c4762a1bSJed Brown   } else {
3299566063dSJacob Faibussowitsch     PetscCall(TSSolve(ts,U));
330c4762a1bSJed Brown   }
3319566063dSJacob Faibussowitsch   PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD));
3329566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
3339566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
336c4762a1bSJed Brown      Adjoint model starts here
337c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
338c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
3399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lambda[0],&y_ptr));
340c4762a1bSJed Brown   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lambda[0],&y_ptr));
342c4762a1bSJed Brown 
3439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(mu[0],&x_ptr));
344c4762a1bSJed Brown   x_ptr[0] = -1.0;
3459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(mu[0],&x_ptr));
346c4762a1bSJed Brown 
347c4762a1bSJed Brown   /*   Set RHS JacobianP */
3489566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&ctx));
349c4762a1bSJed Brown 
3509566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
351c4762a1bSJed Brown 
3529566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n"));
3539566063dSJacob Faibussowitsch   PetscCall(VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD));
3549566063dSJacob Faibussowitsch   PetscCall(VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD));
3559566063dSJacob Faibussowitsch   PetscCall(TSGetCostIntegral(ts,&q));
3569566063dSJacob Faibussowitsch   PetscCall(VecView(q,PETSC_VIEWER_STDOUT_WORLD));
3579566063dSJacob Faibussowitsch   PetscCall(VecGetArray(q,&x_ptr));
3589566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm)));
3599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(q,&x_ptr));
360c4762a1bSJed Brown 
3619566063dSJacob Faibussowitsch   PetscCall(ComputeSensiP(lambda[0],mu[0],&ctx));
362c4762a1bSJed Brown 
363c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
364c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
365c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3669566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
3679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jacp));
3689566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDU));
3699566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&DRDP));
3709566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
3719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lambda[0]));
3729566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mu[0]));
3739566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3749566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
375b122ec5aSJacob Faibussowitsch   return 0;
376c4762a1bSJed Brown }
377c4762a1bSJed Brown 
378c4762a1bSJed Brown /*TEST
379c4762a1bSJed Brown 
380c4762a1bSJed Brown    build:
381c4762a1bSJed Brown       requires: !complex
382c4762a1bSJed Brown 
383c4762a1bSJed Brown    test:
384c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none
385c4762a1bSJed Brown 
386c4762a1bSJed Brown TEST*/
387