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