xref: /petsc/src/ts/tutorials/autodiff/adr_ex1.cxx (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for a nonlinear reaction problem from chemistry.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    REQUIRES configuration of PETSc with option --download-adolc.
5c4762a1bSJed Brown 
6c4762a1bSJed Brown    For documentation on ADOL-C, see
7c4762a1bSJed Brown      $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf
8c4762a1bSJed Brown */
9c4762a1bSJed Brown /* ------------------------------------------------------------------------
10c4762a1bSJed Brown   See ../advection-diffusion-reaction/ex1 for a description of the problem
11c4762a1bSJed Brown   ------------------------------------------------------------------------- */
12c4762a1bSJed Brown #include <petscts.h>
13c4762a1bSJed Brown #include "adolc-utils/drivers.cxx"
14c4762a1bSJed Brown #include <adolc/adolc.h>
15c4762a1bSJed Brown 
16c4762a1bSJed Brown typedef struct {
17c4762a1bSJed Brown   PetscScalar k;
18c4762a1bSJed Brown   Vec         initialsolution;
19c4762a1bSJed Brown   AdolcCtx   *adctx; /* Automatic differentiation support */
20c4762a1bSJed Brown } AppCtx;
21c4762a1bSJed Brown 
22*9371c9d4SSatish Balay PetscErrorCode IFunctionView(AppCtx *ctx, PetscViewer v) {
23c4762a1bSJed Brown   PetscFunctionBegin;
249566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(v, &ctx->k, 1, PETSC_SCALAR));
25c4762a1bSJed Brown   PetscFunctionReturn(0);
26c4762a1bSJed Brown }
27c4762a1bSJed Brown 
28*9371c9d4SSatish Balay PetscErrorCode IFunctionLoad(AppCtx **ctx, PetscViewer v) {
29c4762a1bSJed Brown   PetscFunctionBegin;
309566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
319566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(v, &(*ctx)->k, 1, NULL, PETSC_SCALAR));
32c4762a1bSJed Brown   PetscFunctionReturn(0);
33c4762a1bSJed Brown }
34c4762a1bSJed Brown 
35c4762a1bSJed Brown /*
36c4762a1bSJed Brown   Defines the ODE passed to the ODE solver
37c4762a1bSJed Brown */
38*9371c9d4SSatish Balay PetscErrorCode IFunctionPassive(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) {
39c4762a1bSJed Brown   PetscScalar       *f;
40c4762a1bSJed Brown   const PetscScalar *u, *udot;
41c4762a1bSJed Brown 
42c4762a1bSJed Brown   PetscFunctionBegin;
43c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
47c4762a1bSJed Brown   f[0] = udot[0] + ctx->k * u[0] * u[1];
48c4762a1bSJed Brown   f[1] = udot[1] + ctx->k * u[0] * u[1];
49c4762a1bSJed Brown   f[2] = udot[2] - ctx->k * u[0] * u[1];
509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
53c4762a1bSJed Brown   PetscFunctionReturn(0);
54c4762a1bSJed Brown }
55c4762a1bSJed Brown 
56c4762a1bSJed Brown /*
57c4762a1bSJed Brown   'Active' ADOL-C annotated version, marking dependence upon u.
58c4762a1bSJed Brown */
59*9371c9d4SSatish Balay PetscErrorCode IFunctionActive1(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) {
60c4762a1bSJed Brown   PetscScalar       *f;
61c4762a1bSJed Brown   const PetscScalar *u, *udot;
62c4762a1bSJed Brown 
63c4762a1bSJed Brown   adouble f_a[3]; /* 'active' double for dependent variables */
64c4762a1bSJed Brown   adouble u_a[3]; /* 'active' double for independent variables */
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   PetscFunctionBegin;
67c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
689566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
709566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
71c4762a1bSJed Brown 
72c4762a1bSJed Brown   /* Start of active section */
73c4762a1bSJed Brown   trace_on(1);
74*9371c9d4SSatish Balay   u_a[0] <<= u[0];
75*9371c9d4SSatish Balay   u_a[1] <<= u[1];
76*9371c9d4SSatish Balay   u_a[2] <<= u[2]; /* Mark independence */
77c4762a1bSJed Brown   f_a[0] = udot[0] + ctx->k * u_a[0] * u_a[1];
78c4762a1bSJed Brown   f_a[1] = udot[1] + ctx->k * u_a[0] * u_a[1];
79c4762a1bSJed Brown   f_a[2] = udot[2] - ctx->k * u_a[0] * u_a[1];
80*9371c9d4SSatish Balay   f_a[0] >>= f[0];
81*9371c9d4SSatish Balay   f_a[1] >>= f[1];
82*9371c9d4SSatish Balay   f_a[2] >>= f[2]; /* Mark dependence */
83c4762a1bSJed Brown   trace_off();
84c4762a1bSJed Brown   /* End of active section */
85c4762a1bSJed Brown 
869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
89c4762a1bSJed Brown   PetscFunctionReturn(0);
90c4762a1bSJed Brown }
91c4762a1bSJed Brown 
92c4762a1bSJed Brown /*
93c4762a1bSJed Brown   'Active' ADOL-C annotated version, marking dependence upon udot.
94c4762a1bSJed Brown */
95*9371c9d4SSatish Balay PetscErrorCode IFunctionActive2(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) {
96c4762a1bSJed Brown   PetscScalar       *f;
97c4762a1bSJed Brown   const PetscScalar *u, *udot;
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   adouble f_a[3];    /* 'active' double for dependent variables */
100c4762a1bSJed Brown   adouble udot_a[3]; /* 'active' double for independent variables */
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBegin;
103c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot, &udot));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   /* Start of active section */
109c4762a1bSJed Brown   trace_on(2);
110*9371c9d4SSatish Balay   udot_a[0] <<= udot[0];
111*9371c9d4SSatish Balay   udot_a[1] <<= udot[1];
112*9371c9d4SSatish Balay   udot_a[2] <<= udot[2]; /* Mark independence */
113c4762a1bSJed Brown   f_a[0] = udot_a[0] + ctx->k * u[0] * u[1];
114c4762a1bSJed Brown   f_a[1] = udot_a[1] + ctx->k * u[0] * u[1];
115c4762a1bSJed Brown   f_a[2] = udot_a[2] - ctx->k * u[0] * u[1];
116*9371c9d4SSatish Balay   f_a[0] >>= f[0];
117*9371c9d4SSatish Balay   f_a[1] >>= f[1];
118*9371c9d4SSatish Balay   f_a[2] >>= f[2]; /* Mark dependence */
119c4762a1bSJed Brown   trace_off();
120c4762a1bSJed Brown   /* End of active section */
121c4762a1bSJed Brown 
1229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot, &udot));
1249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
125c4762a1bSJed Brown   PetscFunctionReturn(0);
126c4762a1bSJed Brown }
127c4762a1bSJed Brown 
128c4762a1bSJed Brown /*
129c4762a1bSJed Brown  Defines the Jacobian of the ODE passed to the ODE solver, using the PETSc-ADOL-C driver for
130c4762a1bSJed Brown  implicit TS.
131c4762a1bSJed Brown */
132*9371c9d4SSatish Balay PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) {
133c4762a1bSJed Brown   AppCtx            *appctx = (AppCtx *)ctx;
134410585f6SHong Zhang   const PetscScalar *u;
135c4762a1bSJed Brown 
136c4762a1bSJed Brown   PetscFunctionBegin;
1379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U, &u));
1389566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeIJacobian(1, 2, A, u, a, appctx->adctx));
1399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U, &u));
140c4762a1bSJed Brown   PetscFunctionReturn(0);
141c4762a1bSJed Brown }
142c4762a1bSJed Brown 
143c4762a1bSJed Brown /*
144c4762a1bSJed Brown      Defines the exact (analytic) solution to the ODE
145c4762a1bSJed Brown */
146*9371c9d4SSatish Balay static PetscErrorCode Solution(TS ts, PetscReal t, Vec U, AppCtx *ctx) {
147c4762a1bSJed Brown   const PetscScalar *uinit;
148c4762a1bSJed Brown   PetscScalar       *u, d0, q;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   PetscFunctionBegin;
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->initialsolution, &uinit));
1529566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
153c4762a1bSJed Brown   d0 = uinit[0] - uinit[1];
154c4762a1bSJed Brown   if (d0 == 0.0) q = ctx->k * t;
155c4762a1bSJed Brown   else q = (1.0 - PetscExpScalar(-ctx->k * t * d0)) / d0;
156c4762a1bSJed Brown   u[0] = uinit[0] / (1.0 + uinit[1] * q);
157c4762a1bSJed Brown   u[1] = u[0] - d0;
158c4762a1bSJed Brown   u[2] = uinit[1] + uinit[2] - u[1];
1599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
1609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->initialsolution, &uinit));
161c4762a1bSJed Brown   PetscFunctionReturn(0);
162c4762a1bSJed Brown }
163c4762a1bSJed Brown 
164*9371c9d4SSatish Balay int main(int argc, char **argv) {
165c4762a1bSJed Brown   TS                ts;         /* ODE integrator */
166c4762a1bSJed Brown   Vec               U, Udot, R; /* solution, derivative, residual */
167c4762a1bSJed Brown   Mat               A;          /* Jacobian matrix */
168c4762a1bSJed Brown   PetscMPIInt       size;
169c4762a1bSJed Brown   PetscInt          n = 3;
170c4762a1bSJed Brown   AppCtx            ctx;
171c4762a1bSJed Brown   AdolcCtx         *adctx;
172c4762a1bSJed Brown   PetscScalar      *u;
173c4762a1bSJed Brown   const char *const names[] = {"U1", "U2", "U3", NULL};
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176c4762a1bSJed Brown      Initialize program
177c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178327415f7SBarry Smith   PetscFunctionBeginUser;
1799566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1809566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
18108401ef6SPierre Jolivet   PetscCheck(size <= 1, PETSC_COMM_WORLD, PETSC_ERR_SUP, "Only for sequential runs");
1829566063dSJacob Faibussowitsch   PetscCall(PetscNew(&adctx));
183*9371c9d4SSatish Balay   adctx->m  = n;
184*9371c9d4SSatish Balay   adctx->n  = n;
185*9371c9d4SSatish Balay   adctx->p  = n;
186c4762a1bSJed Brown   ctx.adctx = adctx;
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189c4762a1bSJed Brown     Create necessary matrix and vectors
190c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1919566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1929566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
1939566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1949566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
195c4762a1bSJed Brown 
1969566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199c4762a1bSJed Brown     Set runtime options
200c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201c4762a1bSJed Brown   ctx.k = .9;
2029566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL, NULL, "-k", &ctx.k, NULL));
2039566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &ctx.initialsolution));
2049566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx.initialsolution, &u));
205c4762a1bSJed Brown   u[0] = 1;
206c4762a1bSJed Brown   u[1] = .7;
207c4762a1bSJed Brown   u[2] = 0;
2089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx.initialsolution, &u));
2099566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetVec(NULL, NULL, "-initial", ctx.initialsolution, NULL));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212c4762a1bSJed Brown      Create timestepping solver context
213c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2149566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2159566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
2169566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW));
2179566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunctionPassive, &ctx));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220c4762a1bSJed Brown      Set initial conditions
221c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2229566063dSJacob Faibussowitsch   PetscCall(Solution(ts, 0, U, &ctx));
2239566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
224c4762a1bSJed Brown 
225c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226c4762a1bSJed Brown      Trace just once for each tape
227c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2289566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &Udot));
2299566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U, &R));
2309566063dSJacob Faibussowitsch   PetscCall(IFunctionActive1(ts, 0., U, Udot, R, &ctx));
2319566063dSJacob Faibussowitsch   PetscCall(IFunctionActive2(ts, 0., U, Udot, R, &ctx));
2329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
2339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Udot));
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236c4762a1bSJed Brown      Set Jacobian
237c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2389566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx));
2399566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(ts, (TSSolutionFunction)Solution, &ctx));
240c4762a1bSJed Brown 
241c4762a1bSJed Brown   {
242c4762a1bSJed Brown     DM    dm;
243c4762a1bSJed Brown     void *ptr;
2449566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts, &dm));
2459566063dSJacob Faibussowitsch     PetscCall(PetscDLSym(NULL, "IFunctionView", &ptr));
2469566063dSJacob Faibussowitsch     PetscCall(PetscDLSym(NULL, "IFunctionLoad", &ptr));
2479566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunctionSerialize(dm, (PetscErrorCode(*)(void *, PetscViewer))IFunctionView, (PetscErrorCode(*)(void **, PetscViewer))IFunctionLoad));
2489566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobianSerialize(dm, (PetscErrorCode(*)(void *, PetscViewer))IFunctionView, (PetscErrorCode(*)(void **, PetscViewer))IFunctionLoad));
249c4762a1bSJed Brown   }
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252c4762a1bSJed Brown      Set solver options
253c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2549566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .001));
2559566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 1000));
2569566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 20.0));
2579566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
2589566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2599566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGSetVariableNames(ts, names));
260c4762a1bSJed Brown 
261c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
262c4762a1bSJed Brown      Solve nonlinear system
263c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2649566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
265c4762a1bSJed Brown 
2669566063dSJacob Faibussowitsch   PetscCall(TSView(ts, PETSC_VIEWER_BINARY_WORLD));
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
269c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
270c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.initialsolution));
2729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2739566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2749566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2759566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
2769566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
277b122ec5aSJacob Faibussowitsch   return 0;
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280c4762a1bSJed Brown /*TEST
281c4762a1bSJed Brown 
282c4762a1bSJed Brown   build:
283c4762a1bSJed Brown     requires: double !complex adolc
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   test:
286c4762a1bSJed Brown     suffix: 1
287c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
288c4762a1bSJed Brown     output_file: output/adr_ex1_1.out
289c4762a1bSJed Brown 
290c4762a1bSJed Brown   test:
291c4762a1bSJed Brown     suffix: 2
292c4762a1bSJed Brown     args: -ts_max_steps 1 -snes_test_jacobian
293c4762a1bSJed Brown     output_file: output/adr_ex1_2.out
294c4762a1bSJed Brown 
295c4762a1bSJed Brown TEST*/
296