Actual source code: ex16opt_ic.c
petsc-3.8.4 2018-03-24
1: static char help[] = "Solves an ODE-constrained optimization problem -- finding the optimal initial conditions 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 optimal values for initial conditions.
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;
26: PetscInt steps;
27: PetscReal ftime,x_ob[2];
28: Mat A; /* Jacobian matrix */
29: Vec x,lambda[2]; /* adjoint variables */
30: };
32: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
34: /*
35: * User-defined routines
36: */
37: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
38: {
39: PetscErrorCode ierr;
40: User user = (User)ctx;
41: PetscScalar *f;
42: const PetscScalar *x;
45: VecGetArrayRead(X,&x);
46: VecGetArray(F,&f);
47: f[0] = x[1];
48: f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
49: VecRestoreArrayRead(X,&x);
50: VecRestoreArray(F,&f);
51: return(0);
52: }
54: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
55: {
56: PetscErrorCode ierr;
57: User user = (User)ctx;
58: PetscReal mu = user->mu;
59: PetscInt rowcol[] = {0,1};
60: PetscScalar J[2][2];
61: const PetscScalar *x;
64: VecGetArrayRead(X,&x);
65: J[0][0] = 0;
66: J[1][0] = -2.*mu*x[1]*x[0]-1;
67: J[0][1] = 1.0;
68: J[1][1] = mu*(1.0-x[0]*x[0]);
69: MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
70: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
71: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
72: if (A != B) {
73: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
74: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
75: }
76: VecRestoreArrayRead(X,&x);
77: return(0);
78: }
80: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
81: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
82: {
83: PetscErrorCode ierr;
84: const PetscScalar *x;
85: PetscReal tfinal, dt, tprev;
86: User user = (User)ctx;
89: TSGetTimeStep(ts,&dt);
90: TSGetMaxTime(ts,&tfinal);
91: TSGetPrevTime(ts,&tprev);
92: VecGetArrayRead(X,&x);
93: 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]));
94: PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
95: VecGetArrayRead(X,&x);
96: return(0);
97: }
99: int main(int argc,char **argv)
100: {
101: TS ts; /* nonlinear solver */
102: Vec ic;
103: PetscBool monitor = PETSC_FALSE;
104: PetscScalar *x_ptr;
105: PetscMPIInt size;
106: struct _n_User user;
107: PetscErrorCode ierr;
108: Tao tao;
109: KSP ksp;
110: PC pc;
112: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
113: Initialize program
114: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
115: PetscInitialize(&argc,&argv,NULL,help);
116: MPI_Comm_size(PETSC_COMM_WORLD,&size);
117: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
119: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
120: Set runtime options
121: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
122: user.mu = 1.0;
123: user.next_output = 0.0;
124: user.steps = 0;
125: user.ftime = 0.5;
127: PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
128: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
130: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131: Create necessary matrix and vectors, solve same ODE on every process
132: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: MatCreate(PETSC_COMM_WORLD,&user.A);
134: MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
135: MatSetFromOptions(user.A);
136: MatSetUp(user.A);
137: MatCreateVecs(user.A,&user.x,NULL);
139: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140: Create timestepping solver context
141: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: TSCreate(PETSC_COMM_WORLD,&ts);
143: TSSetType(ts,TSRK);
144: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
145: TSSetMaxTime(ts,user.ftime);
146: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
147: if (monitor) {
148: TSMonitorSet(ts,Monitor,&user,NULL);
149: }
151: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152: Set initial conditions
153: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154: VecGetArray(user.x,&x_ptr);
155: x_ptr[0] = 2.0; x_ptr[1] = 0.66666654321;
156: VecRestoreArray(user.x,&x_ptr);
157: TSSetTime(ts,0.0);
158: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime));
160: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161: Save trajectory of solution so that TSAdjointSolve() may be used
162: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163: TSSetSaveTrajectory(ts);
165: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166: Set runtime options
167: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
168: TSSetFromOptions(ts);
170: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171: Solve nonlinear system
172: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173: TSSolve(ts,user.x);
174: TSGetSolveTime(ts,&(user.ftime));
175: TSGetStepNumber(ts,&user.steps);
176: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);
178: VecGetArray(user.x,&x_ptr);
179: user.x_ob[0] = x_ptr[0];
180: user.x_ob[1] = x_ptr[1];
182: MatCreateVecs(user.A,&user.lambda[0],NULL);
184: /* Create TAO solver and set desired solution method */
185: TaoCreate(PETSC_COMM_WORLD,&tao);
186: TaoSetType(tao,TAOCG);
188: /* Set initial solution guess */
189: MatCreateVecs(user.A,&ic,NULL);
190: VecGetArray(ic,&x_ptr);
191: x_ptr[0] = 2.1;
192: x_ptr[1] = 0.7;
193: VecRestoreArray(ic,&x_ptr);
195: TaoSetInitialVector(tao,ic);
197: /* Set routine for function and gradient evaluation */
198: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);
200: /* Check for any TAO command line options */
201: TaoSetFromOptions(tao);
202: TaoGetKSP(tao,&ksp);
203: if (ksp) {
204: KSPGetPC(ksp,&pc);
205: PCSetType(pc,PCNONE);
206: }
208: TaoSetTolerances(tao,1e-10,PETSC_DEFAULT,PETSC_DEFAULT);
210: /* SOLVE THE APPLICATION */
211: TaoSolve(tao);
213: /* Free TAO data structures */
214: TaoDestroy(&tao);
216: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217: Free work space. All PETSc objects should be destroyed when they
218: are no longer needed.
219: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220: MatDestroy(&user.A);
221: VecDestroy(&user.x);
222: VecDestroy(&user.lambda[0]);
223: TSDestroy(&ts);
225: VecDestroy(&ic);
226: PetscFinalize();
227: return ierr;
228: }
230: /* ------------------------------------------------------------------ */
231: /*
232: FormFunctionGradient - Evaluates the function and corresponding gradient.
234: Input Parameters:
235: tao - the Tao context
236: X - the input vector
237: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
239: Output Parameters:
240: f - the newly evaluated function
241: G - the newly evaluated gradient
242: */
243: PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
244: {
245: User user = (User)ctx;
246: TS ts;
247: PetscScalar *x_ptr,*y_ptr;
248: PetscErrorCode ierr;
249: PetscScalar *ic_ptr;
251: VecCopy(IC,user->x);
253: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
254: Create timestepping solver context
255: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
256: TSCreate(PETSC_COMM_WORLD,&ts);
257: TSSetType(ts,TSRK);
258: TSSetRHSFunction(ts,NULL,RHSFunction,user);
260: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
261: Set time
262: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
263: TSSetTime(ts,0.0);
264: TSSetTimeStep(ts,.001);
265: TSSetMaxTime(ts,0.5);
266: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
268: TSSetTolerances(ts,1e-7,NULL,1e-7,NULL);
270: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
271: Save trajectory of solution so that TSAdjointSolve() may be used
272: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
273: TSSetSaveTrajectory(ts);
275: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276: Set runtime options
277: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278: TSSetFromOptions(ts);
280: TSSolve(ts,user->x);
281: TSGetSolveTime(ts,&user->ftime);
282: TSGetStepNumber(ts,&user->steps);
283: PetscPrintf(PETSC_COMM_WORLD,"mu %.6f, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime);
285: VecGetArray(IC,&ic_ptr);
286: VecGetArray(user->x,&x_ptr);
287: *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]);
288: PetscPrintf(PETSC_COMM_WORLD,"Observed value y_ob=[%f; %f], ODE solution y=[%f;%f], Cost function f=%f\n",(double)user->x_ob[0],(double)user->x_ob[1],(double)x_ptr[0],(double)x_ptr[1],(double)(*f));
290: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
291: Adjoint model starts here
292: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
293: /* Redet initial conditions for the adjoint integration */
294: VecGetArray(user->lambda[0],&y_ptr);
295: y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
296: y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
297: VecRestoreArray(user->lambda[0],&y_ptr);
298: TSSetCostGradients(ts,1,user->lambda,NULL);
300: /* Set RHS Jacobian for the adjoint integration */
301: TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user);
303: TSAdjointSolve(ts);
305: VecCopy(user->lambda[0],G);
307: TSDestroy(&ts);
308: return(0);
309: }