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