xref: /petsc/src/ts/tutorials/power_grid/ex9adj.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n";
3*c4762a1bSJed Brown 
4*c4762a1bSJed Brown /*F
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown \begin{eqnarray}
7*c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8*c4762a1bSJed Brown                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9*c4762a1bSJed Brown \end{eqnarray}
10*c4762a1bSJed Brown 
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown   Ensemble of initial conditions
14*c4762a1bSJed Brown    ./ex9 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown   Fault at .1 seconds
17*c4762a1bSJed Brown    ./ex9 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
18*c4762a1bSJed Brown 
19*c4762a1bSJed Brown   Initial conditions same as when fault is ended
20*c4762a1bSJed Brown    ./ex9 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rk -pc_type lu -ksp_type preonly
21*c4762a1bSJed Brown 
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown F*/
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown /*
26*c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
27*c4762a1bSJed Brown    file automatically includes:
28*c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
29*c4762a1bSJed Brown      petscmat.h - matrices
30*c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
31*c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
32*c4762a1bSJed Brown      petscksp.h   - linear solvers
33*c4762a1bSJed Brown */
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown #include <petscts.h>
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown typedef struct {
38*c4762a1bSJed Brown   PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c;
39*c4762a1bSJed Brown   PetscInt    beta;
40*c4762a1bSJed Brown   PetscReal   tf,tcl;
41*c4762a1bSJed Brown } AppCtx;
42*c4762a1bSJed Brown 
43*c4762a1bSJed Brown PetscErrorCode PostStepFunction(TS ts)
44*c4762a1bSJed Brown {
45*c4762a1bSJed Brown   PetscErrorCode    ierr;
46*c4762a1bSJed Brown   Vec               U;
47*c4762a1bSJed Brown   PetscReal         t;
48*c4762a1bSJed Brown   const PetscScalar *u;
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   PetscFunctionBegin;
51*c4762a1bSJed Brown   ierr = TSGetTime(ts,&t);CHKERRQ(ierr);
52*c4762a1bSJed Brown   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
53*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
54*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",(double)t,(double)u[0]);CHKERRQ(ierr);
55*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
56*c4762a1bSJed Brown   PetscFunctionReturn(0);
57*c4762a1bSJed Brown }
58*c4762a1bSJed Brown 
59*c4762a1bSJed Brown /*
60*c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
61*c4762a1bSJed Brown */
62*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
63*c4762a1bSJed Brown {
64*c4762a1bSJed Brown   PetscErrorCode    ierr;
65*c4762a1bSJed Brown   PetscScalar       *f,Pmax;
66*c4762a1bSJed Brown   const PetscScalar *u;
67*c4762a1bSJed Brown 
68*c4762a1bSJed Brown   PetscFunctionBegin;
69*c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
70*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
71*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
72*c4762a1bSJed 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 */
73*c4762a1bSJed Brown   else Pmax = ctx->Pmax;
74*c4762a1bSJed Brown 
75*c4762a1bSJed Brown   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
76*c4762a1bSJed Brown   f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
79*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
80*c4762a1bSJed Brown   PetscFunctionReturn(0);
81*c4762a1bSJed Brown }
82*c4762a1bSJed Brown 
83*c4762a1bSJed Brown /*
84*c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
85*c4762a1bSJed Brown */
86*c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
87*c4762a1bSJed Brown {
88*c4762a1bSJed Brown   PetscErrorCode    ierr;
89*c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
90*c4762a1bSJed Brown   PetscScalar       J[2][2],Pmax;
91*c4762a1bSJed Brown   const PetscScalar *u;
92*c4762a1bSJed Brown 
93*c4762a1bSJed Brown   PetscFunctionBegin;
94*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
95*c4762a1bSJed 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 */
96*c4762a1bSJed Brown   else Pmax = ctx->Pmax;
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   J[0][0] = 0;                                  J[0][1] = ctx->omega_b;
99*c4762a1bSJed Brown   J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H);  J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H);
100*c4762a1bSJed Brown 
101*c4762a1bSJed Brown   ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
102*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
106*c4762a1bSJed Brown   if (A != B) {
107*c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
108*c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
109*c4762a1bSJed Brown   }
110*c4762a1bSJed Brown   PetscFunctionReturn(0);
111*c4762a1bSJed Brown }
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
114*c4762a1bSJed Brown {
115*c4762a1bSJed Brown   PetscErrorCode ierr;
116*c4762a1bSJed Brown   PetscInt       row[] = {0,1},col[]={0};
117*c4762a1bSJed Brown   PetscScalar    J[2][1];
118*c4762a1bSJed Brown   AppCtx         *ctx=(AppCtx*)ctx0;
119*c4762a1bSJed Brown 
120*c4762a1bSJed Brown   PetscFunctionBeginUser;
121*c4762a1bSJed Brown   J[0][0] = 0;
122*c4762a1bSJed Brown   J[1][0] = ctx->omega_s/(2.0*ctx->H);
123*c4762a1bSJed Brown   ierr    = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
124*c4762a1bSJed Brown   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
126*c4762a1bSJed Brown   PetscFunctionReturn(0);
127*c4762a1bSJed Brown }
128*c4762a1bSJed Brown 
129*c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
130*c4762a1bSJed Brown {
131*c4762a1bSJed Brown   PetscErrorCode    ierr;
132*c4762a1bSJed Brown   PetscScalar       *r;
133*c4762a1bSJed Brown   const PetscScalar *u;
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown   PetscFunctionBegin;
136*c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
138*c4762a1bSJed Brown   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
141*c4762a1bSJed Brown   PetscFunctionReturn(0);
142*c4762a1bSJed Brown }
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx)
145*c4762a1bSJed Brown {
146*c4762a1bSJed Brown   PetscErrorCode    ierr;
147*c4762a1bSJed Brown   PetscScalar       ru[1];
148*c4762a1bSJed Brown   const PetscScalar *u;
149*c4762a1bSJed Brown   PetscInt          row[] = {0},col[] = {0};
150*c4762a1bSJed Brown 
151*c4762a1bSJed Brown   PetscFunctionBegin;
152*c4762a1bSJed Brown   ierr  = VecGetArrayRead(U,&u);CHKERRQ(ierr);
153*c4762a1bSJed Brown   ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr  = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
155*c4762a1bSJed Brown   ierr  = MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES);CHKERRQ(ierr);
156*c4762a1bSJed Brown   ierr  = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157*c4762a1bSJed Brown   ierr  = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158*c4762a1bSJed Brown   PetscFunctionReturn(0);
159*c4762a1bSJed Brown }
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx)
162*c4762a1bSJed Brown {
163*c4762a1bSJed Brown   PetscErrorCode    ierr;
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown   PetscFunctionBegin;
166*c4762a1bSJed Brown   ierr = MatZeroEntries(DRDP);CHKERRQ(ierr);
167*c4762a1bSJed Brown   ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
168*c4762a1bSJed Brown   ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
169*c4762a1bSJed Brown   PetscFunctionReturn(0);
170*c4762a1bSJed Brown }
171*c4762a1bSJed Brown 
172*c4762a1bSJed Brown PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
173*c4762a1bSJed Brown {
174*c4762a1bSJed Brown   PetscErrorCode    ierr;
175*c4762a1bSJed Brown   PetscScalar       sensip;
176*c4762a1bSJed Brown   const PetscScalar *x,*y;
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown   PetscFunctionBegin;
179*c4762a1bSJed Brown   ierr = VecGetArrayRead(lambda,&x);CHKERRQ(ierr);
180*c4762a1bSJed Brown   ierr = VecGetArrayRead(mu,&y);CHKERRQ(ierr);
181*c4762a1bSJed Brown   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
182*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);CHKERRQ(ierr);
183*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(lambda,&x);CHKERRQ(ierr);
184*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(mu,&y);CHKERRQ(ierr);
185*c4762a1bSJed Brown   PetscFunctionReturn(0);
186*c4762a1bSJed Brown }
187*c4762a1bSJed Brown 
188*c4762a1bSJed Brown int main(int argc,char **argv)
189*c4762a1bSJed Brown {
190*c4762a1bSJed Brown   TS             ts,quadts;     /* ODE integrator */
191*c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
192*c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
193*c4762a1bSJed Brown   Mat            Jacp;          /* Jacobian matrix */
194*c4762a1bSJed Brown   Mat            DRDU,DRDP;
195*c4762a1bSJed Brown   PetscErrorCode ierr;
196*c4762a1bSJed Brown   PetscMPIInt    size;
197*c4762a1bSJed Brown   PetscInt       n = 2;
198*c4762a1bSJed Brown   AppCtx         ctx;
199*c4762a1bSJed Brown   PetscScalar    *u;
200*c4762a1bSJed Brown   PetscReal      du[2] = {0.0,0.0};
201*c4762a1bSJed Brown   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
202*c4762a1bSJed Brown   PetscReal      ftime;
203*c4762a1bSJed Brown   PetscInt       steps;
204*c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr;
205*c4762a1bSJed Brown   Vec            lambda[1],q,mu[1];
206*c4762a1bSJed Brown 
207*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208*c4762a1bSJed Brown      Initialize program
209*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
211*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
212*c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215*c4762a1bSJed Brown     Create necessary matrix and vectors
216*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
222*c4762a1bSJed Brown 
223*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
224*c4762a1bSJed Brown 
225*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr);
226*c4762a1bSJed Brown   ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
227*c4762a1bSJed Brown   ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr);
228*c4762a1bSJed Brown   ierr = MatSetUp(Jacp);CHKERRQ(ierr);
229*c4762a1bSJed Brown 
230*c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP);CHKERRQ(ierr);
231*c4762a1bSJed Brown   ierr = MatSetUp(DRDP);CHKERRQ(ierr);
232*c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU);CHKERRQ(ierr);
233*c4762a1bSJed Brown   ierr = MatSetUp(DRDU);CHKERRQ(ierr);
234*c4762a1bSJed Brown 
235*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236*c4762a1bSJed Brown     Set runtime options
237*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
239*c4762a1bSJed Brown   {
240*c4762a1bSJed Brown     ctx.beta    = 2;
241*c4762a1bSJed Brown     ctx.c       = 10000.0;
242*c4762a1bSJed Brown     ctx.u_s     = 1.0;
243*c4762a1bSJed Brown     ctx.omega_s = 1.0;
244*c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
245*c4762a1bSJed Brown     ctx.H       = 5.0;
246*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
247*c4762a1bSJed Brown     ctx.D       = 5.0;
248*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
249*c4762a1bSJed Brown     ctx.E       = 1.1378;
250*c4762a1bSJed Brown     ctx.V       = 1.0;
251*c4762a1bSJed Brown     ctx.X       = 0.545;
252*c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
253*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
254*c4762a1bSJed Brown     ctx.Pm      = 1.1;
255*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
256*c4762a1bSJed Brown     ctx.tf      = 0.1;
257*c4762a1bSJed Brown     ctx.tcl     = 0.2;
258*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
259*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
260*c4762a1bSJed Brown     ierr        = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr);
261*c4762a1bSJed Brown     if (ensemble) {
262*c4762a1bSJed Brown       ctx.tf      = -1;
263*c4762a1bSJed Brown       ctx.tcl     = -1;
264*c4762a1bSJed Brown     }
265*c4762a1bSJed Brown 
266*c4762a1bSJed Brown     ierr = VecGetArray(U,&u);CHKERRQ(ierr);
267*c4762a1bSJed Brown     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
268*c4762a1bSJed Brown     u[1] = 1.0;
269*c4762a1bSJed Brown     ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr);
270*c4762a1bSJed Brown     n    = 2;
271*c4762a1bSJed Brown     ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr);
272*c4762a1bSJed Brown     u[0] += du[0];
273*c4762a1bSJed Brown     u[1] += du[1];
274*c4762a1bSJed Brown     ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
275*c4762a1bSJed Brown     if (flg1 || flg2) {
276*c4762a1bSJed Brown       ctx.tf      = -1;
277*c4762a1bSJed Brown       ctx.tcl     = -1;
278*c4762a1bSJed Brown     }
279*c4762a1bSJed Brown   }
280*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
283*c4762a1bSJed Brown      Create timestepping solver context
284*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
285*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
286*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
287*c4762a1bSJed Brown   ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
288*c4762a1bSJed Brown   ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
289*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr);
290*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr);
291*c4762a1bSJed Brown   ierr = TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts);CHKERRQ(ierr);
292*c4762a1bSJed Brown   ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr);
293*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr);
294*c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx);CHKERRQ(ierr);
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
297*c4762a1bSJed Brown      Set initial conditions
298*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
299*c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
302*c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
303*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
304*c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
305*c4762a1bSJed Brown 
306*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr);
307*c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
308*c4762a1bSJed Brown   ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr);
309*c4762a1bSJed Brown   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
310*c4762a1bSJed Brown   ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr);
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown   ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr);
313*c4762a1bSJed Brown   ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
314*c4762a1bSJed Brown   x_ptr[0] = -1.0;
315*c4762a1bSJed Brown   ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
316*c4762a1bSJed Brown   ierr = TSSetCostGradients(ts,1,lambda,mu);CHKERRQ(ierr);
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
319*c4762a1bSJed Brown      Set solver options
320*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
321*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,10.0);CHKERRQ(ierr);
322*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
323*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr);
324*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
325*c4762a1bSJed Brown 
326*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
327*c4762a1bSJed Brown      Solve nonlinear system
328*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
329*c4762a1bSJed Brown   if (ensemble) {
330*c4762a1bSJed Brown     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
331*c4762a1bSJed Brown       ierr = VecGetArray(U,&u);CHKERRQ(ierr);
332*c4762a1bSJed Brown       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
333*c4762a1bSJed Brown       u[1] = ctx.omega_s;
334*c4762a1bSJed Brown       u[0] += du[0];
335*c4762a1bSJed Brown       u[1] += du[1];
336*c4762a1bSJed Brown       ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
337*c4762a1bSJed Brown       ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr);
338*c4762a1bSJed Brown       ierr = TSSolve(ts,U);CHKERRQ(ierr);
339*c4762a1bSJed Brown     }
340*c4762a1bSJed Brown   } else {
341*c4762a1bSJed Brown     ierr = TSSolve(ts,U);CHKERRQ(ierr);
342*c4762a1bSJed Brown   }
343*c4762a1bSJed Brown   ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
344*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
345*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
346*c4762a1bSJed Brown 
347*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
348*c4762a1bSJed Brown      Adjoint model starts here
349*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
350*c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
351*c4762a1bSJed Brown   ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr);
352*c4762a1bSJed Brown   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
353*c4762a1bSJed Brown   ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr);
354*c4762a1bSJed Brown 
355*c4762a1bSJed Brown   ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
356*c4762a1bSJed Brown   x_ptr[0] = -1.0;
357*c4762a1bSJed Brown   ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
358*c4762a1bSJed Brown 
359*c4762a1bSJed Brown   /*   Set RHS JacobianP */
360*c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr);
361*c4762a1bSJed Brown 
362*c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
363*c4762a1bSJed Brown 
364*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");CHKERRQ(ierr);
365*c4762a1bSJed Brown   ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
366*c4762a1bSJed Brown   ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
367*c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
368*c4762a1bSJed Brown   ierr = VecView(q,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
369*c4762a1bSJed Brown   ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr);
370*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));CHKERRQ(ierr);
371*c4762a1bSJed Brown   ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr);
372*c4762a1bSJed Brown 
373*c4762a1bSJed Brown   ierr = ComputeSensiP(lambda[0],mu[0],&ctx);CHKERRQ(ierr);
374*c4762a1bSJed Brown 
375*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
376*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
377*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
378*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
379*c4762a1bSJed Brown   ierr = MatDestroy(&Jacp);CHKERRQ(ierr);
380*c4762a1bSJed Brown   ierr = MatDestroy(&DRDU);CHKERRQ(ierr);
381*c4762a1bSJed Brown   ierr = MatDestroy(&DRDP);CHKERRQ(ierr);
382*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
383*c4762a1bSJed Brown   ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
384*c4762a1bSJed Brown   ierr = VecDestroy(&mu[0]);CHKERRQ(ierr);
385*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
386*c4762a1bSJed Brown   ierr = PetscFinalize();
387*c4762a1bSJed Brown   return ierr;
388*c4762a1bSJed Brown }
389*c4762a1bSJed Brown 
390*c4762a1bSJed Brown 
391*c4762a1bSJed Brown /*TEST
392*c4762a1bSJed Brown 
393*c4762a1bSJed Brown    build:
394*c4762a1bSJed Brown       requires: !complex
395*c4762a1bSJed Brown 
396*c4762a1bSJed Brown    test:
397*c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none
398*c4762a1bSJed Brown 
399*c4762a1bSJed Brown TEST*/
400