1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\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 Ensemble of initial conditions 12c4762a1bSJed Brown ./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 13c4762a1bSJed Brown 14c4762a1bSJed Brown Fault at .1 seconds 15c4762a1bSJed Brown ./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 16c4762a1bSJed Brown 17c4762a1bSJed Brown Initial conditions same as when fault is ended 18c4762a1bSJed Brown ./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 19c4762a1bSJed Brown 20c4762a1bSJed Brown F*/ 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* 23c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 24c4762a1bSJed Brown file automatically includes: 25c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 26c4762a1bSJed Brown petscmat.h - matrices 27c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 28c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 29c4762a1bSJed Brown petscksp.h - linear solvers 30c4762a1bSJed Brown */ 31c4762a1bSJed Brown 32c4762a1bSJed Brown #include <petscts.h> 33c4762a1bSJed Brown 34c4762a1bSJed Brown typedef struct { 35c4762a1bSJed Brown PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c; 36c4762a1bSJed Brown PetscInt beta; 37c4762a1bSJed Brown PetscReal tf, tcl; 38c4762a1bSJed Brown } AppCtx; 39c4762a1bSJed Brown 40d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStepFunction(TS ts) 41d71ae5a4SJacob Faibussowitsch { 42c4762a1bSJed Brown Vec U; 43c4762a1bSJed Brown PetscReal t; 44c4762a1bSJed Brown const PetscScalar *u; 45c4762a1bSJed Brown 46c4762a1bSJed Brown PetscFunctionBegin; 479566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 489566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 509566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, "delta(%3.2f) = %8.7f\n", (double)t, (double)u[0])); 519566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 52*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 53c4762a1bSJed Brown } 54c4762a1bSJed Brown 55c4762a1bSJed Brown /* 56c4762a1bSJed Brown Defines the ODE passed to the ODE solver 57c4762a1bSJed Brown */ 58d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 59d71ae5a4SJacob Faibussowitsch { 60c4762a1bSJed Brown PetscScalar *f, Pmax; 61c4762a1bSJed Brown const PetscScalar *u; 62c4762a1bSJed Brown 63c4762a1bSJed Brown PetscFunctionBegin; 64c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 669566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 67c4762a1bSJed Brown 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 */ 68c4762a1bSJed Brown else Pmax = ctx->Pmax; 69c4762a1bSJed Brown 70c4762a1bSJed Brown f[0] = ctx->omega_b * (u[1] - ctx->omega_s); 71c4762a1bSJed Brown f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H); 72c4762a1bSJed Brown 739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 749566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 75*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76c4762a1bSJed Brown } 77c4762a1bSJed Brown 78c4762a1bSJed Brown /* 79c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 80c4762a1bSJed Brown */ 81d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) 82d71ae5a4SJacob Faibussowitsch { 83c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 84c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 85c4762a1bSJed Brown const PetscScalar *u; 86c4762a1bSJed Brown 87c4762a1bSJed Brown PetscFunctionBegin; 889566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 89c4762a1bSJed Brown 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 */ 90c4762a1bSJed Brown else Pmax = ctx->Pmax; 91c4762a1bSJed Brown 929371c9d4SSatish Balay J[0][0] = 0; 939371c9d4SSatish Balay J[0][1] = ctx->omega_b; 949371c9d4SSatish Balay J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H); 959371c9d4SSatish Balay J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H); 96c4762a1bSJed Brown 979566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 989566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 99c4762a1bSJed Brown 1009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 102c4762a1bSJed Brown if (A != B) { 1039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 105c4762a1bSJed Brown } 106*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 107c4762a1bSJed Brown } 108c4762a1bSJed Brown 109d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0) 110d71ae5a4SJacob Faibussowitsch { 111c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 112c4762a1bSJed Brown PetscScalar J[2][1]; 113c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 114c4762a1bSJed Brown 115c4762a1bSJed Brown PetscFunctionBeginUser; 116c4762a1bSJed Brown J[0][0] = 0; 117c4762a1bSJed Brown J[1][0] = ctx->omega_s / (2.0 * ctx->H); 1189566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 1199566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1209566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 121*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122c4762a1bSJed Brown } 123c4762a1bSJed Brown 124d71ae5a4SJacob Faibussowitsch static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) 125d71ae5a4SJacob Faibussowitsch { 126c4762a1bSJed Brown PetscScalar *r; 127c4762a1bSJed Brown const PetscScalar *u; 128c4762a1bSJed Brown 129c4762a1bSJed Brown PetscFunctionBegin; 1309566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1319566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 1322f613bf5SBarry Smith r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta); 1339566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 1349566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 135*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 136c4762a1bSJed Brown } 137c4762a1bSJed Brown 138d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) 139d71ae5a4SJacob Faibussowitsch { 140c4762a1bSJed Brown PetscScalar ru[1]; 141c4762a1bSJed Brown const PetscScalar *u; 142c4762a1bSJed Brown PetscInt row[] = {0}, col[] = {0}; 143c4762a1bSJed Brown 144c4762a1bSJed Brown PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1462f613bf5SBarry Smith ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1); 1479566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1489566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES)); 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY)); 1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY)); 151*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx) 155d71ae5a4SJacob Faibussowitsch { 156c4762a1bSJed Brown PetscFunctionBegin; 1579566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 1589566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY)); 1599566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY)); 160*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 161c4762a1bSJed Brown } 162c4762a1bSJed Brown 163d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) 164d71ae5a4SJacob Faibussowitsch { 165c4762a1bSJed Brown PetscScalar sensip; 166c4762a1bSJed Brown const PetscScalar *x, *y; 167c4762a1bSJed Brown 168c4762a1bSJed Brown PetscFunctionBegin; 1699566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(lambda, &x)); 1709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(mu, &y)); 171c4762a1bSJed Brown sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0]; 1729566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt parameter pm: %.7f \n", (double)sensip)); 1739566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(lambda, &x)); 1749566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(mu, &y)); 175*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176c4762a1bSJed Brown } 177c4762a1bSJed Brown 178d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 179d71ae5a4SJacob Faibussowitsch { 180c4762a1bSJed Brown TS ts, quadts; /* ODE integrator */ 181c4762a1bSJed Brown Vec U; /* solution will be stored here */ 182c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 183c4762a1bSJed Brown Mat Jacp; /* Jacobian matrix */ 184c4762a1bSJed Brown Mat DRDU, DRDP; 185c4762a1bSJed Brown PetscMPIInt size; 186c4762a1bSJed Brown PetscInt n = 2; 187c4762a1bSJed Brown AppCtx ctx; 188c4762a1bSJed Brown PetscScalar *u; 189c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0}; 190c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2; 191c4762a1bSJed Brown PetscReal ftime; 192c4762a1bSJed Brown PetscInt steps; 193c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr; 194c4762a1bSJed Brown Vec lambda[1], q, mu[1]; 195c4762a1bSJed Brown 196c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 197c4762a1bSJed Brown Initialize program 198c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 199327415f7SBarry Smith PetscFunctionBeginUser; 2009566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 2019566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2023c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 203c4762a1bSJed Brown 204c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 205c4762a1bSJed Brown Create necessary matrix and vectors 206c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2079566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2089566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 2099566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 2109566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2119566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 212c4762a1bSJed Brown 2139566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 214c4762a1bSJed Brown 2159566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp)); 2169566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1)); 2179566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp)); 2189566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp)); 219c4762a1bSJed Brown 2209566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP)); 2219566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDP)); 2229566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU)); 2239566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDU)); 224c4762a1bSJed Brown 225c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 226c4762a1bSJed Brown Set runtime options 227c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 228d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 229c4762a1bSJed Brown { 230c4762a1bSJed Brown ctx.beta = 2; 231c4762a1bSJed Brown ctx.c = 10000.0; 232c4762a1bSJed Brown ctx.u_s = 1.0; 233c4762a1bSJed Brown ctx.omega_s = 1.0; 234c4762a1bSJed Brown ctx.omega_b = 120.0 * PETSC_PI; 235c4762a1bSJed Brown ctx.H = 5.0; 2369566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 237c4762a1bSJed Brown ctx.D = 5.0; 2389566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 239c4762a1bSJed Brown ctx.E = 1.1378; 240c4762a1bSJed Brown ctx.V = 1.0; 241c4762a1bSJed Brown ctx.X = 0.545; 242c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 2439566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 244c4762a1bSJed Brown ctx.Pm = 1.1; 2459566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 246c4762a1bSJed Brown ctx.tf = 0.1; 247c4762a1bSJed Brown ctx.tcl = 0.2; 2489566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 2499566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 2509566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL)); 251c4762a1bSJed Brown if (ensemble) { 252c4762a1bSJed Brown ctx.tf = -1; 253c4762a1bSJed Brown ctx.tcl = -1; 254c4762a1bSJed Brown } 255c4762a1bSJed Brown 2569566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 257c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 258c4762a1bSJed Brown u[1] = 1.0; 2599566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1)); 260c4762a1bSJed Brown n = 2; 2619566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2)); 262c4762a1bSJed Brown u[0] += du[0]; 263c4762a1bSJed Brown u[1] += du[1]; 2649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 265c4762a1bSJed Brown if (flg1 || flg2) { 266c4762a1bSJed Brown ctx.tf = -1; 267c4762a1bSJed Brown ctx.tcl = -1; 268c4762a1bSJed Brown } 269c4762a1bSJed Brown } 270d0609cedSBarry Smith PetscOptionsEnd(); 271c4762a1bSJed Brown 272c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 273c4762a1bSJed Brown Create timestepping solver context 274c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2759566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 2769566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 2779566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 2789566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK)); 2799566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunction)RHSFunction, &ctx)); 2809566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx)); 2819566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts)); 2829566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, &ctx)); 2839566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobian)DRDUJacobianTranspose, &ctx)); 2849566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianP)DRDPJacobianTranspose, &ctx)); 285c4762a1bSJed Brown 286c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 287c4762a1bSJed Brown Set initial conditions 288c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2899566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 292c4762a1bSJed Brown Save trajectory of solution so that TSAdjointSolve() may be used 293c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2949566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 295c4762a1bSJed Brown 2969566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[0], NULL)); 297c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 2989566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &y_ptr)); 2999371c9d4SSatish Balay y_ptr[0] = 0.0; 3009371c9d4SSatish Balay y_ptr[1] = 0.0; 3019566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &y_ptr)); 302c4762a1bSJed Brown 3039566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp, &mu[0], NULL)); 3049566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr)); 305c4762a1bSJed Brown x_ptr[0] = -1.0; 3069566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr)); 3079566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, lambda, mu)); 308c4762a1bSJed Brown 309c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 310c4762a1bSJed Brown Set solver options 311c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3129566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 10.0)); 3139566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 3149566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 3159566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 316c4762a1bSJed Brown 317c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 318c4762a1bSJed Brown Solve nonlinear system 319c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 320c4762a1bSJed Brown if (ensemble) { 321c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 3229566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 323c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 324c4762a1bSJed Brown u[1] = ctx.omega_s; 325c4762a1bSJed Brown u[0] += du[0]; 326c4762a1bSJed Brown u[1] += du[1]; 3279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 3289566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 3299566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 330c4762a1bSJed Brown } 331c4762a1bSJed Brown } else { 3329566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 333c4762a1bSJed Brown } 3349566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD)); 3359566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 3369566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 337c4762a1bSJed Brown 338c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 339c4762a1bSJed Brown Adjoint model starts here 340c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 341c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 3429566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &y_ptr)); 3439371c9d4SSatish Balay y_ptr[0] = 0.0; 3449371c9d4SSatish Balay y_ptr[1] = 0.0; 3459566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &y_ptr)); 346c4762a1bSJed Brown 3479566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr)); 348c4762a1bSJed Brown x_ptr[0] = -1.0; 3499566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr)); 350c4762a1bSJed Brown 351c4762a1bSJed Brown /* Set RHS JacobianP */ 3529566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &ctx)); 353c4762a1bSJed Brown 3549566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 355c4762a1bSJed Brown 3569566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0] d[Psi(tf)]/d[omega0]\n")); 3579566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD)); 3589566063dSJacob Faibussowitsch PetscCall(VecView(mu[0], PETSC_VIEWER_STDOUT_WORLD)); 3599566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 3609566063dSJacob Faibussowitsch PetscCall(VecView(q, PETSC_VIEWER_STDOUT_WORLD)); 3619566063dSJacob Faibussowitsch PetscCall(VecGetArray(q, &x_ptr)); 3629566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n cost function=%g\n", (double)(x_ptr[0] - ctx.Pm))); 3639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(q, &x_ptr)); 364c4762a1bSJed Brown 3659566063dSJacob Faibussowitsch PetscCall(ComputeSensiP(lambda[0], mu[0], &ctx)); 366c4762a1bSJed Brown 367c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 368c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 369c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3719566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp)); 3729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDU)); 3739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDP)); 3749566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 3769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 3779566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 3789566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 379b122ec5aSJacob Faibussowitsch return 0; 380c4762a1bSJed Brown } 381c4762a1bSJed Brown 382c4762a1bSJed Brown /*TEST 383c4762a1bSJed Brown 384c4762a1bSJed Brown build: 385c4762a1bSJed Brown requires: !complex 386c4762a1bSJed Brown 387c4762a1bSJed Brown test: 388c4762a1bSJed Brown args: -viewer_binary_skip_info -ts_adapt_type none 389c4762a1bSJed Brown 390c4762a1bSJed Brown TEST*/ 391