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