xref: /petsc/src/ts/tutorials/autodiff/ex16opt_ic.cxx (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
55*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
56c4762a1bSJed Brown   f[0] = x[1];
57c4762a1bSJed Brown   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
58*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
59*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
78*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
79*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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 
88*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
89*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
102*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
103*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscAdolcComputeRHSJacobian(1,A,x,user->adctx));
104*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
118*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
119*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetMaxTime(ts,&tfinal));
120*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetPrevTime(ts,&tprev));
121*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
122*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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])));
123*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev));
124*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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   PetscErrorCode     ierr;
138c4762a1bSJed Brown   Tao                tao;
139c4762a1bSJed Brown   KSP                ksp;
140c4762a1bSJed Brown   PC                 pc;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143c4762a1bSJed Brown      Initialize program
144c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
146*5f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1472c71b3e2SJacob Faibussowitsch   PetscCheckFalse(size != 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown     Set runtime options and create AdolcCtx
151c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscNew(&adctx));
153c4762a1bSJed Brown   user.mu          = 1.0;
154c4762a1bSJed Brown   user.next_output = 0.0;
155c4762a1bSJed Brown   user.steps       = 0;
156c4762a1bSJed Brown   user.ftime       = 0.5;
157c4762a1bSJed Brown   adctx->m = 2;adctx->n = 2;adctx->p = 2;
158c4762a1bSJed Brown   user.adctx = adctx;
159c4762a1bSJed Brown 
160*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
161*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
165c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.A));
167*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
168*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.A));
169*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.A));
170*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.x,NULL));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173c4762a1bSJed Brown      Set initial conditions
174c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user.x,&x_ptr));
176c4762a1bSJed Brown   x_ptr[0] = 2.0;   x_ptr[1] = 0.66666654321;
177*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user.x,&x_ptr));
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180c4762a1bSJed Brown      Trace just once on each tape and put zeros on Jacobian diagonal
181c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(user.x,&r));
183*5f80ce2aSJacob Faibussowitsch   CHKERRQ(RHSFunctionActive(ts,0.,user.x,r,&user));
184*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecSet(r,0));
185*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDiagonalSet(user.A,r,INSERT_VALUES));
186*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&r));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189c4762a1bSJed Brown      Create timestepping solver context
190c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
192*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSRK));
193*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user));
194*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user));
195*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,user.ftime));
196*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
197c4762a1bSJed Brown   if (monitor) {
198*5f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL));
199c4762a1bSJed Brown   }
200c4762a1bSJed Brown 
201*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
202*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime)));
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
206c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(ts));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Set runtime options
211c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215c4762a1bSJed Brown      Solve nonlinear system
216c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,user.x));
218*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&(user.ftime)));
219*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&user.steps));
220*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime));
221c4762a1bSJed Brown 
222*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user.x,&x_ptr));
223c4762a1bSJed Brown   user.x_ob[0] = x_ptr[0];
224c4762a1bSJed Brown   user.x_ob[1] = x_ptr[1];
225*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user.x,&x_ptr));
226c4762a1bSJed Brown 
227*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&user.lambda[0],NULL));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
230*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
231*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetType(tao,TAOCG));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* Set initial solution guess */
234*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.A,&ic,NULL));
235*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(ic,&x_ptr));
236c4762a1bSJed Brown   x_ptr[0]  = 2.1;
237c4762a1bSJed Brown   x_ptr[1]  = 0.7;
238*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(ic,&x_ptr));
239c4762a1bSJed Brown 
240*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetSolution(tao,ic));
241c4762a1bSJed Brown 
242c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
243*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user));
244c4762a1bSJed Brown 
245c4762a1bSJed Brown   /* Check for any TAO command line options */
246*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetFromOptions(tao));
247*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoGetKSP(tao,&ksp));
248c4762a1bSJed Brown   if (ksp) {
249*5f80ce2aSJacob Faibussowitsch     CHKERRQ(KSPGetPC(ksp,&pc));
250*5f80ce2aSJacob Faibussowitsch     CHKERRQ(PCSetType(pc,PCNONE));
251c4762a1bSJed Brown   }
252c4762a1bSJed Brown 
253*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSetTolerances(tao,1e-10,PETSC_DEFAULT,PETSC_DEFAULT));
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
256*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoSolve(tao));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   /* Free TAO data structures */
259*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TaoDestroy(&tao));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
263c4762a1bSJed Brown      are no longer needed.
264c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
265*5f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.A));
266*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.x));
267*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.lambda[0]));
268*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
269*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&ic));
270*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscFree(adctx));
271c4762a1bSJed Brown   ierr = PetscFinalize();
272c4762a1bSJed Brown   return ierr;
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
275c4762a1bSJed Brown /* ------------------------------------------------------------------ */
276c4762a1bSJed Brown /*
277c4762a1bSJed Brown    FormFunctionGradient - Evaluates the function and corresponding gradient.
278c4762a1bSJed Brown 
279c4762a1bSJed Brown    Input Parameters:
280c4762a1bSJed Brown    tao - the Tao context
281c4762a1bSJed Brown    X   - the input vector
282a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
283c4762a1bSJed Brown 
284c4762a1bSJed Brown    Output Parameters:
285c4762a1bSJed Brown    f   - the newly evaluated function
286c4762a1bSJed Brown    G   - the newly evaluated gradient
287c4762a1bSJed Brown */
288c4762a1bSJed Brown PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
289c4762a1bSJed Brown {
290c4762a1bSJed Brown   User              user = (User)ctx;
291c4762a1bSJed Brown   TS                ts;
292c4762a1bSJed Brown   PetscScalar       *x_ptr,*y_ptr;
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   PetscFunctionBeginUser;
295*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(IC,user->x));
296c4762a1bSJed Brown 
297c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
298c4762a1bSJed Brown      Create timestepping solver context
299c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
300*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
301*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSRK));
302*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunctionPassive,user));
303c4762a1bSJed Brown   /*   Set RHS Jacobian  for the adjoint integration */
304*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user));
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307c4762a1bSJed Brown      Set time
308c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
309*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTime(ts,0.0));
310*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.001));
311*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,0.5));
312*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
313c4762a1bSJed Brown 
314*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTolerances(ts,1e-7,NULL,1e-7,NULL));
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
317c4762a1bSJed Brown     Save trajectory of solution so that TSAdjointSolve() may be used
318c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
319*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(ts));
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
322c4762a1bSJed Brown      Set runtime options
323c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
324*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
325c4762a1bSJed Brown 
326*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,user->x));
327*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&user->ftime));
328*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&user->steps));
329*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %.6f, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime));
330c4762a1bSJed Brown 
331*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->x,&x_ptr));
332c4762a1bSJed 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]);
333*5f80ce2aSJacob Faibussowitsch   CHKERRQ(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)));
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
336c4762a1bSJed Brown      Adjoint model starts here
337c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
338c4762a1bSJed Brown   /*   Redet initial conditions for the adjoint integration */
339*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user->lambda[0],&y_ptr));
340c4762a1bSJed Brown   y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
341c4762a1bSJed Brown   y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
342*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->lambda[0],&y_ptr));
343*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user->x,&x_ptr));
344*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,user->lambda,NULL));
345c4762a1bSJed Brown 
346*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
347c4762a1bSJed Brown 
348*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(user->lambda[0],G));
349c4762a1bSJed Brown 
350*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
351c4762a1bSJed Brown   PetscFunctionReturn(0);
352c4762a1bSJed Brown }
353c4762a1bSJed Brown 
354c4762a1bSJed Brown /*TEST
355c4762a1bSJed Brown 
356c4762a1bSJed Brown   build:
357c4762a1bSJed Brown     requires: double !complex adolc
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   test:
360c4762a1bSJed Brown     suffix: 1
361c4762a1bSJed Brown     args: -ts_rhs_jacobian_test_mult_transpose FALSE -tao_max_it 2 -ts_rhs_jacobian_test_mult FALSE
362c4762a1bSJed Brown     output_file: output/ex16opt_ic_1.out
363c4762a1bSJed Brown 
364c4762a1bSJed Brown TEST*/
365