static char help[] = "Demonstrates automatic Jacobian generation using ADOL-C for an adjoint sensitivity analysis of the van der Pol equation.\n\ Input parameters include:\n\ -mu : stiffness parameter\n\n"; /* Concepts: TS^time-dependent nonlinear problems Concepts: TS^van der Pol equation Concepts: TS^adjoint sensitivity analysis Concepts: Automatic differentation using ADOL-C Concepts: Automatic differentation w.r.t. a parameter using ADOL-C Processors: 1 */ /* REQUIRES configuration of PETSc with option --download-adolc. For documentation on ADOL-C, see $PETSC_ARCH/externalpackages/ADOL-C-2.6.0/ADOL-C/doc/adolc-manual.pdf */ /* ------------------------------------------------------------------------ See ex16adj for a description of the problem being solved. ------------------------------------------------------------------------- */ #include #include #include "adolc-utils/drivers.cxx" #include typedef struct _n_User *User; struct _n_User { PetscReal mu; PetscReal next_output; PetscReal tprev; /* Automatic differentiation support */ AdolcCtx *adctx; }; /* 'Passive' RHS function, used in residual evaluations during the time integration. */ static PetscErrorCode RHSFunctionPassive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; PetscScalar *f; const PetscScalar *x; PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); f[0] = x[1]; f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0]; ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Trace RHS to mark on tape 1 the dependence of f upon x. This tape is used in generating the Jacobian transform. */ static PetscErrorCode RHSFunctionActive(TS ts,PetscReal t,Vec X,Vec F,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; PetscScalar *f; const PetscScalar *x; adouble f_a[2]; /* 'active' double for dependent variables */ adouble x_a[2]; /* 'active' double for independent variables */ PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); /* Start of active section */ trace_on(1); x_a[0] <<= x[0];x_a[1] <<= x[1]; /* Mark independence */ f_a[0] = x_a[1]; f_a[1] = user->mu*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0]; f_a[0] >>= f[0];f_a[1] >>= f[1]; /* Mark dependence */ trace_off(); /* End of active section */ ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Trace RHS again to mark on tape 2 the dependence of f upon the parameter mu. This tape is used in generating JacobianP. */ static PetscErrorCode RHSFunctionActiveP(TS ts,PetscReal t,Vec X,Vec F,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; PetscScalar *f; const PetscScalar *x; adouble f_a[2]; /* 'active' double for dependent variables */ adouble x_a[2],mu_a; /* 'active' double for independent variables */ PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); ierr = VecGetArray(F,&f);CHKERRQ(ierr); /* Start of active section */ trace_on(3); x_a[0] <<= x[0];x_a[1] <<= x[1];mu_a <<= user->mu; /* Mark independence */ f_a[0] = x_a[1]; f_a[1] = mu_a*(1.-x_a[0]*x_a[0])*x_a[1]-x_a[0]; f_a[0] >>= f[0];f_a[1] >>= f[1]; /* Mark dependence */ trace_off(); /* End of active section */ ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Compute the Jacobian w.r.t. x using PETSc-ADOL-C driver for explicit TS. */ static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; const PetscScalar *x; PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); ierr = PetscAdolcComputeRHSJacobian(1,A,x,user->adctx);CHKERRQ(ierr); ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Compute the Jacobian w.r.t. mu using PETSc-ADOL-C driver for explicit TS. */ static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx) { PetscErrorCode ierr; User user = (User)ctx; const PetscScalar *x; PetscFunctionBeginUser; ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); ierr = PetscAdolcComputeRHSJacobianP(3,A,x,&user->mu,user->adctx);CHKERRQ(ierr); ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); PetscFunctionReturn(0); } /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx) { PetscErrorCode ierr; const PetscScalar *x; PetscReal tfinal, dt, tprev; User user = (User)ctx; PetscFunctionBeginUser; ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr); ierr = TSGetPrevTime(ts,&tprev);CHKERRQ(ierr); ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 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); ierr = PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);CHKERRQ(ierr); ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); PetscFunctionReturn(0); } int main(int argc,char **argv) { TS ts; /* nonlinear solver */ Vec x; /* solution, residual vectors */ Mat A; /* Jacobian matrix */ Mat Jacp; /* JacobianP matrix */ PetscInt steps; PetscReal ftime = 0.5; PetscBool monitor = PETSC_FALSE; PetscScalar *x_ptr; PetscMPIInt size; struct _n_User user; AdolcCtx *adctx; PetscErrorCode ierr; Vec lambda[2],mu[2],r; /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Initialize program - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options and create AdolcCtx - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = PetscNew(&adctx);CHKERRQ(ierr); user.mu = 1; user.next_output = 0.0; adctx->m = 2;adctx->n = 2;adctx->p = 2;adctx->num_params = 1; user.adctx = adctx; ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr); ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create necessary matrix and vectors, solve same ODE on every process - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); ierr = MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); ierr = MatSetFromOptions(A);CHKERRQ(ierr); ierr = MatSetUp(A);CHKERRQ(ierr); ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr); ierr = MatCreate(PETSC_COMM_WORLD,&Jacp);CHKERRQ(ierr); ierr = MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);CHKERRQ(ierr); ierr = MatSetFromOptions(Jacp);CHKERRQ(ierr); ierr = MatSetUp(Jacp);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Create timestepping solver context - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); ierr = TSSetType(ts,TSRK);CHKERRQ(ierr); ierr = TSSetRHSFunction(ts,NULL,RHSFunctionPassive,&user);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set initial conditions - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); x_ptr[0] = 2; x_ptr[1] = 0.66666654321; ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Trace just once on each tape and put zeros on Jacobian diagonal - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = VecDuplicate(x,&r);CHKERRQ(ierr); ierr = RHSFunctionActive(ts,0.,x,r,&user);CHKERRQ(ierr); ierr = RHSFunctionActiveP(ts,0.,x,r,&user);CHKERRQ(ierr); ierr = VecSet(r,0);CHKERRQ(ierr); ierr = MatDiagonalSet(A,r,INSERT_VALUES);CHKERRQ(ierr); ierr = VecDestroy(&r);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set RHS Jacobian for the adjoint integration - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);CHKERRQ(ierr); ierr = TSSetMaxTime(ts,ftime);CHKERRQ(ierr); ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); if (monitor) { ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); } ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr); /* Have the TS save its trajectory so that TSAdjointSolve() may be used */ ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Set runtime options - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSSetFromOptions(ts);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Solve nonlinear system - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = TSSolve(ts,x);CHKERRQ(ierr); ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr); ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr); ierr = PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);CHKERRQ(ierr); ierr = VecView(x,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Start the Adjoint model - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatCreateVecs(A,&lambda[0],NULL);CHKERRQ(ierr); ierr = MatCreateVecs(A,&lambda[1],NULL);CHKERRQ(ierr); /* Reset initial conditions for the adjoint integration */ ierr = VecGetArray(lambda[0],&x_ptr);CHKERRQ(ierr); x_ptr[0] = 1.0; x_ptr[1] = 0.0; ierr = VecRestoreArray(lambda[0],&x_ptr);CHKERRQ(ierr); ierr = VecGetArray(lambda[1],&x_ptr);CHKERRQ(ierr); x_ptr[0] = 0.0; x_ptr[1] = 1.0; ierr = VecRestoreArray(lambda[1],&x_ptr);CHKERRQ(ierr); ierr = MatCreateVecs(Jacp,&mu[0],NULL);CHKERRQ(ierr); ierr = MatCreateVecs(Jacp,&mu[1],NULL);CHKERRQ(ierr); ierr = VecGetArray(mu[0],&x_ptr);CHKERRQ(ierr); x_ptr[0] = 0.0; ierr = VecRestoreArray(mu[0],&x_ptr);CHKERRQ(ierr); ierr = VecGetArray(mu[1],&x_ptr);CHKERRQ(ierr); x_ptr[0] = 0.0; ierr = VecRestoreArray(mu[1],&x_ptr);CHKERRQ(ierr); ierr = TSSetCostGradients(ts,2,lambda,mu);CHKERRQ(ierr); /* Set RHS JacobianP */ ierr = TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user);CHKERRQ(ierr); ierr = TSAdjointSolve(ts);CHKERRQ(ierr); ierr = VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); ierr = VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - Free work space. All PETSc objects should be destroyed when they are no longer needed. - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ ierr = MatDestroy(&A);CHKERRQ(ierr); ierr = MatDestroy(&Jacp);CHKERRQ(ierr); ierr = VecDestroy(&x);CHKERRQ(ierr); ierr = VecDestroy(&lambda[0]);CHKERRQ(ierr); ierr = VecDestroy(&lambda[1]);CHKERRQ(ierr); ierr = VecDestroy(&mu[0]);CHKERRQ(ierr); ierr = VecDestroy(&mu[1]);CHKERRQ(ierr); ierr = TSDestroy(&ts);CHKERRQ(ierr); ierr = PetscFree(adctx);CHKERRQ(ierr); ierr = PetscFinalize(); return ierr; } /*TEST build: requires: double !complex adolc test: suffix: 1 args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor output_file: output/ex16adj_1.out test: suffix: 2 args: -ts_max_steps 10 -ts_monitor -ts_adjoint_monitor -mu 5 output_file: output/ex16adj_2.out test: suffix: 3 args: -ts_max_steps 10 -monitor output_file: output/ex16adj_3.out test: suffix: 4 args: -ts_max_steps 10 -monitor -mu 5 output_file: output/ex16adj_4.out TEST*/