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