1c4762a1bSJed Brown typedef enum {SA_ADJ, SA_TLM} SAMethod; 2c4762a1bSJed Brown static const char *const SAMethods[] = {"ADJ","TLM","SAMethod","SA_",0}; 3c4762a1bSJed Brown 4c4762a1bSJed Brown typedef struct { 5c4762a1bSJed Brown PetscScalar H,D,omega_b,omega_s,Pmax,Pmax_ini,Pm,E,V,X,u_s,c; 6c4762a1bSJed Brown PetscInt beta; 7c4762a1bSJed Brown PetscReal tf,tcl; 8c4762a1bSJed Brown /* Solver context */ 9c4762a1bSJed Brown TS ts,quadts; 10c4762a1bSJed Brown Vec U; /* solution will be stored here */ 11c4762a1bSJed Brown Mat Jac; /* Jacobian matrix */ 12c4762a1bSJed Brown Mat Jacp; /* Jacobianp matrix */ 13c4762a1bSJed Brown Mat DRDU,DRDP; 14c4762a1bSJed Brown SAMethod sa; 15c4762a1bSJed Brown } AppCtx; 16c4762a1bSJed Brown 17c4762a1bSJed Brown /* Event check */ 18c4762a1bSJed Brown PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx) 19c4762a1bSJed Brown { 20c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 21c4762a1bSJed Brown 22c4762a1bSJed Brown PetscFunctionBegin; 23c4762a1bSJed Brown /* Event for fault-on time */ 24c4762a1bSJed Brown fvalue[0] = t - user->tf; 25c4762a1bSJed Brown /* Event for fault-off time */ 26c4762a1bSJed Brown fvalue[1] = t - user->tcl; 27c4762a1bSJed Brown 28c4762a1bSJed Brown PetscFunctionReturn(0); 29c4762a1bSJed Brown } 30c4762a1bSJed Brown 31c4762a1bSJed Brown PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx) 32c4762a1bSJed Brown { 33c4762a1bSJed Brown AppCtx *user=(AppCtx*)ctx; 34c4762a1bSJed Brown 35c4762a1bSJed Brown PetscFunctionBegin; 36c4762a1bSJed Brown if (event_list[0] == 0) { 37c4762a1bSJed Brown if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */ 38c4762a1bSJed Brown else user->Pmax = user->Pmax_ini; /* Going backward, reversal of event */ 39c4762a1bSJed Brown } else if (event_list[0] == 1) { 40c4762a1bSJed Brown if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault - this is done by setting Pmax = Pmax_ini */ 41c4762a1bSJed Brown else user->Pmax = 0.0; /* Going backward, reversal of event */ 42c4762a1bSJed Brown } 43*9566063dSJacob Faibussowitsch PetscCall(TSRestartStep(ts)); /* Must set restart flag to ture, otherwise methods with FSAL will fail */ 44c4762a1bSJed Brown PetscFunctionReturn(0); 45c4762a1bSJed Brown } 46c4762a1bSJed Brown 47c4762a1bSJed Brown /* 48c4762a1bSJed Brown Defines the ODE passed to the ODE solver 49c4762a1bSJed Brown */ 50c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 51c4762a1bSJed Brown { 52c4762a1bSJed Brown PetscScalar *f,Pmax; 53c4762a1bSJed Brown const PetscScalar *u; 54c4762a1bSJed Brown 55c4762a1bSJed Brown PetscFunctionBegin; 56c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 57*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 58*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 59c4762a1bSJed Brown Pmax = ctx->Pmax; 60c4762a1bSJed Brown f[0] = ctx->omega_b*(u[1] - ctx->omega_s); 61c4762a1bSJed Brown f[1] = ctx->omega_s/(2.0*ctx->H)*(ctx->Pm - Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s)); 62c4762a1bSJed Brown 63*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 64*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 65c4762a1bSJed Brown PetscFunctionReturn(0); 66c4762a1bSJed Brown } 67c4762a1bSJed Brown 68c4762a1bSJed Brown /* 69c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of a and the Jacobian. 70c4762a1bSJed Brown */ 71c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx) 72c4762a1bSJed Brown { 73c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 74c4762a1bSJed Brown PetscScalar J[2][2],Pmax; 75c4762a1bSJed Brown const PetscScalar *u; 76c4762a1bSJed Brown 77c4762a1bSJed Brown PetscFunctionBegin; 78*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 79c4762a1bSJed Brown Pmax = ctx->Pmax; 80c4762a1bSJed Brown J[0][0] = 0; 81c4762a1bSJed Brown J[0][1] = ctx->omega_b; 82c4762a1bSJed Brown J[1][0] = -ctx->omega_s/(2.0*ctx->H)*Pmax*PetscCosScalar(u[0]); 83c4762a1bSJed Brown J[1][1] = -ctx->omega_s/(2.0*ctx->H)*ctx->D; 84*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 85*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 86c4762a1bSJed Brown 87*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 88*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 89c4762a1bSJed Brown if (A != B) { 90*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 91*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown PetscFunctionReturn(0); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96c4762a1bSJed Brown /* 97c4762a1bSJed Brown Defines the ODE passed to the ODE solver 98c4762a1bSJed Brown */ 99c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx) 100c4762a1bSJed Brown { 101c4762a1bSJed Brown PetscScalar *f,Pmax; 102c4762a1bSJed Brown const PetscScalar *u,*udot; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBegin; 105c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 106*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 107*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 108*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 109c4762a1bSJed Brown Pmax = ctx->Pmax; 110c4762a1bSJed Brown f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s); 111c4762a1bSJed 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; 112c4762a1bSJed Brown 113*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 114*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 115*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 116c4762a1bSJed Brown PetscFunctionReturn(0); 117c4762a1bSJed Brown } 118c4762a1bSJed Brown 119c4762a1bSJed Brown /* 120c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 121c4762a1bSJed Brown */ 122c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx) 123c4762a1bSJed Brown { 124c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 125c4762a1bSJed Brown PetscScalar J[2][2],Pmax; 126c4762a1bSJed Brown const PetscScalar *u,*udot; 127c4762a1bSJed Brown 128c4762a1bSJed Brown PetscFunctionBegin; 129*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 130*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 131c4762a1bSJed Brown Pmax = ctx->Pmax; 132c4762a1bSJed Brown J[0][0] = a; J[0][1] = -ctx->omega_b; 133c4762a1bSJed Brown J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D; J[1][0] = Pmax*PetscCosScalar(u[0]); 134c4762a1bSJed Brown 135*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 136*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 137*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 138c4762a1bSJed Brown 139*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 140*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 141c4762a1bSJed Brown if (A != B) { 142*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 143*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown PetscFunctionReturn(0); 146c4762a1bSJed Brown } 147c4762a1bSJed Brown 148c4762a1bSJed Brown PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0) 149c4762a1bSJed Brown { 150c4762a1bSJed Brown PetscInt row[] = {0,1},col[] = {0}; 151c4762a1bSJed Brown PetscScalar *x,J[2][1]; 152c4762a1bSJed Brown AppCtx *ctx = (AppCtx*)ctx0; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBeginUser; 155*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(X,&x)); 156c4762a1bSJed Brown J[0][0] = 0; 157c4762a1bSJed Brown J[1][0] = ctx->omega_s/(2.0*ctx->H); 158*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES)); 159c4762a1bSJed Brown 160*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 161*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 162c4762a1bSJed Brown PetscFunctionReturn(0); 163c4762a1bSJed Brown } 164c4762a1bSJed Brown 165c4762a1bSJed Brown PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx) 166c4762a1bSJed Brown { 167c4762a1bSJed Brown PetscScalar *r; 168c4762a1bSJed Brown const PetscScalar *u; 169c4762a1bSJed Brown 170c4762a1bSJed Brown PetscFunctionBegin; 171*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 172*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(R,&r)); 1732f613bf5SBarry Smith r[0] = ctx->c*PetscPowScalarInt(PetscMax(0.,u[0]-ctx->u_s),ctx->beta); 174*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(R,&r)); 175*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 176c4762a1bSJed Brown PetscFunctionReturn(0); 177c4762a1bSJed Brown } 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* Transpose of DRDU */ 180c4762a1bSJed Brown PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx) 181c4762a1bSJed Brown { 182c4762a1bSJed Brown PetscScalar ru[2]; 183c4762a1bSJed Brown PetscInt row[] = {0,1},col[] = {0}; 184c4762a1bSJed Brown const PetscScalar *u; 185c4762a1bSJed Brown 186c4762a1bSJed Brown PetscFunctionBegin; 187*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 1882f613bf5SBarry Smith ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0.,u[0]-ctx->u_s),ctx->beta-1); 189c4762a1bSJed Brown ru[1] = 0; 190*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(DRDU,2,row,1,col,ru,INSERT_VALUES)); 191*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 192*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY)); 193*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY)); 194c4762a1bSJed Brown PetscFunctionReturn(0); 195c4762a1bSJed Brown } 196c4762a1bSJed Brown 197c4762a1bSJed Brown PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,void *ctx) 198c4762a1bSJed Brown { 199c4762a1bSJed Brown PetscFunctionBegin; 200*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(DRDP)); 201*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY)); 202*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY)); 203c4762a1bSJed Brown PetscFunctionReturn(0); 204c4762a1bSJed Brown } 205c4762a1bSJed Brown 206c4762a1bSJed Brown PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx) 207c4762a1bSJed Brown { 208c4762a1bSJed Brown PetscScalar *y,sensip; 209c4762a1bSJed Brown const PetscScalar *x; 210c4762a1bSJed Brown 211c4762a1bSJed Brown PetscFunctionBegin; 212*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(lambda,&x)); 213*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(mu,&y)); 214c4762a1bSJed Brown sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0]; 215c4762a1bSJed Brown y[0] = sensip; 216*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(mu,&y)); 217*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(lambda,&x)); 218c4762a1bSJed Brown PetscFunctionReturn(0); 219c4762a1bSJed Brown } 220