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