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 Concepts: TS^time-dependent nonlinear problems 5c4762a1bSJed Brown Concepts: Automatic differentiation using ADOL-C 6c4762a1bSJed Brown Processors: 1 7c4762a1bSJed Brown */ 8c4762a1bSJed Brown /* 9c4762a1bSJed Brown REQUIRES configuration of PETSc with option --download-adolc. 10c4762a1bSJed Brown 11c4762a1bSJed Brown For documentation on ADOL-C, see 12c4762a1bSJed Brown $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf 13c4762a1bSJed Brown */ 14c4762a1bSJed Brown /* ------------------------------------------------------------------------ 15c4762a1bSJed Brown See ../advection-diffusion-reaction/ex1 for a description of the problem 16c4762a1bSJed Brown ------------------------------------------------------------------------- */ 17c4762a1bSJed Brown #include <petscts.h> 18c4762a1bSJed Brown #include "adolc-utils/drivers.cxx" 19c4762a1bSJed Brown #include <adolc/adolc.h> 20c4762a1bSJed Brown 21c4762a1bSJed Brown typedef struct { 22c4762a1bSJed Brown PetscScalar k; 23c4762a1bSJed Brown Vec initialsolution; 24c4762a1bSJed Brown AdolcCtx *adctx; /* Automatic differentiation support */ 25c4762a1bSJed Brown } AppCtx; 26c4762a1bSJed Brown 27c4762a1bSJed Brown PetscErrorCode IFunctionView(AppCtx *ctx,PetscViewer v) 28c4762a1bSJed Brown { 29c4762a1bSJed Brown PetscFunctionBegin; 30*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryWrite(v,&ctx->k,1,PETSC_SCALAR)); 31c4762a1bSJed Brown PetscFunctionReturn(0); 32c4762a1bSJed Brown } 33c4762a1bSJed Brown 34c4762a1bSJed Brown PetscErrorCode IFunctionLoad(AppCtx **ctx,PetscViewer v) 35c4762a1bSJed Brown { 36c4762a1bSJed Brown PetscFunctionBegin; 37*9566063dSJacob Faibussowitsch PetscCall(PetscNew(ctx)); 38*9566063dSJacob Faibussowitsch PetscCall(PetscViewerBinaryRead(v,&(*ctx)->k,1,NULL,PETSC_SCALAR)); 39c4762a1bSJed Brown PetscFunctionReturn(0); 40c4762a1bSJed Brown } 41c4762a1bSJed Brown 42c4762a1bSJed Brown /* 43c4762a1bSJed Brown Defines the ODE passed to the ODE solver 44c4762a1bSJed Brown */ 45c4762a1bSJed Brown PetscErrorCode IFunctionPassive(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx) 46c4762a1bSJed Brown { 47c4762a1bSJed Brown PetscScalar *f; 48c4762a1bSJed Brown const PetscScalar *u,*udot; 49c4762a1bSJed Brown 50c4762a1bSJed Brown PetscFunctionBegin; 51c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 52*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 53*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 54*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 55c4762a1bSJed Brown f[0] = udot[0] + ctx->k*u[0]*u[1]; 56c4762a1bSJed Brown f[1] = udot[1] + ctx->k*u[0]*u[1]; 57c4762a1bSJed Brown f[2] = udot[2] - ctx->k*u[0]*u[1]; 58*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 59*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 60*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 61c4762a1bSJed Brown PetscFunctionReturn(0); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown 'Active' ADOL-C annotated version, marking dependence upon u. 66c4762a1bSJed Brown */ 67c4762a1bSJed Brown PetscErrorCode IFunctionActive1(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown PetscScalar *f; 70c4762a1bSJed Brown const PetscScalar *u,*udot; 71c4762a1bSJed Brown 72c4762a1bSJed Brown adouble f_a[3]; /* 'active' double for dependent variables */ 73c4762a1bSJed Brown adouble u_a[3]; /* 'active' double for independent variables */ 74c4762a1bSJed Brown 75c4762a1bSJed Brown PetscFunctionBegin; 76c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 77*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 78*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 79*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 80c4762a1bSJed Brown 81c4762a1bSJed Brown /* Start of active section */ 82c4762a1bSJed Brown trace_on(1); 83c4762a1bSJed Brown u_a[0] <<= u[0]; u_a[1] <<= u[1]; u_a[2] <<= u[2]; /* Mark independence */ 84c4762a1bSJed Brown f_a[0] = udot[0] + ctx->k*u_a[0]*u_a[1]; 85c4762a1bSJed Brown f_a[1] = udot[1] + ctx->k*u_a[0]*u_a[1]; 86c4762a1bSJed Brown f_a[2] = udot[2] - ctx->k*u_a[0]*u_a[1]; 87c4762a1bSJed Brown f_a[0] >>= f[0]; f_a[1] >>= f[1]; f_a[2] >>= f[2]; /* Mark dependence */ 88c4762a1bSJed Brown trace_off(); 89c4762a1bSJed Brown /* End of active section */ 90c4762a1bSJed Brown 91*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 92*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 93*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 94c4762a1bSJed Brown PetscFunctionReturn(0); 95c4762a1bSJed Brown } 96c4762a1bSJed Brown 97c4762a1bSJed Brown /* 98c4762a1bSJed Brown 'Active' ADOL-C annotated version, marking dependence upon udot. 99c4762a1bSJed Brown */ 100c4762a1bSJed Brown PetscErrorCode IFunctionActive2(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx) 101c4762a1bSJed Brown { 102c4762a1bSJed Brown PetscScalar *f; 103c4762a1bSJed Brown const PetscScalar *u,*udot; 104c4762a1bSJed Brown 105c4762a1bSJed Brown adouble f_a[3]; /* 'active' double for dependent variables */ 106c4762a1bSJed Brown adouble udot_a[3]; /* 'active' double for independent variables */ 107c4762a1bSJed Brown 108c4762a1bSJed Brown PetscFunctionBegin; 109c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 110*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 111*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot,&udot)); 112*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(F,&f)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Start of active section */ 115c4762a1bSJed Brown trace_on(2); 116c4762a1bSJed Brown udot_a[0] <<= udot[0]; udot_a[1] <<= udot[1]; udot_a[2] <<= udot[2]; /* Mark independence */ 117c4762a1bSJed Brown f_a[0] = udot_a[0] + ctx->k*u[0]*u[1]; 118c4762a1bSJed Brown f_a[1] = udot_a[1] + ctx->k*u[0]*u[1]; 119c4762a1bSJed Brown f_a[2] = udot_a[2] - ctx->k*u[0]*u[1]; 120c4762a1bSJed Brown f_a[0] >>= f[0]; f_a[1] >>= f[1]; f_a[2] >>= f[2]; /* Mark dependence */ 121c4762a1bSJed Brown trace_off(); 122c4762a1bSJed Brown /* End of active section */ 123c4762a1bSJed Brown 124*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F,&f)); 125*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot,&udot)); 126*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 127c4762a1bSJed Brown PetscFunctionReturn(0); 128c4762a1bSJed Brown } 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* 131c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver, using the PETSc-ADOL-C driver for 132c4762a1bSJed Brown implicit TS. 133c4762a1bSJed Brown */ 134c4762a1bSJed Brown PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx) 135c4762a1bSJed Brown { 136c4762a1bSJed Brown AppCtx *appctx = (AppCtx*)ctx; 137410585f6SHong Zhang const PetscScalar *u; 138c4762a1bSJed Brown 139c4762a1bSJed Brown PetscFunctionBegin; 140*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U,&u)); 141*9566063dSJacob Faibussowitsch PetscCall(PetscAdolcComputeIJacobian(1,2,A,u,a,appctx->adctx)); 142*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U,&u)); 143c4762a1bSJed Brown PetscFunctionReturn(0); 144c4762a1bSJed Brown } 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* 147c4762a1bSJed Brown Defines the exact (analytic) solution to the ODE 148c4762a1bSJed Brown */ 149c4762a1bSJed Brown static PetscErrorCode Solution(TS ts,PetscReal t,Vec U,AppCtx *ctx) 150c4762a1bSJed Brown { 151c4762a1bSJed Brown const PetscScalar *uinit; 152c4762a1bSJed Brown PetscScalar *u,d0,q; 153c4762a1bSJed Brown 154c4762a1bSJed Brown PetscFunctionBegin; 155*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(ctx->initialsolution,&uinit)); 156*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(U,&u)); 157c4762a1bSJed Brown d0 = uinit[0] - uinit[1]; 158c4762a1bSJed Brown if (d0 == 0.0) q = ctx->k*t; 159c4762a1bSJed Brown else q = (1.0 - PetscExpScalar(-ctx->k*t*d0))/d0; 160c4762a1bSJed Brown u[0] = uinit[0]/(1.0 + uinit[1]*q); 161c4762a1bSJed Brown u[1] = u[0] - d0; 162c4762a1bSJed Brown u[2] = uinit[1] + uinit[2] - u[1]; 163*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U,&u)); 164*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(ctx->initialsolution,&uinit)); 165c4762a1bSJed Brown PetscFunctionReturn(0); 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 168c4762a1bSJed Brown int main(int argc,char **argv) 169c4762a1bSJed Brown { 170c4762a1bSJed Brown TS ts; /* ODE integrator */ 171c4762a1bSJed Brown Vec U,Udot,R; /* solution, derivative, residual */ 172c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 173c4762a1bSJed Brown PetscMPIInt size; 174c4762a1bSJed Brown PetscInt n = 3; 175c4762a1bSJed Brown AppCtx ctx; 176c4762a1bSJed Brown AdolcCtx *adctx; 177c4762a1bSJed Brown PetscScalar *u; 178c4762a1bSJed Brown const char * const names[] = {"U1","U2","U3",NULL}; 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Initialize program 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 183*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 184*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1852c71b3e2SJacob Faibussowitsch PetscCheckFalse(size > 1,PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs"); 186*9566063dSJacob Faibussowitsch PetscCall(PetscNew(&adctx)); 187c4762a1bSJed Brown adctx->m = n;adctx->n = n;adctx->p = n; 188c4762a1bSJed Brown ctx.adctx = adctx; 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191c4762a1bSJed Brown Create necessary matrix and vectors 192c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 193*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&A)); 194*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 195*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 196*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 197c4762a1bSJed Brown 198*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A,&U,NULL)); 199c4762a1bSJed Brown 200c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 201c4762a1bSJed Brown Set runtime options 202c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 203c4762a1bSJed Brown ctx.k = .9; 204*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetScalar(NULL,NULL,"-k",&ctx.k,NULL)); 205*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U,&ctx.initialsolution)); 206*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(ctx.initialsolution,&u)); 207c4762a1bSJed Brown u[0] = 1; 208c4762a1bSJed Brown u[1] = .7; 209c4762a1bSJed Brown u[2] = 0; 210*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(ctx.initialsolution,&u)); 211*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetVec(NULL,NULL,"-initial",ctx.initialsolution,NULL)); 212c4762a1bSJed Brown 213c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 214c4762a1bSJed Brown Create timestepping solver context 215c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 216*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 217*9566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 218*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSROSW)); 219*9566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunctionPassive,&ctx)); 220c4762a1bSJed Brown 221c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 222c4762a1bSJed Brown Set initial conditions 223c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 224*9566063dSJacob Faibussowitsch PetscCall(Solution(ts,0,U,&ctx)); 225*9566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts,U)); 226c4762a1bSJed Brown 227c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 228c4762a1bSJed Brown Trace just once for each tape 229c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 230*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U,&Udot)); 231*9566063dSJacob Faibussowitsch PetscCall(VecDuplicate(U,&R)); 232*9566063dSJacob Faibussowitsch PetscCall(IFunctionActive1(ts,0.,U,Udot,R,&ctx)); 233*9566063dSJacob Faibussowitsch PetscCall(IFunctionActive2(ts,0.,U,Udot,R,&ctx)); 234*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&R)); 235*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&Udot)); 236c4762a1bSJed Brown 237c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 238c4762a1bSJed Brown Set Jacobian 239c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 240*9566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx)); 241*9566063dSJacob Faibussowitsch PetscCall(TSSetSolutionFunction(ts,(TSSolutionFunction)Solution,&ctx)); 242c4762a1bSJed Brown 243c4762a1bSJed Brown { 244c4762a1bSJed Brown DM dm; 245c4762a1bSJed Brown void *ptr; 246*9566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 247*9566063dSJacob Faibussowitsch PetscCall(PetscDLSym(NULL,"IFunctionView",&ptr)); 248*9566063dSJacob Faibussowitsch PetscCall(PetscDLSym(NULL,"IFunctionLoad",&ptr)); 249*9566063dSJacob Faibussowitsch PetscCall(DMTSSetIFunctionSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad)); 250*9566063dSJacob Faibussowitsch PetscCall(DMTSSetIJacobianSerialize(dm,(PetscErrorCode (*)(void*,PetscViewer))IFunctionView,(PetscErrorCode (*)(void**,PetscViewer))IFunctionLoad)); 251c4762a1bSJed Brown } 252c4762a1bSJed Brown 253c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 254c4762a1bSJed Brown Set solver options 255c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 256*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,.001)); 257*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxSteps(ts,1000)); 258*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,20.0)); 259*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER)); 260*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 261*9566063dSJacob Faibussowitsch PetscCall(TSMonitorLGSetVariableNames(ts,names)); 262c4762a1bSJed Brown 263c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264c4762a1bSJed Brown Solve nonlinear system 265c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 266*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,U)); 267c4762a1bSJed Brown 268*9566063dSJacob Faibussowitsch PetscCall(TSView(ts,PETSC_VIEWER_BINARY_WORLD)); 269c4762a1bSJed Brown 270c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 271c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 272c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 273*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ctx.initialsolution)); 274*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 275*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 276*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 277*9566063dSJacob Faibussowitsch PetscCall(PetscFree(adctx)); 278*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 279b122ec5aSJacob Faibussowitsch return 0; 280c4762a1bSJed Brown } 281c4762a1bSJed Brown 282c4762a1bSJed Brown /*TEST 283c4762a1bSJed Brown 284c4762a1bSJed Brown build: 285c4762a1bSJed Brown requires: double !complex adolc 286c4762a1bSJed Brown 287c4762a1bSJed Brown test: 288c4762a1bSJed Brown suffix: 1 289c4762a1bSJed Brown args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor 290c4762a1bSJed Brown output_file: output/adr_ex1_1.out 291c4762a1bSJed Brown 292c4762a1bSJed Brown test: 293c4762a1bSJed Brown suffix: 2 294c4762a1bSJed Brown args: -ts_max_steps 1 -snes_test_jacobian 295c4762a1bSJed Brown output_file: output/adr_ex1_2.out 296c4762a1bSJed Brown 297c4762a1bSJed Brown TEST*/ 298