xref: /petsc/src/ts/tutorials/power_grid/ex3opt_fd.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown 
2*c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\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 F*/
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown /*
14*c4762a1bSJed Brown   Solve the same optimization problem as in ex3opt.c.
15*c4762a1bSJed Brown   Use finite difference to approximate the gradients.
16*c4762a1bSJed Brown */
17*c4762a1bSJed Brown #include <petsctao.h>
18*c4762a1bSJed Brown #include <petscts.h>
19*c4762a1bSJed Brown #include "ex3.h"
20*c4762a1bSJed Brown 
21*c4762a1bSJed Brown PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);
22*c4762a1bSJed Brown 
23*c4762a1bSJed Brown PetscErrorCode monitor(Tao tao,AppCtx *ctx)
24*c4762a1bSJed Brown {
25*c4762a1bSJed Brown   FILE               *fp;
26*c4762a1bSJed Brown   PetscInt           iterate;
27*c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
28*c4762a1bSJed Brown   Vec                X,G;
29*c4762a1bSJed Brown   const PetscScalar  *x,*g;
30*c4762a1bSJed Brown   TaoConvergedReason reason;
31*c4762a1bSJed Brown   PetscErrorCode     ierr;
32*c4762a1bSJed Brown 
33*c4762a1bSJed Brown   PetscFunctionBeginUser;
34*c4762a1bSJed Brown   ierr = TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason);CHKERRQ(ierr);
35*c4762a1bSJed Brown   ierr = TaoGetSolutionVector(tao,&X);CHKERRQ(ierr);
36*c4762a1bSJed Brown   ierr = TaoGetGradientVector(tao,&G);CHKERRQ(ierr);
37*c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
38*c4762a1bSJed Brown   ierr = VecGetArrayRead(G,&g);CHKERRQ(ierr);
39*c4762a1bSJed Brown   fp = fopen("ex3opt_fd_conv.out","a");
40*c4762a1bSJed Brown   ierr = PetscFPrintf(PETSC_COMM_WORLD,fp,"%d %g %.12lf %.12lf\n",iterate,gnorm,x[0],g[0]);CHKERRQ(ierr);
41*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
42*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(G,&g);CHKERRQ(ierr);
43*c4762a1bSJed Brown   fclose(fp);
44*c4762a1bSJed Brown   PetscFunctionReturn(0);
45*c4762a1bSJed Brown }
46*c4762a1bSJed Brown 
47*c4762a1bSJed Brown int main(int argc,char **argv)
48*c4762a1bSJed Brown {
49*c4762a1bSJed Brown   Vec                p;
50*c4762a1bSJed Brown   PetscScalar        *x_ptr;
51*c4762a1bSJed Brown   PetscErrorCode     ierr;
52*c4762a1bSJed Brown   PetscMPIInt        size;
53*c4762a1bSJed Brown   AppCtx             ctx;
54*c4762a1bSJed Brown   Vec                lowerb,upperb;
55*c4762a1bSJed Brown   Tao                tao;
56*c4762a1bSJed Brown   KSP                ksp;
57*c4762a1bSJed Brown   PC                 pc;
58*c4762a1bSJed Brown   PetscBool          printtofile;
59*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60*c4762a1bSJed Brown      Initialize program
61*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
63*c4762a1bSJed Brown   PetscFunctionBeginUser;
64*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
65*c4762a1bSJed Brown   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
66*c4762a1bSJed Brown 
67*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
68*c4762a1bSJed Brown     Set runtime options
69*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
71*c4762a1bSJed Brown   {
72*c4762a1bSJed Brown     ctx.beta    = 2;
73*c4762a1bSJed Brown     ctx.c       = 10000.0;
74*c4762a1bSJed Brown     ctx.u_s     = 1.0;
75*c4762a1bSJed Brown     ctx.omega_s = 1.0;
76*c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
77*c4762a1bSJed Brown     ctx.H       = 5.0;
78*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr);
79*c4762a1bSJed Brown     ctx.D       = 5.0;
80*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr);
81*c4762a1bSJed Brown     ctx.E       = 1.1378;
82*c4762a1bSJed Brown     ctx.V       = 1.0;
83*c4762a1bSJed Brown     ctx.X       = 0.545;
84*c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
85*c4762a1bSJed Brown     ctx.Pmax_ini = ctx.Pmax;
86*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr);
87*c4762a1bSJed Brown     ctx.Pm      = 1.06;
88*c4762a1bSJed Brown     ierr        = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr);
89*c4762a1bSJed Brown     ctx.tf      = 0.1;
90*c4762a1bSJed Brown     ctx.tcl     = 0.2;
91*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr);
92*c4762a1bSJed Brown     ierr        = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr);
93*c4762a1bSJed Brown     printtofile = PETSC_FALSE;
94*c4762a1bSJed Brown     ierr        = PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL);CHKERRQ(ierr);
95*c4762a1bSJed Brown   }
96*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
99*c4762a1bSJed Brown   ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr);
100*c4762a1bSJed Brown   ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr);
101*c4762a1bSJed Brown   if(printtofile) {
102*c4762a1bSJed Brown     ierr = TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL);CHKERRQ(ierr);
103*c4762a1bSJed Brown   }
104*c4762a1bSJed Brown   ierr = TaoSetMaximumIterations(tao,30);CHKERRQ(ierr);
105*c4762a1bSJed Brown   /*
106*c4762a1bSJed Brown      Optimization starts
107*c4762a1bSJed Brown   */
108*c4762a1bSJed Brown   /* Set initial solution guess */
109*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_WORLD,1,&p);CHKERRQ(ierr);
110*c4762a1bSJed Brown   ierr = VecGetArray(p,&x_ptr);CHKERRQ(ierr);
111*c4762a1bSJed Brown   x_ptr[0]   = ctx.Pm;
112*c4762a1bSJed Brown   ierr = VecRestoreArray(p,&x_ptr);CHKERRQ(ierr);
113*c4762a1bSJed Brown 
114*c4762a1bSJed Brown   ierr = TaoSetInitialVector(tao,p);CHKERRQ(ierr);
115*c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
116*c4762a1bSJed Brown   ierr = TaoSetObjectiveRoutine(tao,FormFunction,(void *)&ctx);CHKERRQ(ierr);
117*c4762a1bSJed Brown   ierr = TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,(void *)&ctx);CHKERRQ(ierr);
118*c4762a1bSJed Brown 
119*c4762a1bSJed Brown   /* Set bounds for the optimization */
120*c4762a1bSJed Brown   ierr = VecDuplicate(p,&lowerb);CHKERRQ(ierr);
121*c4762a1bSJed Brown   ierr = VecDuplicate(p,&upperb);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = VecGetArray(lowerb,&x_ptr);CHKERRQ(ierr);
123*c4762a1bSJed Brown   x_ptr[0] = 0.;
124*c4762a1bSJed Brown   ierr = VecRestoreArray(lowerb,&x_ptr);CHKERRQ(ierr);
125*c4762a1bSJed Brown   ierr = VecGetArray(upperb,&x_ptr);CHKERRQ(ierr);
126*c4762a1bSJed Brown   x_ptr[0] = 1.1;
127*c4762a1bSJed Brown   ierr = VecRestoreArray(upperb,&x_ptr);CHKERRQ(ierr);
128*c4762a1bSJed Brown   ierr = TaoSetVariableBounds(tao,lowerb,upperb);CHKERRQ(ierr);
129*c4762a1bSJed Brown 
130*c4762a1bSJed Brown   /* Check for any TAO command line options */
131*c4762a1bSJed Brown   ierr = TaoSetFromOptions(tao);CHKERRQ(ierr);
132*c4762a1bSJed Brown   ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr);
133*c4762a1bSJed Brown   if (ksp) {
134*c4762a1bSJed Brown     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
135*c4762a1bSJed Brown     ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
136*c4762a1bSJed Brown   }
137*c4762a1bSJed Brown 
138*c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
139*c4762a1bSJed Brown   ierr = TaoSolve(tao); CHKERRQ(ierr);
140*c4762a1bSJed Brown 
141*c4762a1bSJed Brown   ierr = VecView(p,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
142*c4762a1bSJed Brown 
143*c4762a1bSJed Brown   /* Free TAO data structures */
144*c4762a1bSJed Brown   ierr = TaoDestroy(&tao);CHKERRQ(ierr);
145*c4762a1bSJed Brown   ierr = VecDestroy(&p);CHKERRQ(ierr);
146*c4762a1bSJed Brown   ierr = VecDestroy(&lowerb);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = VecDestroy(&upperb);CHKERRQ(ierr);
148*c4762a1bSJed Brown   ierr = PetscFinalize();
149*c4762a1bSJed Brown   return ierr;
150*c4762a1bSJed Brown }
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown /* ------------------------------------------------------------------ */
153*c4762a1bSJed Brown /*
154*c4762a1bSJed Brown    FormFunction - Evaluates the function and corresponding gradient.
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown    Input Parameters:
157*c4762a1bSJed Brown    tao - the Tao context
158*c4762a1bSJed Brown    X   - the input vector
159*c4762a1bSJed Brown    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown    Output Parameters:
162*c4762a1bSJed Brown    f   - the newly evaluated function
163*c4762a1bSJed Brown */
164*c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
165*c4762a1bSJed Brown {
166*c4762a1bSJed Brown   AppCtx            *ctx = (AppCtx*)ctx0;
167*c4762a1bSJed Brown   TS                ts,quadts;
168*c4762a1bSJed Brown   Vec               U;             /* solution will be stored here */
169*c4762a1bSJed Brown   Mat               A;             /* Jacobian matrix */
170*c4762a1bSJed Brown   PetscErrorCode    ierr;
171*c4762a1bSJed Brown   PetscInt          n = 2;
172*c4762a1bSJed Brown   PetscReal         ftime;
173*c4762a1bSJed Brown   PetscInt          steps;
174*c4762a1bSJed Brown   PetscScalar       *u;
175*c4762a1bSJed Brown   const PetscScalar *x_ptr,*qx_ptr;
176*c4762a1bSJed Brown   Vec               q;
177*c4762a1bSJed Brown   PetscInt          direction[2];
178*c4762a1bSJed Brown   PetscBool         terminate[2];
179*c4762a1bSJed Brown 
180*c4762a1bSJed Brown   ierr = VecGetArrayRead(P,&x_ptr);CHKERRQ(ierr);
181*c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
182*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(P,&x_ptr);CHKERRQ(ierr);
183*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
184*c4762a1bSJed Brown     Create necessary matrix and vectors
185*c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
186*c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
187*c4762a1bSJed Brown   ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
188*c4762a1bSJed Brown   ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
190*c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
191*c4762a1bSJed Brown 
192*c4762a1bSJed Brown   ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr);
193*c4762a1bSJed Brown 
194*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195*c4762a1bSJed Brown      Create timestepping solver context
196*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
198*c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
199*c4762a1bSJed Brown   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
200*c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);CHKERRQ(ierr);
201*c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx);CHKERRQ(ierr);
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204*c4762a1bSJed Brown      Set initial conditions
205*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206*c4762a1bSJed Brown   ierr = VecGetArray(U,&u);CHKERRQ(ierr);
207*c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
208*c4762a1bSJed Brown   u[1] = 1.0;
209*c4762a1bSJed Brown   ierr = VecRestoreArray(U,&u);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = TSSetSolution(ts,U);CHKERRQ(ierr);
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213*c4762a1bSJed Brown      Set solver options
214*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
215*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,1.0);CHKERRQ(ierr);
216*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
217*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,0.03125);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts);CHKERRQ(ierr);
219*c4762a1bSJed Brown   ierr = TSGetSolution(quadts,&q);CHKERRQ(ierr);
220*c4762a1bSJed Brown   ierr = VecSet(q,0.0);CHKERRQ(ierr);
221*c4762a1bSJed Brown   ierr = TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
223*c4762a1bSJed Brown 
224*c4762a1bSJed Brown   direction[0] = direction[1] = 1;
225*c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
226*c4762a1bSJed Brown 
227*c4762a1bSJed Brown   ierr = TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)ctx);CHKERRQ(ierr);
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230*c4762a1bSJed Brown      Solve nonlinear system
231*c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
232*c4762a1bSJed Brown   ierr = TSSolve(ts,U);CHKERRQ(ierr);
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = VecGetArrayRead(q,&qx_ptr);CHKERRQ(ierr);
237*c4762a1bSJed Brown   *f   = -ctx->Pm + qx_ptr[0];
238*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(q,&qx_ptr);CHKERRQ(ierr);
239*c4762a1bSJed Brown 
240*c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
241*c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
242*c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
243*c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = VecDestroy(&U);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
246*c4762a1bSJed Brown   return 0;
247*c4762a1bSJed Brown }
248*c4762a1bSJed Brown 
249*c4762a1bSJed Brown /*TEST
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown    build:
252*c4762a1bSJed Brown       requires: !complex !single
253*c4762a1bSJed Brown 
254*c4762a1bSJed Brown    test:
255*c4762a1bSJed Brown       args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3
256*c4762a1bSJed Brown 
257*c4762a1bSJed Brown TEST*/
258