xref: /petsc/src/ts/tutorials/power_grid/ex3opt_fd.c (revision d0609ced746bc51b019815ca91d747429db24893)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\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 F*/
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown   Solve the same optimization problem as in ex3opt.c.
15c4762a1bSJed Brown   Use finite difference to approximate the gradients.
16c4762a1bSJed Brown */
17c4762a1bSJed Brown #include <petsctao.h>
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown #include "ex3.h"
20c4762a1bSJed Brown 
21c4762a1bSJed Brown PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);
22c4762a1bSJed Brown 
23c4762a1bSJed Brown PetscErrorCode monitor(Tao tao,AppCtx *ctx)
24c4762a1bSJed Brown {
25c4762a1bSJed Brown   FILE               *fp;
26c4762a1bSJed Brown   PetscInt           iterate;
27c4762a1bSJed Brown   PetscReal          f,gnorm,cnorm,xdiff;
28c4762a1bSJed Brown   Vec                X,G;
29c4762a1bSJed Brown   const PetscScalar  *x,*g;
30c4762a1bSJed Brown   TaoConvergedReason reason;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   PetscFunctionBeginUser;
339566063dSJacob Faibussowitsch   PetscCall(TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason));
349566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao,&X));
359566063dSJacob Faibussowitsch   PetscCall(TaoGetGradient(tao,&G,NULL,NULL));
369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(G,&g));
38c4762a1bSJed Brown   fp = fopen("ex3opt_fd_conv.out","a");
399566063dSJacob Faibussowitsch   PetscCall(PetscFPrintf(PETSC_COMM_WORLD,fp,"%d %g %.12lf %.12lf\n",iterate,gnorm,x[0],g[0]));
409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(G,&g));
42c4762a1bSJed Brown   fclose(fp);
43c4762a1bSJed Brown   PetscFunctionReturn(0);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46c4762a1bSJed Brown int main(int argc,char **argv)
47c4762a1bSJed Brown {
48c4762a1bSJed Brown   Vec                p;
49c4762a1bSJed Brown   PetscScalar        *x_ptr;
50c4762a1bSJed Brown   PetscMPIInt        size;
51c4762a1bSJed Brown   AppCtx             ctx;
52c4762a1bSJed Brown   Vec                lowerb,upperb;
53c4762a1bSJed Brown   Tao                tao;
54c4762a1bSJed Brown   KSP                ksp;
55c4762a1bSJed Brown   PC                 pc;
56c4762a1bSJed Brown   PetscBool          printtofile;
57c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Initialize program
59c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
609566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
61c4762a1bSJed Brown   PetscFunctionBeginUser;
629566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
633c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
64c4762a1bSJed Brown 
65c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
66c4762a1bSJed Brown     Set runtime options
67c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
68*d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
69c4762a1bSJed Brown   {
70c4762a1bSJed Brown     ctx.beta    = 2;
71c4762a1bSJed Brown     ctx.c       = 10000.0;
72c4762a1bSJed Brown     ctx.u_s     = 1.0;
73c4762a1bSJed Brown     ctx.omega_s = 1.0;
74c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
75c4762a1bSJed Brown     ctx.H       = 5.0;
769566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL));
77c4762a1bSJed Brown     ctx.D       = 5.0;
789566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL));
79c4762a1bSJed Brown     ctx.E       = 1.1378;
80c4762a1bSJed Brown     ctx.V       = 1.0;
81c4762a1bSJed Brown     ctx.X       = 0.545;
82c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
83c4762a1bSJed Brown     ctx.Pmax_ini = ctx.Pmax;
849566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL));
85c4762a1bSJed Brown     ctx.Pm      = 1.06;
869566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL));
87c4762a1bSJed Brown     ctx.tf      = 0.1;
88c4762a1bSJed Brown     ctx.tcl     = 0.2;
899566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL));
909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL));
91c4762a1bSJed Brown     printtofile = PETSC_FALSE;
929566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL));
93c4762a1bSJed Brown   }
94*d0609cedSBarry Smith   PetscOptionsEnd();
95c4762a1bSJed Brown 
96c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
979566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
989566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOBLMVM));
99c4762a1bSJed Brown   if (printtofile) {
1009566063dSJacob Faibussowitsch     PetscCall(TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL));
101c4762a1bSJed Brown   }
1029566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(tao,30));
103c4762a1bSJed Brown   /*
104c4762a1bSJed Brown      Optimization starts
105c4762a1bSJed Brown   */
106c4762a1bSJed Brown   /* Set initial solution guess */
1079566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,1,&p));
1089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(p,&x_ptr));
109c4762a1bSJed Brown   x_ptr[0]   = ctx.Pm;
1109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(p,&x_ptr));
111c4762a1bSJed Brown 
1129566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,p));
113c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
1149566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(tao,FormFunction,(void *)&ctx));
1159566063dSJacob Faibussowitsch   PetscCall(TaoSetGradient(tao,NULL,TaoDefaultComputeGradient,(void *)&ctx));
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   /* Set bounds for the optimization */
1189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p,&lowerb));
1199566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p,&upperb));
1209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lowerb,&x_ptr));
121c4762a1bSJed Brown   x_ptr[0] = 0.;
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lowerb,&x_ptr));
1239566063dSJacob Faibussowitsch   PetscCall(VecGetArray(upperb,&x_ptr));
124c4762a1bSJed Brown   x_ptr[0] = 1.1;
1259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(upperb,&x_ptr));
1269566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao,lowerb,upperb));
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   /* Check for any TAO command line options */
1299566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
1309566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao,&ksp));
131c4762a1bSJed Brown   if (ksp) {
1329566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp,&pc));
1339566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc,PCNONE));
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1379566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
138c4762a1bSJed Brown 
1399566063dSJacob Faibussowitsch   PetscCall(VecView(p,PETSC_VIEWER_STDOUT_WORLD));
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* Free TAO data structures */
1429566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&p));
1449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lowerb));
1459566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&upperb));
1469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
147b122ec5aSJacob Faibussowitsch   return 0;
148c4762a1bSJed Brown }
149c4762a1bSJed Brown 
150c4762a1bSJed Brown /* ------------------------------------------------------------------ */
151c4762a1bSJed Brown /*
152c4762a1bSJed Brown    FormFunction - Evaluates the function and corresponding gradient.
153c4762a1bSJed Brown 
154c4762a1bSJed Brown    Input Parameters:
155c4762a1bSJed Brown    tao - the Tao context
156c4762a1bSJed Brown    X   - the input vector
157a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
158c4762a1bSJed Brown 
159c4762a1bSJed Brown    Output Parameters:
160c4762a1bSJed Brown    f   - the newly evaluated function
161c4762a1bSJed Brown */
162c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
163c4762a1bSJed Brown {
164c4762a1bSJed Brown   AppCtx            *ctx = (AppCtx*)ctx0;
165c4762a1bSJed Brown   TS                ts,quadts;
166c4762a1bSJed Brown   Vec               U;             /* solution will be stored here */
167c4762a1bSJed Brown   Mat               A;             /* Jacobian matrix */
168c4762a1bSJed Brown   PetscInt          n = 2;
169c4762a1bSJed Brown   PetscReal         ftime;
170c4762a1bSJed Brown   PetscInt          steps;
171c4762a1bSJed Brown   PetscScalar       *u;
172c4762a1bSJed Brown   const PetscScalar *x_ptr,*qx_ptr;
173c4762a1bSJed Brown   Vec               q;
174c4762a1bSJed Brown   PetscInt          direction[2];
175c4762a1bSJed Brown   PetscBool         terminate[2];
176c4762a1bSJed Brown 
1779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P,&x_ptr));
178c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
1799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P,&x_ptr));
180c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown     Create necessary matrix and vectors
182c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1839566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
1849566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
1859566063dSJacob Faibussowitsch   PetscCall(MatSetType(A,MATDENSE));
1869566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1879566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
188c4762a1bSJed Brown 
1899566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&U,NULL));
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192c4762a1bSJed Brown      Create timestepping solver context
193c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1949566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1959566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
1969566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSCN));
1979566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx));
1989566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown      Set initial conditions
202c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
204c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
205c4762a1bSJed Brown   u[1] = 1.0;
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
2079566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,U));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Set solver options
211c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2129566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,1.0));
2139566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
2149566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,0.03125));
2159566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts));
2169566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(quadts,&q));
2179566063dSJacob Faibussowitsch   PetscCall(VecSet(q,0.0));
2189566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx));
2199566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   direction[0] = direction[1] = 1;
222c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
223c4762a1bSJed Brown 
2249566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)ctx));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227c4762a1bSJed Brown      Solve nonlinear system
228c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2299566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
2329566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(q,&qx_ptr));
234c4762a1bSJed Brown   *f   = -ctx->Pm + qx_ptr[0];
2359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(q,&qx_ptr));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
239c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2429566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
243c4762a1bSJed Brown   return 0;
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown /*TEST
247c4762a1bSJed Brown 
248c4762a1bSJed Brown    build:
249c4762a1bSJed Brown       requires: !complex !single
250c4762a1bSJed Brown 
251c4762a1bSJed Brown    test:
252c4762a1bSJed Brown       args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3
253c4762a1bSJed Brown 
254c4762a1bSJed Brown TEST*/
255