1c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*F 4c4762a1bSJed Brown 5c4762a1bSJed Brown \begin{eqnarray} 6c4762a1bSJed Brown \frac{d \theta}{dt} = \omega_b (\omega - \omega_s) 7c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\ 8c4762a1bSJed Brown \end{eqnarray} 9c4762a1bSJed Brown 10c4762a1bSJed Brown Ensemble of initial conditions 11c4762a1bSJed 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 12c4762a1bSJed Brown 13c4762a1bSJed Brown Fault at .1 seconds 14c4762a1bSJed 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 15c4762a1bSJed Brown 16c4762a1bSJed Brown Initial conditions same as when fault is ended 17c4762a1bSJed 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 18c4762a1bSJed Brown 19c4762a1bSJed Brown F*/ 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* 22c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 23c4762a1bSJed Brown file automatically includes: 24c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 25c4762a1bSJed Brown petscmat.h - matrices 26c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 27c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 28c4762a1bSJed Brown petscksp.h - linear solvers 29c4762a1bSJed Brown */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown #include <petsctao.h> 32c4762a1bSJed Brown #include <petscts.h> 33c4762a1bSJed Brown 34c4762a1bSJed Brown typedef struct { 35c4762a1bSJed Brown TS ts; 36c4762a1bSJed Brown PetscScalar H, D, omega_b, omega_s, Pmax, Pm, E, V, X, u_s, c; 37c4762a1bSJed Brown PetscInt beta; 38c4762a1bSJed Brown PetscReal tf, tcl, dt; 39c4762a1bSJed Brown } AppCtx; 40c4762a1bSJed Brown 41c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 42c4762a1bSJed Brown PetscErrorCode FormGradient(Tao, Vec, Vec, void *); 43c4762a1bSJed Brown 44c4762a1bSJed Brown /* 45c4762a1bSJed Brown Defines the ODE passed to the ODE solver 46c4762a1bSJed Brown */ 47d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 48d71ae5a4SJacob Faibussowitsch { 49c4762a1bSJed Brown PetscScalar *f, Pmax; 50c4762a1bSJed Brown const PetscScalar *u; 51c4762a1bSJed Brown 52c4762a1bSJed Brown PetscFunctionBegin; 53c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 549566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 559566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 56c4762a1bSJed 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 */ 57c4762a1bSJed Brown else Pmax = ctx->Pmax; 58c4762a1bSJed Brown 59c4762a1bSJed Brown f[0] = ctx->omega_b * (u[1] - ctx->omega_s); 60c4762a1bSJed Brown f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H); 61c4762a1bSJed Brown 629566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 639566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* 68c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 69c4762a1bSJed Brown */ 70d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) 71d71ae5a4SJacob Faibussowitsch { 72c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 73c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 74c4762a1bSJed Brown const PetscScalar *u; 75c4762a1bSJed Brown 76c4762a1bSJed Brown PetscFunctionBegin; 779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 78c4762a1bSJed 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 */ 79c4762a1bSJed Brown else Pmax = ctx->Pmax; 80c4762a1bSJed Brown 819371c9d4SSatish Balay J[0][0] = 0; 829371c9d4SSatish Balay J[0][1] = ctx->omega_b; 839371c9d4SSatish Balay J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H); 849371c9d4SSatish Balay J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H); 85c4762a1bSJed Brown 869566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 91c4762a1bSJed Brown if (A != B) { 929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 94c4762a1bSJed Brown } 953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx0) 99d71ae5a4SJacob Faibussowitsch { 100c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 101c4762a1bSJed Brown PetscScalar J[2][1]; 102c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBeginUser; 105c4762a1bSJed Brown J[0][0] = 0; 106c4762a1bSJed Brown J[1][0] = ctx->omega_s / (2.0 * ctx->H); 1079566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 1089566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1099566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 111c4762a1bSJed Brown } 112c4762a1bSJed Brown 113d71ae5a4SJacob Faibussowitsch static PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) 114d71ae5a4SJacob Faibussowitsch { 115c4762a1bSJed Brown PetscScalar *r; 116c4762a1bSJed Brown const PetscScalar *u; 117c4762a1bSJed Brown 118c4762a1bSJed Brown PetscFunctionBegin; 1199566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1209566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 1212f613bf5SBarry Smith r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta); 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 1239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 125c4762a1bSJed Brown } 126c4762a1bSJed Brown 127d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) 128d71ae5a4SJacob Faibussowitsch { 129c4762a1bSJed Brown PetscScalar ru[1]; 130c4762a1bSJed Brown const PetscScalar *u; 131c4762a1bSJed Brown PetscInt row[] = {0}, col[] = {0}; 132c4762a1bSJed Brown 133c4762a1bSJed Brown PetscFunctionBegin; 1349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1352f613bf5SBarry Smith ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1); 1369566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1379566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU, 1, row, 1, col, ru, INSERT_VALUES)); 1389566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY)); 1399566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY)); 1403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 141c4762a1bSJed Brown } 142c4762a1bSJed Brown 143d71ae5a4SJacob Faibussowitsch static PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, AppCtx *ctx) 144d71ae5a4SJacob Faibussowitsch { 145c4762a1bSJed Brown PetscFunctionBegin; 1469566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 1479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY)); 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY)); 1493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 152d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) 153d71ae5a4SJacob Faibussowitsch { 154c4762a1bSJed Brown PetscScalar *y, sensip; 155c4762a1bSJed Brown const PetscScalar *x; 156c4762a1bSJed Brown 157c4762a1bSJed Brown PetscFunctionBegin; 1589566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(lambda, &x)); 1599566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu, &y)); 160c4762a1bSJed Brown sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0]; 161c4762a1bSJed Brown y[0] = sensip; 1629566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu, &y)); 1639566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(lambda, &x)); 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 165c4762a1bSJed Brown } 166c4762a1bSJed Brown 167d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 168d71ae5a4SJacob Faibussowitsch { 169c4762a1bSJed Brown Vec p; 170c4762a1bSJed Brown PetscScalar *x_ptr; 171c4762a1bSJed Brown PetscMPIInt size; 172c4762a1bSJed Brown AppCtx ctx; 173c4762a1bSJed Brown Vec lowerb, upperb; 174c4762a1bSJed Brown Tao tao; 175c4762a1bSJed Brown KSP ksp; 176c4762a1bSJed Brown PC pc; 177c4762a1bSJed Brown Vec U, lambda[1], mu[1]; 178c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 179c4762a1bSJed Brown Mat Jacp; /* Jacobian matrix */ 180c4762a1bSJed Brown Mat DRDU, DRDP; 181c4762a1bSJed Brown PetscInt n = 2; 182c4762a1bSJed Brown TS quadts; 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185c4762a1bSJed Brown Initialize program 186c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 187327415f7SBarry Smith PetscFunctionBeginUser; 1889566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 1899566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1903c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 191c4762a1bSJed Brown 192c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 193c4762a1bSJed Brown Set runtime options 194c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 195d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 196c4762a1bSJed Brown { 197c4762a1bSJed Brown ctx.beta = 2; 198c4762a1bSJed Brown ctx.c = PetscRealConstant(10000.0); 199c4762a1bSJed Brown ctx.u_s = PetscRealConstant(1.0); 200c4762a1bSJed Brown ctx.omega_s = PetscRealConstant(1.0); 201c4762a1bSJed Brown ctx.omega_b = PetscRealConstant(120.0) * PETSC_PI; 202c4762a1bSJed Brown ctx.H = PetscRealConstant(5.0); 2039566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 204c4762a1bSJed Brown ctx.D = PetscRealConstant(5.0); 2059566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 206c4762a1bSJed Brown ctx.E = PetscRealConstant(1.1378); 207c4762a1bSJed Brown ctx.V = PetscRealConstant(1.0); 208c4762a1bSJed Brown ctx.X = PetscRealConstant(0.545); 209c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 2109566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 211c4762a1bSJed Brown ctx.Pm = PetscRealConstant(1.0194); 2129566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 213c4762a1bSJed Brown ctx.tf = PetscRealConstant(0.1); 214c4762a1bSJed Brown ctx.tcl = PetscRealConstant(0.2); 2159566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 2169566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 217c4762a1bSJed Brown } 218d0609cedSBarry Smith PetscOptionsEnd(); 219c4762a1bSJed Brown 220c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 221c4762a1bSJed Brown Create necessary matrix and vectors 222c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2239566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 2249566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 2259566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 2269566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 2279566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 228c4762a1bSJed Brown 2299566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 230c4762a1bSJed Brown 2319566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp)); 2329566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 1)); 2339566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp)); 2349566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp)); 2359566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 1, NULL, &DRDP)); 2369566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDP)); 2379566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 1, 2, NULL, &DRDU)); 2389566063dSJacob Faibussowitsch PetscCall(MatSetUp(DRDU)); 239c4762a1bSJed Brown 240c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 241c4762a1bSJed Brown Create timestepping solver context 242c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2439566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ctx.ts)); 2449566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ctx.ts, TS_NONLINEAR)); 2459566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ctx.ts, TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 2469566063dSJacob Faibussowitsch PetscCall(TSSetType(ctx.ts, TSRK)); 2478434afd1SBarry Smith PetscCall(TSSetRHSFunction(ctx.ts, NULL, (TSRHSFunctionFn *)RHSFunction, &ctx)); 2488434afd1SBarry Smith PetscCall(TSSetRHSJacobian(ctx.ts, A, A, (TSRHSJacobianFn *)RHSJacobian, &ctx)); 2499566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ctx.ts, TS_EXACTFINALTIME_MATCHSTEP)); 250c4762a1bSJed Brown 2519566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &lambda[0], NULL)); 2529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(Jacp, &mu[0], NULL)); 2539566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ctx.ts, 1, lambda, mu)); 2549566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ctx.ts, Jacp, RHSJacobianP, &ctx)); 255c4762a1bSJed Brown 256c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 257c4762a1bSJed Brown Set solver options 258c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2599566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ctx.ts, PetscRealConstant(1.0))); 2609566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ctx.ts, PetscRealConstant(0.01))); 2619566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ctx.ts)); 262c4762a1bSJed Brown 2639566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ctx.ts, &ctx.dt)); /* save the stepsize */ 264c4762a1bSJed Brown 2659566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ctx.ts, PETSC_TRUE, &quadts)); 2668434afd1SBarry Smith PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, &ctx)); 2678434afd1SBarry Smith PetscCall(TSSetRHSJacobian(quadts, DRDU, DRDU, (TSRHSJacobianFn *)DRDUJacobianTranspose, &ctx)); 2688434afd1SBarry Smith PetscCall(TSSetRHSJacobianP(quadts, DRDP, (TSRHSJacobianPFn *)DRDPJacobianTranspose, &ctx)); 2699566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ctx.ts, U)); 270c4762a1bSJed Brown 271c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 2729566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 2739566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM)); 274c4762a1bSJed Brown 275c4762a1bSJed Brown /* 276c4762a1bSJed Brown Optimization starts 277c4762a1bSJed Brown */ 278c4762a1bSJed Brown /* Set initial solution guess */ 2799566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p)); 2809566063dSJacob Faibussowitsch PetscCall(VecGetArray(p, &x_ptr)); 281c4762a1bSJed Brown x_ptr[0] = ctx.Pm; 2829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(p, &x_ptr)); 283c4762a1bSJed Brown 2849566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, p)); 285c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 2869566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx)); 2879566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, FormGradient, (void *)&ctx)); 288c4762a1bSJed Brown 289c4762a1bSJed Brown /* Set bounds for the optimization */ 2909566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &lowerb)); 2919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &upperb)); 2929566063dSJacob Faibussowitsch PetscCall(VecGetArray(lowerb, &x_ptr)); 293c4762a1bSJed Brown x_ptr[0] = 0.; 2949566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lowerb, &x_ptr)); 2959566063dSJacob Faibussowitsch PetscCall(VecGetArray(upperb, &x_ptr)); 296c4762a1bSJed Brown x_ptr[0] = PetscRealConstant(1.1); 2979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(upperb, &x_ptr)); 2989566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, lowerb, upperb)); 299c4762a1bSJed Brown 300c4762a1bSJed Brown /* Check for any TAO command line options */ 3019566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 3029566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 303c4762a1bSJed Brown if (ksp) { 3049566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 3059566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 306c4762a1bSJed Brown } 307c4762a1bSJed Brown 308c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 3099566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 310c4762a1bSJed Brown 3119566063dSJacob Faibussowitsch PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD)); 312c4762a1bSJed Brown /* Free TAO data structures */ 3139566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 3149566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p)); 3159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lowerb)); 3169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&upperb)); 317c4762a1bSJed Brown 3189566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ctx.ts)); 3199566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 3209566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 3219566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp)); 3229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDU)); 3239566063dSJacob Faibussowitsch PetscCall(MatDestroy(&DRDP)); 3249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lambda[0])); 3259566063dSJacob Faibussowitsch PetscCall(VecDestroy(&mu[0])); 3269566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 327b122ec5aSJacob Faibussowitsch return 0; 328c4762a1bSJed Brown } 329c4762a1bSJed Brown 330c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 331c4762a1bSJed Brown /* 332c4762a1bSJed Brown FormFunction - Evaluates the function 333c4762a1bSJed Brown 334c4762a1bSJed Brown Input Parameters: 335c4762a1bSJed Brown tao - the Tao context 336c4762a1bSJed Brown X - the input vector 337a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 338c4762a1bSJed Brown 339c4762a1bSJed Brown Output Parameters: 340c4762a1bSJed Brown f - the newly evaluated function 341c4762a1bSJed Brown */ 342*2a8381b2SBarry Smith PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, PetscCtx ctx0) 343d71ae5a4SJacob Faibussowitsch { 344c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 345c4762a1bSJed Brown TS ts = ctx->ts; 346c4762a1bSJed Brown Vec U; /* solution will be stored here */ 347c4762a1bSJed Brown PetscScalar *u; 348c4762a1bSJed Brown PetscScalar *x_ptr; 349c4762a1bSJed Brown Vec q; 350c4762a1bSJed Brown 3513ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 3529566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr)); 353c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 3549566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr)); 355c4762a1bSJed Brown 356c4762a1bSJed Brown /* reset time */ 3579566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 358c4762a1bSJed Brown /* reset step counter, this is critical for adjoint solver */ 3599566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 360c4762a1bSJed Brown /* reset step size, the step size becomes negative after TSAdjointSolve */ 3619566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx->dt)); 362c4762a1bSJed Brown /* reinitialize the integral value */ 3639566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 3649566063dSJacob Faibussowitsch PetscCall(VecSet(q, 0.0)); 365c4762a1bSJed Brown 366c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 367c4762a1bSJed Brown Set initial conditions 368c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3699566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 3709566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 371c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax); 372c4762a1bSJed Brown u[1] = PetscRealConstant(1.0); 3739566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 374c4762a1bSJed Brown 375c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 376c4762a1bSJed Brown Solve nonlinear system 377c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3789566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 3799566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 3809566063dSJacob Faibussowitsch PetscCall(VecGetArray(q, &x_ptr)); 381c4762a1bSJed Brown *f = -ctx->Pm + x_ptr[0]; 3829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(q, &x_ptr)); 3833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 384c4762a1bSJed Brown } 385c4762a1bSJed Brown 386*2a8381b2SBarry Smith PetscErrorCode FormGradient(Tao tao, Vec P, Vec G, PetscCtx ctx0) 387d71ae5a4SJacob Faibussowitsch { 388c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 389c4762a1bSJed Brown TS ts = ctx->ts; 390c4762a1bSJed Brown Vec U; /* solution will be stored here */ 391c4762a1bSJed Brown PetscReal ftime; 392c4762a1bSJed Brown PetscInt steps; 393c4762a1bSJed Brown PetscScalar *u; 394c4762a1bSJed Brown PetscScalar *x_ptr, *y_ptr; 395c4762a1bSJed Brown Vec *lambda, q, *mu; 396c4762a1bSJed Brown 3973ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 3989566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, (const PetscScalar **)&x_ptr)); 399c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 4009566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, (const PetscScalar **)&x_ptr)); 401c4762a1bSJed Brown 402c4762a1bSJed Brown /* reset time */ 4039566063dSJacob Faibussowitsch PetscCall(TSSetTime(ts, 0.0)); 404c4762a1bSJed Brown /* reset step counter, this is critical for adjoint solver */ 4059566063dSJacob Faibussowitsch PetscCall(TSSetStepNumber(ts, 0)); 406c4762a1bSJed Brown /* reset step size, the step size becomes negative after TSAdjointSolve */ 4079566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, ctx->dt)); 408c4762a1bSJed Brown /* reinitialize the integral value */ 4099566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 4109566063dSJacob Faibussowitsch PetscCall(VecSet(q, 0.0)); 411c4762a1bSJed Brown 412c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 413c4762a1bSJed Brown Set initial conditions 414c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4159566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &U)); 4169566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 417c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax); 418c4762a1bSJed Brown u[1] = PetscRealConstant(1.0); 4199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 420c4762a1bSJed Brown 421f32d6360SSatish Balay /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */ 4229566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 4239566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 424c4762a1bSJed Brown 425c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 426c4762a1bSJed Brown Solve nonlinear system 427c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4289566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 429c4762a1bSJed Brown 4309566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 4319566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 432c4762a1bSJed Brown 433c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 434c4762a1bSJed Brown Adjoint model starts here 435c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4369566063dSJacob Faibussowitsch PetscCall(TSGetCostGradients(ts, NULL, &lambda, &mu)); 437c4762a1bSJed Brown /* Set initial conditions for the adjoint integration */ 4389566063dSJacob Faibussowitsch PetscCall(VecGetArray(lambda[0], &y_ptr)); 4399371c9d4SSatish Balay y_ptr[0] = 0.0; 4409371c9d4SSatish Balay y_ptr[1] = 0.0; 4419566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lambda[0], &y_ptr)); 4429566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu[0], &x_ptr)); 443c4762a1bSJed Brown x_ptr[0] = PetscRealConstant(-1.0); 4449566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu[0], &x_ptr)); 445c4762a1bSJed Brown 4469566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 4479566063dSJacob Faibussowitsch PetscCall(TSGetCostIntegral(ts, &q)); 4489566063dSJacob Faibussowitsch PetscCall(ComputeSensiP(lambda[0], mu[0], ctx)); 4499566063dSJacob Faibussowitsch PetscCall(VecCopy(mu[0], G)); 4503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 451c4762a1bSJed Brown } 452c4762a1bSJed Brown 453c4762a1bSJed Brown /*TEST 454c4762a1bSJed Brown 455c4762a1bSJed Brown build: 4569d5502f9SJunchao Zhang requires: !complex !single 457c4762a1bSJed Brown 458c4762a1bSJed Brown test: 459c4762a1bSJed Brown args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason 460c4762a1bSJed Brown 461c4762a1bSJed Brown test: 462c4762a1bSJed Brown suffix: 2 463c4762a1bSJed 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 464c4762a1bSJed Brown 465c4762a1bSJed Brown TEST*/ 466