1 2 static char help[] = "Basic equation for generator stability analysis.\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 Ensemble of initial conditions 12 ./ex9 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly 13 14 Fault at .1 seconds 15 ./ex9 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly 16 17 Initial conditions same as when fault is ended 18 ./ex9 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly 19 20 F*/ 21 22 /* 23 Include "petscts.h" so that we can use TS solvers. Note that this 24 file automatically includes: 25 petscsys.h - base PETSc routines petscvec.h - vectors 26 petscmat.h - matrices 27 petscis.h - index sets petscksp.h - Krylov subspace methods 28 petscviewer.h - viewers petscpc.h - preconditioners 29 petscksp.h - linear solvers 30 */ 31 32 #include <petscts.h> 33 34 typedef struct { 35 PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c; 36 PetscInt beta; 37 PetscReal tf,tcl; 38 } AppCtx; 39 40 PetscErrorCode PostStepFunction(TS ts) 41 { 42 PetscErrorCode ierr; 43 Vec U; 44 PetscReal t; 45 const PetscScalar *u; 46 47 PetscFunctionBegin; 48 ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 49 ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 50 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 51 ierr = PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",(double)t,(double)u[0]);CHKERRQ(ierr); 52 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 53 PetscFunctionReturn(0); 54 } 55 56 /* 57 Defines the ODE passed to the ODE solver 58 */ 59 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 60 { 61 PetscErrorCode ierr; 62 PetscScalar *f,Pmax; 63 const PetscScalar *u; 64 65 PetscFunctionBegin; 66 /* The next three lines allow us to access the entries of the vectors directly */ 67 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 68 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 69 if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */ 70 else Pmax = ctx->Pmax; 71 72 f[0] = ctx->omega_b*(u[1] - ctx->omega_s); 73 f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H); 74 75 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 76 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 77 PetscFunctionReturn(0); 78 } 79 80 /* 81 Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 82 */ 83 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx) 84 { 85 PetscErrorCode ierr; 86 PetscInt rowcol[] = {0,1}; 87 PetscScalar J[2][2],Pmax; 88 const PetscScalar *u; 89 90 PetscFunctionBegin; 91 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 92 if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */ 93 else Pmax = ctx->Pmax; 94 95 J[0][0] = 0; J[0][1] = ctx->omega_b; 96 J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H); J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H); 97 98 ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 99 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 100 101 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 102 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 103 if (A != B) { 104 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 105 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 106 } 107 PetscFunctionReturn(0); 108 } 109 110 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0) 111 { 112 PetscErrorCode ierr; 113 PetscInt row[] = {0,1},col[]={0}; 114 PetscScalar J[2][1]; 115 AppCtx *ctx=(AppCtx*)ctx0; 116 117 PetscFunctionBeginUser; 118 J[0][0] = 0; 119 J[1][0] = ctx->omega_s/(2.0*ctx->H); 120 ierr = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 121 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 122 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 123 PetscFunctionReturn(0); 124 } 125 126 static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx) 127 { 128 PetscErrorCode ierr; 129 PetscScalar *r; 130 const PetscScalar *u; 131 132 PetscFunctionBegin; 133 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 134 ierr = VecGetArray(R,&r);CHKERRQ(ierr); 135 r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);CHKERRQ(ierr); 136 ierr = VecRestoreArray(R,&r);CHKERRQ(ierr); 137 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 138 PetscFunctionReturn(0); 139 } 140 141 static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx) 142 { 143 PetscErrorCode ierr; 144 PetscScalar ru[1]; 145 const PetscScalar *u; 146 PetscInt row[] = {0},col[] = {0}; 147 148 PetscFunctionBegin; 149 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 150 ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);CHKERRQ(ierr); 151 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 152 ierr = MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES);CHKERRQ(ierr); 153 ierr = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 154 ierr = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 155 PetscFunctionReturn(0); 156 } 157 158 static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx) 159 { 160 PetscErrorCode ierr; 161 162 PetscFunctionBegin; 163 ierr = MatZeroEntries(DRDP);CHKERRQ(ierr); 164 ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 165 ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166 PetscFunctionReturn(0); 167 } 168 169 PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx) 170 { 171 PetscErrorCode ierr; 172 PetscScalar sensip; 173 const PetscScalar *x,*y; 174 175 PetscFunctionBegin; 176 ierr = VecGetArrayRead(lambda,&x);CHKERRQ(ierr); 177 ierr = VecGetArrayRead(mu,&y);CHKERRQ(ierr); 178 sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0]; 179 ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);CHKERRQ(ierr); 180 ierr = VecRestoreArrayRead(lambda,&x);CHKERRQ(ierr); 181 ierr = VecRestoreArrayRead(mu,&y);CHKERRQ(ierr); 182 PetscFunctionReturn(0); 183 } 184 185 int main(int argc,char **argv) 186 { 187 TS ts,quadts; /* ODE integrator */ 188 Vec U; /* solution will be stored here */ 189 Mat A; /* Jacobian matrix */ 190 Mat Jacp; /* Jacobian matrix */ 191 Mat DRDU,DRDP; 192 PetscErrorCode ierr; 193 PetscMPIInt size; 194 PetscInt n = 2; 195 AppCtx ctx; 196 PetscScalar *u; 197 PetscReal du[2] = {0.0,0.0}; 198 PetscBool ensemble = PETSC_FALSE,flg1,flg2; 199 PetscReal ftime; 200 PetscInt steps; 201 PetscScalar *x_ptr,*y_ptr; 202 Vec lambda[1],q,mu[1]; 203 204 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 205 Initialize program 206 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 207 ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 208 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 209 if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 210 211 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 212 Create necessary matrix and vectors 213 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 214 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 215 ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 216 ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 217 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 218 ierr = MatSetUp(A);CHKERRQ(ierr); 219 220 ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 221 222 ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr); 223 ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr); 224 ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr); 225 ierr = MatSetUp(Jacp);CHKERRQ(ierr); 226 227 ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP);CHKERRQ(ierr); 228 ierr = MatSetUp(DRDP);CHKERRQ(ierr); 229 ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU);CHKERRQ(ierr); 230 ierr = MatSetUp(DRDU);CHKERRQ(ierr); 231 232 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 233 Set runtime options 234 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 235 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr); 236 { 237 ctx.beta = 2; 238 ctx.c = 10000.0; 239 ctx.u_s = 1.0; 240 ctx.omega_s = 1.0; 241 ctx.omega_b = 120.0*PETSC_PI; 242 ctx.H = 5.0; 243 ierr = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr); 244 ctx.D = 5.0; 245 ierr = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr); 246 ctx.E = 1.1378; 247 ctx.V = 1.0; 248 ctx.X = 0.545; 249 ctx.Pmax = ctx.E*ctx.V/ctx.X; 250 ierr = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr); 251 ctx.Pm = 1.1; 252 ierr = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr); 253 ctx.tf = 0.1; 254 ctx.tcl = 0.2; 255 ierr = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr); 256 ierr = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr); 257 ierr = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr); 258 if (ensemble) { 259 ctx.tf = -1; 260 ctx.tcl = -1; 261 } 262 263 ierr = VecGetArray(U,&u);CHKERRQ(ierr); 264 u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 265 u[1] = 1.0; 266 ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr); 267 n = 2; 268 ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr); 269 u[0] += du[0]; 270 u[1] += du[1]; 271 ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 272 if (flg1 || flg2) { 273 ctx.tf = -1; 274 ctx.tcl = -1; 275 } 276 } 277 ierr = PetscOptionsEnd();CHKERRQ(ierr); 278 279 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 280 Create timestepping solver context 281 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 282 ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 283 ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 284 ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 285 ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); 286 ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr); 287 ierr = TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr); 288 ierr = TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts);CHKERRQ(ierr); 289 ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr); 290 ierr = TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr); 291 ierr = TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx);CHKERRQ(ierr); 292 293 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 294 Set initial conditions 295 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 296 ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 297 298 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 299 Save trajectory of solution so that TSAdjointSolve() may be used 300 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 301 ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 302 303 ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr); 304 /* Set initial conditions for the adjoint integration */ 305 ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr); 306 y_ptr[0] = 0.0; y_ptr[1] = 0.0; 307 ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr); 308 309 ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr); 310 ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 311 x_ptr[0] = -1.0; 312 ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 313 ierr = TSSetCostGradients(ts,1,lambda,mu);CHKERRQ(ierr); 314 315 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 316 Set solver options 317 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 318 ierr = TSSetMaxTime(ts,10.0);CHKERRQ(ierr); 319 ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 320 ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr); 321 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 322 323 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 324 Solve nonlinear system 325 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 326 if (ensemble) { 327 for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 328 ierr = VecGetArray(U,&u);CHKERRQ(ierr); 329 u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 330 u[1] = ctx.omega_s; 331 u[0] += du[0]; 332 u[1] += du[1]; 333 ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 334 ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr); 335 ierr = TSSolve(ts,U);CHKERRQ(ierr); 336 } 337 } else { 338 ierr = TSSolve(ts,U);CHKERRQ(ierr); 339 } 340 ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 341 ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 342 ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 343 344 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 345 Adjoint model starts here 346 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 347 /* Set initial conditions for the adjoint integration */ 348 ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr); 349 y_ptr[0] = 0.0; y_ptr[1] = 0.0; 350 ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr); 351 352 ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 353 x_ptr[0] = -1.0; 354 ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 355 356 /* Set RHS JacobianP */ 357 ierr = TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr); 358 359 ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 360 361 ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0] d[Psi(tf)]/d[omega0]\n");CHKERRQ(ierr); 362 ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 363 ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 364 ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 365 ierr = VecView(q,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 366 ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr); 367 ierr = PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));CHKERRQ(ierr); 368 ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr); 369 370 ierr = ComputeSensiP(lambda[0],mu[0],&ctx);CHKERRQ(ierr); 371 372 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 373 Free work space. All PETSc objects should be destroyed when they are no longer needed. 374 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 375 ierr = MatDestroy(&A);CHKERRQ(ierr); 376 ierr = MatDestroy(&Jacp);CHKERRQ(ierr); 377 ierr = MatDestroy(&DRDU);CHKERRQ(ierr); 378 ierr = MatDestroy(&DRDP);CHKERRQ(ierr); 379 ierr = VecDestroy(&U);CHKERRQ(ierr); 380 ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 381 ierr = VecDestroy(&mu[0]);CHKERRQ(ierr); 382 ierr = TSDestroy(&ts);CHKERRQ(ierr); 383 ierr = PetscFinalize(); 384 return ierr; 385 } 386 387 /*TEST 388 389 build: 390 requires: !complex 391 392 test: 393 args: -viewer_binary_skip_info -ts_adapt_type none 394 395 TEST*/ 396