1a4963045SJacob Faibussowitsch #pragma once 23ea99036SJacob Faibussowitsch 39371c9d4SSatish Balay typedef enum { 49371c9d4SSatish Balay SA_ADJ, 59371c9d4SSatish Balay SA_TLM 69371c9d4SSatish Balay } SAMethod; 7c4762a1bSJed Brown static const char *const SAMethods[] = {"ADJ", "TLM", "SAMethod", "SA_", 0}; 8c4762a1bSJed Brown 9c4762a1bSJed Brown typedef struct { 10c4762a1bSJed Brown PetscScalar H, D, omega_b, omega_s, Pmax, Pmax_ini, Pm, E, V, X, u_s, c; 11c4762a1bSJed Brown PetscInt beta; 12c4762a1bSJed Brown PetscReal tf, tcl; 13c4762a1bSJed Brown /* Solver context */ 14c4762a1bSJed Brown TS ts, quadts; 15c4762a1bSJed Brown Vec U; /* solution will be stored here */ 16c4762a1bSJed Brown Mat Jac; /* Jacobian matrix */ 17c4762a1bSJed Brown Mat Jacp; /* Jacobianp matrix */ 18c4762a1bSJed Brown Mat DRDU, DRDP; 19c4762a1bSJed Brown SAMethod sa; 20c4762a1bSJed Brown } AppCtx; 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* Event check */ 23*2a8381b2SBarry Smith PetscErrorCode EventFunction(TS ts, PetscReal t, Vec X, PetscReal *fvalue, PetscCtx ctx) 24d71ae5a4SJacob Faibussowitsch { 25c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 26c4762a1bSJed Brown 27c4762a1bSJed Brown PetscFunctionBegin; 28c4762a1bSJed Brown /* Event for fault-on time */ 29c4762a1bSJed Brown fvalue[0] = t - user->tf; 30c4762a1bSJed Brown /* Event for fault-off time */ 31c4762a1bSJed Brown fvalue[1] = t - user->tcl; 323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 33c4762a1bSJed Brown } 34c4762a1bSJed Brown 35*2a8381b2SBarry Smith PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec X, PetscBool forwardsolve, PetscCtx ctx) 36d71ae5a4SJacob Faibussowitsch { 37c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 38c4762a1bSJed Brown 39c4762a1bSJed Brown PetscFunctionBegin; 40c4762a1bSJed Brown if (event_list[0] == 0) { 41c4762a1bSJed Brown if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */ 42c4762a1bSJed Brown else user->Pmax = user->Pmax_ini; /* Going backward, reversal of event */ 43c4762a1bSJed Brown } else if (event_list[0] == 1) { 44c4762a1bSJed Brown if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault - this is done by setting Pmax = Pmax_ini */ 45c4762a1bSJed Brown else user->Pmax = 0.0; /* Going backward, reversal of event */ 46c4762a1bSJed Brown } 4769d47153SPierre Jolivet PetscCall(TSRestartStep(ts)); /* Must set restart flag to true, otherwise methods with FSAL will fail */ 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51c4762a1bSJed Brown /* 52c4762a1bSJed Brown Defines the ODE passed to the ODE solver 53c4762a1bSJed Brown */ 54d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 55d71ae5a4SJacob Faibussowitsch { 56c4762a1bSJed Brown PetscScalar *f, Pmax; 57c4762a1bSJed Brown const PetscScalar *u; 58c4762a1bSJed Brown 59c4762a1bSJed Brown PetscFunctionBegin; 60c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 619566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 629566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 63c4762a1bSJed Brown Pmax = ctx->Pmax; 64c4762a1bSJed Brown f[0] = ctx->omega_b * (u[1] - ctx->omega_s); 65c4762a1bSJed Brown f[1] = ctx->omega_s / (2.0 * ctx->H) * (ctx->Pm - Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s)); 66c4762a1bSJed Brown 679566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 689566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 70c4762a1bSJed Brown } 71c4762a1bSJed Brown 72c4762a1bSJed Brown /* 73c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of a and the Jacobian. 74c4762a1bSJed Brown */ 75d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) 76d71ae5a4SJacob Faibussowitsch { 77c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 78c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 79c4762a1bSJed Brown const PetscScalar *u; 80c4762a1bSJed Brown 81c4762a1bSJed Brown PetscFunctionBegin; 829566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 83c4762a1bSJed Brown Pmax = ctx->Pmax; 84c4762a1bSJed Brown J[0][0] = 0; 85c4762a1bSJed Brown J[0][1] = ctx->omega_b; 86c4762a1bSJed Brown J[1][0] = -ctx->omega_s / (2.0 * ctx->H) * Pmax * PetscCosScalar(u[0]); 87c4762a1bSJed Brown J[1][1] = -ctx->omega_s / (2.0 * ctx->H) * ctx->D; 889566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 899566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 90c4762a1bSJed Brown 919566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 929566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 93c4762a1bSJed Brown if (A != B) { 949566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 959566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 96c4762a1bSJed Brown } 973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 98c4762a1bSJed Brown } 99c4762a1bSJed Brown 100c4762a1bSJed Brown /* 101c4762a1bSJed Brown Defines the ODE passed to the ODE solver 102c4762a1bSJed Brown */ 103d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) 104d71ae5a4SJacob Faibussowitsch { 105c4762a1bSJed Brown PetscScalar *f, Pmax; 106c4762a1bSJed Brown const PetscScalar *u, *udot; 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscFunctionBegin; 109c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1119566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 1129566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 113c4762a1bSJed Brown Pmax = ctx->Pmax; 114c4762a1bSJed Brown f[0] = udot[0] - ctx->omega_b * (u[1] - ctx->omega_s); 115c4762a1bSJed Brown f[1] = 2.0 * ctx->H / ctx->omega_s * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - ctx->omega_s) - ctx->Pm; 116c4762a1bSJed Brown 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1189566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 121c4762a1bSJed Brown } 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* 124c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 125c4762a1bSJed Brown */ 126d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) 127d71ae5a4SJacob Faibussowitsch { 128c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 129c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 130c4762a1bSJed Brown const PetscScalar *u, *udot; 131c4762a1bSJed Brown 132c4762a1bSJed Brown PetscFunctionBegin; 1339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 135c4762a1bSJed Brown Pmax = ctx->Pmax; 1369371c9d4SSatish Balay J[0][0] = a; 1379371c9d4SSatish Balay J[0][1] = -ctx->omega_b; 1389371c9d4SSatish Balay J[1][1] = 2.0 * ctx->H / ctx->omega_s * a + ctx->D; 1399371c9d4SSatish Balay J[1][0] = Pmax * PetscCosScalar(u[0]); 140c4762a1bSJed Brown 1419566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 1429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 144c4762a1bSJed Brown 1459566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1469566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 147c4762a1bSJed Brown if (A != B) { 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1499566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 150c4762a1bSJed Brown } 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown 154*2a8381b2SBarry Smith PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx0) 155d71ae5a4SJacob Faibussowitsch { 156c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 157c4762a1bSJed Brown PetscScalar *x, J[2][1]; 158c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 159c4762a1bSJed Brown 160c4762a1bSJed Brown PetscFunctionBeginUser; 1619566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 162c4762a1bSJed Brown J[0][0] = 0; 163c4762a1bSJed Brown J[1][0] = ctx->omega_s / (2.0 * ctx->H); 1649566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 165c4762a1bSJed Brown 1669566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1679566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 169c4762a1bSJed Brown } 170c4762a1bSJed Brown 171d71ae5a4SJacob Faibussowitsch PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) 172d71ae5a4SJacob Faibussowitsch { 173c4762a1bSJed Brown PetscScalar *r; 174c4762a1bSJed Brown const PetscScalar *u; 175c4762a1bSJed Brown 176c4762a1bSJed Brown PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1789566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 1792f613bf5SBarry Smith r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta); 1809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 1819566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 183c4762a1bSJed Brown } 184c4762a1bSJed Brown 185c4762a1bSJed Brown /* Transpose of DRDU */ 186d71ae5a4SJacob Faibussowitsch PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) 187d71ae5a4SJacob Faibussowitsch { 188c4762a1bSJed Brown PetscScalar ru[2]; 189c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 190c4762a1bSJed Brown const PetscScalar *u; 191c4762a1bSJed Brown 192c4762a1bSJed Brown PetscFunctionBegin; 1939566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1942f613bf5SBarry Smith ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1); 195c4762a1bSJed Brown ru[1] = 0; 1969566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU, 2, row, 1, col, ru, INSERT_VALUES)); 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1989566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY)); 1999566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY)); 2003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown 203*2a8381b2SBarry Smith PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, PetscCtx ctx) 204d71ae5a4SJacob Faibussowitsch { 205c4762a1bSJed Brown PetscFunctionBegin; 2069566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 2079566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY)); 2089566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY)); 2093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 210c4762a1bSJed Brown } 211c4762a1bSJed Brown 212d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) 213d71ae5a4SJacob Faibussowitsch { 214c4762a1bSJed Brown PetscScalar *y, sensip; 215c4762a1bSJed Brown const PetscScalar *x; 216c4762a1bSJed Brown 217c4762a1bSJed Brown PetscFunctionBegin; 2189566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(lambda, &x)); 2199566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu, &y)); 220c4762a1bSJed Brown sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0]; 221c4762a1bSJed Brown y[0] = sensip; 2229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu, &y)); 2239566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(lambda, &x)); 2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 225c4762a1bSJed Brown } 226