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