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