xref: /petsc/src/ts/tutorials/power_grid/ex9opt.c (revision f32d636044d38655283e0a09d0e1a46bb2053733)
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 
12c4762a1bSJed Brown 
13c4762a1bSJed Brown   Ensemble of initial conditions
14c4762a1bSJed Brown    ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3      -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
15c4762a1bSJed Brown 
16c4762a1bSJed Brown   Fault at .1 seconds
17c4762a1bSJed Brown    ./ex2           -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   Initial conditions same as when fault is ended
20c4762a1bSJed Brown    ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05  -ts_adapt_dt_max .01  -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly
21c4762a1bSJed Brown 
22c4762a1bSJed Brown 
23c4762a1bSJed Brown F*/
24c4762a1bSJed Brown 
25c4762a1bSJed Brown /*
26c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
27c4762a1bSJed Brown    file automatically includes:
28c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
29c4762a1bSJed Brown      petscmat.h - matrices
30c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
31c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
32c4762a1bSJed Brown      petscksp.h   - linear solvers
33c4762a1bSJed Brown */
34c4762a1bSJed Brown 
35c4762a1bSJed Brown #include <petsctao.h>
36c4762a1bSJed Brown #include <petscts.h>
37c4762a1bSJed Brown 
38c4762a1bSJed Brown typedef struct {
39c4762a1bSJed Brown   TS          ts;
40c4762a1bSJed Brown   PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c;
41c4762a1bSJed Brown   PetscInt    beta;
42c4762a1bSJed Brown   PetscReal   tf,tcl,dt;
43c4762a1bSJed Brown } AppCtx;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);
46c4762a1bSJed Brown PetscErrorCode FormGradient(Tao,Vec,Vec,void*);
47c4762a1bSJed Brown 
48c4762a1bSJed Brown /*
49c4762a1bSJed Brown      Defines the ODE passed to the ODE solver
50c4762a1bSJed Brown */
51c4762a1bSJed Brown static 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   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 */
62c4762a1bSJed Brown   else Pmax = ctx->Pmax;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown   f[0] = ctx->omega_b*(u[1] - ctx->omega_s);
65c4762a1bSJed Brown   f[1] = (-Pmax*PetscSinScalar(u[0]) - ctx->D*(u[1] - ctx->omega_s) + ctx->Pm)*ctx->omega_s/(2.0*ctx->H);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
69c4762a1bSJed Brown   PetscFunctionReturn(0);
70c4762a1bSJed Brown }
71c4762a1bSJed Brown 
72c4762a1bSJed Brown /*
73c4762a1bSJed Brown      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
74c4762a1bSJed Brown */
75c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
76c4762a1bSJed Brown {
77c4762a1bSJed Brown   PetscErrorCode    ierr;
78c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
79c4762a1bSJed Brown   PetscScalar       J[2][2],Pmax;
80c4762a1bSJed Brown   const PetscScalar *u;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   PetscFunctionBegin;
83c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
84c4762a1bSJed 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 */
85c4762a1bSJed Brown   else Pmax = ctx->Pmax;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   J[0][0] = 0;                                  J[0][1] = ctx->omega_b;
88c4762a1bSJed 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);
89c4762a1bSJed Brown 
90c4762a1bSJed Brown   ierr    = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
91c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
94c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
95c4762a1bSJed Brown   if (A != B) {
96c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
97c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
98c4762a1bSJed Brown   }
99c4762a1bSJed Brown   PetscFunctionReturn(0);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
103c4762a1bSJed Brown {
104c4762a1bSJed Brown   PetscErrorCode ierr;
105c4762a1bSJed Brown   PetscInt       row[] = {0,1},col[]={0};
106c4762a1bSJed Brown   PetscScalar    J[2][1];
107c4762a1bSJed Brown   AppCtx         *ctx=(AppCtx*)ctx0;
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   PetscFunctionBeginUser;
110c4762a1bSJed Brown   J[0][0] = 0;
111c4762a1bSJed Brown   J[1][0] = ctx->omega_s/(2.0*ctx->H);
112c4762a1bSJed Brown   ierr    = MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
113c4762a1bSJed Brown   ierr    = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
114c4762a1bSJed Brown   ierr    = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
115c4762a1bSJed Brown   PetscFunctionReturn(0);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
118c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
119c4762a1bSJed Brown {
120c4762a1bSJed Brown   PetscErrorCode    ierr;
121c4762a1bSJed Brown   PetscScalar       *r;
122c4762a1bSJed Brown   const PetscScalar *u;
123c4762a1bSJed Brown 
124c4762a1bSJed Brown   PetscFunctionBegin;
125c4762a1bSJed Brown   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
126c4762a1bSJed Brown   ierr = VecGetArray(R,&r);CHKERRQ(ierr);
127c4762a1bSJed Brown   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);CHKERRQ(ierr);
128c4762a1bSJed Brown   ierr = VecRestoreArray(R,&r);CHKERRQ(ierr);
129c4762a1bSJed Brown   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
130c4762a1bSJed Brown   PetscFunctionReturn(0);
131c4762a1bSJed Brown }
132c4762a1bSJed Brown 
133c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,AppCtx *ctx)
134c4762a1bSJed Brown {
135c4762a1bSJed Brown   PetscErrorCode    ierr;
136c4762a1bSJed Brown   PetscScalar       ru[1];
137c4762a1bSJed Brown   const PetscScalar *u;
138c4762a1bSJed Brown   PetscInt          row[] = {0},col[] = {0};
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   PetscFunctionBegin;
141c4762a1bSJed Brown   ierr  = VecGetArrayRead(U,&u);CHKERRQ(ierr);
142c4762a1bSJed Brown   ru[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);CHKERRQ(ierr);
143c4762a1bSJed Brown   ierr  = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
144c4762a1bSJed Brown   ierr  = MatSetValues(DRDU,1,row,1,col,ru,INSERT_VALUES);CHKERRQ(ierr);
145c4762a1bSJed Brown   ierr  = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
146c4762a1bSJed Brown   ierr  = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
147c4762a1bSJed Brown   PetscFunctionReturn(0);
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,AppCtx *ctx)
151c4762a1bSJed Brown {
152c4762a1bSJed Brown   PetscErrorCode ierr;
153c4762a1bSJed Brown 
154c4762a1bSJed Brown   PetscFunctionBegin;
155c4762a1bSJed Brown   ierr = MatZeroEntries(DRDP);CHKERRQ(ierr);
156c4762a1bSJed Brown   ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
157c4762a1bSJed Brown   ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
158c4762a1bSJed Brown   PetscFunctionReturn(0);
159c4762a1bSJed Brown }
160c4762a1bSJed Brown 
161c4762a1bSJed Brown PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
162c4762a1bSJed Brown {
163c4762a1bSJed Brown   PetscErrorCode    ierr;
164c4762a1bSJed Brown   PetscScalar       *y,sensip;
165c4762a1bSJed Brown   const PetscScalar *x;
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   PetscFunctionBegin;
168c4762a1bSJed Brown   ierr = VecGetArrayRead(lambda,&x);CHKERRQ(ierr);
169c4762a1bSJed Brown   ierr = VecGetArray(mu,&y);CHKERRQ(ierr);
170c4762a1bSJed Brown   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
171c4762a1bSJed Brown   y[0] = sensip;
172c4762a1bSJed Brown   ierr = VecRestoreArray(mu,&y);CHKERRQ(ierr);
173c4762a1bSJed Brown   ierr = VecRestoreArrayRead(lambda,&x);CHKERRQ(ierr);
174c4762a1bSJed Brown   PetscFunctionReturn(0);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
177c4762a1bSJed Brown int main(int argc,char **argv)
178c4762a1bSJed Brown {
179c4762a1bSJed Brown   Vec            p;
180c4762a1bSJed Brown   PetscScalar    *x_ptr;
181c4762a1bSJed Brown   PetscErrorCode ierr;
182c4762a1bSJed Brown   PetscMPIInt    size;
183c4762a1bSJed Brown   AppCtx         ctx;
184c4762a1bSJed Brown   Vec            lowerb,upperb;
185c4762a1bSJed Brown   Tao            tao;
186c4762a1bSJed Brown   KSP            ksp;
187c4762a1bSJed Brown   PC             pc;
188c4762a1bSJed Brown   Vec            U,lambda[1],mu[1];
189c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
190c4762a1bSJed Brown   Mat            Jacp;          /* Jacobian matrix */
191c4762a1bSJed Brown   Mat            DRDU,DRDP;
192c4762a1bSJed Brown   PetscInt       n = 2;
193c4762a1bSJed Brown   TS             quadts;
194c4762a1bSJed Brown 
195c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196c4762a1bSJed Brown      Initialize program
197c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
199c4762a1bSJed Brown   PetscFunctionBeginUser;
200c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
201c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown     Set runtime options
205c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
207c4762a1bSJed Brown   {
208c4762a1bSJed Brown     ctx.beta    = 2;
209c4762a1bSJed Brown     ctx.c       = PetscRealConstant(10000.0);
210c4762a1bSJed Brown     ctx.u_s     = PetscRealConstant(1.0);
211c4762a1bSJed Brown     ctx.omega_s = PetscRealConstant(1.0);
212c4762a1bSJed Brown     ctx.omega_b = PetscRealConstant(120.0)*PETSC_PI;
213c4762a1bSJed Brown     ctx.H       = PetscRealConstant(5.0);
214c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
215c4762a1bSJed Brown     ctx.D       = PetscRealConstant(5.0);
216c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
217c4762a1bSJed Brown     ctx.E       = PetscRealConstant(1.1378);
218c4762a1bSJed Brown     ctx.V       = PetscRealConstant(1.0);
219c4762a1bSJed Brown     ctx.X       = PetscRealConstant(0.545);
220c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
221c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
222c4762a1bSJed Brown     ctx.Pm      = PetscRealConstant(1.0194);
223c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
224c4762a1bSJed Brown     ctx.tf      = PetscRealConstant(0.1);
225c4762a1bSJed Brown     ctx.tcl     = PetscRealConstant(0.2);
226c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
227c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   }
230c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233c4762a1bSJed Brown     Create necessary matrix and vectors
234c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
235c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
236c4762a1bSJed Brown   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
238c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
239c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr);
244c4762a1bSJed Brown   ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr);
245c4762a1bSJed Brown   ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr);
246c4762a1bSJed Brown   ierr = MatSetUp(Jacp);CHKERRQ(ierr);
247c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&DRDP);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr = MatSetUp(DRDP);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,2,NULL,&DRDU);CHKERRQ(ierr);
250c4762a1bSJed Brown   ierr = MatSetUp(DRDU);CHKERRQ(ierr);
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253c4762a1bSJed Brown      Create timestepping solver context
254c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
255c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ctx.ts);CHKERRQ(ierr);
256c4762a1bSJed Brown   ierr = TSSetProblemType(ctx.ts,TS_NONLINEAR);CHKERRQ(ierr);
257c4762a1bSJed Brown   ierr = TSSetEquationType(ctx.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */
258c4762a1bSJed Brown   ierr = TSSetType(ctx.ts,TSRK);CHKERRQ(ierr);
259c4762a1bSJed Brown   ierr = TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx);CHKERRQ(ierr);
260c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ctx.ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx);CHKERRQ(ierr);
261c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr);
264c4762a1bSJed Brown   ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr);
265c4762a1bSJed Brown   ierr = TSSetCostGradients(ctx.ts,1,lambda,mu);CHKERRQ(ierr);
266c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(ctx.ts,Jacp,RHSJacobianP,&ctx);CHKERRQ(ierr);
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269c4762a1bSJed Brown      Set solver options
270c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
271c4762a1bSJed Brown   ierr = TSSetMaxTime(ctx.ts,PetscRealConstant(1.0));CHKERRQ(ierr);
272c4762a1bSJed Brown   ierr = TSSetTimeStep(ctx.ts,PetscRealConstant(0.01));CHKERRQ(ierr);
273c4762a1bSJed Brown   ierr = TSSetFromOptions(ctx.ts);CHKERRQ(ierr);
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   ierr = TSGetTimeStep(ctx.ts,&ctx.dt);CHKERRQ(ierr); /* save the stepsize */
276c4762a1bSJed Brown 
277c4762a1bSJed Brown   ierr = TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&quadts);CHKERRQ(ierr);
278c4762a1bSJed Brown   ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx);CHKERRQ(ierr);
279c4762a1bSJed Brown   ierr = TSSetRHSJacobian(quadts,DRDU,DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx);CHKERRQ(ierr);
280c4762a1bSJed Brown   ierr = TSSetRHSJacobianP(quadts,DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&ctx);CHKERRQ(ierr);
281c4762a1bSJed Brown   ierr = TSSetSolution(ctx.ts,U);CHKERRQ(ierr);
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
284c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
285c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
286c4762a1bSJed Brown 
287c4762a1bSJed Brown   /*
288c4762a1bSJed Brown      Optimization starts
289c4762a1bSJed Brown   */
290c4762a1bSJed Brown   /* Set initial solution guess */
291c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_WORLD,1,&p);CHKERRQ(ierr);
292c4762a1bSJed Brown   ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr);
293c4762a1bSJed Brown   x_ptr[0]   = ctx.Pm;
294c4762a1bSJed Brown   ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr);
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,p);CHKERRQ(ierr);
297c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
298c4762a1bSJed Brown   ierr = TaoSetObjectiveRoutine(tao,FormFunction,(void *)&ctx);CHKERRQ(ierr);
299c4762a1bSJed Brown   ierr = TaoSetGradientRoutine(tao,FormGradient,(void *)&ctx);CHKERRQ(ierr);
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   /* Set bounds for the optimization */
302c4762a1bSJed Brown   ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr);
303c4762a1bSJed Brown   ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr);
304c4762a1bSJed Brown   ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr);
305c4762a1bSJed Brown   x_ptr[0] = 0.;
306c4762a1bSJed Brown   ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr);
307c4762a1bSJed Brown   ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr);
308c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(1.1);
309c4762a1bSJed Brown   ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr);
310c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr);
311c4762a1bSJed Brown 
312c4762a1bSJed Brown   /* Check for any TAO command line options */
313c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
314c4762a1bSJed Brown   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
315c4762a1bSJed Brown   if (ksp) {
316c4762a1bSJed Brown     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
317c4762a1bSJed Brown     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
318c4762a1bSJed Brown   }
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
321c4762a1bSJed Brown   ierr = TaoSolve(tao); CHKERRQ(ierr);
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
324c4762a1bSJed Brown   /* Free TAO data structures */
325c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
326c4762a1bSJed Brown   ierr = VecDestroy(&p);CHKERRQ(ierr);
327c4762a1bSJed Brown   ierr = VecDestroy(&lowerb);CHKERRQ(ierr);
328c4762a1bSJed Brown   ierr = VecDestroy(&upperb);CHKERRQ(ierr);
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   ierr = TSDestroy(&ctx.ts);CHKERRQ(ierr);
331c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
332c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
333c4762a1bSJed Brown   ierr = MatDestroy(&Jacp);CHKERRQ(ierr);
334c4762a1bSJed Brown   ierr = MatDestroy(&DRDU);CHKERRQ(ierr);
335c4762a1bSJed Brown   ierr = MatDestroy(&DRDP);CHKERRQ(ierr);
336c4762a1bSJed Brown   ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr);
337c4762a1bSJed Brown   ierr = VecDestroy(&mu[0]);CHKERRQ(ierr);
338c4762a1bSJed Brown   ierr = PetscFinalize();
339c4762a1bSJed Brown   return ierr;
340c4762a1bSJed Brown }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown /* ------------------------------------------------------------------ */
343c4762a1bSJed Brown /*
344c4762a1bSJed Brown    FormFunction - Evaluates the function
345c4762a1bSJed Brown 
346c4762a1bSJed Brown    Input Parameters:
347c4762a1bSJed Brown    tao - the Tao context
348c4762a1bSJed Brown    X   - the input vector
349c4762a1bSJed Brown    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
350c4762a1bSJed Brown 
351c4762a1bSJed Brown    Output Parameters:
352c4762a1bSJed Brown    f   - the newly evaluated function
353c4762a1bSJed Brown */
354c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
355c4762a1bSJed Brown {
356c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)ctx0;
357c4762a1bSJed Brown   TS             ts = ctx->ts;
358c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
359c4762a1bSJed Brown   PetscErrorCode ierr;
360c4762a1bSJed Brown   PetscScalar    *u;
361c4762a1bSJed Brown   PetscScalar    *x_ptr;
362c4762a1bSJed Brown   Vec            q;
363c4762a1bSJed Brown 
364c4762a1bSJed Brown   ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
365c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
366c4762a1bSJed Brown   ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
367c4762a1bSJed Brown 
368c4762a1bSJed Brown   /* reset time */
369c4762a1bSJed Brown   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
370c4762a1bSJed Brown   /* reset step counter, this is critical for adjoint solver */
371c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
372c4762a1bSJed Brown   /* reset step size, the step size becomes negative after TSAdjointSolve */
373c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,ctx->dt);CHKERRQ(ierr);
374c4762a1bSJed Brown   /* reinitialize the integral value */
375c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
376c4762a1bSJed Brown   ierr = VecSet(q,0.0);CHKERRQ(ierr);
377c4762a1bSJed Brown 
378c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
379c4762a1bSJed Brown      Set initial conditions
380c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
381c4762a1bSJed Brown   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
382c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
383c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
384c4762a1bSJed Brown   u[1] = PetscRealConstant(1.0);
385c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
388c4762a1bSJed Brown      Solve nonlinear system
389c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
390c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
391c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
392c4762a1bSJed Brown   ierr = VecGetArray(q,&x_ptr);CHKERRQ(ierr);
393c4762a1bSJed Brown   *f   = -ctx->Pm + x_ptr[0];
394c4762a1bSJed Brown   ierr = VecRestoreArray(q,&x_ptr);CHKERRQ(ierr);
395c4762a1bSJed Brown   return 0;
396c4762a1bSJed Brown }
397c4762a1bSJed Brown 
398c4762a1bSJed Brown PetscErrorCode FormGradient(Tao tao,Vec P,Vec G,void *ctx0)
399c4762a1bSJed Brown {
400c4762a1bSJed Brown   AppCtx         *ctx = (AppCtx*)ctx0;
401c4762a1bSJed Brown   TS             ts = ctx->ts;
402c4762a1bSJed Brown   Vec            U;             /* solution will be stored here */
403c4762a1bSJed Brown   PetscErrorCode ierr;
404c4762a1bSJed Brown   PetscReal      ftime;
405c4762a1bSJed Brown   PetscInt       steps;
406c4762a1bSJed Brown   PetscScalar    *u;
407c4762a1bSJed Brown   PetscScalar    *x_ptr,*y_ptr;
408c4762a1bSJed Brown   Vec            *lambda,q,*mu;
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   ierr = VecGetArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
411c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
412c4762a1bSJed Brown   ierr = VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);CHKERRQ(ierr);
413c4762a1bSJed Brown 
414c4762a1bSJed Brown   /* reset time */
415c4762a1bSJed Brown   ierr = TSSetTime(ts,0.0);CHKERRQ(ierr);
416c4762a1bSJed Brown   /* reset step counter, this is critical for adjoint solver */
417c4762a1bSJed Brown   ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr);
418c4762a1bSJed Brown   /* reset step size, the step size becomes negative after TSAdjointSolve */
419c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,ctx->dt);CHKERRQ(ierr);
420c4762a1bSJed Brown   /* reinitialize the integral value */
421c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
422c4762a1bSJed Brown   ierr = VecSet(q,0.0);CHKERRQ(ierr);
423c4762a1bSJed Brown 
424c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
425c4762a1bSJed Brown      Set initial conditions
426c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
427c4762a1bSJed Brown   ierr = TSGetSolution(ts,&U);CHKERRQ(ierr);
428c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
429c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
430c4762a1bSJed Brown   u[1] = PetscRealConstant(1.0);
431c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
432c4762a1bSJed Brown 
433*f32d6360SSatish Balay   /* Set up to save trajectory before TSSetFromOptions() so that TSTrajectory options can be captured */
434c4762a1bSJed Brown   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
435*f32d6360SSatish Balay   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
436c4762a1bSJed Brown 
437c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
438c4762a1bSJed Brown      Solve nonlinear system
439c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
440c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
441c4762a1bSJed Brown 
442c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
443c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
444c4762a1bSJed Brown 
445c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
446c4762a1bSJed Brown      Adjoint model starts here
447c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
448c4762a1bSJed Brown   ierr = TSGetCostGradients(ts,NULL,&lambda,&mu);CHKERRQ(ierr);
449c4762a1bSJed Brown   /*   Set initial conditions for the adjoint integration */
450c4762a1bSJed Brown   ierr = VecGetArray(lambda[0],&y_ptr);CHKERRQ(ierr);
451c4762a1bSJed Brown   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
452c4762a1bSJed Brown   ierr = VecRestoreArray(lambda[0],&y_ptr);CHKERRQ(ierr);
453c4762a1bSJed Brown   ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr);
454c4762a1bSJed Brown   x_ptr[0] = PetscRealConstant(-1.0);
455c4762a1bSJed Brown   ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr);
456c4762a1bSJed Brown 
457c4762a1bSJed Brown   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
458c4762a1bSJed Brown   ierr = TSGetCostIntegral(ts,&q);CHKERRQ(ierr);
459c4762a1bSJed Brown   ierr = ComputeSensiP(lambda[0],mu[0],ctx);CHKERRQ(ierr);
460c4762a1bSJed Brown   ierr = VecCopy(mu[0],G);CHKERRQ(ierr);
461c4762a1bSJed Brown   return 0;
462c4762a1bSJed Brown }
463c4762a1bSJed Brown 
464c4762a1bSJed Brown 
465c4762a1bSJed Brown /*TEST
466c4762a1bSJed Brown 
467c4762a1bSJed Brown    build:
468c4762a1bSJed Brown       requires: !complex
469c4762a1bSJed Brown 
470c4762a1bSJed Brown    test:
471c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason
472c4762a1bSJed Brown 
473c4762a1bSJed Brown    test:
474c4762a1bSJed Brown       suffix: 2
475c4762a1bSJed Brown       args: -viewer_binary_skip_info -ts_adapt_type none -tao_monitor -tao_gatol 0.0 -tao_grtol 1.e-3 -tao_converged_reason -tao_test_gradient
476c4762a1bSJed Brown 
477c4762a1bSJed Brown TEST*/
478