xref: /petsc/src/ts/tutorials/autodiff/ex16opt_ic.cxx (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
1c4762a1bSJed Brown static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for an ODE-constrained optimization problem.\n\
2c4762a1bSJed Brown Input parameters include:\n\
3c4762a1bSJed Brown       -mu : stiffness parameter\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
7c4762a1bSJed Brown    Concepts: TS^van der Pol equation
8c4762a1bSJed Brown    Concepts: Optimization using adjoint sensitivities
9c4762a1bSJed Brown    Concepts: Automatic differentation using ADOL-C
10c4762a1bSJed Brown    Processors: 1
11c4762a1bSJed Brown */
12c4762a1bSJed Brown /*
13c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
14c4762a1bSJed Brown 
15c4762a1bSJed Brown    For documentation on ADOL-C, see
16c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
17c4762a1bSJed Brown */
18c4762a1bSJed Brown /* ------------------------------------------------------------------------
19c4762a1bSJed Brown   See ex16opt_ic for a description of the problem being solved.
20c4762a1bSJed Brown   ------------------------------------------------------------------------- */
21c4762a1bSJed Brown #include <petsctao.h>
22c4762a1bSJed Brown #include <petscts.h>
23c4762a1bSJed Brown #include <petscmat.h>
24c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
25c4762a1bSJed Brown #include <adolc/adolc.h>
26c4762a1bSJed Brown 
27c4762a1bSJed Brown typedef struct _n_User *User;
28c4762a1bSJed Brown struct _n_User {
29c4762a1bSJed Brown   PetscReal mu;
30c4762a1bSJed Brown   PetscReal next_output;
31c4762a1bSJed Brown   PetscInt  steps;
32c4762a1bSJed Brown 
33c4762a1bSJed Brown   /* Sensitivity analysis support */
34c4762a1bSJed Brown   PetscReal ftime,x_ob[2];
35c4762a1bSJed Brown   Mat       A;             /* Jacobian matrix */
36c4762a1bSJed Brown   Vec       x,lambda[2];   /* adjoint variables */
37c4762a1bSJed Brown 
38c4762a1bSJed Brown   /* Automatic differentiation support */
39c4762a1bSJed Brown   AdolcCtx  *adctx;
40c4762a1bSJed Brown };
41c4762a1bSJed Brown 
42c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
43c4762a1bSJed Brown 
44c4762a1bSJed Brown /*
45c4762a1bSJed Brown   'Passive' RHS function, used in residual evaluations during the time integration.
46c4762a1bSJed Brown */
47c4762a1bSJed Brown static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
48c4762a1bSJed Brown {
49c4762a1bSJed Brown   User              user = (User)ctx;
50c4762a1bSJed Brown   PetscScalar       *f;
51c4762a1bSJed Brown   const PetscScalar *x;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   PetscFunctionBeginUser;
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
56c4762a1bSJed Brown   f[0] = x[1];
57c4762a1bSJed Brown   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
60c4762a1bSJed Brown   PetscFunctionReturn(0);
61c4762a1bSJed Brown }
62c4762a1bSJed Brown 
63c4762a1bSJed Brown /*
64c4762a1bSJed Brown   Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the
65c4762a1bSJed Brown   Jacobian transform.
66c4762a1bSJed Brown */
67c4762a1bSJed Brown static PetscErrorCode RHSFunctionActive(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
68c4762a1bSJed Brown {
69c4762a1bSJed Brown   User              user = (User)ctx;
70c4762a1bSJed Brown   PetscReal         mu   = user->mu;
71c4762a1bSJed Brown   PetscScalar       *f;
72c4762a1bSJed Brown   const PetscScalar *x;
73c4762a1bSJed Brown 
74c4762a1bSJed Brown   adouble           f_a[2];                     /* adouble for dependent variables */
75c4762a1bSJed Brown   adouble           x_a[2];                     /* adouble for independent variables */
76c4762a1bSJed Brown 
77c4762a1bSJed Brown   PetscFunctionBeginUser;
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   trace_on(1);                                  /* Start of active section */
82c4762a1bSJed Brown   x_a[0] <<= x[0]; x_a[1] <<= x[1];             /* Mark as independent */
83c4762a1bSJed Brown   f_a[0] = x_a[1];
84c4762a1bSJed Brown   f_a[1] = mu*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0];
85c4762a1bSJed Brown   f_a[0] >>= f[0]; f_a[1] >>= f[1];             /* Mark as dependent */
86c4762a1bSJed Brown   trace_off(1);                                 /* End of active section */
87c4762a1bSJed Brown 
889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
90c4762a1bSJed Brown   PetscFunctionReturn(0);
91c4762a1bSJed Brown }
92c4762a1bSJed Brown 
93c4762a1bSJed Brown /*
94c4762a1bSJed Brown   Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver.
95c4762a1bSJed Brown */
96c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
97c4762a1bSJed Brown {
98c4762a1bSJed Brown   User              user=(User)ctx;
99a8c08197SHong Zhang   const PetscScalar *x;
100c4762a1bSJed Brown 
101c4762a1bSJed Brown   PetscFunctionBeginUser;
1029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
1039566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeRHSJacobian(1,A,x,user->adctx));
1049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
105c4762a1bSJed Brown   PetscFunctionReturn(0);
106c4762a1bSJed Brown }
107c4762a1bSJed Brown 
108c4762a1bSJed Brown /*
109c4762a1bSJed Brown   Monitor timesteps and use interpolation to output at integer multiples of 0.1
110c4762a1bSJed Brown */
111c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
112c4762a1bSJed Brown {
113c4762a1bSJed Brown   const PetscScalar *x;
114c4762a1bSJed Brown   PetscReal         tfinal, dt, tprev;
115c4762a1bSJed Brown   User              user = (User)ctx;
116c4762a1bSJed Brown 
117c4762a1bSJed Brown   PetscFunctionBeginUser;
1189566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts,&dt));
1199566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts,&tfinal));
1209566063dSJacob Faibussowitsch   PetscCall(TSGetPrevTime(ts,&tprev));
1219566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
1229566063dSJacob Faibussowitsch   PetscCall(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])));
1239566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev));
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
125c4762a1bSJed Brown   PetscFunctionReturn(0);
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown int main(int argc,char **argv)
129c4762a1bSJed Brown {
130410585f6SHong Zhang   TS                 ts = NULL;          /* nonlinear solver */
131c4762a1bSJed Brown   Vec                ic,r;
132c4762a1bSJed Brown   PetscBool          monitor = PETSC_FALSE;
133c4762a1bSJed Brown   PetscScalar        *x_ptr;
134c4762a1bSJed Brown   PetscMPIInt        size;
135c4762a1bSJed Brown   struct _n_User     user;
136c4762a1bSJed Brown   AdolcCtx           *adctx;
137c4762a1bSJed Brown   Tao                tao;
138c4762a1bSJed Brown   KSP                ksp;
139c4762a1bSJed Brown   PC                 pc;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown      Initialize program
143c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1449566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
1459566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
146*08401ef6SPierre Jolivet   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149c4762a1bSJed Brown     Set runtime options and create AdolcCtx
150c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1519566063dSJacob Faibussowitsch   PetscCall(PetscNew(&adctx));
152c4762a1bSJed Brown   user.mu          = 1.0;
153c4762a1bSJed Brown   user.next_output = 0.0;
154c4762a1bSJed Brown   user.steps       = 0;
155c4762a1bSJed Brown   user.ftime       = 0.5;
156c4762a1bSJed Brown   adctx->m = 2;adctx->n = 2;adctx->p = 2;
157c4762a1bSJed Brown   user.adctx = adctx;
158c4762a1bSJed Brown 
1599566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
1609566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
164c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1659566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.A));
1669566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
1679566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.A));
1689566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.A));
1699566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.x,NULL));
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172c4762a1bSJed Brown      Set initial conditions
173c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.x,&x_ptr));
175c4762a1bSJed Brown   x_ptr[0] = 2.0;   x_ptr[1] = 0.66666654321;
1769566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.x,&x_ptr));
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179c4762a1bSJed Brown      Trace just once on each tape and put zeros on Jacobian diagonal
180c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1819566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(user.x,&r));
1829566063dSJacob Faibussowitsch   PetscCall(RHSFunctionActive(ts,0.,user.x,r,&user));
1839566063dSJacob Faibussowitsch   PetscCall(VecSet(r,0));
1849566063dSJacob Faibussowitsch   PetscCall(MatDiagonalSet(user.A,r,INSERT_VALUES));
1859566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown      Create timestepping solver context
189c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1909566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1919566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSRK));
1929566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user));
1939566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user));
1949566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,user.ftime));
1959566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
196c4762a1bSJed Brown   if (monitor) {
1979566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
198c4762a1bSJed Brown   }
199c4762a1bSJed Brown 
2009566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0.0));
2019566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime)));
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
205c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2069566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209c4762a1bSJed Brown      Set runtime options
210c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2119566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
212c4762a1bSJed Brown 
213c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214c4762a1bSJed Brown      Solve nonlinear system
215c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2169566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,user.x));
2179566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&(user.ftime)));
2189566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&user.steps));
2199566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime));
220c4762a1bSJed Brown 
2219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.x,&x_ptr));
222c4762a1bSJed Brown   user.x_ob[0] = x_ptr[0];
223c4762a1bSJed Brown   user.x_ob[1] = x_ptr[1];
2249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.x,&x_ptr));
225c4762a1bSJed Brown 
2269566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.lambda[0],NULL));
227c4762a1bSJed Brown 
228c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
2299566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD,&tao));
2309566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao,TAOCG));
231c4762a1bSJed Brown 
232c4762a1bSJed Brown   /* Set initial solution guess */
2339566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&ic,NULL));
2349566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ic,&x_ptr));
235c4762a1bSJed Brown   x_ptr[0]  = 2.1;
236c4762a1bSJed Brown   x_ptr[1]  = 0.7;
2379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ic,&x_ptr));
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao,ic));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
2429566063dSJacob Faibussowitsch   PetscCall(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /* Check for any TAO command line options */
2459566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
2469566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao,&ksp));
247c4762a1bSJed Brown   if (ksp) {
2489566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp,&pc));
2499566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc,PCNONE));
250c4762a1bSJed Brown   }
251c4762a1bSJed Brown 
2529566063dSJacob Faibussowitsch   PetscCall(TaoSetTolerances(tao,1e-10,PETSC_DEFAULT,PETSC_DEFAULT));
253c4762a1bSJed Brown 
254c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
2559566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
256c4762a1bSJed Brown 
257c4762a1bSJed Brown   /* Free TAO data structures */
2589566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
261c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
262c4762a1bSJed Brown      are no longer needed.
263c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2649566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
2659566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.x));
2669566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.lambda[0]));
2679566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2689566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ic));
2699566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
2709566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
271b122ec5aSJacob Faibussowitsch   return 0;
272c4762a1bSJed Brown }
273c4762a1bSJed Brown 
274c4762a1bSJed Brown /* ------------------------------------------------------------------ */
275c4762a1bSJed Brown /*
276c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
277c4762a1bSJed Brown 
278c4762a1bSJed Brown    Input Parameters:
279c4762a1bSJed Brown    tao - the Tao context
280c4762a1bSJed Brown    X   - the input vector
281a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
282c4762a1bSJed Brown 
283c4762a1bSJed Brown    Output Parameters:
284c4762a1bSJed Brown    f   - the newly evaluated function
285c4762a1bSJed Brown    G   - the newly evaluated gradient
286c4762a1bSJed Brown */
287c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
288c4762a1bSJed Brown {
289c4762a1bSJed Brown   User              user = (User)ctx;
290c4762a1bSJed Brown   TS                ts;
291c4762a1bSJed Brown   PetscScalar       *x_ptr,*y_ptr;
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   PetscFunctionBeginUser;
2949566063dSJacob Faibussowitsch   PetscCall(VecCopy(IC,user->x));
295c4762a1bSJed Brown 
296c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
297c4762a1bSJed Brown      Create timestepping solver context
298c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2999566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
3009566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSRK));
3019566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,user));
302c4762a1bSJed Brown   /*   Set RHS Jacobian  for the adjoint integration */
3039566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user));
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
306c4762a1bSJed Brown      Set time
307c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3089566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts,0.0));
3099566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.001));
3109566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,0.5));
3119566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
312c4762a1bSJed Brown 
3139566063dSJacob Faibussowitsch   PetscCall(TSSetTolerances(ts,1e-7,NULL,1e-7,NULL));
314c4762a1bSJed Brown 
315c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
316c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
317c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3189566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
319c4762a1bSJed Brown 
320c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
321c4762a1bSJed Brown      Set runtime options
322c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3239566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
324c4762a1bSJed Brown 
3259566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,user->x));
3269566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&user->ftime));
3279566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&user->steps));
3289566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %.6f, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime));
329c4762a1bSJed Brown 
3309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->x,&x_ptr));
331c4762a1bSJed Brown   *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]);
3329566063dSJacob Faibussowitsch   PetscCall(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)));
333c4762a1bSJed Brown 
334c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
335c4762a1bSJed Brown      Adjoint model starts here
336c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
337c4762a1bSJed Brown   /*   Redet initial conditions for the adjoint integration */
3389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user->lambda[0],&y_ptr));
339c4762a1bSJed Brown   y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
340c4762a1bSJed Brown   y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->lambda[0],&y_ptr));
3429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user->x,&x_ptr));
3439566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,user->lambda,NULL));
344c4762a1bSJed Brown 
3459566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
346c4762a1bSJed Brown 
3479566063dSJacob Faibussowitsch   PetscCall(VecCopy(user->lambda[0],G));
348c4762a1bSJed Brown 
3499566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
350c4762a1bSJed Brown   PetscFunctionReturn(0);
351c4762a1bSJed Brown }
352c4762a1bSJed Brown 
353c4762a1bSJed Brown /*TEST
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   build:
356c4762a1bSJed Brown     requires: double !complex adolc
357c4762a1bSJed Brown 
358c4762a1bSJed Brown   test:
359c4762a1bSJed Brown     suffix: 1
360c4762a1bSJed Brown     args: -ts_rhs_jacobian_test_mult_transpose FALSE -tao_max_it 2 -ts_rhs_jacobian_test_mult FALSE
361c4762a1bSJed Brown     output_file: output/ex16opt_ic_1.out
362c4762a1bSJed Brown 
363c4762a1bSJed Brown TEST*/
364