xref: /petsc/src/ts/tutorials/power_grid/ex3opt_fd.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
335f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetSolution(tao,&X));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetGradient(tao,&G,NULL,NULL));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(G,&g));
38c4762a1bSJed Brown   fp = fopen("ex3opt_fd_conv.out","a");
395f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFPrintf(PETSC_COMM_WORLD,fp,"%d %g %.12lf %.12lf\n",iterate,gnorm,x[0],g[0]));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode     ierr;
51c4762a1bSJed Brown   PetscMPIInt        size;
52c4762a1bSJed Brown   AppCtx             ctx;
53c4762a1bSJed Brown   Vec                lowerb,upperb;
54c4762a1bSJed Brown   Tao                tao;
55c4762a1bSJed Brown   KSP                ksp;
56c4762a1bSJed Brown   PC                 pc;
57c4762a1bSJed Brown   PetscBool          printtofile;
58c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59c4762a1bSJed Brown      Initialize program
60c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
62c4762a1bSJed Brown   PetscFunctionBeginUser;
635f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
643c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4762a1bSJed Brown     Set runtime options
68c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
70c4762a1bSJed Brown   {
71c4762a1bSJed Brown     ctx.beta    = 2;
72c4762a1bSJed Brown     ctx.c       = 10000.0;
73c4762a1bSJed Brown     ctx.u_s     = 1.0;
74c4762a1bSJed Brown     ctx.omega_s = 1.0;
75c4762a1bSJed Brown     ctx.omega_b = 120.0*PETSC_PI;
76c4762a1bSJed Brown     ctx.H       = 5.0;
775f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL));
78c4762a1bSJed Brown     ctx.D       = 5.0;
795f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL));
80c4762a1bSJed Brown     ctx.E       = 1.1378;
81c4762a1bSJed Brown     ctx.V       = 1.0;
82c4762a1bSJed Brown     ctx.X       = 0.545;
83c4762a1bSJed Brown     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
84c4762a1bSJed Brown     ctx.Pmax_ini = ctx.Pmax;
855f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL));
86c4762a1bSJed Brown     ctx.Pm      = 1.06;
875f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL));
88c4762a1bSJed Brown     ctx.tf      = 0.1;
89c4762a1bSJed Brown     ctx.tcl     = 0.2;
905f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL));
915f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL));
92c4762a1bSJed Brown     printtofile = PETSC_FALSE;
935f80ce2aSJacob Faibussowitsch     CHKERRQ(PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL));
94c4762a1bSJed Brown   }
95c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
985f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
995f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOBLMVM));
100c4762a1bSJed Brown   if (printtofile) {
1015f80ce2aSJacob Faibussowitsch     CHKERRQ(TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL));
102c4762a1bSJed Brown   }
1035f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetMaximumIterations(tao,30));
104c4762a1bSJed Brown   /*
105c4762a1bSJed Brown      Optimization starts
106c4762a1bSJed Brown   */
107c4762a1bSJed Brown   /* Set initial solution guess */
1085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD,1,&p));
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(p,&x_ptr));
110c4762a1bSJed Brown   x_ptr[0]   = ctx.Pm;
1115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(p,&x_ptr));
112c4762a1bSJed Brown 
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,p));
114c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjective(tao,FormFunction,(void *)&ctx));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetGradient(tao,NULL,TaoDefaultComputeGradient,(void *)&ctx));
117c4762a1bSJed Brown 
118c4762a1bSJed Brown   /* Set bounds for the optimization */
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(p,&lowerb));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(p,&upperb));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(lowerb,&x_ptr));
122c4762a1bSJed Brown   x_ptr[0] = 0.;
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(lowerb,&x_ptr));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(upperb,&x_ptr));
125c4762a1bSJed Brown   x_ptr[0] = 1.1;
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(upperb,&x_ptr));
1275f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetVariableBounds(tao,lowerb,upperb));
128c4762a1bSJed Brown 
129c4762a1bSJed Brown   /* Check for any TAO command line options */
1305f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
132c4762a1bSJed Brown   if (ksp) {
1335f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
1345f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCNONE));
135c4762a1bSJed Brown   }
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
139c4762a1bSJed Brown 
1405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(p,PETSC_VIEWER_STDOUT_WORLD));
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* Free TAO data structures */
1435f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&p));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&lowerb));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&upperb));
147*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
148*b122ec5aSJacob Faibussowitsch   return 0;
149c4762a1bSJed Brown }
150c4762a1bSJed Brown 
151c4762a1bSJed Brown /* ------------------------------------------------------------------ */
152c4762a1bSJed Brown /*
153c4762a1bSJed Brown    FormFunction - Evaluates the function and corresponding gradient.
154c4762a1bSJed Brown 
155c4762a1bSJed Brown    Input Parameters:
156c4762a1bSJed Brown    tao - the Tao context
157c4762a1bSJed Brown    X   - the input vector
158a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
159c4762a1bSJed Brown 
160c4762a1bSJed Brown    Output Parameters:
161c4762a1bSJed Brown    f   - the newly evaluated function
162c4762a1bSJed Brown */
163c4762a1bSJed Brown PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
164c4762a1bSJed Brown {
165c4762a1bSJed Brown   AppCtx            *ctx = (AppCtx*)ctx0;
166c4762a1bSJed Brown   TS                ts,quadts;
167c4762a1bSJed Brown   Vec               U;             /* solution will be stored here */
168c4762a1bSJed Brown   Mat               A;             /* Jacobian matrix */
169c4762a1bSJed Brown   PetscInt          n = 2;
170c4762a1bSJed Brown   PetscReal         ftime;
171c4762a1bSJed Brown   PetscInt          steps;
172c4762a1bSJed Brown   PetscScalar       *u;
173c4762a1bSJed Brown   const PetscScalar *x_ptr,*qx_ptr;
174c4762a1bSJed Brown   Vec               q;
175c4762a1bSJed Brown   PetscInt          direction[2];
176c4762a1bSJed Brown   PetscBool         terminate[2];
177c4762a1bSJed Brown 
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(P,&x_ptr));
179c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(P,&x_ptr));
181c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182c4762a1bSJed Brown     Create necessary matrix and vectors
183c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1845f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetType(A,MATDENSE));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
1885f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
189c4762a1bSJed Brown 
1905f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&U,NULL));
191c4762a1bSJed Brown 
192c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193c4762a1bSJed Brown      Create timestepping solver context
194c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1955f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSCN));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx));
1995f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx));
200c4762a1bSJed Brown 
201c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
202c4762a1bSJed Brown      Set initial conditions
203c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2045f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(U,&u));
205c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
206c4762a1bSJed Brown   u[1] = 1.0;
2075f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(U,&u));
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,U));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
211c4762a1bSJed Brown      Set solver options
212c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2135f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,1.0));
2145f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
2155f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,0.03125));
2165f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreateQuadratureTS(ts,PETSC_TRUE,&quadts));
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolution(quadts,&q));
2185f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(q,0.0));
2195f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(quadts,NULL,(TSRHSFunction)CostIntegrand,ctx));
2205f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
221c4762a1bSJed Brown 
222c4762a1bSJed Brown   direction[0] = direction[1] = 1;
223c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
224c4762a1bSJed Brown 
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)ctx));
226c4762a1bSJed Brown 
227c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
228c4762a1bSJed Brown      Solve nonlinear system
229c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2305f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,U));
231c4762a1bSJed Brown 
2325f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
2335f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
2345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(q,&qx_ptr));
235c4762a1bSJed Brown   *f   = -ctx->Pm + qx_ptr[0];
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(q,&qx_ptr));
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
240c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2415f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2425f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&U));
2435f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
244c4762a1bSJed Brown   return 0;
245c4762a1bSJed Brown }
246c4762a1bSJed Brown 
247c4762a1bSJed Brown /*TEST
248c4762a1bSJed Brown 
249c4762a1bSJed Brown    build:
250c4762a1bSJed Brown       requires: !complex !single
251c4762a1bSJed Brown 
252c4762a1bSJed Brown    test:
253c4762a1bSJed Brown       args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3
254c4762a1bSJed Brown 
255c4762a1bSJed Brown TEST*/
256