1c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*F 4c4762a1bSJed Brown 5c4762a1bSJed Brown \begin{eqnarray} 6c4762a1bSJed Brown \frac{d \theta}{dt} = \omega_b (\omega - \omega_s) 7c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\ 8c4762a1bSJed Brown \end{eqnarray} 9c4762a1bSJed Brown 10c4762a1bSJed Brown F*/ 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown This code demonstrates how to solve a ODE-constrained optimization problem with TAO, TSEvent, TSAdjoint and TS. 14c4762a1bSJed Brown The problem features discontinuities and a cost function in integral form. 15c4762a1bSJed Brown The gradient is computed with the discrete adjoint of an implicit theta method, see ex3adj.c for details. 16c4762a1bSJed Brown */ 17c4762a1bSJed Brown 18c4762a1bSJed Brown #include <petsctao.h> 19c4762a1bSJed Brown #include <petscts.h> 20c4762a1bSJed Brown #include "ex3.h" 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao, Vec, PetscReal *, Vec, void *); 23c4762a1bSJed Brown 24d71ae5a4SJacob Faibussowitsch PetscErrorCode monitor(Tao tao, AppCtx *ctx) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown FILE *fp; 27c4762a1bSJed Brown PetscInt iterate; 28c4762a1bSJed Brown PetscReal f, gnorm, cnorm, xdiff; 29c4762a1bSJed Brown TaoConvergedReason reason; 30c4762a1bSJed Brown 31c4762a1bSJed Brown PetscFunctionBeginUser; 329566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao, &iterate, &f, &gnorm, &cnorm, &xdiff, &reason)); 33c4762a1bSJed Brown 34c4762a1bSJed Brown fp = fopen("ex3opt_conv.out", "a"); 3563a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fp, "%" PetscInt_FMT " %g\n", iterate, (double)gnorm)); 36c4762a1bSJed Brown fclose(fp); 373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 41d71ae5a4SJacob Faibussowitsch { 42c4762a1bSJed Brown Vec p; 43c4762a1bSJed Brown PetscScalar *x_ptr; 44c4762a1bSJed Brown PetscMPIInt size; 45c4762a1bSJed Brown AppCtx ctx; 46c4762a1bSJed Brown Tao tao; 47c4762a1bSJed Brown KSP ksp; 48c4762a1bSJed Brown PC pc; 49c4762a1bSJed Brown Vec lambda[1], mu[1], lowerb, upperb; 50c4762a1bSJed Brown PetscBool printtofile; 51c4762a1bSJed Brown PetscInt direction[2]; 52c4762a1bSJed Brown PetscBool terminate[2]; 53c4762a1bSJed Brown Mat qgrad; /* Forward sesivitiy */ 54c4762a1bSJed Brown Mat sp; /* Forward sensitivity matrix */ 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 57c4762a1bSJed Brown Initialize program 58c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 59327415f7SBarry Smith PetscFunctionBeginUser; 609566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 623c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65c4762a1bSJed Brown Set runtime options 66c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 68c4762a1bSJed Brown { 69c4762a1bSJed Brown ctx.beta = 2; 70c4762a1bSJed Brown ctx.c = 10000.0; 71c4762a1bSJed Brown ctx.u_s = 1.0; 72c4762a1bSJed Brown ctx.omega_s = 1.0; 73c4762a1bSJed Brown ctx.omega_b = 120.0 * PETSC_PI; 74c4762a1bSJed Brown ctx.H = 5.0; 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 76c4762a1bSJed Brown ctx.D = 5.0; 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 78c4762a1bSJed Brown ctx.E = 1.1378; 79c4762a1bSJed Brown ctx.V = 1.0; 80c4762a1bSJed Brown ctx.X = 0.545; 81c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 82c4762a1bSJed Brown ctx.Pmax_ini = ctx.Pmax; 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 84c4762a1bSJed Brown ctx.Pm = 1.06; 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 86c4762a1bSJed Brown ctx.tf = 0.1; 87c4762a1bSJed Brown ctx.tcl = 0.2; 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 90c4762a1bSJed Brown printtofile = PETSC_FALSE; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-printtofile", "Print convergence results to file", "", printtofile, &printtofile, NULL)); 92c4762a1bSJed Brown ctx.sa = SA_ADJ; 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sa_method", "Sensitivity analysis method (adj or tlm)", "", SAMethods, (PetscEnum)ctx.sa, (PetscEnum *)&ctx.sa, NULL)); 94c4762a1bSJed Brown } 95d0609cedSBarry Smith PetscOptionsEnd(); 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 98c4762a1bSJed Brown Create necessary matrix and vectors 99c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1009566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &ctx.Jac)); 1019566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx.Jac, 2, 2, PETSC_DETERMINE, PETSC_DETERMINE)); 1029566063dSJacob Faibussowitsch PetscCall(MatSetType(ctx.Jac, MATDENSE)); 1039566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ctx.Jac)); 1049566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.Jac)); 1059566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &ctx.Jacp)); 1069566063dSJacob Faibussowitsch PetscCall(MatSetSizes(ctx.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1)); 1079566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(ctx.Jacp)); 1089566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.Jacp)); 1099566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx.Jac, &ctx.U, NULL)); 1109566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &ctx.DRDP)); 1119566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.DRDP)); 1129566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, &ctx.DRDU)); 1139566063dSJacob Faibussowitsch PetscCall(MatSetUp(ctx.DRDU)); 114c4762a1bSJed Brown 115c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 116c4762a1bSJed Brown Create timestepping solver context 117c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1189566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ctx.ts)); 1199566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ctx.ts, TS_NONLINEAR)); 1209566063dSJacob Faibussowitsch PetscCall(TSSetType(ctx.ts, TSCN)); 1218434afd1SBarry Smith PetscCall(TSSetRHSFunction(ctx.ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx)); 1228434afd1SBarry Smith PetscCall(TSSetRHSJacobian(ctx.ts, ctx.Jac, ctx.Jac, (TSRHSJacobianFn *)RHSJacobian, &ctx)); 1239566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ctx.ts, ctx.Jacp, RHSJacobianP, &ctx)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown if (ctx.sa == SA_ADJ) { 1269566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx.Jac, &lambda[0], NULL)); 1279566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(ctx.Jacp, &mu[0], NULL)); 1289566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ctx.ts)); 1299566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ctx.ts, 1, lambda, mu)); 1309566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ctx.ts, PETSC_FALSE, &ctx.quadts)); 1318434afd1SBarry Smith PetscCall(TSSetRHSFunction(ctx.quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &ctx)); 1328434afd1SBarry Smith PetscCall(TSSetRHSJacobian(ctx.quadts, ctx.DRDU, ctx.DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, &ctx)); 1339566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ctx.quadts, ctx.DRDP, DRDPJacobianTranspose, &ctx)); 134c4762a1bSJed Brown } 135c4762a1bSJed Brown if (ctx.sa == SA_TLM) { 1369566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &qgrad)); 1379566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 1, NULL, &sp)); 1389566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ctx.ts, 1, sp)); 1399566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ctx.ts, PETSC_TRUE, &ctx.quadts)); 1409566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ctx.quadts, 1, qgrad)); 1418434afd1SBarry Smith PetscCall(TSSetRHSFunction(ctx.quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &ctx)); 1428434afd1SBarry Smith PetscCall(TSSetRHSJacobian(ctx.quadts, ctx.DRDU, ctx.DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, &ctx)); 1439566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ctx.quadts, ctx.DRDP, DRDPJacobianTranspose, &ctx)); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Set solver options 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1499566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ctx.ts, 1.0)); 1509566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ctx.ts, TS_EXACTFINALTIME_MATCHSTEP)); 1519566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ctx.ts, 0.03125)); 1529566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ctx.ts)); 153c4762a1bSJed Brown 154c4762a1bSJed Brown direction[0] = direction[1] = 1; 155c4762a1bSJed Brown terminate[0] = terminate[1] = PETSC_FALSE; 1569566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ctx.ts, 2, direction, terminate, EventFunction, PostEventFunction, &ctx)); 157c4762a1bSJed Brown 158c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 1599566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 1609566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM)); 16110978b7dSBarry Smith if (printtofile) PetscCall(TaoMonitorSet(tao, (PetscErrorCode (*)(Tao, void *))monitor, (void *)&ctx, NULL)); 162c4762a1bSJed Brown /* 163c4762a1bSJed Brown Optimization starts 164c4762a1bSJed Brown */ 165c4762a1bSJed Brown /* Set initial solution guess */ 1669566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p)); 1679566063dSJacob Faibussowitsch PetscCall(VecGetArray(p, &x_ptr)); 168c4762a1bSJed Brown x_ptr[0] = ctx.Pm; 1699566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(p, &x_ptr)); 170c4762a1bSJed Brown 1719566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, p)); 172c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 1739566063dSJacob Faibussowitsch PetscCall(TaoSetObjectiveAndGradient(tao, NULL, FormFunctionGradient, (void *)&ctx)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* Set bounds for the optimization */ 1769566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &lowerb)); 1779566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &upperb)); 1789566063dSJacob Faibussowitsch PetscCall(VecGetArray(lowerb, &x_ptr)); 179c4762a1bSJed Brown x_ptr[0] = 0.; 1809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lowerb, &x_ptr)); 1819566063dSJacob Faibussowitsch PetscCall(VecGetArray(upperb, &x_ptr)); 182c4762a1bSJed Brown x_ptr[0] = 1.1; 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(upperb, &x_ptr)); 1849566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, lowerb, upperb)); 185c4762a1bSJed Brown 186c4762a1bSJed Brown /* Check for any TAO command line options */ 1879566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 1889566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 189c4762a1bSJed Brown if (ksp) { 1909566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1919566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 192c4762a1bSJed Brown } 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1959566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 196c4762a1bSJed Brown 1979566063dSJacob Faibussowitsch PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD)); 198c4762a1bSJed Brown 199c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 200c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 201c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2029566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.Jac)); 2039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.Jacp)); 2049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.DRDU)); 2059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ctx.DRDP)); 2069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.U)); 207c4762a1bSJed Brown if (ctx.sa == SA_ADJ) { 2089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 2099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown if (ctx.sa == SA_TLM) { 2129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&qgrad)); 2139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sp)); 214c4762a1bSJed Brown } 2159566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ctx.ts)); 2169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p)); 2179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lowerb)); 2189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&upperb)); 2199566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 2209566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 221b122ec5aSJacob Faibussowitsch return 0; 222c4762a1bSJed Brown } 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 225c4762a1bSJed Brown /* 226c4762a1bSJed Brown FormFunctionGradient - Evaluates the function and corresponding gradient. 227c4762a1bSJed Brown 228c4762a1bSJed Brown Input Parameters: 229c4762a1bSJed Brown tao - the Tao context 230c4762a1bSJed Brown X - the input vector 231a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 232c4762a1bSJed Brown 233c4762a1bSJed Brown Output Parameters: 234c4762a1bSJed Brown f - the newly evaluated function 235c4762a1bSJed Brown G - the newly evaluated gradient 236c4762a1bSJed Brown */ 237*2a8381b2SBarry Smith PetscErrorCode FormFunctionGradient(Tao tao, Vec P, PetscReal *f, Vec G, PetscCtx ctx0) 238d71ae5a4SJacob Faibussowitsch { 239c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 240c4762a1bSJed Brown PetscInt nadj; 241c4762a1bSJed Brown PetscReal ftime; 242c4762a1bSJed Brown PetscInt steps; 243c4762a1bSJed Brown PetscScalar *u; 244c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr; 245c4762a1bSJed Brown Vec q; 246c4762a1bSJed Brown Mat qgrad; 247c4762a1bSJed Brown 2483ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 2499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr)); 250c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 2519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr)); 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* reinitialize the solution vector */ 2549566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx->U, &u)); 255c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax); 256c4762a1bSJed Brown u[1] = 1.0; 2579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx->U, &u)); 2589566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ctx->ts, ctx->U)); 259c4762a1bSJed Brown 260c4762a1bSJed Brown /* reset time */ 2619566063dSJacob Faibussowitsch PetscCall(TSSetTime(ctx->ts, 0.0)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* reset step counter, this is critical for adjoint solver */ 2649566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ctx->ts, 0)); 265c4762a1bSJed Brown 266c4762a1bSJed Brown /* reset step size, the step size becomes negative after TSAdjointSolve */ 2679566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ctx->ts, 0.03125)); 268c4762a1bSJed Brown 269c4762a1bSJed Brown /* reinitialize the integral value */ 2709566063dSJacob Faibussowitsch PetscCall(TSGetQuadratureTS(ctx->ts, NULL, &ctx->quadts)); 2719566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ctx->quadts, &q)); 2729566063dSJacob Faibussowitsch PetscCall(VecSet(q, 0.0)); 273c4762a1bSJed Brown 274c4762a1bSJed Brown if (ctx->sa == SA_TLM) { /* reset the forward sensitivities */ 275c4762a1bSJed Brown TS quadts; 276c4762a1bSJed Brown Mat sp; 277c4762a1bSJed Brown PetscScalar val[2]; 278c4762a1bSJed Brown const PetscInt row[] = {0, 1}, col[] = {0}; 279c4762a1bSJed Brown 2809566063dSJacob Faibussowitsch PetscCall(TSGetQuadratureTS(ctx->ts, NULL, &quadts)); 2819566063dSJacob Faibussowitsch PetscCall(TSForwardGetSensitivities(quadts, NULL, &qgrad)); 2829566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(qgrad)); 2839566063dSJacob Faibussowitsch PetscCall(TSForwardGetSensitivities(ctx->ts, NULL, &sp)); 284c4762a1bSJed Brown val[0] = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax; 285c4762a1bSJed Brown val[1] = 0.0; 2869566063dSJacob Faibussowitsch PetscCall(MatSetValues(sp, 2, row, 1, col, val, INSERT_VALUES)); 2879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(sp, MAT_FINAL_ASSEMBLY)); 2889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(sp, MAT_FINAL_ASSEMBLY)); 289c4762a1bSJed Brown } 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* solve the ODE */ 2929566063dSJacob Faibussowitsch PetscCall(TSSolve(ctx->ts, ctx->U)); 2939566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ctx->ts, &ftime)); 2949566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ctx->ts, &steps)); 295c4762a1bSJed Brown 296c4762a1bSJed Brown if (ctx->sa == SA_ADJ) { 297c4762a1bSJed Brown Vec *lambda, *mu; 298c4762a1bSJed Brown /* reset the terminal condition for adjoint */ 2999566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(ctx->ts, &nadj, &lambda, &mu)); 3009566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &y_ptr)); 3019371c9d4SSatish Balay y_ptr[0] = 0.0; 3029371c9d4SSatish Balay y_ptr[1] = 0.0; 3039566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &y_ptr)); 3049566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr)); 305c4762a1bSJed Brown x_ptr[0] = -1.0; 3069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr)); 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* solve the adjont */ 3099566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ctx->ts)); 310c4762a1bSJed Brown 3119566063dSJacob Faibussowitsch PetscCall(ComputeSensiP(lambda[0], mu[0], ctx)); 3129566063dSJacob Faibussowitsch PetscCall(VecCopy(mu[0], G)); 313c4762a1bSJed Brown } 314c4762a1bSJed Brown 315c4762a1bSJed Brown if (ctx->sa == SA_TLM) { 3169566063dSJacob Faibussowitsch PetscCall(VecGetArray(G, &x_ptr)); 3179566063dSJacob Faibussowitsch PetscCall(MatDenseGetArray(qgrad, &y_ptr)); 318c4762a1bSJed Brown x_ptr[0] = y_ptr[0] - 1.; 3199566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(qgrad, &y_ptr)); 3209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(G, &x_ptr)); 321c4762a1bSJed Brown } 322c4762a1bSJed Brown 3239566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ctx->quadts, &q)); 3249566063dSJacob Faibussowitsch PetscCall(VecGetArray(q, &x_ptr)); 325c4762a1bSJed Brown *f = -ctx->Pm + x_ptr[0]; 3269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(q, &x_ptr)); 3273ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 328c4762a1bSJed Brown } 329c4762a1bSJed Brown 330c4762a1bSJed Brown /*TEST 331c4762a1bSJed Brown 332c4762a1bSJed Brown build: 333c4762a1bSJed Brown requires: !complex !single 334c4762a1bSJed Brown 335c4762a1bSJed Brown test: 336c4762a1bSJed Brown args: -viewer_binary_skip_info -ts_type cn -pc_type lu -tao_monitor 337c4762a1bSJed Brown 338c4762a1bSJed Brown test: 339c4762a1bSJed Brown suffix: 2 340c4762a1bSJed Brown output_file: output/ex3opt_1.out 341c4762a1bSJed Brown args: -sa_method tlm -ts_type cn -pc_type lu -tao_monitor 342c4762a1bSJed Brown TEST*/ 343