xref: /petsc/src/ts/tutorials/power_grid/ex3opt.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
335f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason));
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   fp = fopen("ex3opt_conv.out","a");
365f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode     ierr;
46c4762a1bSJed Brown   PetscMPIInt        size;
47c4762a1bSJed Brown   AppCtx             ctx;
48c4762a1bSJed Brown   Tao                tao;
49c4762a1bSJed Brown   KSP                ksp;
50c4762a1bSJed Brown   PC                 pc;
51c4762a1bSJed Brown   Vec                lambda[1],mu[1],lowerb,upperb;
52c4762a1bSJed Brown   PetscBool          printtofile;
53c4762a1bSJed Brown   PetscInt           direction[2];
54c4762a1bSJed Brown   PetscBool          terminate[2];
55c4762a1bSJed Brown   Mat                qgrad;         /* Forward sesivitiy */
56c4762a1bSJed Brown   Mat                sp;            /* Forward sensitivity matrix */
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59c4762a1bSJed Brown      Initialize program
60c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
62c4762a1bSJed Brown   PetscFunctionBeginUser;
635f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
643c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4762a1bSJed Brown     Set runtime options
68c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
70c4762a1bSJed Brown   {
71c4762a1bSJed Brown     ctx.beta    = 2;
72c4762a1bSJed Brown     ctx.c       = 10000.0;
73c4762a1bSJed Brown     ctx.u_s     = 1.0;
74c4762a1bSJed Brown     ctx.omega_s = 1.0;
75c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
76c4762a1bSJed Brown     ctx.H       = 5.0;
775f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL));
78c4762a1bSJed Brown     ctx.D       = 5.0;
795f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL));
80c4762a1bSJed Brown     ctx.E       = 1.1378;
81c4762a1bSJed Brown     ctx.V       = 1.0;
82c4762a1bSJed Brown     ctx.X       = 0.545;
83c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
84c4762a1bSJed Brown     ctx.Pmax_ini = ctx.Pmax;
855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL));
86c4762a1bSJed Brown     ctx.Pm      = 1.06;
875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL));
88c4762a1bSJed Brown     ctx.tf      = 0.1;
89c4762a1bSJed Brown     ctx.tcl     = 0.2;
905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL));
915f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL));
92c4762a1bSJed Brown     printtofile = PETSC_FALSE;
935f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL));
94c4762a1bSJed Brown     ctx.sa      = SA_ADJ;
955f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsEnum("-sa_method","Sensitivity analysis method (adj or tlm)","",SAMethods,(PetscEnum)ctx.sa,(PetscEnum*)&ctx.sa,NULL));
96c4762a1bSJed Brown   }
97c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100c4762a1bSJed Brown     Create necessary matrix and vectors
101c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1025f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&ctx.Jac));
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(ctx.Jac,2,2,PETSC_DETERMINE,PETSC_DETERMINE));
1045f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(ctx.Jac,MATDENSE));
1055f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(ctx.Jac));
1065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(ctx.Jac));
1075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&ctx.Jacp));
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(ctx.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(ctx.Jacp));
1105f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(ctx.Jacp));
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(ctx.Jac,&ctx.U,NULL));
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&ctx.DRDP));
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(ctx.DRDP));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&ctx.DRDU));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(ctx.DRDU));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118c4762a1bSJed Brown      Create timestepping solver context
119c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ctx.ts));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ctx.ts,TS_NONLINEAR));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ctx.ts,TSCN));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ctx.ts,ctx.Jac,ctx.Jac,(TSRHSJacobian)RHSJacobian,&ctx));
1255f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(ctx.ts,ctx.Jacp,RHSJacobianP,&ctx));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   if (ctx.sa == SA_ADJ) {
1285f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(ctx.Jac,&lambda[0],NULL));
1295f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateVecs(ctx.Jacp,&mu[0],NULL));
1305f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetSaveTrajectory(ctx.ts));
1315f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetCostGradients(ctx.ts,1,lambda,mu));
1325f80ce2aSJacob Faibussowitsch     CHKERRQ(TSCreateQuadratureTS(ctx.ts,PETSC_FALSE,&ctx.quadts));
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx));
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx));
1355f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx));
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown   if (ctx.sa == SA_TLM) {
1385f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&qgrad));
1395f80ce2aSJacob Faibussowitsch     CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&sp));
1405f80ce2aSJacob Faibussowitsch     CHKERRQ(TSForwardSetSensitivities(ctx.ts,1,sp));
1415f80ce2aSJacob Faibussowitsch     CHKERRQ(TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&ctx.quadts));
1425f80ce2aSJacob Faibussowitsch     CHKERRQ(TSForwardSetSensitivities(ctx.quadts,1,qgrad));
1435f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx));
1445f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx));
1455f80ce2aSJacob Faibussowitsch     CHKERRQ(TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx));
146c4762a1bSJed Brown   }
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown      Set solver options
150c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1515f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ctx.ts,1.0));
1525f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP));
1535f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ctx.ts,0.03125));
1545f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ctx.ts));
155c4762a1bSJed Brown 
156c4762a1bSJed Brown   direction[0] = direction[1] = 1;
157c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetEventHandler(ctx.ts,2,direction,terminate,EventFunction,PostEventFunction,&ctx));
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
1625f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBLMVM));
163c4762a1bSJed Brown   if (printtofile) {
1645f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL));
165c4762a1bSJed Brown   }
166c4762a1bSJed Brown   /*
167c4762a1bSJed Brown      Optimization starts
168c4762a1bSJed Brown   */
169c4762a1bSJed Brown   /* Set initial solution guess */
1705f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD,1,&p));
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(p,&x_ptr));
172c4762a1bSJed Brown   x_ptr[0] = ctx.Pm;
1735f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(p,&x_ptr));
174c4762a1bSJed Brown 
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,p));
176c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&ctx));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* Set bounds for the optimization */
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(p,&lowerb));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(p,&upperb));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(lowerb,&x_ptr));
183c4762a1bSJed Brown   x_ptr[0] = 0.;
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(lowerb,&x_ptr));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(upperb,&x_ptr));
186c4762a1bSJed Brown   x_ptr[0] = 1.1;
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(upperb,&x_ptr));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetVariableBounds(tao,lowerb,upperb));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* Check for any TAO command line options */
1915f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
193c4762a1bSJed Brown   if (ksp) {
1945f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
1955f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCNONE));
196c4762a1bSJed Brown   }
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
200c4762a1bSJed Brown 
2015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(p,PETSC_VIEWER_STDOUT_WORLD));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
205c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2065f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&ctx.Jac));
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&ctx.Jacp));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&ctx.DRDU));
2095f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&ctx.DRDP));
2105f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ctx.U));
211c4762a1bSJed Brown   if (ctx.sa == SA_ADJ) {
2125f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&lambda[0]));
2135f80ce2aSJacob Faibussowitsch     CHKERRQ(VecDestroy(&mu[0]));
214c4762a1bSJed Brown   }
215c4762a1bSJed Brown   if (ctx.sa == SA_TLM) {
2165f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&qgrad));
2175f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDestroy(&sp));
218c4762a1bSJed Brown   }
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ctx.ts));
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&p));
2215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&lowerb));
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&upperb));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
224*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
225*b122ec5aSJacob Faibussowitsch   return 0;
226c4762a1bSJed Brown }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown /* ------------------------------------------------------------------ */
229c4762a1bSJed Brown /*
230c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
231c4762a1bSJed Brown 
232c4762a1bSJed Brown    Input Parameters:
233c4762a1bSJed Brown    tao - the Tao context
234c4762a1bSJed Brown    X   - the input vector
235a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
236c4762a1bSJed Brown 
237c4762a1bSJed Brown    Output Parameters:
238c4762a1bSJed Brown    f   - the newly evaluated function
239c4762a1bSJed Brown    G   - the newly evaluated gradient
240c4762a1bSJed Brown */
241c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
242c4762a1bSJed Brown {
243c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)ctx0;
244c4762a1bSJed Brown   PetscInt       nadj;
245c4762a1bSJed Brown   PetscReal      ftime;
246c4762a1bSJed Brown   PetscInt       steps;
247c4762a1bSJed Brown   PetscScalar    *u;
248c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr;
249c4762a1bSJed Brown   Vec            q;
250c4762a1bSJed Brown   Mat            qgrad;
251c4762a1bSJed Brown 
2525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,(const PetscScalar**)&x_ptr));
253c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
2545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr));
255c4762a1bSJed Brown 
256c4762a1bSJed Brown   /* reinitialize the solution vector */
2575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ctx->U,&u));
258c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
259c4762a1bSJed Brown   u[1] = 1.0;
2605f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ctx->U,&u));
2615f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ctx->ts,ctx->U));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /* reset time */
2645f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ctx->ts,0.0));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* reset step counter, this is critical for adjoint solver */
2675f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetStepNumber(ctx->ts,0));
268c4762a1bSJed Brown 
269c4762a1bSJed Brown   /* reset step size, the step size becomes negative after TSAdjointSolve */
2705f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ctx->ts,0.03125));
271c4762a1bSJed Brown 
272c4762a1bSJed Brown   /* reinitialize the integral value */
2735f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetQuadratureTS(ctx->ts,NULL,&ctx->quadts));
2745f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ctx->quadts,&q));
2755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(q,0.0));
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   if (ctx->sa == SA_TLM) { /* reset the forward sensitivities */
278c4762a1bSJed Brown     TS             quadts;
279c4762a1bSJed Brown     Mat            sp;
280c4762a1bSJed Brown     PetscScalar    val[2];
281c4762a1bSJed Brown     const PetscInt row[]={0,1},col[]={0};
282c4762a1bSJed Brown 
2835f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetQuadratureTS(ctx->ts,NULL,&quadts));
2845f80ce2aSJacob Faibussowitsch     CHKERRQ(TSForwardGetSensitivities(quadts,NULL,&qgrad));
2855f80ce2aSJacob Faibussowitsch     CHKERRQ(MatZeroEntries(qgrad));
2865f80ce2aSJacob Faibussowitsch     CHKERRQ(TSForwardGetSensitivities(ctx->ts,NULL,&sp));
287c4762a1bSJed Brown     val[0] = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax;
288c4762a1bSJed Brown     val[1] = 0.0;
2895f80ce2aSJacob Faibussowitsch     CHKERRQ(MatSetValues(sp,2,row,1,col,val,INSERT_VALUES));
2905f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(sp,MAT_FINAL_ASSEMBLY));
2915f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(sp,MAT_FINAL_ASSEMBLY));
292c4762a1bSJed Brown   }
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   /* solve the ODE */
2955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ctx->ts,ctx->U));
2965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ctx->ts,&ftime));
2975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ctx->ts,&steps));
298c4762a1bSJed Brown 
299c4762a1bSJed Brown   if (ctx->sa == SA_ADJ) {
300c4762a1bSJed Brown     Vec *lambda,*mu;
301c4762a1bSJed Brown     /* reset the terminal condition for adjoint */
3025f80ce2aSJacob Faibussowitsch     CHKERRQ(TSGetCostGradients(ctx->ts,&nadj,&lambda,&mu));
3035f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(lambda[0],&y_ptr));
304c4762a1bSJed Brown     y_ptr[0] = 0.0; y_ptr[1] = 0.0;
3055f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(lambda[0],&y_ptr));
3065f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(mu[0],&x_ptr));
307c4762a1bSJed Brown     x_ptr[0] = -1.0;
3085f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(mu[0],&x_ptr));
309c4762a1bSJed Brown 
310c4762a1bSJed Brown     /* solve the adjont */
3115f80ce2aSJacob Faibussowitsch     CHKERRQ(TSAdjointSolve(ctx->ts));
312c4762a1bSJed Brown 
3135f80ce2aSJacob Faibussowitsch     CHKERRQ(ComputeSensiP(lambda[0],mu[0],ctx));
3145f80ce2aSJacob Faibussowitsch     CHKERRQ(VecCopy(mu[0],G));
315c4762a1bSJed Brown   }
316c4762a1bSJed Brown 
317c4762a1bSJed Brown   if (ctx->sa == SA_TLM) {
3185f80ce2aSJacob Faibussowitsch     CHKERRQ(VecGetArray(G,&x_ptr));
3195f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseGetArray(qgrad,&y_ptr));
320c4762a1bSJed Brown     x_ptr[0] = y_ptr[0]-1.;
3215f80ce2aSJacob Faibussowitsch     CHKERRQ(MatDenseRestoreArray(qgrad,&y_ptr));
3225f80ce2aSJacob Faibussowitsch     CHKERRQ(VecRestoreArray(G,&x_ptr));
323c4762a1bSJed Brown   }
324c4762a1bSJed Brown 
3255f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(ctx->quadts,&q));
3265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(q,&x_ptr));
327c4762a1bSJed Brown   *f   = -ctx->Pm + x_ptr[0];
3285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(q,&x_ptr));
329c4762a1bSJed Brown   return 0;
330c4762a1bSJed Brown }
331c4762a1bSJed Brown 
332c4762a1bSJed Brown /*TEST
333c4762a1bSJed Brown 
334c4762a1bSJed Brown    build:
335c4762a1bSJed Brown       requires: !complex !single
336c4762a1bSJed Brown 
337c4762a1bSJed Brown    test:
338c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_type cn -pc_type lu -tao_monitor
339c4762a1bSJed Brown 
340c4762a1bSJed Brown    test:
341c4762a1bSJed Brown       suffix: 2
342c4762a1bSJed Brown       output_file: output/ex3opt_1.out
343c4762a1bSJed Brown       args: -sa_method tlm -ts_type cn -pc_type lu -tao_monitor
344c4762a1bSJed Brown TEST*/
345