1*3ea99036SJacob Faibussowitsch #ifndef EX3_H 2*3ea99036SJacob Faibussowitsch #define EX3_H 3*3ea99036SJacob Faibussowitsch 49371c9d4SSatish Balay typedef enum { 59371c9d4SSatish Balay SA_ADJ, 69371c9d4SSatish Balay SA_TLM 79371c9d4SSatish Balay } SAMethod; 8c4762a1bSJed Brown static const char *const SAMethods[] = {"ADJ", "TLM", "SAMethod", "SA_", 0}; 9c4762a1bSJed Brown 10c4762a1bSJed Brown typedef struct { 11c4762a1bSJed Brown PetscScalar H, D, omega_b, omega_s, Pmax, Pmax_ini, Pm, E, V, X, u_s, c; 12c4762a1bSJed Brown PetscInt beta; 13c4762a1bSJed Brown PetscReal tf, tcl; 14c4762a1bSJed Brown /* Solver context */ 15c4762a1bSJed Brown TS ts, quadts; 16c4762a1bSJed Brown Vec U; /* solution will be stored here */ 17c4762a1bSJed Brown Mat Jac; /* Jacobian matrix */ 18c4762a1bSJed Brown Mat Jacp; /* Jacobianp matrix */ 19c4762a1bSJed Brown Mat DRDU, DRDP; 20c4762a1bSJed Brown SAMethod sa; 21c4762a1bSJed Brown } AppCtx; 22c4762a1bSJed Brown 23c4762a1bSJed Brown /* Event check */ 24d71ae5a4SJacob Faibussowitsch PetscErrorCode EventFunction(TS ts, PetscReal t, Vec X, PetscScalar *fvalue, void *ctx) 25d71ae5a4SJacob Faibussowitsch { 26c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 27c4762a1bSJed Brown 28c4762a1bSJed Brown PetscFunctionBegin; 29c4762a1bSJed Brown /* Event for fault-on time */ 30c4762a1bSJed Brown fvalue[0] = t - user->tf; 31c4762a1bSJed Brown /* Event for fault-off time */ 32c4762a1bSJed Brown fvalue[1] = t - user->tcl; 33c4762a1bSJed Brown 343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 35c4762a1bSJed Brown } 36c4762a1bSJed Brown 37d71ae5a4SJacob Faibussowitsch PetscErrorCode PostEventFunction(TS ts, PetscInt nevents, PetscInt event_list[], PetscReal t, Vec X, PetscBool forwardsolve, void *ctx) 38d71ae5a4SJacob Faibussowitsch { 39c4762a1bSJed Brown AppCtx *user = (AppCtx *)ctx; 40c4762a1bSJed Brown 41c4762a1bSJed Brown PetscFunctionBegin; 42c4762a1bSJed Brown if (event_list[0] == 0) { 43c4762a1bSJed Brown if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */ 44c4762a1bSJed Brown else user->Pmax = user->Pmax_ini; /* Going backward, reversal of event */ 45c4762a1bSJed Brown } else if (event_list[0] == 1) { 46c4762a1bSJed Brown if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault - this is done by setting Pmax = Pmax_ini */ 47c4762a1bSJed Brown else user->Pmax = 0.0; /* Going backward, reversal of event */ 48c4762a1bSJed Brown } 4969d47153SPierre Jolivet PetscCall(TSRestartStep(ts)); /* Must set restart flag to true, otherwise methods with FSAL will fail */ 503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 51c4762a1bSJed Brown } 52c4762a1bSJed Brown 53c4762a1bSJed Brown /* 54c4762a1bSJed Brown Defines the ODE passed to the ODE solver 55c4762a1bSJed Brown */ 56d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 57d71ae5a4SJacob Faibussowitsch { 58c4762a1bSJed Brown PetscScalar *f, Pmax; 59c4762a1bSJed Brown const PetscScalar *u; 60c4762a1bSJed Brown 61c4762a1bSJed Brown PetscFunctionBegin; 62c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 639566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 649566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 65c4762a1bSJed Brown Pmax = ctx->Pmax; 66c4762a1bSJed Brown f[0] = ctx->omega_b * (u[1] - ctx->omega_s); 67c4762a1bSJed Brown f[1] = ctx->omega_s / (2.0 * ctx->H) * (ctx->Pm - Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s)); 68c4762a1bSJed Brown 699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72c4762a1bSJed Brown } 73c4762a1bSJed Brown 74c4762a1bSJed Brown /* 75c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of a and the Jacobian. 76c4762a1bSJed Brown */ 77d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) 78d71ae5a4SJacob Faibussowitsch { 79c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 80c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 81c4762a1bSJed Brown const PetscScalar *u; 82c4762a1bSJed Brown 83c4762a1bSJed Brown PetscFunctionBegin; 849566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 85c4762a1bSJed Brown Pmax = ctx->Pmax; 86c4762a1bSJed Brown J[0][0] = 0; 87c4762a1bSJed Brown J[0][1] = ctx->omega_b; 88c4762a1bSJed Brown J[1][0] = -ctx->omega_s / (2.0 * ctx->H) * Pmax * PetscCosScalar(u[0]); 89c4762a1bSJed Brown J[1][1] = -ctx->omega_s / (2.0 * ctx->H) * ctx->D; 909566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 919566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 92c4762a1bSJed Brown 939566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 949566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 95c4762a1bSJed Brown if (A != B) { 969566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 979566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 98c4762a1bSJed Brown } 993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 100c4762a1bSJed Brown } 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* 103c4762a1bSJed Brown Defines the ODE passed to the ODE solver 104c4762a1bSJed Brown */ 105d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) 106d71ae5a4SJacob Faibussowitsch { 107c4762a1bSJed Brown PetscScalar *f, Pmax; 108c4762a1bSJed Brown const PetscScalar *u, *udot; 109c4762a1bSJed Brown 110c4762a1bSJed Brown PetscFunctionBegin; 111c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 1129566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1139566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 1149566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 115c4762a1bSJed Brown Pmax = ctx->Pmax; 116c4762a1bSJed Brown f[0] = udot[0] - ctx->omega_b * (u[1] - ctx->omega_s); 117c4762a1bSJed 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; 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1209566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 1219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123c4762a1bSJed Brown } 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* 126c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 127c4762a1bSJed Brown */ 128d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) 129d71ae5a4SJacob Faibussowitsch { 130c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 131c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 132c4762a1bSJed Brown const PetscScalar *u, *udot; 133c4762a1bSJed Brown 134c4762a1bSJed Brown PetscFunctionBegin; 1359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 137c4762a1bSJed Brown Pmax = ctx->Pmax; 1389371c9d4SSatish Balay J[0][0] = a; 1399371c9d4SSatish Balay J[0][1] = -ctx->omega_b; 1409371c9d4SSatish Balay J[1][1] = 2.0 * ctx->H / ctx->omega_s * a + ctx->D; 1419371c9d4SSatish Balay J[1][0] = Pmax * PetscCosScalar(u[0]); 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 146c4762a1bSJed Brown 1479566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1489566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 149c4762a1bSJed Brown if (A != B) { 1509566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 1519566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 152c4762a1bSJed Brown } 1533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx0) 157d71ae5a4SJacob Faibussowitsch { 158c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 159c4762a1bSJed Brown PetscScalar *x, J[2][1]; 160c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 161c4762a1bSJed Brown 162c4762a1bSJed Brown PetscFunctionBeginUser; 1639566063dSJacob Faibussowitsch PetscCall(VecGetArray(X, &x)); 164c4762a1bSJed Brown J[0][0] = 0; 165c4762a1bSJed Brown J[1][0] = ctx->omega_s / (2.0 * ctx->H); 1669566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES)); 167c4762a1bSJed Brown 1689566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1699566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173d71ae5a4SJacob Faibussowitsch PetscErrorCode CostIntegrand(TS ts, PetscReal t, Vec U, Vec R, AppCtx *ctx) 174d71ae5a4SJacob Faibussowitsch { 175c4762a1bSJed Brown PetscScalar *r; 176c4762a1bSJed Brown const PetscScalar *u; 177c4762a1bSJed Brown 178c4762a1bSJed Brown PetscFunctionBegin; 1799566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1809566063dSJacob Faibussowitsch PetscCall(VecGetArray(R, &r)); 1812f613bf5SBarry Smith r[0] = ctx->c * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta); 1829566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R, &r)); 1839566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 185c4762a1bSJed Brown } 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* Transpose of DRDU */ 188d71ae5a4SJacob Faibussowitsch PetscErrorCode DRDUJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDU, Mat B, AppCtx *ctx) 189d71ae5a4SJacob Faibussowitsch { 190c4762a1bSJed Brown PetscScalar ru[2]; 191c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {0}; 192c4762a1bSJed Brown const PetscScalar *u; 193c4762a1bSJed Brown 194c4762a1bSJed Brown PetscFunctionBegin; 1959566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1962f613bf5SBarry Smith ru[0] = ctx->c * ctx->beta * PetscPowScalarInt(PetscMax(0., u[0] - ctx->u_s), ctx->beta - 1); 197c4762a1bSJed Brown ru[1] = 0; 1989566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU, 2, row, 1, col, ru, INSERT_VALUES)); 1999566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 2009566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU, MAT_FINAL_ASSEMBLY)); 2019566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU, MAT_FINAL_ASSEMBLY)); 2023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown 205d71ae5a4SJacob Faibussowitsch PetscErrorCode DRDPJacobianTranspose(TS ts, PetscReal t, Vec U, Mat DRDP, void *ctx) 206d71ae5a4SJacob Faibussowitsch { 207c4762a1bSJed Brown PetscFunctionBegin; 2089566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 2099566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDP, MAT_FINAL_ASSEMBLY)); 2109566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDP, MAT_FINAL_ASSEMBLY)); 2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 212c4762a1bSJed Brown } 213c4762a1bSJed Brown 214d71ae5a4SJacob Faibussowitsch PetscErrorCode ComputeSensiP(Vec lambda, Vec mu, AppCtx *ctx) 215d71ae5a4SJacob Faibussowitsch { 216c4762a1bSJed Brown PetscScalar *y, sensip; 217c4762a1bSJed Brown const PetscScalar *x; 218c4762a1bSJed Brown 219c4762a1bSJed Brown PetscFunctionBegin; 2209566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(lambda, &x)); 2219566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu, &y)); 222c4762a1bSJed Brown sensip = 1. / PetscSqrtScalar(1. - (ctx->Pm / ctx->Pmax) * (ctx->Pm / ctx->Pmax)) / ctx->Pmax * x[0] + y[0]; 223c4762a1bSJed Brown y[0] = sensip; 2249566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu, &y)); 2259566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(lambda, &x)); 2263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 227c4762a1bSJed Brown } 228*3ea99036SJacob Faibussowitsch #endif 229