Actual source code: ex16opt_p.c
petsc-3.8.4 2018-03-24
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: */
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: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
81: {
82: PetscErrorCode ierr;
83: PetscInt row[] = {0,1},col[]={0};
84: PetscScalar J[2][1];
85: const PetscScalar *x;
88: VecGetArrayRead(X,&x);
89: J[0][0] = 0;
90: J[1][0] = (1.-x[0]*x[0])*x[1];
91: MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
92: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
93: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
94: VecRestoreArrayRead(X,&x);
95: return(0);
96: }
98: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
99: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
100: {
101: PetscErrorCode ierr;
102: const PetscScalar *x;
103: PetscReal tfinal, dt, tprev;
104: User user = (User)ctx;
107: TSGetTimeStep(ts,&dt);
108: TSGetMaxTime(ts,&tfinal);
109: TSGetPrevTime(ts,&tprev);
110: VecGetArrayRead(X,&x);
111: 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]));
112: PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
113: VecRestoreArrayRead(X,&x);
114: return(0);
115: }
117: int main(int argc,char **argv)
118: {
119: TS ts; /* nonlinear solver */
120: Vec p;
121: PetscBool monitor = PETSC_FALSE;
122: PetscScalar *x_ptr;
123: PetscMPIInt size;
124: struct _n_User user;
125: PetscErrorCode ierr;
126: Tao tao;
127: Vec lowerb,upperb;
128: KSP ksp;
129: PC pc;
131: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132: Initialize program
133: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134: PetscInitialize(&argc,&argv,NULL,help);
135: MPI_Comm_size(PETSC_COMM_WORLD,&size);
136: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
138: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: Set runtime options
140: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141: user.mu = 1.0;
142: user.next_output = 0.0;
143: user.steps = 0;
144: user.ftime = 0.5;
146: PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
147: PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Create necessary matrix and vectors, solve same ODE on every process
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152: MatCreate(PETSC_COMM_WORLD,&user.A);
153: MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
154: MatSetFromOptions(user.A);
155: MatSetUp(user.A);
156: MatCreateVecs(user.A,&user.x,NULL);
158: MatCreate(PETSC_COMM_WORLD,&user.Jacp);
159: MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
160: MatSetFromOptions(user.Jacp);
161: MatSetUp(user.Jacp);
163: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164: Create timestepping solver context
165: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166: TSCreate(PETSC_COMM_WORLD,&ts);
167: TSSetType(ts,TSRK);
168: TSSetRHSFunction(ts,NULL,RHSFunction,&user);
169: TSSetMaxTime(ts,user.ftime);
170: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
171: if (monitor) {
172: TSMonitorSet(ts,Monitor,&user,NULL);
173: }
175: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176: Set initial conditions
177: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178: VecGetArray(user.x,&x_ptr);
179: x_ptr[0] = 2.0; x_ptr[1] = 0.66666654321;
180: VecRestoreArray(user.x,&x_ptr);
181: TSSetTime(ts,0.0);
182: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime));
184: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185: Save trajectory of solution so that TSAdjointSolve() may be used
186: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187: TSSetSaveTrajectory(ts);
189: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190: Set runtime options
191: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192: TSSetFromOptions(ts);
194: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195: Solve nonlinear system
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197: TSSolve(ts,user.x);
198: TSGetSolveTime(ts,&(user.ftime));
199: TSGetStepNumber(ts,&user.steps);
200: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);
202: VecGetArray(user.x,&x_ptr);
203: user.x_ob[0] = x_ptr[0];
204: user.x_ob[1] = x_ptr[1];
206: MatCreateVecs(user.A,&user.lambda[0],NULL);
207: MatCreateVecs(user.A,&user.lambda[1],NULL);
208: MatCreateVecs(user.Jacp,&user.mup[0],NULL);
209: MatCreateVecs(user.Jacp,&user.mup[1],NULL);
211: /* Create TAO solver and set desired solution method */
212: TaoCreate(PETSC_COMM_WORLD,&tao);
213: TaoSetType(tao,TAOCG);
215: /* Set initial solution guess */
216: MatCreateVecs(user.Jacp,&p,NULL);
217: VecGetArray(p,&x_ptr);
218: x_ptr[0] = 6.0;
219: VecRestoreArray(p,&x_ptr);
221: TaoSetInitialVector(tao,p);
223: /* Set routine for function and gradient evaluation */
224: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);
226: VecDuplicate(p,&lowerb);
227: VecDuplicate(p,&upperb);
228: VecGetArray(lowerb,&x_ptr);
229: x_ptr[0] = 0.0;
230: VecRestoreArray(lowerb,&x_ptr);
231: VecGetArray(upperb,&x_ptr);
232: x_ptr[0] = 100.0;
233: VecRestoreArray(upperb,&x_ptr);
235: TaoSetVariableBounds(tao,lowerb,upperb);
237: /* Check for any TAO command line options */
238: TaoSetFromOptions(tao);
239: TaoGetKSP(tao,&ksp);
240: if (ksp) {
241: KSPGetPC(ksp,&pc);
242: PCSetType(pc,PCNONE);
243: }
245: TaoSetTolerances(tao,1e-13,PETSC_DEFAULT,PETSC_DEFAULT);
247: /* SOLVE THE APPLICATION */
248: TaoSolve(tao);
250: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
251: /* Free TAO data structures */
252: TaoDestroy(&tao);
254: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
255: Free work space. All PETSc objects should be destroyed when they
256: are no longer needed.
257: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258: MatDestroy(&user.A);
259: MatDestroy(&user.Jacp);
260: VecDestroy(&user.x);
261: VecDestroy(&user.lambda[0]);
262: VecDestroy(&user.lambda[1]);
263: VecDestroy(&user.mup[0]);
264: VecDestroy(&user.mup[1]);
265: TSDestroy(&ts);
267: VecDestroy(&lowerb);
268: VecDestroy(&upperb);
269: VecDestroy(&p);
270: PetscFinalize();
271: return ierr;
272: }
274: /* ------------------------------------------------------------------ */
275: /*
276: FormFunctionGradient - Evaluates the function and corresponding gradient.
278: Input Parameters:
279: tao - the Tao context
280: X - the input vector
281: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
283: Output Parameters:
284: f - the newly evaluated function
285: G - the newly evaluated gradient
286: */
287: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
288: {
289: User user = (User)ctx;
290: TS ts;
291: PetscScalar *x_ptr,*y_ptr;
292: PetscErrorCode ierr;
294: VecGetArray(P,&x_ptr);
295: user->mu = x_ptr[0];
297: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
298: Create timestepping solver context
299: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
300: TSCreate(PETSC_COMM_WORLD,&ts);
301: TSSetType(ts,TSRK);
302: TSSetRHSFunction(ts,NULL,RHSFunction,user);
304: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
305: Set initial conditions
306: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
307: VecGetArray(user->x,&x_ptr);
308: x_ptr[0] = 2; x_ptr[1] = 0.66666654321;
309: VecRestoreArray(user->x,&x_ptr);
310: TSSetTime(ts,0.0);
311: TSSetMaxTime(ts,0.5);
312: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
314: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
315: Save trajectory of solution so that TSAdjointSolve() may be used
316: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
317: TSSetSaveTrajectory(ts);
319: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
320: Set runtime options
321: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
322: TSSetFromOptions(ts);
324: TSSolve(ts,user->x);
325: TSGetSolveTime(ts,&user->ftime);
326: TSGetStepNumber(ts,&user->steps);
327: PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime);
329: VecGetArray(user->x,&x_ptr);
330: *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]);
331: 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));
333: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
334: Adjoint model starts here
335: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
336: /* Redet initial conditions for the adjoint integration */
337: VecGetArray(user->lambda[0],&y_ptr);
338: y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
339: y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
340: VecRestoreArray(user->lambda[0],&y_ptr);
341: VecGetArray(user->lambda[1],&x_ptr);
342: x_ptr[0] = 0.0; x_ptr[1] = 1.0;
343: VecRestoreArray(user->lambda[1],&x_ptr);
345: VecGetArray(user->mup[0],&x_ptr);
346: x_ptr[0] = 0.0;
347: VecRestoreArray(user->mup[0],&x_ptr);
348: VecGetArray(user->mup[1],&x_ptr);
349: x_ptr[0] = 0.0;
350: VecRestoreArray(user->mup[1],&x_ptr);
351: TSSetCostGradients(ts,1,user->lambda,user->mup);
353: TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user);
354: TSAdjointSetRHSJacobian(ts,user->Jacp,RHSJacobianP,user);
356: TSAdjointSolve(ts);
358: VecCopy(user->mup[0],G);
360: TSDestroy(&ts);
361: return(0);
362: }