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 ./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 13c4762a1bSJed Brown 14c4762a1bSJed Brown Fault at .1 seconds 15c4762a1bSJed Brown ./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 16c4762a1bSJed Brown 17c4762a1bSJed Brown Initial conditions same as when fault is ended 18c4762a1bSJed Brown ./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 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 <petsctao.h> 33c4762a1bSJed Brown #include <petscts.h> 34c4762a1bSJed Brown 35c4762a1bSJed Brown typedef struct { 36c4762a1bSJed Brown TS ts; 37c4762a1bSJed Brown PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c; 38c4762a1bSJed Brown PetscInt beta; 39c4762a1bSJed Brown PetscReal tf, tcl, dt; 40c4762a1bSJed Brown } AppCtx; 41c4762a1bSJed Brown 42c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 43c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 44c4762a1bSJed Brown 45c4762a1bSJed Brown /* 46c4762a1bSJed Brown Defines the ODE passed to the ODE solver 47c4762a1bSJed Brown */ 48d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 49d71ae5a4SJacob Faibussowitsch { 50c4762a1bSJed Brown PetscScalar *f, Pmax; 51c4762a1bSJed Brown const PetscScalar *u; 52c4762a1bSJed Brown 53c4762a1bSJed Brown PetscFunctionBegin; 54c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 559566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 569566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 57c4762a1bSJed 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 */ 58c4762a1bSJed Brown else Pmax = ctx->Pmax; 59c4762a1bSJed Brown 60c4762a1bSJed Brown f[0] = ctx->omega_b * (u[1] - ctx->omega_s); 61c4762a1bSJed Brown f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H); 62c4762a1bSJed Brown 639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 649566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 65*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* 69c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 70c4762a1bSJed Brown */ 71d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) 72d71ae5a4SJacob Faibussowitsch { 73c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 74c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 75c4762a1bSJed Brown const PetscScalar *u; 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscFunctionBegin; 789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 79c4762a1bSJed 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 */ 80c4762a1bSJed Brown else Pmax = ctx->Pmax; 81c4762a1bSJed Brown 829371c9d4SSatish Balay J[0][0] = 0; 839371c9d4SSatish Balay J[0][1] = ctx->omega_b; 849371c9d4SSatish Balay J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H); 859371c9d4SSatish Balay J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H); 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 89c4762a1bSJed Brown 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 92c4762a1bSJed Brown if (A != B) { 939566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 949566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 95c4762a1bSJed Brown } 96*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 97c4762a1bSJed Brown } 98c4762a1bSJed Brown 99d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0) 100d71ae5a4SJacob Faibussowitsch { 101c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 102c4762a1bSJed Brown PetscScalar J[2][1]; 103c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 104c4762a1bSJed Brown 105c4762a1bSJed Brown PetscFunctionBeginUser; 106c4762a1bSJed Brown J[0][0] = 0; 107c4762a1bSJed Brown J[1][0] = ctx->omega_s / (2.0 * ctx->H); 1089566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 1099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 111*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114d71ae5a4SJacob Faibussowitsch static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) 115d71ae5a4SJacob Faibussowitsch { 116c4762a1bSJed Brown PetscScalar *r; 117c4762a1bSJed Brown const PetscScalar *u; 118c4762a1bSJed Brown 119c4762a1bSJed Brown PetscFunctionBegin; 1209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1219566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 1222f613bf5SBarry Smith r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta); 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 1249566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 125*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 126c4762a1bSJed Brown } 127c4762a1bSJed Brown 128d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) 129d71ae5a4SJacob Faibussowitsch { 130c4762a1bSJed Brown PetscScalar ru[1]; 131c4762a1bSJed Brown const PetscScalar *u; 132c4762a1bSJed Brown PetscInt row[] = {0}, col[] = {0}; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBegin; 1359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1362f613bf5SBarry Smith ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1); 1379566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1389566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES)); 1399566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY)); 1409566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY)); 141*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 142c4762a1bSJed Brown } 143c4762a1bSJed Brown 144d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx) 145d71ae5a4SJacob Faibussowitsch { 146c4762a1bSJed Brown PetscFunctionBegin; 1479566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY)); 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY)); 150*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151c4762a1bSJed Brown } 152c4762a1bSJed Brown 153d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) 154d71ae5a4SJacob Faibussowitsch { 155c4762a1bSJed Brown PetscScalar *y, sensip; 156c4762a1bSJed Brown const PetscScalar *x; 157c4762a1bSJed Brown 158c4762a1bSJed Brown PetscFunctionBegin; 1599566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(lambda, &x)); 1609566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu, &y)); 161c4762a1bSJed Brown sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0]; 162c4762a1bSJed Brown y[0] = sensip; 1639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu, &y)); 1649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(lambda, &x)); 165*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 169d71ae5a4SJacob Faibussowitsch { 170c4762a1bSJed Brown Vec p; 171c4762a1bSJed Brown PetscScalar *x_ptr; 172c4762a1bSJed Brown PetscMPIInt size; 173c4762a1bSJed Brown AppCtx ctx; 174c4762a1bSJed Brown Vec lowerb, upperb; 175c4762a1bSJed Brown Tao tao; 176c4762a1bSJed Brown KSP ksp; 177c4762a1bSJed Brown PC pc; 178c4762a1bSJed Brown Vec U, lambda[1], mu[1]; 179c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 180c4762a1bSJed Brown Mat Jacp; /* Jacobian matrix */ 181c4762a1bSJed Brown Mat DRDU, DRDP; 182c4762a1bSJed Brown PetscInt n = 2; 183c4762a1bSJed Brown TS quadts; 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 186c4762a1bSJed Brown Initialize program 187c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 188327415f7SBarry Smith PetscFunctionBeginUser; 1899566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 190c4762a1bSJed Brown PetscFunctionBeginUser; 1919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1923c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 193c4762a1bSJed Brown 194c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 195c4762a1bSJed Brown Set runtime options 196c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 197d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 198c4762a1bSJed Brown { 199c4762a1bSJed Brown ctx.beta = 2; 200c4762a1bSJed Brown ctx.c = PetscRealConstant(10000.0); 201c4762a1bSJed Brown ctx.u_s = PetscRealConstant(1.0); 202c4762a1bSJed Brown ctx.omega_s = PetscRealConstant(1.0); 203c4762a1bSJed Brown ctx.omega_b = PetscRealConstant(120.0) * PETSC_PI; 204c4762a1bSJed Brown ctx.H = PetscRealConstant(5.0); 2059566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 206c4762a1bSJed Brown ctx.D = PetscRealConstant(5.0); 2079566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 208c4762a1bSJed Brown ctx.E = PetscRealConstant(1.1378); 209c4762a1bSJed Brown ctx.V = PetscRealConstant(1.0); 210c4762a1bSJed Brown ctx.X = PetscRealConstant(0.545); 211c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 213c4762a1bSJed Brown ctx.Pm = PetscRealConstant(1.0194); 2149566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 215c4762a1bSJed Brown ctx.tf = PetscRealConstant(0.1); 216c4762a1bSJed Brown ctx.tcl = PetscRealConstant(0.2); 2179566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 2189566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 219c4762a1bSJed Brown } 220d0609cedSBarry Smith PetscOptionsEnd(); 221c4762a1bSJed Brown 222c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 223c4762a1bSJed Brown Create necessary matrix and vectors 224c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2259566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2269566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 2279566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 2289566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2299566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 230c4762a1bSJed Brown 2319566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 232c4762a1bSJed Brown 2339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp)); 2349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1)); 2359566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp)); 2369566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp)); 2379566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDP)); 2399566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU)); 2409566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDU)); 241c4762a1bSJed Brown 242c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 243c4762a1bSJed Brown Create timestepping solver context 244c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2459566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ctx.ts)); 2469566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ctx.ts, TS_NONLINEAR)); 2479566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ctx.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 2489566063dSJacob Faibussowitsch PetscCall(TSSetType(ctx.ts, TSRK)); 2499566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ctx.ts, NULL, (TSRHSFunction)RHSFunction, &ctx)); 2509566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ctx.ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx)); 2519566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ctx.ts, TS_EXACTFINALTIME_MATCHSTEP)); 252c4762a1bSJed Brown 2539566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[0], NULL)); 2549566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp, &mu[0], NULL)); 2559566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ctx.ts, 1, lambda, mu)); 2569566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ctx.ts, Jacp, RHSJacobianP, &ctx)); 257c4762a1bSJed Brown 258c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 259c4762a1bSJed Brown Set solver options 260c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2619566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ctx.ts, PetscRealConstant(1.0))); 2629566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ctx.ts, PetscRealConstant(0.01))); 2639566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ctx.ts)); 264c4762a1bSJed Brown 2659566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ctx.ts, &ctx.dt)); /* save the stepsize */ 266c4762a1bSJed Brown 2679566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ctx.ts, PETSC_TRUE, &quadts)); 2689566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, &ctx)); 2699566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobian)DRDUJacobianTranspose, &ctx)); 2709566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianP)DRDPJacobianTranspose, &ctx)); 2719566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ctx.ts, U)); 272c4762a1bSJed Brown 273c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2749566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2759566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM)); 276c4762a1bSJed Brown 277c4762a1bSJed Brown /* 278c4762a1bSJed Brown Optimization starts 279c4762a1bSJed Brown */ 280c4762a1bSJed Brown /* Set initial solution guess */ 2819566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p)); 2829566063dSJacob Faibussowitsch PetscCall(VecGetArray(p, &x_ptr)); 283c4762a1bSJed Brown x_ptr[0] = ctx.Pm; 2849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(p, &x_ptr)); 285c4762a1bSJed Brown 2869566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, p)); 287c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2889566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx)); 2899566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, (void *)&ctx)); 290c4762a1bSJed Brown 291c4762a1bSJed Brown /* Set bounds for the optimization */ 2929566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &lowerb)); 2939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &upperb)); 2949566063dSJacob Faibussowitsch PetscCall(VecGetArray(lowerb, &x_ptr)); 295c4762a1bSJed Brown x_ptr[0] = 0.; 2969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lowerb, &x_ptr)); 2979566063dSJacob Faibussowitsch PetscCall(VecGetArray(upperb, &x_ptr)); 298c4762a1bSJed Brown x_ptr[0] = PetscRealConstant(1.1); 2999566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(upperb, &x_ptr)); 3009566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, lowerb, upperb)); 301c4762a1bSJed Brown 302c4762a1bSJed Brown /* Check for any TAO command line options */ 3039566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 3049566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 305c4762a1bSJed Brown if (ksp) { 3069566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 3079566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 308c4762a1bSJed Brown } 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 3119566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 312c4762a1bSJed Brown 3139566063dSJacob Faibussowitsch PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD)); 314c4762a1bSJed Brown /* Free TAO data structures */ 3159566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 3169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p)); 3179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lowerb)); 3189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&upperb)); 319c4762a1bSJed Brown 3209566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ctx.ts)); 3219566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp)); 3249566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDU)); 3259566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDP)); 3269566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 3279566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 3289566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 329b122ec5aSJacob Faibussowitsch return 0; 330c4762a1bSJed Brown } 331c4762a1bSJed Brown 332c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 333c4762a1bSJed Brown /* 334c4762a1bSJed Brown FormFunction - Evaluates the function 335c4762a1bSJed Brown 336c4762a1bSJed Brown Input Parameters: 337c4762a1bSJed Brown tao - the Tao context 338c4762a1bSJed Brown X - the input vector 339a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 340c4762a1bSJed Brown 341c4762a1bSJed Brown Output Parameters: 342c4762a1bSJed Brown f - the newly evaluated function 343c4762a1bSJed Brown */ 344d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0) 345d71ae5a4SJacob Faibussowitsch { 346c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 347c4762a1bSJed Brown TS ts = ctx->ts; 348c4762a1bSJed Brown Vec U; /* solution will be stored here */ 349c4762a1bSJed Brown PetscScalar *u; 350c4762a1bSJed Brown PetscScalar *x_ptr; 351c4762a1bSJed Brown Vec q; 352c4762a1bSJed Brown 353*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 3549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr)); 355c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 3569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr)); 357c4762a1bSJed Brown 358c4762a1bSJed Brown /* reset time */ 3599566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 360c4762a1bSJed Brown /* reset step counter, this is critical for adjoint solver */ 3619566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 362c4762a1bSJed Brown /* reset step size, the step size becomes negative after TSAdjointSolve */ 3639566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx->dt)); 364c4762a1bSJed Brown /* reinitialize the integral value */ 3659566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 3669566063dSJacob Faibussowitsch PetscCall(VecSet(q, 0.0)); 367c4762a1bSJed Brown 368c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 369c4762a1bSJed Brown Set initial conditions 370c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3719566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 3729566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 373c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax); 374c4762a1bSJed Brown u[1] = PetscRealConstant(1.0); 3759566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 376c4762a1bSJed Brown 377c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 378c4762a1bSJed Brown Solve nonlinear system 379c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3809566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 3819566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 3829566063dSJacob Faibussowitsch PetscCall(VecGetArray(q, &x_ptr)); 383c4762a1bSJed Brown *f = -ctx->Pm + x_ptr[0]; 3849566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(q, &x_ptr)); 385*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 386c4762a1bSJed Brown } 387c4762a1bSJed Brown 388d71ae5a4SJacob Faibussowitsch PetscErrorCode FormGradient(Tao tao, Vec P, Vec G, void *ctx0) 389d71ae5a4SJacob Faibussowitsch { 390c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 391c4762a1bSJed Brown TS ts = ctx->ts; 392c4762a1bSJed Brown Vec U; /* solution will be stored here */ 393c4762a1bSJed Brown PetscReal ftime; 394c4762a1bSJed Brown PetscInt steps; 395c4762a1bSJed Brown PetscScalar *u; 396c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr; 397c4762a1bSJed Brown Vec *lambda, q, *mu; 398c4762a1bSJed Brown 399*3ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 4009566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr)); 401c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 4029566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr)); 403c4762a1bSJed Brown 404c4762a1bSJed Brown /* reset time */ 4059566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 406c4762a1bSJed Brown /* reset step counter, this is critical for adjoint solver */ 4079566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 408c4762a1bSJed Brown /* reset step size, the step size becomes negative after TSAdjointSolve */ 4099566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx->dt)); 410c4762a1bSJed Brown /* reinitialize the integral value */ 4119566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 4129566063dSJacob Faibussowitsch PetscCall(VecSet(q, 0.0)); 413c4762a1bSJed Brown 414c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 415c4762a1bSJed Brown Set initial conditions 416c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4179566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 4189566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 419c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax); 420c4762a1bSJed Brown u[1] = PetscRealConstant(1.0); 4219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 422c4762a1bSJed Brown 423f32d6360SSatish Balay /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */ 4249566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 4259566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 426c4762a1bSJed Brown 427c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 428c4762a1bSJed Brown Solve nonlinear system 429c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4309566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 431c4762a1bSJed Brown 4329566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 4339566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 436c4762a1bSJed Brown Adjoint model starts here 437c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4389566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(ts, NULL, &lambda, &mu)); 439c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 4409566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &y_ptr)); 4419371c9d4SSatish Balay y_ptr[0] = 0.0; 4429371c9d4SSatish Balay y_ptr[1] = 0.0; 4439566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &y_ptr)); 4449566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr)); 445c4762a1bSJed Brown x_ptr[0] = PetscRealConstant(-1.0); 4469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr)); 447c4762a1bSJed Brown 4489566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 4499566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 4509566063dSJacob Faibussowitsch PetscCall(ComputeSensiP(lambda[0], mu[0], ctx)); 4519566063dSJacob Faibussowitsch PetscCall(VecCopy(mu[0], G)); 452*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 453c4762a1bSJed Brown } 454c4762a1bSJed Brown 455c4762a1bSJed Brown /*TEST 456c4762a1bSJed Brown 457c4762a1bSJed Brown build: 458c4762a1bSJed Brown requires: !complex 459c4762a1bSJed Brown 460c4762a1bSJed Brown test: 461c4762a1bSJed Brown args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason 462c4762a1bSJed Brown 463c4762a1bSJed Brown test: 464c4762a1bSJed Brown suffix: 2 465c4762a1bSJed Brown 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 466c4762a1bSJed Brown 467c4762a1bSJed Brown TEST*/ 468