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