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 12 13 Ensemble of initial conditions 14 ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 15 16 Fault at .1 seconds 17 ./ex2 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 18 19 Initial conditions same as when fault is ended 20 ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 21 22 23 F*/ 24 25 /* 26 Include "petscts.h" so that we can use TS solvers. Note that this 27 file automatically includes: 28 petscsys.h - base PETSc routines petscvec.h - vectors 29 petscmat.h - matrices 30 petscis.h - index sets petscksp.h - Krylov subspace methods 31 petscviewer.h - viewers petscpc.h - preconditioners 32 petscksp.h - linear solvers 33 */ 34 35 #include <petsctao.h> 36 #include <petscts.h> 37 38 typedef struct { 39 TS ts; 40 PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c; 41 PetscInt beta; 42 PetscReal tf,tcl,dt; 43 } AppCtx; 44 45 PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*); 46 PetscErrorCode FormGradient(Tao,Vec,Vec,void*); 47 48 /* 49 Defines the ODE passed to the ODE solver 50 */ 51 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 52 { 53 PetscErrorCode ierr; 54 PetscScalar *f,Pmax; 55 const PetscScalar *u; 56 57 PetscFunctionBegin; 58 /* The next three lines allow us to access the entries of the vectors directly */ 59 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 60 ierr = VecGetArray(F,&f);CHKERRQ(ierr); 61 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 */ 62 else Pmax = ctx->Pmax; 63 64 f[0] = ctx->omega_b*(u[1] - ctx->omega_s); 65 f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H); 66 67 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 68 ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 69 PetscFunctionReturn(0); 70 } 71 72 /* 73 Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 74 */ 75 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx) 76 { 77 PetscErrorCode ierr; 78 PetscInt rowcol[] = {0,1}; 79 PetscScalar J[2][2],Pmax; 80 const PetscScalar *u; 81 82 PetscFunctionBegin; 83 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 84 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 */ 85 else Pmax = ctx->Pmax; 86 87 J[0][0] = 0; J[0][1] = ctx->omega_b; 88 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); 89 90 ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 91 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 92 93 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 95 if (A != B) { 96 ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 97 ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 98 } 99 PetscFunctionReturn(0); 100 } 101 102 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0) 103 { 104 PetscErrorCode ierr; 105 PetscInt row[] = {0,1},col[]={0}; 106 PetscScalar J[2][1]; 107 AppCtx *ctx=(AppCtx*)ctx0; 108 109 PetscFunctionBeginUser; 110 J[0][0] = 0; 111 J[1][0] = ctx->omega_s/(2.0*ctx->H); 112 ierr = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 113 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 114 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 115 PetscFunctionReturn(0); 116 } 117 118 static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx) 119 { 120 PetscErrorCode ierr; 121 PetscScalar *r; 122 const PetscScalar *u; 123 124 PetscFunctionBegin; 125 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 126 ierr = VecGetArray(R,&r);CHKERRQ(ierr); 127 r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);CHKERRQ(ierr); 128 ierr = VecRestoreArray(R,&r);CHKERRQ(ierr); 129 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 130 PetscFunctionReturn(0); 131 } 132 133 static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx) 134 { 135 PetscErrorCode ierr; 136 PetscScalar ru[1]; 137 const PetscScalar *u; 138 PetscInt row[] = {0},col[] = {0}; 139 140 PetscFunctionBegin; 141 ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 142 ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);CHKERRQ(ierr); 143 ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 144 ierr = MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES);CHKERRQ(ierr); 145 ierr = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 146 ierr = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 147 PetscFunctionReturn(0); 148 } 149 150 static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx) 151 { 152 PetscErrorCode ierr; 153 154 PetscFunctionBegin; 155 ierr = MatZeroEntries(DRDP);CHKERRQ(ierr); 156 ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157 ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 158 PetscFunctionReturn(0); 159 } 160 161 PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx) 162 { 163 PetscErrorCode ierr; 164 PetscScalar *y,sensip; 165 const PetscScalar *x; 166 167 PetscFunctionBegin; 168 ierr = VecGetArrayRead(lambda,&x);CHKERRQ(ierr); 169 ierr = VecGetArray(mu,&y);CHKERRQ(ierr); 170 sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0]; 171 y[0] = sensip; 172 ierr = VecRestoreArray(mu,&y);CHKERRQ(ierr); 173 ierr = VecRestoreArrayRead(lambda,&x);CHKERRQ(ierr); 174 PetscFunctionReturn(0); 175 } 176 177 int main(int argc,char **argv) 178 { 179 Vec p; 180 PetscScalar *x_ptr; 181 PetscErrorCode ierr; 182 PetscMPIInt size; 183 AppCtx ctx; 184 Vec lowerb,upperb; 185 Tao tao; 186 KSP ksp; 187 PC pc; 188 Vec U,lambda[1],mu[1]; 189 Mat A; /* Jacobian matrix */ 190 Mat Jacp; /* Jacobian matrix */ 191 Mat DRDU,DRDP; 192 PetscInt n = 2; 193 TS quadts; 194 195 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 196 Initialize program 197 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 198 ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 199 PetscFunctionBeginUser; 200 ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 201 if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 202 203 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204 Set runtime options 205 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206 ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr); 207 { 208 ctx.beta = 2; 209 ctx.c = PetscRealConstant(10000.0); 210 ctx.u_s = PetscRealConstant(1.0); 211 ctx.omega_s = PetscRealConstant(1.0); 212 ctx.omega_b = PetscRealConstant(120.0)*PETSC_PI; 213 ctx.H = PetscRealConstant(5.0); 214 ierr = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr); 215 ctx.D = PetscRealConstant(5.0); 216 ierr = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr); 217 ctx.E = PetscRealConstant(1.1378); 218 ctx.V = PetscRealConstant(1.0); 219 ctx.X = PetscRealConstant(0.545); 220 ctx.Pmax = ctx.E*ctx.V/ctx.X; 221 ierr = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr); 222 ctx.Pm = PetscRealConstant(1.0194); 223 ierr = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr); 224 ctx.tf = PetscRealConstant(0.1); 225 ctx.tcl = PetscRealConstant(0.2); 226 ierr = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr); 227 ierr = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr); 228 229 } 230 ierr = PetscOptionsEnd();CHKERRQ(ierr); 231 232 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 233 Create necessary matrix and vectors 234 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 235 ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 236 ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 237 ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 238 ierr = MatSetFromOptions(A);CHKERRQ(ierr); 239 ierr = MatSetUp(A);CHKERRQ(ierr); 240 241 ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 242 243 ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr); 244 ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr); 245 ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr); 246 ierr = MatSetUp(Jacp);CHKERRQ(ierr); 247 ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP);CHKERRQ(ierr); 248 ierr = MatSetUp(DRDP);CHKERRQ(ierr); 249 ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU);CHKERRQ(ierr); 250 ierr = MatSetUp(DRDU);CHKERRQ(ierr); 251 252 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 253 Create timestepping solver context 254 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 255 ierr = TSCreate(PETSC_COMM_WORLD,&ctx.ts);CHKERRQ(ierr); 256 ierr = TSSetProblemType(ctx.ts,TS_NONLINEAR);CHKERRQ(ierr); 257 ierr = TSSetEquationType(ctx.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 258 ierr = TSSetType(ctx.ts,TSRK);CHKERRQ(ierr); 259 ierr = TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr); 260 ierr = TSSetRHSJacobian(ctx.ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr); 261 ierr = TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 262 263 ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr); 264 ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr); 265 ierr = TSSetCostGradients(ctx.ts,1,lambda,mu);CHKERRQ(ierr); 266 ierr = TSSetRHSJacobianP(ctx.ts,Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr); 267 268 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 269 Set solver options 270 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 271 ierr = TSSetMaxTime(ctx.ts,PetscRealConstant(1.0));CHKERRQ(ierr); 272 ierr = TSSetTimeStep(ctx.ts,PetscRealConstant(0.01));CHKERRQ(ierr); 273 ierr = TSSetFromOptions(ctx.ts);CHKERRQ(ierr); 274 275 ierr = TSGetTimeStep(ctx.ts,&ctx.dt);CHKERRQ(ierr); /* save the stepsize */ 276 277 ierr = TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&quadts);CHKERRQ(ierr); 278 ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr); 279 ierr = TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr); 280 ierr = TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx);CHKERRQ(ierr); 281 ierr = TSSetSolution(ctx.ts,U);CHKERRQ(ierr); 282 283 /* Create TAO solver and set desired solution method */ 284 ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 285 ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr); 286 287 /* 288 Optimization starts 289 */ 290 /* Set initial solution guess */ 291 ierr = VecCreateSeq(PETSC_COMM_WORLD,1,&p);CHKERRQ(ierr); 292 ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr); 293 x_ptr[0] = ctx.Pm; 294 ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr); 295 296 ierr = TaoSetInitialVector(tao,p);CHKERRQ(ierr); 297 /* Set routine for function and gradient evaluation */ 298 ierr = TaoSetObjectiveRoutine(tao,FormFunction,(void *)&ctx);CHKERRQ(ierr); 299 ierr = TaoSetGradientRoutine(tao,FormGradient,(void *)&ctx);CHKERRQ(ierr); 300 301 /* Set bounds for the optimization */ 302 ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr); 303 ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr); 304 ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr); 305 x_ptr[0] = 0.; 306 ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr); 307 ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr); 308 x_ptr[0] = PetscRealConstant(1.1); 309 ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr); 310 ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr); 311 312 /* Check for any TAO command line options */ 313 ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 314 ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 315 if (ksp) { 316 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 317 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 318 } 319 320 /* SOLVE THE APPLICATION */ 321 ierr = TaoSolve(tao);CHKERRQ(ierr); 322 323 ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 324 /* Free TAO data structures */ 325 ierr = TaoDestroy(&tao);CHKERRQ(ierr); 326 ierr = VecDestroy(&p);CHKERRQ(ierr); 327 ierr = VecDestroy(&lowerb);CHKERRQ(ierr); 328 ierr = VecDestroy(&upperb);CHKERRQ(ierr); 329 330 ierr = TSDestroy(&ctx.ts);CHKERRQ(ierr); 331 ierr = VecDestroy(&U);CHKERRQ(ierr); 332 ierr = MatDestroy(&A);CHKERRQ(ierr); 333 ierr = MatDestroy(&Jacp);CHKERRQ(ierr); 334 ierr = MatDestroy(&DRDU);CHKERRQ(ierr); 335 ierr = MatDestroy(&DRDP);CHKERRQ(ierr); 336 ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); 337 ierr = VecDestroy(&mu[0]);CHKERRQ(ierr); 338 ierr = PetscFinalize(); 339 return ierr; 340 } 341 342 /* ------------------------------------------------------------------ */ 343 /* 344 FormFunction - Evaluates the function 345 346 Input Parameters: 347 tao - the Tao context 348 X - the input vector 349 ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine() 350 351 Output Parameters: 352 f - the newly evaluated function 353 */ 354 PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0) 355 { 356 AppCtx *ctx = (AppCtx*)ctx0; 357 TS ts = ctx->ts; 358 Vec U; /* solution will be stored here */ 359 PetscErrorCode ierr; 360 PetscScalar *u; 361 PetscScalar *x_ptr; 362 Vec q; 363 364 ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr); 365 ctx->Pm = x_ptr[0]; 366 ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr); 367 368 /* reset time */ 369 ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 370 /* reset step counter, this is critical for adjoint solver */ 371 ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 372 /* reset step size, the step size becomes negative after TSAdjointSolve */ 373 ierr = TSSetTimeStep(ts,ctx->dt);CHKERRQ(ierr); 374 /* reinitialize the integral value */ 375 ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 376 ierr = VecSet(q,0.0);CHKERRQ(ierr); 377 378 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 379 Set initial conditions 380 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 381 ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 382 ierr = VecGetArray(U,&u);CHKERRQ(ierr); 383 u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax); 384 u[1] = PetscRealConstant(1.0); 385 ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 386 387 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 388 Solve nonlinear system 389 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 390 ierr = TSSolve(ts,U);CHKERRQ(ierr); 391 ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 392 ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr); 393 *f = -ctx->Pm + x_ptr[0]; 394 ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr); 395 return 0; 396 } 397 398 PetscErrorCode FormGradient(Tao tao,Vec P,Vec G,void *ctx0) 399 { 400 AppCtx *ctx = (AppCtx*)ctx0; 401 TS ts = ctx->ts; 402 Vec U; /* solution will be stored here */ 403 PetscErrorCode ierr; 404 PetscReal ftime; 405 PetscInt steps; 406 PetscScalar *u; 407 PetscScalar *x_ptr,*y_ptr; 408 Vec *lambda,q,*mu; 409 410 ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr); 411 ctx->Pm = x_ptr[0]; 412 ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr); 413 414 /* reset time */ 415 ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 416 /* reset step counter, this is critical for adjoint solver */ 417 ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 418 /* reset step size, the step size becomes negative after TSAdjointSolve */ 419 ierr = TSSetTimeStep(ts,ctx->dt);CHKERRQ(ierr); 420 /* reinitialize the integral value */ 421 ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 422 ierr = VecSet(q,0.0);CHKERRQ(ierr); 423 424 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 425 Set initial conditions 426 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 427 ierr = TSGetSolution(ts,&U);CHKERRQ(ierr); 428 ierr = VecGetArray(U,&u);CHKERRQ(ierr); 429 u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax); 430 u[1] = PetscRealConstant(1.0); 431 ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 432 433 /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */ 434 ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 435 ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 436 437 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 438 Solve nonlinear system 439 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 440 ierr = TSSolve(ts,U);CHKERRQ(ierr); 441 442 ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); 443 ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); 444 445 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 446 Adjoint model starts here 447 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 448 ierr = TSGetCostGradients(ts,NULL,&lambda,&mu);CHKERRQ(ierr); 449 /* Set initial conditions for the adjoint integration */ 450 ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr); 451 y_ptr[0] = 0.0; y_ptr[1] = 0.0; 452 ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr); 453 ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); 454 x_ptr[0] = PetscRealConstant(-1.0); 455 ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); 456 457 ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 458 ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr); 459 ierr = ComputeSensiP(lambda[0],mu[0],ctx);CHKERRQ(ierr); 460 ierr = VecCopy(mu[0],G);CHKERRQ(ierr); 461 return 0; 462 } 463 464 465 /*TEST 466 467 build: 468 requires: !complex 469 470 test: 471 args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason 472 473 test: 474 suffix: 2 475 args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason -tao_test_gradient 476 477 TEST*/ 478