xref: /petsc/src/ts/tutorials/power_grid/ex3opt.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\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 F*/
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown   This code demonstrates how to solve a ODE-constrained optimization problem with TAO, TSEvent, TSAdjoint and TS.
15c4762a1bSJed Brown   The problem features discontinuities and a cost function in integral form.
16c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit theta method, see ex3adj.c for details.
17c4762a1bSJed Brown */
18c4762a1bSJed Brown 
19c4762a1bSJed Brown #include <petsctao.h>
20c4762a1bSJed Brown #include <petscts.h>
21c4762a1bSJed Brown #include "ex3.h"
22c4762a1bSJed Brown 
23c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
24c4762a1bSJed Brown 
25c4762a1bSJed Brown PetscErrorCode monitor(Tao tao,AppCtx *ctx)
26c4762a1bSJed Brown {
27c4762a1bSJed Brown   FILE               *fp;
28c4762a1bSJed Brown   PetscInt           iterate;
29c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
30c4762a1bSJed Brown   TaoConvergedReason reason;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   PetscFunctionBeginUser;
339566063dSJacob Faibussowitsch   PetscCall(TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   fp = fopen("ex3opt_conv.out","a");
369566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fp,"%D %g\n",iterate,(double)gnorm));
37c4762a1bSJed Brown   fclose(fp);
38c4762a1bSJed Brown   PetscFunctionReturn(0);
39c4762a1bSJed Brown }
40c4762a1bSJed Brown 
41c4762a1bSJed Brown int main(int argc,char **argv)
42c4762a1bSJed Brown {
43c4762a1bSJed Brown   Vec                p;
44c4762a1bSJed Brown   PetscScalar        *x_ptr;
45c4762a1bSJed Brown   PetscMPIInt        size;
46c4762a1bSJed Brown   AppCtx             ctx;
47c4762a1bSJed Brown   Tao                tao;
48c4762a1bSJed Brown   KSP                ksp;
49c4762a1bSJed Brown   PC                 pc;
50c4762a1bSJed Brown   Vec                lambda[1],mu[1],lowerb,upperb;
51c4762a1bSJed Brown   PetscBool          printtofile;
52c4762a1bSJed Brown   PetscInt           direction[2];
53c4762a1bSJed Brown   PetscBool          terminate[2];
54c4762a1bSJed Brown   Mat                qgrad;         /* Forward sesivitiy */
55c4762a1bSJed Brown   Mat                sp;            /* Forward sensitivity matrix */
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Initialize program
59c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
609566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
61c4762a1bSJed Brown   PetscFunctionBeginUser;
629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
633c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4762a1bSJed Brown     Set runtime options
67c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
69c4762a1bSJed Brown   {
70c4762a1bSJed Brown     ctx.beta    = 2;
71c4762a1bSJed Brown     ctx.c       = 10000.0;
72c4762a1bSJed Brown     ctx.u_s     = 1.0;
73c4762a1bSJed Brown     ctx.omega_s = 1.0;
74c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
75c4762a1bSJed Brown     ctx.H       = 5.0;
769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL));
77c4762a1bSJed Brown     ctx.D       = 5.0;
789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL));
79c4762a1bSJed Brown     ctx.E       = 1.1378;
80c4762a1bSJed Brown     ctx.V       = 1.0;
81c4762a1bSJed Brown     ctx.X       = 0.545;
82c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
83c4762a1bSJed Brown     ctx.Pmax_ini = ctx.Pmax;
849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL));
85c4762a1bSJed Brown     ctx.Pm      = 1.06;
869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL));
87c4762a1bSJed Brown     ctx.tf      = 0.1;
88c4762a1bSJed Brown     ctx.tcl     = 0.2;
899566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL));
909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL));
91c4762a1bSJed Brown     printtofile = PETSC_FALSE;
929566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL));
93c4762a1bSJed Brown     ctx.sa      = SA_ADJ;
949566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-sa_method","Sensitivity analysis method (adj or tlm)","",SAMethods,(PetscEnum)ctx.sa,(PetscEnum*)&ctx.sa,NULL));
95c4762a1bSJed Brown   }
96*d0609cedSBarry Smith   PetscOptionsEnd();
97c4762a1bSJed Brown 
98c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
99c4762a1bSJed Brown     Create necessary matrix and vectors
100c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1019566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx.Jac));
1029566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(ctx.Jac,2,2,PETSC_DETERMINE,PETSC_DETERMINE));
1039566063dSJacob Faibussowitsch   PetscCall(MatSetType(ctx.Jac,MATDENSE));
1049566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(ctx.Jac));
1059566063dSJacob Faibussowitsch   PetscCall(MatSetUp(ctx.Jac));
1069566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&ctx.Jacp));
1079566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(ctx.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
1089566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(ctx.Jacp));
1099566063dSJacob Faibussowitsch   PetscCall(MatSetUp(ctx.Jacp));
1109566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(ctx.Jac,&ctx.U,NULL));
1119566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&ctx.DRDP));
1129566063dSJacob Faibussowitsch   PetscCall(MatSetUp(ctx.DRDP));
1139566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&ctx.DRDU));
1149566063dSJacob Faibussowitsch   PetscCall(MatSetUp(ctx.DRDU));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117c4762a1bSJed Brown      Create timestepping solver context
118c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1199566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ctx.ts));
1209566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ctx.ts,TS_NONLINEAR));
1219566063dSJacob Faibussowitsch   PetscCall(TSSetType(ctx.ts,TSCN));
1229566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx));
1239566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ctx.ts,ctx.Jac,ctx.Jac,(TSRHSJacobian)RHSJacobian,&ctx));
1249566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ctx.ts,ctx.Jacp,RHSJacobianP,&ctx));
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   if (ctx.sa == SA_ADJ) {
1279566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(ctx.Jac,&lambda[0],NULL));
1289566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(ctx.Jacp,&mu[0],NULL));
1299566063dSJacob Faibussowitsch     PetscCall(TSSetSaveTrajectory(ctx.ts));
1309566063dSJacob Faibussowitsch     PetscCall(TSSetCostGradients(ctx.ts,1,lambda,mu));
1319566063dSJacob Faibussowitsch     PetscCall(TSCreateQuadratureTS(ctx.ts,PETSC_FALSE,&ctx.quadts));
1329566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx));
1339566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx));
1349566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx));
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown   if (ctx.sa == SA_TLM) {
1379566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&qgrad));
1389566063dSJacob Faibussowitsch     PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&sp));
1399566063dSJacob Faibussowitsch     PetscCall(TSForwardSetSensitivities(ctx.ts,1,sp));
1409566063dSJacob Faibussowitsch     PetscCall(TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&ctx.quadts));
1419566063dSJacob Faibussowitsch     PetscCall(TSForwardSetSensitivities(ctx.quadts,1,qgrad));
1429566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx));
1439566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx));
1449566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx));
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148c4762a1bSJed Brown      Set solver options
149c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1509566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ctx.ts,1.0));
1519566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP));
1529566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ctx.ts,0.03125));
1539566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ctx.ts));
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   direction[0] = direction[1] = 1;
156c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
1579566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ctx.ts,2,direction,terminate,EventFunction,PostEventFunction,&ctx));
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1609566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
1619566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBLMVM));
162c4762a1bSJed Brown   if (printtofile) {
1639566063dSJacob Faibussowitsch     PetscCall(TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL));
164c4762a1bSJed Brown   }
165c4762a1bSJed Brown   /*
166c4762a1bSJed Brown      Optimization starts
167c4762a1bSJed Brown   */
168c4762a1bSJed Brown   /* Set initial solution guess */
1699566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,1,&p));
1709566063dSJacob Faibussowitsch   PetscCall(VecGetArray(p,&x_ptr));
171c4762a1bSJed Brown   x_ptr[0] = ctx.Pm;
1729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(p,&x_ptr));
173c4762a1bSJed Brown 
1749566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,p));
175c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
1769566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&ctx));
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* Set bounds for the optimization */
1799566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p,&lowerb));
1809566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p,&upperb));
1819566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lowerb,&x_ptr));
182c4762a1bSJed Brown   x_ptr[0] = 0.;
1839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lowerb,&x_ptr));
1849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(upperb,&x_ptr));
185c4762a1bSJed Brown   x_ptr[0] = 1.1;
1869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(upperb,&x_ptr));
1879566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao,lowerb,upperb));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* Check for any TAO command line options */
1909566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
1919566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao,&ksp));
192c4762a1bSJed Brown   if (ksp) {
1939566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp,&pc));
1949566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc,PCNONE));
195c4762a1bSJed Brown   }
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1989566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch   PetscCall(VecView(p,PETSC_VIEWER_STDOUT_WORLD));
201c4762a1bSJed Brown 
202c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
204c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2059566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx.Jac));
2069566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx.Jacp));
2079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx.DRDU));
2089566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&ctx.DRDP));
2099566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.U));
210c4762a1bSJed Brown   if (ctx.sa == SA_ADJ) {
2119566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&lambda[0]));
2129566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&mu[0]));
213c4762a1bSJed Brown   }
214c4762a1bSJed Brown   if (ctx.sa == SA_TLM) {
2159566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&qgrad));
2169566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&sp));
217c4762a1bSJed Brown   }
2189566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ctx.ts));
2199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&p));
2209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lowerb));
2219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&upperb));
2229566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
2239566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
224b122ec5aSJacob Faibussowitsch   return 0;
225c4762a1bSJed Brown }
226c4762a1bSJed Brown 
227c4762a1bSJed Brown /* ------------------------------------------------------------------ */
228c4762a1bSJed Brown /*
229c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
230c4762a1bSJed Brown 
231c4762a1bSJed Brown    Input Parameters:
232c4762a1bSJed Brown    tao - the Tao context
233c4762a1bSJed Brown    X   - the input vector
234a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
235c4762a1bSJed Brown 
236c4762a1bSJed Brown    Output Parameters:
237c4762a1bSJed Brown    f   - the newly evaluated function
238c4762a1bSJed Brown    G   - the newly evaluated gradient
239c4762a1bSJed Brown */
240c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
241c4762a1bSJed Brown {
242c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)ctx0;
243c4762a1bSJed Brown   PetscInt       nadj;
244c4762a1bSJed Brown   PetscReal      ftime;
245c4762a1bSJed Brown   PetscInt       steps;
246c4762a1bSJed Brown   PetscScalar    *u;
247c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr;
248c4762a1bSJed Brown   Vec            q;
249c4762a1bSJed Brown   Mat            qgrad;
250c4762a1bSJed Brown 
2519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P,(const PetscScalar**)&x_ptr));
252c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
2539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr));
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   /* reinitialize the solution vector */
2569566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx->U,&u));
257c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
258c4762a1bSJed Brown   u[1] = 1.0;
2599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx->U,&u));
2609566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ctx->ts,ctx->U));
261c4762a1bSJed Brown 
262c4762a1bSJed Brown   /* reset time */
2639566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ctx->ts,0.0));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* reset step counter, this is critical for adjoint solver */
2669566063dSJacob Faibussowitsch   PetscCall(TSSetStepNumber(ctx->ts,0));
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   /* reset step size, the step size becomes negative after TSAdjointSolve */
2699566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ctx->ts,0.03125));
270c4762a1bSJed Brown 
271c4762a1bSJed Brown   /* reinitialize the integral value */
2729566063dSJacob Faibussowitsch   PetscCall(TSGetQuadratureTS(ctx->ts,NULL,&ctx->quadts));
2739566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ctx->quadts,&q));
2749566063dSJacob Faibussowitsch   PetscCall(VecSet(q,0.0));
275c4762a1bSJed Brown 
276c4762a1bSJed Brown   if (ctx->sa == SA_TLM) { /* reset the forward sensitivities */
277c4762a1bSJed Brown     TS             quadts;
278c4762a1bSJed Brown     Mat            sp;
279c4762a1bSJed Brown     PetscScalar    val[2];
280c4762a1bSJed Brown     const PetscInt row[]={0,1},col[]={0};
281c4762a1bSJed Brown 
2829566063dSJacob Faibussowitsch     PetscCall(TSGetQuadratureTS(ctx->ts,NULL,&quadts));
2839566063dSJacob Faibussowitsch     PetscCall(TSForwardGetSensitivities(quadts,NULL,&qgrad));
2849566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(qgrad));
2859566063dSJacob Faibussowitsch     PetscCall(TSForwardGetSensitivities(ctx->ts,NULL,&sp));
286c4762a1bSJed Brown     val[0] = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax;
287c4762a1bSJed Brown     val[1] = 0.0;
2889566063dSJacob Faibussowitsch     PetscCall(MatSetValues(sp,2,row,1,col,val,INSERT_VALUES));
2899566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(sp,MAT_FINAL_ASSEMBLY));
2909566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(sp,MAT_FINAL_ASSEMBLY));
291c4762a1bSJed Brown   }
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   /* solve the ODE */
2949566063dSJacob Faibussowitsch   PetscCall(TSSolve(ctx->ts,ctx->U));
2959566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ctx->ts,&ftime));
2969566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ctx->ts,&steps));
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   if (ctx->sa == SA_ADJ) {
299c4762a1bSJed Brown     Vec *lambda,*mu;
300c4762a1bSJed Brown     /* reset the terminal condition for adjoint */
3019566063dSJacob Faibussowitsch     PetscCall(TSGetCostGradients(ctx->ts,&nadj,&lambda,&mu));
3029566063dSJacob Faibussowitsch     PetscCall(VecGetArray(lambda[0],&y_ptr));
303c4762a1bSJed Brown     y_ptr[0] = 0.0; y_ptr[1] = 0.0;
3049566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(lambda[0],&y_ptr));
3059566063dSJacob Faibussowitsch     PetscCall(VecGetArray(mu[0],&x_ptr));
306c4762a1bSJed Brown     x_ptr[0] = -1.0;
3079566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(mu[0],&x_ptr));
308c4762a1bSJed Brown 
309c4762a1bSJed Brown     /* solve the adjont */
3109566063dSJacob Faibussowitsch     PetscCall(TSAdjointSolve(ctx->ts));
311c4762a1bSJed Brown 
3129566063dSJacob Faibussowitsch     PetscCall(ComputeSensiP(lambda[0],mu[0],ctx));
3139566063dSJacob Faibussowitsch     PetscCall(VecCopy(mu[0],G));
314c4762a1bSJed Brown   }
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   if (ctx->sa == SA_TLM) {
3179566063dSJacob Faibussowitsch     PetscCall(VecGetArray(G,&x_ptr));
3189566063dSJacob Faibussowitsch     PetscCall(MatDenseGetArray(qgrad,&y_ptr));
319c4762a1bSJed Brown     x_ptr[0] = y_ptr[0]-1.;
3209566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(qgrad,&y_ptr));
3219566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(G,&x_ptr));
322c4762a1bSJed Brown   }
323c4762a1bSJed Brown 
3249566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(ctx->quadts,&q));
3259566063dSJacob Faibussowitsch   PetscCall(VecGetArray(q,&x_ptr));
326c4762a1bSJed Brown   *f   = -ctx->Pm + x_ptr[0];
3279566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(q,&x_ptr));
328c4762a1bSJed Brown   return 0;
329c4762a1bSJed Brown }
330c4762a1bSJed Brown 
331c4762a1bSJed Brown /*TEST
332c4762a1bSJed Brown 
333c4762a1bSJed Brown    build:
334c4762a1bSJed Brown       requires: !complex !single
335c4762a1bSJed Brown 
336c4762a1bSJed Brown    test:
337c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_type cn -pc_type lu -tao_monitor
338c4762a1bSJed Brown 
339c4762a1bSJed Brown    test:
340c4762a1bSJed Brown       suffix: 2
341c4762a1bSJed Brown       output_file: output/ex3opt_1.out
342c4762a1bSJed Brown       args: -sa_method tlm -ts_type cn -pc_type lu -tao_monitor
343c4762a1bSJed Brown TEST*/
344