xref: /petsc/src/ts/tutorials/autodiff/adr_ex1.cxx (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 
22c4762a1bSJed Brown PetscErrorCode IFunctionView(AppCtx *ctx,PetscViewer v)
23c4762a1bSJed Brown {
24c4762a1bSJed Brown   PetscFunctionBegin;
259566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryWrite(v,&ctx->k,1,PETSC_SCALAR));
26c4762a1bSJed Brown   PetscFunctionReturn(0);
27c4762a1bSJed Brown }
28c4762a1bSJed Brown 
29c4762a1bSJed Brown PetscErrorCode IFunctionLoad(AppCtx **ctx,PetscViewer v)
30c4762a1bSJed Brown {
31c4762a1bSJed Brown   PetscFunctionBegin;
329566063dSJacob Faibussowitsch   PetscCall(PetscNew(ctx));
339566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryRead(v,&(*ctx)->k,1,NULL,PETSC_SCALAR));
34c4762a1bSJed Brown   PetscFunctionReturn(0);
35c4762a1bSJed Brown }
36c4762a1bSJed Brown 
37c4762a1bSJed Brown /*
38c4762a1bSJed Brown   Defines the ODE passed to the ODE solver
39c4762a1bSJed Brown */
40c4762a1bSJed Brown PetscErrorCode IFunctionPassive(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
41c4762a1bSJed Brown {
42c4762a1bSJed Brown   PetscScalar       *f;
43c4762a1bSJed Brown   const PetscScalar *u,*udot;
44c4762a1bSJed Brown 
45c4762a1bSJed Brown   PetscFunctionBegin;
46c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
50c4762a1bSJed Brown   f[0] = udot[0] + ctx->k*u[0]*u[1];
51c4762a1bSJed Brown   f[1] = udot[1] + ctx->k*u[0]*u[1];
52c4762a1bSJed Brown   f[2] = udot[2] - ctx->k*u[0]*u[1];
539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
56c4762a1bSJed Brown   PetscFunctionReturn(0);
57c4762a1bSJed Brown }
58c4762a1bSJed Brown 
59c4762a1bSJed Brown /*
60c4762a1bSJed Brown   'Active' ADOL-C annotated version, marking dependence upon u.
61c4762a1bSJed Brown */
62c4762a1bSJed Brown PetscErrorCode IFunctionActive1(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
63c4762a1bSJed Brown {
64c4762a1bSJed Brown   PetscScalar       *f;
65c4762a1bSJed Brown   const PetscScalar *u,*udot;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   adouble           f_a[3]; /* 'active' double for dependent variables */
68c4762a1bSJed Brown   adouble           u_a[3]; /* 'active' double for independent variables */
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBegin;
71c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
729566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
739566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
749566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   /* Start of active section */
77c4762a1bSJed Brown   trace_on(1);
78c4762a1bSJed Brown   u_a[0] <<= u[0]; u_a[1] <<= u[1]; u_a[2] <<= u[2]; /* Mark independence */
79c4762a1bSJed Brown   f_a[0] = udot[0] + ctx->k*u_a[0]*u_a[1];
80c4762a1bSJed Brown   f_a[1] = udot[1] + ctx->k*u_a[0]*u_a[1];
81c4762a1bSJed Brown   f_a[2] = udot[2] - ctx->k*u_a[0]*u_a[1];
82c4762a1bSJed Brown   f_a[0] >>= f[0]; f_a[1] >>= f[1]; 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 */
95c4762a1bSJed Brown PetscErrorCode IFunctionActive2(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
96c4762a1bSJed Brown {
97c4762a1bSJed Brown   PetscScalar       *f;
98c4762a1bSJed Brown   const PetscScalar *u,*udot;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   adouble           f_a[3];    /* 'active' double for dependent variables */
101c4762a1bSJed Brown   adouble           udot_a[3]; /* 'active' double for independent variables */
102c4762a1bSJed Brown 
103c4762a1bSJed Brown   PetscFunctionBegin;
104c4762a1bSJed Brown   /*  The next three lines allow us to access the entries of the vectors directly */
1059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Udot,&udot));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
108c4762a1bSJed Brown 
109c4762a1bSJed Brown   /* Start of active section */
110c4762a1bSJed Brown   trace_on(2);
111c4762a1bSJed Brown   udot_a[0] <<= udot[0]; udot_a[1] <<= udot[1]; udot_a[2] <<= udot[2]; /* Mark independence */
112c4762a1bSJed Brown   f_a[0] = udot_a[0] + ctx->k*u[0]*u[1];
113c4762a1bSJed Brown   f_a[1] = udot_a[1] + ctx->k*u[0]*u[1];
114c4762a1bSJed Brown   f_a[2] = udot_a[2] - ctx->k*u[0]*u[1];
115c4762a1bSJed Brown   f_a[0] >>= f[0]; f_a[1] >>= f[1]; f_a[2] >>= f[2];                   /* Mark dependence */
116c4762a1bSJed Brown   trace_off();
117c4762a1bSJed Brown   /* End of active section */
118c4762a1bSJed Brown 
1199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
1209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Udot,&udot));
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
122c4762a1bSJed Brown   PetscFunctionReturn(0);
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown /*
126c4762a1bSJed Brown  Defines the Jacobian of the ODE passed to the ODE solver, using the PETSc-ADOL-C driver for
127c4762a1bSJed Brown  implicit TS.
128c4762a1bSJed Brown */
129c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
130c4762a1bSJed Brown {
131c4762a1bSJed Brown   AppCtx            *appctx = (AppCtx*)ctx;
132410585f6SHong Zhang   const PetscScalar *u;
133c4762a1bSJed Brown 
134c4762a1bSJed Brown   PetscFunctionBegin;
1359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
1369566063dSJacob Faibussowitsch   PetscCall(PetscAdolcComputeIJacobian(1,2,A,u,a,appctx->adctx));
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
138c4762a1bSJed Brown   PetscFunctionReturn(0);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown /*
142c4762a1bSJed Brown      Defines the exact (analytic) solution to the ODE
143c4762a1bSJed Brown */
144c4762a1bSJed Brown static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx)
145c4762a1bSJed Brown {
146c4762a1bSJed Brown   const PetscScalar *uinit;
147c4762a1bSJed Brown   PetscScalar       *u,d0,q;
148c4762a1bSJed Brown 
149c4762a1bSJed Brown   PetscFunctionBegin;
1509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(ctx->initialsolution,&uinit));
1519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U,&u));
152c4762a1bSJed Brown   d0   = uinit[0] - uinit[1];
153c4762a1bSJed Brown   if (d0 == 0.0) q = ctx->k*t;
154c4762a1bSJed Brown   else q = (1.0 - PetscExpScalar(-ctx->k*t*d0))/d0;
155c4762a1bSJed Brown   u[0] = uinit[0]/(1.0 + uinit[1]*q);
156c4762a1bSJed Brown   u[1] = u[0] - d0;
157c4762a1bSJed Brown   u[2] = uinit[1] + uinit[2] - u[1];
1589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U,&u));
1599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(ctx->initialsolution,&uinit));
160c4762a1bSJed Brown   PetscFunctionReturn(0);
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown int main(int argc,char **argv)
164c4762a1bSJed Brown {
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      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178*327415f7SBarry 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));
183c4762a1bSJed Brown   adctx->m = n;adctx->n = n;adctx->p = n;
184c4762a1bSJed Brown   ctx.adctx = adctx;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown     Create necessary matrix and vectors
188c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1899566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
1909566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE));
1919566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1929566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
193c4762a1bSJed Brown 
1949566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&U,NULL));
195c4762a1bSJed Brown 
196c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
197c4762a1bSJed Brown     Set runtime options
198c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
199c4762a1bSJed Brown   ctx.k = .9;
2009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalar(NULL,NULL,"-k",&ctx.k,NULL));
2019566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U,&ctx.initialsolution));
2029566063dSJacob Faibussowitsch   PetscCall(VecGetArray(ctx.initialsolution,&u));
203c4762a1bSJed Brown   u[0]  = 1;
204c4762a1bSJed Brown   u[1]  = .7;
205c4762a1bSJed Brown   u[2]  = 0;
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(ctx.initialsolution,&u));
2079566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetVec(NULL,NULL,"-initial",ctx.initialsolution,NULL));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Create timestepping solver context
211c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2129566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
2139566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
2149566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSROSW));
2159566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunctionPassive,&ctx));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218c4762a1bSJed Brown      Set initial conditions
219c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2209566063dSJacob Faibussowitsch   PetscCall(Solution(ts,0,U,&ctx));
2219566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,U));
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224c4762a1bSJed Brown      Trace just once for each tape
225c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2269566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U,&Udot));
2279566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(U,&R));
2289566063dSJacob Faibussowitsch   PetscCall(IFunctionActive1(ts,0.,U,Udot,R,&ctx));
2299566063dSJacob Faibussowitsch   PetscCall(IFunctionActive2(ts,0.,U,Udot,R,&ctx));
2309566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&R));
2319566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Udot));
232c4762a1bSJed Brown 
233c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234c4762a1bSJed Brown      Set Jacobian
235c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2369566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx));
2379566063dSJacob Faibussowitsch   PetscCall(TSSetSolutionFunction(ts,(TSSolutionFunction)Solution,&ctx));
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   {
240c4762a1bSJed Brown     DM   dm;
241c4762a1bSJed Brown     void *ptr;
2429566063dSJacob Faibussowitsch     PetscCall(TSGetDM(ts,&dm));
2439566063dSJacob Faibussowitsch     PetscCall(PetscDLSym(NULL,"IFunctionView",&ptr));
2449566063dSJacob Faibussowitsch     PetscCall(PetscDLSym(NULL,"IFunctionLoad",&ptr));
2459566063dSJacob Faibussowitsch     PetscCall(DMTSSetIFunctionSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad));
2469566063dSJacob Faibussowitsch     PetscCall(DMTSSetIJacobianSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad));
247c4762a1bSJed Brown   }
248c4762a1bSJed Brown 
249c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250c4762a1bSJed Brown      Set solver options
251c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2529566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.001));
2539566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts,1000));
2549566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,20.0));
2559566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
2569566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
2579566063dSJacob Faibussowitsch   PetscCall(TSMonitorLGSetVariableNames(ts,names));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260c4762a1bSJed Brown      Solve nonlinear system
261c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2629566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,U));
263c4762a1bSJed Brown 
2649566063dSJacob Faibussowitsch   PetscCall(TSView(ts,PETSC_VIEWER_BINARY_WORLD));
265c4762a1bSJed Brown 
266c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
267c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
268c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&ctx.initialsolution));
2709566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2719566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2729566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
2739566063dSJacob Faibussowitsch   PetscCall(PetscFree(adctx));
2749566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
275b122ec5aSJacob Faibussowitsch   return 0;
276c4762a1bSJed Brown }
277c4762a1bSJed Brown 
278c4762a1bSJed Brown /*TEST
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   build:
281c4762a1bSJed Brown     requires: double !complex adolc
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   test:
284c4762a1bSJed Brown     suffix: 1
285c4762a1bSJed Brown     args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor
286c4762a1bSJed Brown     output_file: output/adr_ex1_1.out
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   test:
289c4762a1bSJed Brown     suffix: 2
290c4762a1bSJed Brown     args: -ts_max_steps 1 -snes_test_jacobian
291c4762a1bSJed Brown     output_file: output/adr_ex1_2.out
292c4762a1bSJed Brown 
293c4762a1bSJed Brown TEST*/
294