Actual source code: ex16opt_p.c
petsc-3.7.7 2017-09-25
1: static char help[] = "Solves an ODE-constrained optimization problem -- finding the optimal stiffness parameter for the van der Pol equation.\n\
2: Input parameters include:\n\
3: -mu : stiffness parameter\n\n";
5: /*
6: Concepts: TS^time-dependent nonlinear problems
7: Concepts: TS^van der Pol equation
8: Concepts: Optimization using adjoint sensitivities
9: Processors: 1
10: */
11: /* ------------------------------------------------------------------------
13: Notes:
14: This code demonstrates how to solve an ODE-constrained optimization problem with TAO, TSAdjoint and TS.
15: The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu.
16: The gradient is computed with the discrete adjoint of an explicit Runge-Kutta method, see ex16adj.c for details.
17: ------------------------------------------------------------------------- */
18: #include <petsctao.h>
19: #include <petscts.h>
20: #include <petscmat.h>
21: typedef struct _n_User *User;
22: struct _n_User {
23: PetscReal mu;
24: PetscReal next_output;
25: PetscInt steps;
26: PetscReal ftime,x_ob[2];
27: Mat A; /* Jacobian matrix */
28: Mat Jacp; /* JacobianP matrix */
29: Vec x,lambda[2],mup[2]; /* adjoint variables */
30: };
32: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
34: /*
35: * User-defined routines
36: */
39: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
40: {
41: PetscErrorCode ierr;
42: User user = (User)ctx;
43: PetscScalar *f;
44: const PetscScalar *x;
47: VecGetArrayRead(X,&x);
48: VecGetArray(F,&f);
49: f[0] = x[1];
50: f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
51: VecRestoreArrayRead(X,&x);
52: VecRestoreArray(F,&f);
53: return(0);
54: }
58: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
59: {
60: PetscErrorCode ierr;
61: User user = (User)ctx;
62: PetscReal mu = user->mu;
63: PetscInt rowcol[] = {0,1};
64: PetscScalar J[2][2];
65: const PetscScalar *x;
68: VecGetArrayRead(X,&x);
69: J[0][0] = 0;
70: J[1][0] = -2.*mu*x[1]*x[0]-1.;
71: J[0][1] = 1.0;
72: J[1][1] = mu*(1.0-x[0]*x[0]);
73: MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
74: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
75: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
76: if (A != B) {
77: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
78: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
79: }
80: VecRestoreArrayRead(X,&x);
81: return(0);
82: }
86: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
87: {
88: PetscErrorCode ierr;
89: PetscInt row[] = {0,1},col[]={0};
90: PetscScalar J[2][1];
91: const PetscScalar *x;
94: VecGetArrayRead(X,&x);
95: J[0][0] = 0;
96: J[1][0] = (1.-x[0]*x[0])*x[1];
97: MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
98: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
99: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
100: VecRestoreArrayRead(X,&x);
101: return(0);
102: }
106: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
107: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
108: {
109: PetscErrorCode ierr;
110: const PetscScalar *x;
111: PetscReal tfinal, dt, tprev;
112: User user = (User)ctx;
115: TSGetTimeStep(ts,&dt);
116: TSGetDuration(ts,NULL,&tfinal);
117: TSGetPrevTime(ts,&tprev);
118: VecGetArrayRead(X,&x);
119: PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",(double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(x[0]),(double)PetscRealPart(x[1]));
120: PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
121: VecRestoreArrayRead(X,&x);
122: return(0);
123: }
127: int main(int argc,char **argv)
128: {
129: TS ts; /* nonlinear solver */
130: Vec p;
131: PetscBool monitor = PETSC_FALSE;
132: PetscScalar *x_ptr;
133: PetscMPIInt size;
134: struct _n_User user;
135: PetscErrorCode ierr;
136: Tao tao;
137: Vec lowerb,upperb;
138: KSP ksp;
139: PC pc;
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Initialize program
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
144: PetscInitialize(&argc,&argv,NULL,help);
146: MPI_Comm_size(PETSC_COMM_WORLD,&size);
147: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Set runtime options
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152: user.mu = 1.0;
153: user.next_output = 0.0;
154: user.steps = 0;
155: user.ftime = 0.5;
157: PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
158: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Create necessary matrix and vectors, solve same ODE on every process
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: MatCreate(PETSC_COMM_WORLD,&user.A);
164: MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
165: MatSetFromOptions(user.A);
166: MatSetUp(user.A);
167: MatCreateVecs(user.A,&user.x,NULL);
169: MatCreate(PETSC_COMM_WORLD,&user.Jacp);
170: MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
171: MatSetFromOptions(user.Jacp);
172: MatSetUp(user.Jacp);
174: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175: Create timestepping solver context
176: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177: TSCreate(PETSC_COMM_WORLD,&ts);
178: TSSetType(ts,TSRK);
179: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
180: TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
181: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
182: if (monitor) {
183: TSMonitorSet(ts,Monitor,&user,NULL);
184: }
186: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187: Set initial conditions
188: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189: VecGetArray(user.x,&x_ptr);
190: x_ptr[0] = 2.0; x_ptr[1] = 0.66666654321;
191: VecRestoreArray(user.x,&x_ptr);
192: TSSetTime(ts,0.0);
193: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime));
195: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
196: Save trajectory of solution so that TSAdjointSolve() may be used
197: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
198: TSSetSaveTrajectory(ts);
200: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201: Set runtime options
202: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203: TSSetFromOptions(ts);
205: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: Solve nonlinear system
207: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
208: TSSolve(ts,user.x);
209: TSGetSolveTime(ts,&(user.ftime));
210: TSGetTimeStepNumber(ts,&user.steps);
211: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);
213: VecGetArray(user.x,&x_ptr);
214: user.x_ob[0] = x_ptr[0];
215: user.x_ob[1] = x_ptr[1];
217: MatCreateVecs(user.A,&user.lambda[0],NULL);
218: MatCreateVecs(user.A,&user.lambda[1],NULL);
219: MatCreateVecs(user.Jacp,&user.mup[0],NULL);
220: MatCreateVecs(user.Jacp,&user.mup[1],NULL);
222: /* Create TAO solver and set desired solution method */
223: TaoCreate(PETSC_COMM_WORLD,&tao);
224: TaoSetType(tao,TAOCG);
226: /* Set initial solution guess */
227: MatCreateVecs(user.Jacp,&p,NULL);
228: VecGetArray(p,&x_ptr);
229: x_ptr[0] = 6.0;
230: VecRestoreArray(p,&x_ptr);
232: TaoSetInitialVector(tao,p);
234: /* Set routine for function and gradient evaluation */
235: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);
237: VecDuplicate(p,&lowerb);
238: VecDuplicate(p,&upperb);
239: VecGetArray(lowerb,&x_ptr);
240: x_ptr[0] = 0.0;
241: VecRestoreArray(lowerb,&x_ptr);
242: VecGetArray(upperb,&x_ptr);
243: x_ptr[0] = 100.0;
244: VecRestoreArray(upperb,&x_ptr);
246: TaoSetVariableBounds(tao,lowerb,upperb);
248: /* Check for any TAO command line options */
249: TaoSetFromOptions(tao);
250: TaoGetKSP(tao,&ksp);
251: if (ksp) {
252: KSPGetPC(ksp,&pc);
253: PCSetType(pc,PCNONE);
254: }
256: TaoSetTolerances(tao,1e-13,PETSC_DEFAULT,PETSC_DEFAULT);
258: /* SOLVE THE APPLICATION */
259: TaoSolve(tao);
261: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
262: /* Free TAO data structures */
263: TaoDestroy(&tao);
265: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
266: Free work space. All PETSc objects should be destroyed when they
267: are no longer needed.
268: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
269: MatDestroy(&user.A);
270: MatDestroy(&user.Jacp);
271: VecDestroy(&user.x);
272: VecDestroy(&user.lambda[0]);
273: VecDestroy(&user.lambda[1]);
274: VecDestroy(&user.mup[0]);
275: VecDestroy(&user.mup[1]);
276: TSDestroy(&ts);
278: VecDestroy(&lowerb);
279: VecDestroy(&upperb);
280: VecDestroy(&p);
281: PetscFinalize();
282: return(0);
283: }
285: /* ------------------------------------------------------------------ */
288: /*
289: FormFunctionGradient - Evaluates the function and corresponding gradient.
291: Input Parameters:
292: tao - the Tao context
293: X - the input vector
294: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
296: Output Parameters:
297: f - the newly evaluated function
298: G - the newly evaluated gradient
299: */
300: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
301: {
302: User user = (User)ctx;
303: TS ts;
304: PetscScalar *x_ptr,*y_ptr;
305: PetscErrorCode ierr;
307: VecGetArray(P,&x_ptr);
308: user->mu = x_ptr[0];
310: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
311: Create timestepping solver context
312: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
313: TSCreate(PETSC_COMM_WORLD,&ts);
314: TSSetType(ts,TSRK);
315: TSSetRHSFunction(ts,NULL,RHSFunction,user);
317: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318: Set initial conditions
319: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
320: VecGetArray(user->x,&x_ptr);
321: x_ptr[0] = 2; x_ptr[1] = 0.66666654321;
322: VecRestoreArray(user->x,&x_ptr);
323: TSSetTime(ts,0.0);
324: TSSetDuration(ts,PETSC_DEFAULT,0.5);
325: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
327: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
328: Save trajectory of solution so that TSAdjointSolve() may be used
329: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
330: TSSetSaveTrajectory(ts);
332: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
333: Set runtime options
334: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
335: TSSetFromOptions(ts);
337: TSSolve(ts,user->x);
338: TSGetSolveTime(ts,&user->ftime);
339: TSGetTimeStepNumber(ts,&user->steps);
340: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime);
342: VecGetArray(user->x,&x_ptr);
343: *f = (x_ptr[0]-user->x_ob[0])*(x_ptr[0]-user->x_ob[0])+(x_ptr[1]-user->x_ob[1])*(x_ptr[1]-user->x_ob[1]);
344: PetscPrintf(PETSC_COMM_WORLD,"Observed value y_ob=[%f; %f], ODE solution y=%f, Cost function f=%f\n",(double)user->x_ob[0],(double)user->x_ob[1],(double)x_ptr[0],(double)(*f));
346: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
347: Adjoint model starts here
348: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
349: /* Redet initial conditions for the adjoint integration */
350: VecGetArray(user->lambda[0],&y_ptr);
351: y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
352: y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
353: VecRestoreArray(user->lambda[0],&y_ptr);
354: VecGetArray(user->lambda[1],&x_ptr);
355: x_ptr[0] = 0.0; x_ptr[1] = 1.0;
356: VecRestoreArray(user->lambda[1],&x_ptr);
358: VecGetArray(user->mup[0],&x_ptr);
359: x_ptr[0] = 0.0;
360: VecRestoreArray(user->mup[0],&x_ptr);
361: VecGetArray(user->mup[1],&x_ptr);
362: x_ptr[0] = 0.0;
363: VecRestoreArray(user->mup[1],&x_ptr);
364: TSSetCostGradients(ts,1,user->lambda,user->mup);
366: TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user);
367: TSAdjointSetRHSJacobian(ts,user->Jacp,RHSJacobianP,user);
369: TSAdjointSolve(ts);
371: VecCopy(user->mup[0],G);
373: TSDestroy(&ts);
374: return(0);
375: }