xref: /petsc/src/ts/tutorials/power_grid/ex3.h (revision 2f613bf53f46f9356e00a2ca2bd69453be72fc31)
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   PetscErrorCode    ierr;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   PetscFunctionBegin;
37c4762a1bSJed Brown   if (event_list[0] == 0) {
38c4762a1bSJed Brown     if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
39c4762a1bSJed Brown     else user->Pmax = user->Pmax_ini; /* Going backward, reversal of event */
40c4762a1bSJed Brown   } else if (event_list[0] == 1) {
41c4762a1bSJed Brown     if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault  - this is done by setting Pmax = Pmax_ini */
42c4762a1bSJed Brown     else user->Pmax = 0.0; /* Going backward, reversal of event */
43c4762a1bSJed Brown   }
44c4762a1bSJed Brown   ierr = TSRestartStep(ts);CHKERRQ(ierr); /* Must set restart flag to ture, otherwise methods with FSAL will fail */
45c4762a1bSJed Brown   PetscFunctionReturn(0);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*
49c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
50c4762a1bSJed Brown */
51c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
52c4762a1bSJed Brown {
53c4762a1bSJed Brown   PetscErrorCode    ierr;
54c4762a1bSJed Brown   PetscScalar       *f,Pmax;
55c4762a1bSJed Brown   const PetscScalar *u;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   PetscFunctionBegin;
58c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
59c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
60c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
61c4762a1bSJed Brown   Pmax = ctx->Pmax;
62c4762a1bSJed Brown   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
63c4762a1bSJed Brown   f[1] = ctx->omega_s/(2.0*ctx->H)*(ctx->Pm - Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s));
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
66c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
67c4762a1bSJed Brown   PetscFunctionReturn(0);
68c4762a1bSJed Brown }
69c4762a1bSJed Brown 
70c4762a1bSJed Brown /*
71c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetRHSJacobian() for the meaning of a and the Jacobian.
72c4762a1bSJed Brown */
73c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
74c4762a1bSJed Brown {
75c4762a1bSJed Brown   PetscErrorCode    ierr;
76c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
77c4762a1bSJed Brown   PetscScalar       J[2][2],Pmax;
78c4762a1bSJed Brown   const PetscScalar *u;
79c4762a1bSJed Brown 
80c4762a1bSJed Brown   PetscFunctionBegin;
81c4762a1bSJed Brown   ierr    = VecGetArrayRead(U,&u);CHKERRQ(ierr);
82c4762a1bSJed Brown   Pmax    = ctx->Pmax;
83c4762a1bSJed Brown   J[0][0] = 0;
84c4762a1bSJed Brown   J[0][1] = ctx->omega_b;
85c4762a1bSJed Brown   J[1][0] = -ctx->omega_s/(2.0*ctx->H)*Pmax*PetscCosScalar(u[0]);
86c4762a1bSJed Brown   J[1][1] = -ctx->omega_s/(2.0*ctx->H)*ctx->D;
87c4762a1bSJed Brown   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
88c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
92c4762a1bSJed Brown   if (A != B) {
93c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95c4762a1bSJed Brown   }
96c4762a1bSJed Brown   PetscFunctionReturn(0);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown 
99c4762a1bSJed Brown /*
100c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
101c4762a1bSJed Brown */
102c4762a1bSJed Brown PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
103c4762a1bSJed Brown {
104c4762a1bSJed Brown   PetscErrorCode    ierr;
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 */
110c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
111c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
112c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
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 
117c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
118c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
119c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
120c4762a1bSJed Brown   PetscFunctionReturn(0);
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 */
126c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
127c4762a1bSJed Brown {
128c4762a1bSJed Brown   PetscErrorCode    ierr;
129c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
130c4762a1bSJed Brown   PetscScalar       J[2][2],Pmax;
131c4762a1bSJed Brown   const PetscScalar *u,*udot;
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   PetscFunctionBegin;
134c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
135c4762a1bSJed Brown   ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr);
136c4762a1bSJed Brown   Pmax = ctx->Pmax;
137c4762a1bSJed Brown   J[0][0] = a;                       J[0][1] = -ctx->omega_b;
138c4762a1bSJed Brown   J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   ierr    = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
141c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
142c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr);
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146c4762a1bSJed Brown   if (A != B) {
147c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
148c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
149c4762a1bSJed Brown   }
150c4762a1bSJed Brown   PetscFunctionReturn(0);
151c4762a1bSJed Brown }
152c4762a1bSJed Brown 
153c4762a1bSJed Brown PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
154c4762a1bSJed Brown {
155c4762a1bSJed Brown   PetscErrorCode ierr;
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;
161c4762a1bSJed Brown   ierr    = VecGetArray(X,&x);CHKERRQ(ierr);
162c4762a1bSJed Brown   J[0][0] = 0;
163c4762a1bSJed Brown   J[1][0] = ctx->omega_s/(2.0*ctx->H);
164c4762a1bSJed Brown   ierr    = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
165c4762a1bSJed Brown 
166c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
167c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168c4762a1bSJed Brown   PetscFunctionReturn(0);
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
172c4762a1bSJed Brown {
173c4762a1bSJed Brown   PetscErrorCode    ierr;
174c4762a1bSJed Brown   PetscScalar       *r;
175c4762a1bSJed Brown   const PetscScalar *u;
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   PetscFunctionBegin;
178c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
179c4762a1bSJed Brown   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
180*2f613bf5SBarry Smith   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0.,u[0]-ctx->u_s),ctx->beta);
181c4762a1bSJed Brown   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
182c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
183c4762a1bSJed Brown   PetscFunctionReturn(0);
184c4762a1bSJed Brown }
185c4762a1bSJed Brown 
186c4762a1bSJed Brown /* Transpose of DRDU */
187c4762a1bSJed Brown PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx)
188c4762a1bSJed Brown {
189c4762a1bSJed Brown   PetscScalar       ru[2];
190c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[] = {0};
191c4762a1bSJed Brown   const PetscScalar *u;
192c4762a1bSJed Brown   PetscErrorCode    ierr;
193c4762a1bSJed Brown 
194c4762a1bSJed Brown   PetscFunctionBegin;
195c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
196*2f613bf5SBarry Smith   ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0.,u[0]-ctx->u_s),ctx->beta-1);
197c4762a1bSJed Brown   ru[1] = 0;
198c4762a1bSJed Brown   ierr = MatSetValues(DRDU,2,row,1,col,ru,INSERT_VALUES);CHKERRQ(ierr);
199c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
200c4762a1bSJed Brown   ierr = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
201c4762a1bSJed Brown   ierr = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
202c4762a1bSJed Brown   PetscFunctionReturn(0);
203c4762a1bSJed Brown }
204c4762a1bSJed Brown 
205c4762a1bSJed Brown PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,void *ctx)
206c4762a1bSJed Brown {
207c4762a1bSJed Brown   PetscErrorCode ierr;
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   PetscFunctionBegin;
210c4762a1bSJed Brown   ierr = MatZeroEntries(DRDP);CHKERRQ(ierr);
211c4762a1bSJed Brown   ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
212c4762a1bSJed Brown   ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
213c4762a1bSJed Brown   PetscFunctionReturn(0);
214c4762a1bSJed Brown }
215c4762a1bSJed Brown 
216c4762a1bSJed Brown PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
217c4762a1bSJed Brown {
218c4762a1bSJed Brown   PetscErrorCode    ierr;
219c4762a1bSJed Brown   PetscScalar       *y,sensip;
220c4762a1bSJed Brown   const PetscScalar *x;
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   PetscFunctionBegin;
223c4762a1bSJed Brown   ierr = VecGetArrayRead(lambda,&x);CHKERRQ(ierr);
224c4762a1bSJed Brown   ierr = VecGetArray(mu,&y);CHKERRQ(ierr);
225c4762a1bSJed Brown   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
226c4762a1bSJed Brown   y[0] = sensip;
227c4762a1bSJed Brown   ierr = VecRestoreArray(mu,&y);CHKERRQ(ierr);
228c4762a1bSJed Brown   ierr = VecRestoreArrayRead(lambda,&x);CHKERRQ(ierr);
229c4762a1bSJed Brown   PetscFunctionReturn(0);
230c4762a1bSJed Brown }
231