xref: /petsc/src/ts/tutorials/power_grid/ex3opt.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*F
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown \begin{eqnarray}
7*c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8*c4762a1bSJed Brown                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9*c4762a1bSJed Brown \end{eqnarray}
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown F*/
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown /*
14*c4762a1bSJed Brown   This code demonstrates how to solve a ODE-constrained optimization problem with TAO, TSEvent, TSAdjoint and TS.
15*c4762a1bSJed Brown   The problem features discontinuities and a cost function in integral form.
16*c4762a1bSJed Brown   The gradient is computed with the discrete adjoint of an implicit theta method, see ex3adj.c for details.
17*c4762a1bSJed Brown */
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown #include <petsctao.h>
20*c4762a1bSJed Brown #include <petscts.h>
21*c4762a1bSJed Brown #include "ex3.h"
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown PetscErrorCode monitor(Tao tao,AppCtx *ctx)
26*c4762a1bSJed Brown {
27*c4762a1bSJed Brown   FILE               *fp;
28*c4762a1bSJed Brown   PetscInt           iterate;
29*c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
30*c4762a1bSJed Brown   TaoConvergedReason reason;
31*c4762a1bSJed Brown   PetscErrorCode     ierr;
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   PetscFunctionBeginUser;
34*c4762a1bSJed Brown   ierr = TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason);CHKERRQ(ierr);
35*c4762a1bSJed Brown 
36*c4762a1bSJed Brown   fp = fopen("ex3opt_conv.out","a");
37*c4762a1bSJed Brown   ierr = PetscFPrintf(PETSC_COMM_WORLD,fp,"%D %g\n",iterate,(double)gnorm);CHKERRQ(ierr);
38*c4762a1bSJed Brown   fclose(fp);
39*c4762a1bSJed Brown   PetscFunctionReturn(0);
40*c4762a1bSJed Brown }
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown int main(int argc,char **argv)
43*c4762a1bSJed Brown {
44*c4762a1bSJed Brown   Vec                p;
45*c4762a1bSJed Brown   PetscScalar        *x_ptr;
46*c4762a1bSJed Brown   PetscErrorCode     ierr;
47*c4762a1bSJed Brown   PetscMPIInt        size;
48*c4762a1bSJed Brown   AppCtx             ctx;
49*c4762a1bSJed Brown   Tao                tao;
50*c4762a1bSJed Brown   KSP                ksp;
51*c4762a1bSJed Brown   PC                 pc;
52*c4762a1bSJed Brown   Vec                lambda[1],mu[1],lowerb,upperb;
53*c4762a1bSJed Brown   PetscBool          printtofile;
54*c4762a1bSJed Brown   PetscInt           direction[2];
55*c4762a1bSJed Brown   PetscBool          terminate[2];
56*c4762a1bSJed Brown   Mat                qgrad;         /* Forward sesivitiy */
57*c4762a1bSJed Brown   Mat                sp;            /* Forward sensitivity matrix */
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60*c4762a1bSJed Brown      Initialize program
61*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
63*c4762a1bSJed Brown   PetscFunctionBeginUser;
64*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
65*c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68*c4762a1bSJed Brown     Set runtime options
69*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
71*c4762a1bSJed Brown   {
72*c4762a1bSJed Brown     ctx.beta    = 2;
73*c4762a1bSJed Brown     ctx.c       = 10000.0;
74*c4762a1bSJed Brown     ctx.u_s     = 1.0;
75*c4762a1bSJed Brown     ctx.omega_s = 1.0;
76*c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
77*c4762a1bSJed Brown     ctx.H       = 5.0;
78*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
79*c4762a1bSJed Brown     ctx.D       = 5.0;
80*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
81*c4762a1bSJed Brown     ctx.E       = 1.1378;
82*c4762a1bSJed Brown     ctx.V       = 1.0;
83*c4762a1bSJed Brown     ctx.X       = 0.545;
84*c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
85*c4762a1bSJed Brown     ctx.Pmax_ini = ctx.Pmax;
86*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
87*c4762a1bSJed Brown     ctx.Pm      = 1.06;
88*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
89*c4762a1bSJed Brown     ctx.tf      = 0.1;
90*c4762a1bSJed Brown     ctx.tcl     = 0.2;
91*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
92*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
93*c4762a1bSJed Brown     printtofile = PETSC_FALSE;
94*c4762a1bSJed Brown     ierr        = PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL);CHKERRQ(ierr);
95*c4762a1bSJed Brown     ctx.sa      = SA_ADJ;
96*c4762a1bSJed Brown     ierr        = PetscOptionsEnum("-sa_method","Sensitivity analysis method (adj or tlm)","",SAMethods,(PetscEnum)ctx.sa,(PetscEnum*)&ctx.sa,NULL);CHKERRQ(ierr);
97*c4762a1bSJed Brown   }
98*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
101*c4762a1bSJed Brown     Create necessary matrix and vectors
102*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&ctx.Jac);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = MatSetSizes(ctx.Jac,2,2,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatSetType(ctx.Jac,MATDENSE);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = MatSetFromOptions(ctx.Jac);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = MatSetUp(ctx.Jac);CHKERRQ(ierr);
108*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&ctx.Jacp);CHKERRQ(ierr);
109*c4762a1bSJed Brown   ierr = MatSetSizes(ctx.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = MatSetFromOptions(ctx.Jacp);CHKERRQ(ierr);
111*c4762a1bSJed Brown   ierr = MatSetUp(ctx.Jacp);CHKERRQ(ierr);
112*c4762a1bSJed Brown   ierr = MatCreateVecs(ctx.Jac,&ctx.U,NULL);CHKERRQ(ierr);
113*c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&ctx.DRDP);CHKERRQ(ierr);
114*c4762a1bSJed Brown   ierr = MatSetUp(ctx.DRDP);CHKERRQ(ierr);
115*c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&ctx.DRDU);CHKERRQ(ierr);
116*c4762a1bSJed Brown   ierr = MatSetUp(ctx.DRDU);CHKERRQ(ierr);
117*c4762a1bSJed Brown 
118*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
119*c4762a1bSJed Brown      Create timestepping solver context
120*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
121*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ctx.ts);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = TSSetProblemType(ctx.ts,TS_NONLINEAR);CHKERRQ(ierr);
123*c4762a1bSJed Brown   ierr = TSSetType(ctx.ts,TSCN);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ctx.ts,ctx.Jac,ctx.Jac,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr);
126*c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(ctx.ts,ctx.Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr);
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown   if (ctx.sa == SA_ADJ) {
129*c4762a1bSJed Brown     ierr = MatCreateVecs(ctx.Jac,&lambda[0],NULL);CHKERRQ(ierr);
130*c4762a1bSJed Brown     ierr = MatCreateVecs(ctx.Jacp,&mu[0],NULL);CHKERRQ(ierr);
131*c4762a1bSJed Brown     ierr = TSSetSaveTrajectory(ctx.ts);CHKERRQ(ierr);
132*c4762a1bSJed Brown     ierr = TSSetCostGradients(ctx.ts,1,lambda,mu);CHKERRQ(ierr);
133*c4762a1bSJed Brown     ierr = TSCreateQuadratureTS(ctx.ts,PETSC_FALSE,&ctx.quadts);CHKERRQ(ierr);
134*c4762a1bSJed Brown     ierr = TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr);
135*c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr);
136*c4762a1bSJed Brown     ierr = TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx);CHKERRQ(ierr);
137*c4762a1bSJed Brown   }
138*c4762a1bSJed Brown   if (ctx.sa == SA_TLM) {
139*c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&qgrad);CHKERRQ(ierr);
140*c4762a1bSJed Brown     ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&sp);CHKERRQ(ierr);
141*c4762a1bSJed Brown     ierr = TSForwardSetSensitivities(ctx.ts,1,sp);CHKERRQ(ierr);
142*c4762a1bSJed Brown     ierr = TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&ctx.quadts);CHKERRQ(ierr);
143*c4762a1bSJed Brown     ierr = TSForwardSetSensitivities(ctx.quadts,1,qgrad);CHKERRQ(ierr);
144*c4762a1bSJed Brown     ierr = TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr);
145*c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr);
146*c4762a1bSJed Brown     ierr = TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx);CHKERRQ(ierr);
147*c4762a1bSJed Brown   }
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150*c4762a1bSJed Brown      Set solver options
151*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152*c4762a1bSJed Brown   ierr = TSSetMaxTime(ctx.ts,1.0);CHKERRQ(ierr);
153*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = TSSetTimeStep(ctx.ts,0.03125);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr = TSSetFromOptions(ctx.ts);CHKERRQ(ierr);
156*c4762a1bSJed Brown 
157*c4762a1bSJed Brown   direction[0] = direction[1] = 1;
158*c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
159*c4762a1bSJed Brown   ierr = TSSetEventHandler(ctx.ts,2,direction,terminate,EventFunction,PostEventFunction,&ctx);CHKERRQ(ierr);
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
162*c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
164*c4762a1bSJed Brown   if(printtofile) {
165*c4762a1bSJed Brown     ierr = TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL);CHKERRQ(ierr);
166*c4762a1bSJed Brown   }
167*c4762a1bSJed Brown   /*
168*c4762a1bSJed Brown      Optimization starts
169*c4762a1bSJed Brown   */
170*c4762a1bSJed Brown   /* Set initial solution guess */
171*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_WORLD,1,&p);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr);
173*c4762a1bSJed Brown   x_ptr[0] = ctx.Pm;
174*c4762a1bSJed Brown   ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr);
175*c4762a1bSJed Brown 
176*c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,p);CHKERRQ(ierr);
177*c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
178*c4762a1bSJed Brown   ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&ctx);CHKERRQ(ierr);
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   /* Set bounds for the optimization */
181*c4762a1bSJed Brown   ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr);
182*c4762a1bSJed Brown   ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr);
184*c4762a1bSJed Brown   x_ptr[0] = 0.;
185*c4762a1bSJed Brown   ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr);
187*c4762a1bSJed Brown   x_ptr[0] = 1.1;
188*c4762a1bSJed Brown   ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr);
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown   /* Check for any TAO command line options */
192*c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
193*c4762a1bSJed Brown   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
194*c4762a1bSJed Brown   if (ksp) {
195*c4762a1bSJed Brown     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
196*c4762a1bSJed Brown     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
197*c4762a1bSJed Brown   }
198*c4762a1bSJed Brown 
199*c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
200*c4762a1bSJed Brown   ierr = TaoSolve(tao);CHKERRQ(ierr);
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown   ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
203*c4762a1bSJed Brown 
204*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
206*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207*c4762a1bSJed Brown   ierr = MatDestroy(&ctx.Jac);CHKERRQ(ierr);
208*c4762a1bSJed Brown   ierr = MatDestroy(&ctx.Jacp);CHKERRQ(ierr);
209*c4762a1bSJed Brown   ierr = MatDestroy(&ctx.DRDU);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = MatDestroy(&ctx.DRDP);CHKERRQ(ierr);
211*c4762a1bSJed Brown   ierr = VecDestroy(&ctx.U);CHKERRQ(ierr);
212*c4762a1bSJed Brown   if (ctx.sa == SA_ADJ) {
213*c4762a1bSJed Brown     ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
214*c4762a1bSJed Brown     ierr = VecDestroy(&mu[0]);CHKERRQ(ierr);
215*c4762a1bSJed Brown   }
216*c4762a1bSJed Brown   if (ctx.sa == SA_TLM) {
217*c4762a1bSJed Brown     ierr = MatDestroy(&qgrad);CHKERRQ(ierr);
218*c4762a1bSJed Brown     ierr = MatDestroy(&sp);CHKERRQ(ierr);
219*c4762a1bSJed Brown   }
220*c4762a1bSJed Brown   ierr = TSDestroy(&ctx.ts);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = VecDestroy(&p);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = VecDestroy(&lowerb);CHKERRQ(ierr);
223*c4762a1bSJed Brown   ierr = VecDestroy(&upperb);CHKERRQ(ierr);
224*c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
225*c4762a1bSJed Brown   ierr = PetscFinalize();
226*c4762a1bSJed Brown   return ierr;
227*c4762a1bSJed Brown }
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown /* ------------------------------------------------------------------ */
230*c4762a1bSJed Brown /*
231*c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
232*c4762a1bSJed Brown 
233*c4762a1bSJed Brown    Input Parameters:
234*c4762a1bSJed Brown    tao - the Tao context
235*c4762a1bSJed Brown    X   - the input vector
236*c4762a1bSJed Brown    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown    Output Parameters:
239*c4762a1bSJed Brown    f   - the newly evaluated function
240*c4762a1bSJed Brown    G   - the newly evaluated gradient
241*c4762a1bSJed Brown */
242*c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
243*c4762a1bSJed Brown {
244*c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)ctx0;
245*c4762a1bSJed Brown   PetscInt       nadj;
246*c4762a1bSJed Brown   PetscReal      ftime;
247*c4762a1bSJed Brown   PetscInt       steps;
248*c4762a1bSJed Brown   PetscScalar    *u;
249*c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr;
250*c4762a1bSJed Brown   Vec            q;
251*c4762a1bSJed Brown   Mat            qgrad;
252*c4762a1bSJed Brown   PetscErrorCode ierr;
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown   ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
255*c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
256*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
257*c4762a1bSJed Brown 
258*c4762a1bSJed Brown   /* reinitialize the solution vector */
259*c4762a1bSJed Brown   ierr = VecGetArray(ctx->U,&u);CHKERRQ(ierr);
260*c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
261*c4762a1bSJed Brown   u[1] = 1.0;
262*c4762a1bSJed Brown   ierr = VecRestoreArray(ctx->U,&u);CHKERRQ(ierr);
263*c4762a1bSJed Brown   ierr = TSSetSolution(ctx->ts,ctx->U);CHKERRQ(ierr);
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown   /* reset time */
266*c4762a1bSJed Brown   ierr = TSSetTime(ctx->ts,0.0);CHKERRQ(ierr);
267*c4762a1bSJed Brown 
268*c4762a1bSJed Brown   /* reset step counter, this is critical for adjoint solver */
269*c4762a1bSJed Brown   ierr = TSSetStepNumber(ctx->ts,0);CHKERRQ(ierr);
270*c4762a1bSJed Brown 
271*c4762a1bSJed Brown   /* reset step size, the step size becomes negative after TSAdjointSolve */
272*c4762a1bSJed Brown   ierr = TSSetTimeStep(ctx->ts,0.03125);CHKERRQ(ierr);
273*c4762a1bSJed Brown 
274*c4762a1bSJed Brown   /* reinitialize the integral value */
275*c4762a1bSJed Brown   ierr = TSGetQuadratureTS(ctx->ts,NULL,&ctx->quadts);CHKERRQ(ierr);
276*c4762a1bSJed Brown   ierr = TSGetSolution(ctx->quadts,&q);CHKERRQ(ierr);
277*c4762a1bSJed Brown   ierr = VecSet(q,0.0);CHKERRQ(ierr);
278*c4762a1bSJed Brown 
279*c4762a1bSJed Brown   if (ctx->sa == SA_TLM) { /* reset the forward sensitivities */
280*c4762a1bSJed Brown     TS             quadts;
281*c4762a1bSJed Brown     Mat            sp;
282*c4762a1bSJed Brown     PetscScalar    val[2];
283*c4762a1bSJed Brown     const PetscInt row[]={0,1},col[]={0};
284*c4762a1bSJed Brown 
285*c4762a1bSJed Brown     ierr = TSGetQuadratureTS(ctx->ts,NULL,&quadts);CHKERRQ(ierr);
286*c4762a1bSJed Brown     ierr = TSForwardGetSensitivities(quadts,NULL,&qgrad);CHKERRQ(ierr);
287*c4762a1bSJed Brown     ierr = MatZeroEntries(qgrad);CHKERRQ(ierr);
288*c4762a1bSJed Brown     ierr = TSForwardGetSensitivities(ctx->ts,NULL,&sp);CHKERRQ(ierr);
289*c4762a1bSJed Brown     val[0] = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax;
290*c4762a1bSJed Brown     val[1] = 0.0;
291*c4762a1bSJed Brown     ierr = MatSetValues(sp,2,row,1,col,val,INSERT_VALUES);CHKERRQ(ierr);
292*c4762a1bSJed Brown     ierr = MatAssemblyBegin(sp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
293*c4762a1bSJed Brown     ierr = MatAssemblyEnd(sp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
294*c4762a1bSJed Brown   }
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown   /* reset the trajectory */
297*c4762a1bSJed Brown   ierr = TSResetTrajectory(ctx->ts);CHKERRQ(ierr);
298*c4762a1bSJed Brown 
299*c4762a1bSJed Brown   /* solve the ODE */
300*c4762a1bSJed Brown   ierr = TSSolve(ctx->ts,ctx->U);CHKERRQ(ierr);
301*c4762a1bSJed Brown   ierr = TSGetSolveTime(ctx->ts,&ftime);CHKERRQ(ierr);
302*c4762a1bSJed Brown   ierr = TSGetStepNumber(ctx->ts,&steps);CHKERRQ(ierr);
303*c4762a1bSJed Brown 
304*c4762a1bSJed Brown   if (ctx->sa == SA_ADJ) {
305*c4762a1bSJed Brown     Vec *lambda,*mu;
306*c4762a1bSJed Brown     /* reset the terminal condition for adjoint */
307*c4762a1bSJed Brown     ierr = TSGetCostGradients(ctx->ts,&nadj,&lambda,&mu);CHKERRQ(ierr);
308*c4762a1bSJed Brown     ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr);
309*c4762a1bSJed Brown     y_ptr[0] = 0.0; y_ptr[1] = 0.0;
310*c4762a1bSJed Brown     ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr);
311*c4762a1bSJed Brown     ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
312*c4762a1bSJed Brown     x_ptr[0] = -1.0;
313*c4762a1bSJed Brown     ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
314*c4762a1bSJed Brown 
315*c4762a1bSJed Brown     /* solve the adjont */
316*c4762a1bSJed Brown     ierr = TSAdjointSolve(ctx->ts);CHKERRQ(ierr);
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown     ierr = ComputeSensiP(lambda[0],mu[0],ctx);CHKERRQ(ierr);
319*c4762a1bSJed Brown     ierr = VecCopy(mu[0],G);CHKERRQ(ierr);
320*c4762a1bSJed Brown   }
321*c4762a1bSJed Brown 
322*c4762a1bSJed Brown   if (ctx->sa == SA_TLM) {
323*c4762a1bSJed Brown     ierr = VecGetArray(G,&x_ptr);CHKERRQ(ierr);
324*c4762a1bSJed Brown     ierr = MatDenseGetArray(qgrad,&y_ptr);CHKERRQ(ierr);
325*c4762a1bSJed Brown     x_ptr[0] = y_ptr[0]-1.;
326*c4762a1bSJed Brown     ierr = MatDenseRestoreArray(qgrad,&y_ptr);CHKERRQ(ierr);
327*c4762a1bSJed Brown     ierr = VecRestoreArray(G,&x_ptr);CHKERRQ(ierr);
328*c4762a1bSJed Brown   }
329*c4762a1bSJed Brown 
330*c4762a1bSJed Brown   ierr = TSGetSolution(ctx->quadts,&q);CHKERRQ(ierr);
331*c4762a1bSJed Brown   ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr);
332*c4762a1bSJed Brown   *f   = -ctx->Pm + x_ptr[0];
333*c4762a1bSJed Brown   ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr);
334*c4762a1bSJed Brown   return 0;
335*c4762a1bSJed Brown }
336*c4762a1bSJed Brown 
337*c4762a1bSJed Brown /*TEST
338*c4762a1bSJed Brown 
339*c4762a1bSJed Brown    build:
340*c4762a1bSJed Brown       requires: !complex !single
341*c4762a1bSJed Brown 
342*c4762a1bSJed Brown    test:
343*c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_type cn -pc_type lu -tao_monitor
344*c4762a1bSJed Brown 
345*c4762a1bSJed Brown    test:
346*c4762a1bSJed Brown       suffix: 2
347*c4762a1bSJed Brown       output_file: output/ex3opt_1.out
348*c4762a1bSJed Brown       args: -sa_method tlm -ts_type cn -pc_type lu -tao_monitor
349*c4762a1bSJed Brown TEST*/
350