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