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