1*c4762a1bSJed Brown static char help[] = "Solves a ODE-constrained optimization problem -- finding the optimal initial conditions for the van der Pol equation.\n"; 2*c4762a1bSJed Brown 3*c4762a1bSJed Brown /** 4*c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 5*c4762a1bSJed Brown Concepts: TS^van der Pol equation DAE equivalent 6*c4762a1bSJed Brown Concepts: Optimization using adjoint sensitivities 7*c4762a1bSJed Brown Processors: 1 8*c4762a1bSJed Brown */ 9*c4762a1bSJed Brown /** 10*c4762a1bSJed Brown Notes: 11*c4762a1bSJed Brown This code demonstrates how to solve an ODE-constrained optimization problem with TAO, TSAdjoint and TS. 12*c4762a1bSJed Brown The nonlinear problem is written in an ODE equivalent form. 13*c4762a1bSJed Brown The objective is to minimize the difference between observation and model prediction by finding optimal values for initial conditions. 14*c4762a1bSJed Brown The gradient is computed with the discrete adjoint of an implicit method or an explicit method, see ex20adj.c for details. 15*c4762a1bSJed Brown */ 16*c4762a1bSJed Brown 17*c4762a1bSJed Brown #include <petsctao.h> 18*c4762a1bSJed Brown #include <petscts.h> 19*c4762a1bSJed Brown 20*c4762a1bSJed Brown typedef struct _n_User *User; 21*c4762a1bSJed Brown struct _n_User { 22*c4762a1bSJed Brown TS ts; 23*c4762a1bSJed Brown PetscReal mu; 24*c4762a1bSJed Brown PetscReal next_output; 25*c4762a1bSJed Brown 26*c4762a1bSJed Brown /* Sensitivity analysis support */ 27*c4762a1bSJed Brown PetscInt steps; 28*c4762a1bSJed Brown PetscReal ftime; 29*c4762a1bSJed Brown Mat A; /* Jacobian matrix for ODE */ 30*c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix for ODE*/ 31*c4762a1bSJed Brown Mat H; /* Hessian matrix for optimization */ 32*c4762a1bSJed Brown Vec U,Lambda[1],Mup[1]; /* first-order adjoint variables */ 33*c4762a1bSJed Brown Vec Lambda2[2]; /* second-order adjoint variables */ 34*c4762a1bSJed Brown Vec Ihp1[1]; /* working space for Hessian evaluations */ 35*c4762a1bSJed Brown Vec Dir; /* direction vector */ 36*c4762a1bSJed Brown PetscReal ob[2]; /* observation used by the cost function */ 37*c4762a1bSJed Brown PetscBool implicitform; /* implicit ODE? */ 38*c4762a1bSJed Brown }; 39*c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec,PetscScalar[],User); 40*c4762a1bSJed Brown 41*c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE -------------------- */ 42*c4762a1bSJed Brown 43*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 44*c4762a1bSJed Brown { 45*c4762a1bSJed Brown PetscErrorCode ierr; 46*c4762a1bSJed Brown User user = (User)ctx; 47*c4762a1bSJed Brown PetscScalar *f; 48*c4762a1bSJed Brown const PetscScalar *u; 49*c4762a1bSJed Brown 50*c4762a1bSJed Brown PetscFunctionBeginUser; 51*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 52*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 53*c4762a1bSJed Brown f[0] = u[1]; 54*c4762a1bSJed Brown f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]); 55*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 56*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 57*c4762a1bSJed Brown PetscFunctionReturn(0); 58*c4762a1bSJed Brown } 59*c4762a1bSJed Brown 60*c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 61*c4762a1bSJed Brown { 62*c4762a1bSJed Brown PetscErrorCode ierr; 63*c4762a1bSJed Brown User user = (User)ctx; 64*c4762a1bSJed Brown PetscReal mu = user->mu; 65*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 66*c4762a1bSJed Brown PetscScalar J[2][2]; 67*c4762a1bSJed Brown const PetscScalar *u; 68*c4762a1bSJed Brown 69*c4762a1bSJed Brown PetscFunctionBeginUser; 70*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 71*c4762a1bSJed Brown J[0][0] = 0; 72*c4762a1bSJed Brown J[1][0] = -mu*(2.0*u[1]*u[0]+1.); 73*c4762a1bSJed Brown J[0][1] = 1.0; 74*c4762a1bSJed Brown J[1][1] = mu*(1.0-u[0]*u[0]); 75*c4762a1bSJed Brown ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 77*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 78*c4762a1bSJed Brown if (A != B) { 79*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 80*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 81*c4762a1bSJed Brown } 82*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 83*c4762a1bSJed Brown PetscFunctionReturn(0); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 87*c4762a1bSJed Brown { 88*c4762a1bSJed Brown const PetscScalar *vl,*vr,*u; 89*c4762a1bSJed Brown PetscScalar *vhv; 90*c4762a1bSJed Brown PetscScalar dJdU[2][2][2]={{{0}}}; 91*c4762a1bSJed Brown PetscInt i,j,k; 92*c4762a1bSJed Brown User user = (User)ctx; 93*c4762a1bSJed Brown PetscErrorCode ierr; 94*c4762a1bSJed Brown 95*c4762a1bSJed Brown PetscFunctionBeginUser; 96*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 97*c4762a1bSJed Brown ierr = VecGetArrayRead(Vl[0],&vl);CHKERRQ(ierr); 98*c4762a1bSJed Brown ierr = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr); 99*c4762a1bSJed Brown ierr = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr); 100*c4762a1bSJed Brown 101*c4762a1bSJed Brown dJdU[1][0][0] = -2.*user->mu*u[1]; 102*c4762a1bSJed Brown dJdU[1][1][0] = -2.*user->mu*u[0]; 103*c4762a1bSJed Brown dJdU[1][0][1] = -2.*user->mu*u[0]; 104*c4762a1bSJed Brown for (j=0;j<2;j++) { 105*c4762a1bSJed Brown vhv[j] = 0; 106*c4762a1bSJed Brown for (k=0;k<2;k++) 107*c4762a1bSJed Brown for (i=0;i<2;i++ ) 108*c4762a1bSJed Brown vhv[j] += vl[i]*dJdU[i][j][k]*vr[k]; 109*c4762a1bSJed Brown } 110*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Vl[0],&vl);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr); 114*c4762a1bSJed Brown PetscFunctionReturn(0); 115*c4762a1bSJed Brown } 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE -------------------- */ 118*c4762a1bSJed Brown 119*c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 120*c4762a1bSJed Brown { 121*c4762a1bSJed Brown PetscErrorCode ierr; 122*c4762a1bSJed Brown User user = (User)ctx; 123*c4762a1bSJed Brown const PetscScalar *u,*udot; 124*c4762a1bSJed Brown PetscScalar *f; 125*c4762a1bSJed Brown 126*c4762a1bSJed Brown PetscFunctionBeginUser; 127*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 128*c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 129*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 130*c4762a1bSJed Brown f[0] = udot[0] - u[1]; 131*c4762a1bSJed Brown f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]) ; 132*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 134*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 135*c4762a1bSJed Brown PetscFunctionReturn(0); 136*c4762a1bSJed Brown } 137*c4762a1bSJed Brown 138*c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 139*c4762a1bSJed Brown { 140*c4762a1bSJed Brown PetscErrorCode ierr; 141*c4762a1bSJed Brown User user = (User)ctx; 142*c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 143*c4762a1bSJed Brown PetscScalar J[2][2]; 144*c4762a1bSJed Brown const PetscScalar *u; 145*c4762a1bSJed Brown 146*c4762a1bSJed Brown PetscFunctionBeginUser; 147*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 148*c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 149*c4762a1bSJed Brown J[1][0] = user->mu*(1.0 + 2.0*u[0]*u[1]); J[1][1] = a - user->mu*(1.0-u[0]*u[0]); 150*c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 154*c4762a1bSJed Brown if (A != B) { 155*c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 156*c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 157*c4762a1bSJed Brown } 158*c4762a1bSJed Brown PetscFunctionReturn(0); 159*c4762a1bSJed Brown } 160*c4762a1bSJed Brown 161*c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 162*c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx) 163*c4762a1bSJed Brown { 164*c4762a1bSJed Brown PetscErrorCode ierr; 165*c4762a1bSJed Brown const PetscScalar *u; 166*c4762a1bSJed Brown PetscReal tfinal, dt; 167*c4762a1bSJed Brown User user = (User)ctx; 168*c4762a1bSJed Brown Vec interpolatedU; 169*c4762a1bSJed Brown 170*c4762a1bSJed Brown PetscFunctionBeginUser; 171*c4762a1bSJed Brown ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 172*c4762a1bSJed Brown ierr = TSGetMaxTime(ts,&tfinal);CHKERRQ(ierr); 173*c4762a1bSJed Brown 174*c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 175*c4762a1bSJed Brown ierr = VecDuplicate(U,&interpolatedU);CHKERRQ(ierr); 176*c4762a1bSJed Brown ierr = TSInterpolate(ts,user->next_output,interpolatedU);CHKERRQ(ierr); 177*c4762a1bSJed Brown ierr = VecGetArrayRead(interpolatedU,&u);CHKERRQ(ierr); 178*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n", 179*c4762a1bSJed Brown (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(u[0]), 180*c4762a1bSJed Brown (double)PetscRealPart(u[1]));CHKERRQ(ierr); 181*c4762a1bSJed Brown ierr = VecRestoreArrayRead(interpolatedU,&u);CHKERRQ(ierr); 182*c4762a1bSJed Brown ierr = VecDestroy(&interpolatedU);CHKERRQ(ierr); 183*c4762a1bSJed Brown user->next_output += 0.1; 184*c4762a1bSJed Brown } 185*c4762a1bSJed Brown PetscFunctionReturn(0); 186*c4762a1bSJed Brown } 187*c4762a1bSJed Brown 188*c4762a1bSJed Brown static PetscErrorCode IHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 189*c4762a1bSJed Brown { 190*c4762a1bSJed Brown const PetscScalar *vl,*vr,*u; 191*c4762a1bSJed Brown PetscScalar *vhv; 192*c4762a1bSJed Brown PetscScalar dJdU[2][2][2]={{{0}}}; 193*c4762a1bSJed Brown PetscInt i,j,k; 194*c4762a1bSJed Brown User user = (User)ctx; 195*c4762a1bSJed Brown PetscErrorCode ierr; 196*c4762a1bSJed Brown 197*c4762a1bSJed Brown PetscFunctionBeginUser; 198*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 199*c4762a1bSJed Brown ierr = VecGetArrayRead(Vl[0],&vl);CHKERRQ(ierr); 200*c4762a1bSJed Brown ierr = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr); 201*c4762a1bSJed Brown ierr = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr); 202*c4762a1bSJed Brown dJdU[1][0][0] = 2.*user->mu*u[1]; 203*c4762a1bSJed Brown dJdU[1][1][0] = 2.*user->mu*u[0]; 204*c4762a1bSJed Brown dJdU[1][0][1] = 2.*user->mu*u[0]; 205*c4762a1bSJed Brown for (j=0;j<2;j++) { 206*c4762a1bSJed Brown vhv[j] = 0; 207*c4762a1bSJed Brown for (k=0;k<2;k++) 208*c4762a1bSJed Brown for (i=0;i<2;i++ ) 209*c4762a1bSJed Brown vhv[j] += vl[i]*dJdU[i][j][k]*vr[k]; 210*c4762a1bSJed Brown } 211*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 212*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Vl[0],&vl);CHKERRQ(ierr); 213*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr); 214*c4762a1bSJed Brown ierr = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr); 215*c4762a1bSJed Brown PetscFunctionReturn(0); 216*c4762a1bSJed Brown } 217*c4762a1bSJed Brown 218*c4762a1bSJed Brown /* ------------------ User-defined routines for TAO -------------------------- */ 219*c4762a1bSJed Brown 220*c4762a1bSJed Brown static PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx) 221*c4762a1bSJed Brown { 222*c4762a1bSJed Brown User user_ptr = (User)ctx; 223*c4762a1bSJed Brown TS ts = user_ptr->ts; 224*c4762a1bSJed Brown const PetscScalar *x_ptr; 225*c4762a1bSJed Brown PetscScalar *y_ptr; 226*c4762a1bSJed Brown PetscErrorCode ierr; 227*c4762a1bSJed Brown 228*c4762a1bSJed Brown PetscFunctionBeginUser; 229*c4762a1bSJed Brown ierr = VecCopy(IC,user_ptr->U);CHKERRQ(ierr); /* set up the initial condition */ 230*c4762a1bSJed Brown 231*c4762a1bSJed Brown ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 232*c4762a1bSJed Brown ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 233*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr); /* can be overwritten by command line options */ 234*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = TSSolve(ts,user_ptr->U);CHKERRQ(ierr); 236*c4762a1bSJed Brown 237*c4762a1bSJed Brown ierr = VecGetArrayRead(user_ptr->U,&x_ptr);CHKERRQ(ierr); 238*c4762a1bSJed Brown ierr = VecGetArray(user_ptr->Lambda[0],&y_ptr);CHKERRQ(ierr); 239*c4762a1bSJed Brown *f = (x_ptr[0]-user_ptr->ob[0])*(x_ptr[0]-user_ptr->ob[0])+(x_ptr[1]-user_ptr->ob[1])*(x_ptr[1]-user_ptr->ob[1]); 240*c4762a1bSJed Brown y_ptr[0] = 2.*(x_ptr[0]-user_ptr->ob[0]); 241*c4762a1bSJed Brown y_ptr[1] = 2.*(x_ptr[1]-user_ptr->ob[1]); 242*c4762a1bSJed Brown ierr = VecRestoreArray(user_ptr->Lambda[0],&y_ptr);CHKERRQ(ierr); 243*c4762a1bSJed Brown ierr = VecRestoreArrayRead(user_ptr->U,&x_ptr);CHKERRQ(ierr); 244*c4762a1bSJed Brown 245*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,user_ptr->Lambda,NULL);CHKERRQ(ierr); 246*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 247*c4762a1bSJed Brown ierr = VecCopy(user_ptr->Lambda[0],G);CHKERRQ(ierr); 248*c4762a1bSJed Brown PetscFunctionReturn(0); 249*c4762a1bSJed Brown } 250*c4762a1bSJed Brown 251*c4762a1bSJed Brown static PetscErrorCode FormHessian(Tao tao,Vec U,Mat H,Mat Hpre,void *ctx) 252*c4762a1bSJed Brown { 253*c4762a1bSJed Brown User user_ptr = (User)ctx; 254*c4762a1bSJed Brown PetscScalar harr[2]; 255*c4762a1bSJed Brown PetscScalar *x_ptr; 256*c4762a1bSJed Brown const PetscInt rows[2] = {0,1}; 257*c4762a1bSJed Brown PetscInt col; 258*c4762a1bSJed Brown PetscErrorCode ierr; 259*c4762a1bSJed Brown 260*c4762a1bSJed Brown PetscFunctionBeginUser; 261*c4762a1bSJed Brown ierr = VecCopy(U,user_ptr->U);CHKERRQ(ierr); 262*c4762a1bSJed Brown ierr = VecGetArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr); 263*c4762a1bSJed Brown x_ptr[0] = 1.; 264*c4762a1bSJed Brown x_ptr[1] = 0.; 265*c4762a1bSJed Brown ierr = VecRestoreArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = Adjoint2(user_ptr->U,harr,user_ptr);CHKERRQ(ierr); 267*c4762a1bSJed Brown col = 0; 268*c4762a1bSJed Brown ierr = MatSetValues(H,2,rows,1,&col,harr,INSERT_VALUES);CHKERRQ(ierr); 269*c4762a1bSJed Brown 270*c4762a1bSJed Brown ierr = VecCopy(U,user_ptr->U);CHKERRQ(ierr); 271*c4762a1bSJed Brown ierr = VecGetArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr); 272*c4762a1bSJed Brown x_ptr[0] = 0.; 273*c4762a1bSJed Brown x_ptr[1] = 1.; 274*c4762a1bSJed Brown ierr = VecRestoreArray(user_ptr->Dir,&x_ptr);CHKERRQ(ierr); 275*c4762a1bSJed Brown ierr = Adjoint2(user_ptr->U,harr,user_ptr);CHKERRQ(ierr); 276*c4762a1bSJed Brown col = 1; 277*c4762a1bSJed Brown ierr = MatSetValues(H,2,rows,1,&col,harr,INSERT_VALUES);CHKERRQ(ierr); 278*c4762a1bSJed Brown 279*c4762a1bSJed Brown ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 281*c4762a1bSJed Brown if (H != Hpre) { 282*c4762a1bSJed Brown ierr = MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 283*c4762a1bSJed Brown ierr = MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 284*c4762a1bSJed Brown } 285*c4762a1bSJed Brown PetscFunctionReturn(0); 286*c4762a1bSJed Brown } 287*c4762a1bSJed Brown 288*c4762a1bSJed Brown static PetscErrorCode MatrixFreeHessian(Tao tao,Vec U,Mat H,Mat Hpre,void *ctx) 289*c4762a1bSJed Brown { 290*c4762a1bSJed Brown User user_ptr = (User)ctx; 291*c4762a1bSJed Brown PetscErrorCode ierr; 292*c4762a1bSJed Brown 293*c4762a1bSJed Brown PetscFunctionBeginUser; 294*c4762a1bSJed Brown ierr = VecCopy(U,user_ptr->U);CHKERRQ(ierr); 295*c4762a1bSJed Brown PetscFunctionReturn(0); 296*c4762a1bSJed Brown } 297*c4762a1bSJed Brown 298*c4762a1bSJed Brown /* ------------ Routines calculating second-order derivatives -------------- */ 299*c4762a1bSJed Brown 300*c4762a1bSJed Brown /** 301*c4762a1bSJed Brown Compute the Hessian-vector product for the cost function using Second-order adjoint 302*c4762a1bSJed Brown */ 303*c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec U,PetscScalar arr[],User ctx) 304*c4762a1bSJed Brown { 305*c4762a1bSJed Brown TS ts = ctx->ts; 306*c4762a1bSJed Brown PetscScalar *x_ptr,*y_ptr; 307*c4762a1bSJed Brown Mat tlmsen; 308*c4762a1bSJed Brown PetscErrorCode ierr; 309*c4762a1bSJed Brown 310*c4762a1bSJed Brown PetscFunctionBeginUser; 311*c4762a1bSJed Brown ierr = TSAdjointReset(ts);CHKERRQ(ierr); 312*c4762a1bSJed Brown 313*c4762a1bSJed Brown ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 314*c4762a1bSJed Brown ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 315*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr); 316*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 317*c4762a1bSJed Brown ierr = TSSetCostHessianProducts(ts,1,ctx->Lambda2,NULL,ctx->Dir);CHKERRQ(ierr); 318*c4762a1bSJed Brown ierr = TSAdjointSetForward(ts,NULL);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 320*c4762a1bSJed Brown 321*c4762a1bSJed Brown /* Set terminal conditions for first- and second-order adjonts */ 322*c4762a1bSJed Brown ierr = VecGetArray(U,&x_ptr);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = VecGetArray(ctx->Lambda[0],&y_ptr);CHKERRQ(ierr); 324*c4762a1bSJed Brown y_ptr[0] = 2.*(x_ptr[0]-ctx->ob[0]); 325*c4762a1bSJed Brown y_ptr[1] = 2.*(x_ptr[1]-ctx->ob[1]); 326*c4762a1bSJed Brown ierr = VecRestoreArray(ctx->Lambda[0],&y_ptr);CHKERRQ(ierr); 327*c4762a1bSJed Brown ierr = VecRestoreArray(U,&x_ptr);CHKERRQ(ierr); 328*c4762a1bSJed Brown ierr = TSForwardGetSensitivities(ts,NULL,&tlmsen);CHKERRQ(ierr); 329*c4762a1bSJed Brown ierr = MatDenseGetColumn(tlmsen,0,&x_ptr);CHKERRQ(ierr); 330*c4762a1bSJed Brown ierr = VecGetArray(ctx->Lambda2[0],&y_ptr);CHKERRQ(ierr); 331*c4762a1bSJed Brown y_ptr[0] = 2.*x_ptr[0]; 332*c4762a1bSJed Brown y_ptr[1] = 2.*x_ptr[1]; 333*c4762a1bSJed Brown ierr = VecRestoreArray(ctx->Lambda2[0],&y_ptr);CHKERRQ(ierr); 334*c4762a1bSJed Brown ierr = MatDenseRestoreColumn(tlmsen,&x_ptr);CHKERRQ(ierr); 335*c4762a1bSJed Brown 336*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,ctx->Lambda,NULL);CHKERRQ(ierr); 337*c4762a1bSJed Brown if (ctx->implicitform) { 338*c4762a1bSJed Brown ierr = TSSetIHessianProduct(ts,ctx->Ihp1,IHessianProductUU,NULL,NULL,NULL,NULL,NULL,NULL,ctx);CHKERRQ(ierr); 339*c4762a1bSJed Brown } else { 340*c4762a1bSJed Brown ierr = TSSetRHSHessianProduct(ts,ctx->Ihp1,RHSHessianProductUU,NULL,NULL,NULL,NULL,NULL,NULL,ctx);CHKERRQ(ierr); 341*c4762a1bSJed Brown } 342*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 343*c4762a1bSJed Brown 344*c4762a1bSJed Brown ierr = VecGetArray(ctx->Lambda2[0],&x_ptr);CHKERRQ(ierr); 345*c4762a1bSJed Brown arr[0] = x_ptr[0]; 346*c4762a1bSJed Brown arr[1] = x_ptr[1]; 347*c4762a1bSJed Brown ierr = VecRestoreArray(ctx->Lambda2[0],&x_ptr);CHKERRQ(ierr); 348*c4762a1bSJed Brown 349*c4762a1bSJed Brown ierr = TSAdjointReset(ts);CHKERRQ(ierr); 350*c4762a1bSJed Brown ierr = TSAdjointResetForward(ts);CHKERRQ(ierr); 351*c4762a1bSJed Brown PetscFunctionReturn(0); 352*c4762a1bSJed Brown } 353*c4762a1bSJed Brown 354*c4762a1bSJed Brown PetscErrorCode FiniteDiff(Vec U,PetscScalar arr[],User ctx) 355*c4762a1bSJed Brown { 356*c4762a1bSJed Brown Vec Up,G,Gp; 357*c4762a1bSJed Brown const PetscScalar eps = PetscRealConstant(1e-7); 358*c4762a1bSJed Brown PetscScalar *u; 359*c4762a1bSJed Brown Tao tao = NULL; 360*c4762a1bSJed Brown PetscReal f; 361*c4762a1bSJed Brown PetscErrorCode ierr; 362*c4762a1bSJed Brown 363*c4762a1bSJed Brown PetscFunctionBeginUser; 364*c4762a1bSJed Brown ierr = VecDuplicate(U,&Up);CHKERRQ(ierr); 365*c4762a1bSJed Brown ierr = VecDuplicate(U,&G);CHKERRQ(ierr); 366*c4762a1bSJed Brown ierr = VecDuplicate(U,&Gp);CHKERRQ(ierr); 367*c4762a1bSJed Brown 368*c4762a1bSJed Brown ierr = FormFunctionGradient(tao,U,&f,G,ctx);CHKERRQ(ierr); 369*c4762a1bSJed Brown 370*c4762a1bSJed Brown ierr = VecCopy(U,Up);CHKERRQ(ierr); 371*c4762a1bSJed Brown ierr = VecGetArray(Up,&u);CHKERRQ(ierr); 372*c4762a1bSJed Brown u[0] += eps; 373*c4762a1bSJed Brown ierr = VecRestoreArray(Up,&u);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = FormFunctionGradient(tao,Up,&f,Gp,ctx);CHKERRQ(ierr); 375*c4762a1bSJed Brown ierr = VecAXPY(Gp,-1,G);CHKERRQ(ierr); 376*c4762a1bSJed Brown ierr = VecScale(Gp,1./eps);CHKERRQ(ierr); 377*c4762a1bSJed Brown ierr = VecGetArray(Gp,&u);CHKERRQ(ierr); 378*c4762a1bSJed Brown arr[0] = u[0]; 379*c4762a1bSJed Brown arr[1] = u[1]; 380*c4762a1bSJed Brown ierr = VecRestoreArray(Gp,&u);CHKERRQ(ierr); 381*c4762a1bSJed Brown 382*c4762a1bSJed Brown ierr = VecCopy(U,Up);CHKERRQ(ierr); 383*c4762a1bSJed Brown ierr = VecGetArray(Up,&u);CHKERRQ(ierr); 384*c4762a1bSJed Brown u[1] += eps; 385*c4762a1bSJed Brown ierr = VecRestoreArray(Up,&u);CHKERRQ(ierr); 386*c4762a1bSJed Brown ierr = FormFunctionGradient(tao,Up,&f,Gp,ctx);CHKERRQ(ierr); 387*c4762a1bSJed Brown ierr = VecAXPY(Gp,-1,G);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = VecScale(Gp,1./eps);CHKERRQ(ierr); 389*c4762a1bSJed Brown ierr = VecGetArray(Gp,&u);CHKERRQ(ierr); 390*c4762a1bSJed Brown arr[2] = u[0]; 391*c4762a1bSJed Brown arr[3] = u[1]; 392*c4762a1bSJed Brown ierr = VecRestoreArray(Gp,&u);CHKERRQ(ierr); 393*c4762a1bSJed Brown 394*c4762a1bSJed Brown ierr = VecDestroy(&G);CHKERRQ(ierr); 395*c4762a1bSJed Brown ierr = VecDestroy(&Gp);CHKERRQ(ierr); 396*c4762a1bSJed Brown ierr = VecDestroy(&Up);CHKERRQ(ierr); 397*c4762a1bSJed Brown PetscFunctionReturn(0); 398*c4762a1bSJed Brown } 399*c4762a1bSJed Brown 400*c4762a1bSJed Brown static PetscErrorCode HessianProductMat(Mat mat,Vec svec,Vec y) 401*c4762a1bSJed Brown { 402*c4762a1bSJed Brown User user_ptr; 403*c4762a1bSJed Brown PetscScalar *y_ptr; 404*c4762a1bSJed Brown PetscErrorCode ierr; 405*c4762a1bSJed Brown 406*c4762a1bSJed Brown PetscFunctionBeginUser; 407*c4762a1bSJed Brown ierr = MatShellGetContext(mat,(void*)&user_ptr);CHKERRQ(ierr); 408*c4762a1bSJed Brown ierr = VecCopy(svec,user_ptr->Dir);CHKERRQ(ierr); 409*c4762a1bSJed Brown ierr = VecGetArray(y,&y_ptr);CHKERRQ(ierr); 410*c4762a1bSJed Brown ierr = Adjoint2(user_ptr->U,y_ptr,user_ptr);CHKERRQ(ierr); 411*c4762a1bSJed Brown ierr = VecRestoreArray(y,&y_ptr);CHKERRQ(ierr); 412*c4762a1bSJed Brown PetscFunctionReturn(0); 413*c4762a1bSJed Brown } 414*c4762a1bSJed Brown 415*c4762a1bSJed Brown int main(int argc,char **argv) 416*c4762a1bSJed Brown { 417*c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE,mf = PETSC_TRUE; 418*c4762a1bSJed Brown PetscInt mode = 0; 419*c4762a1bSJed Brown PetscMPIInt size; 420*c4762a1bSJed Brown struct _n_User user; 421*c4762a1bSJed Brown Vec x; /* working vector for TAO */ 422*c4762a1bSJed Brown PetscScalar *x_ptr,arr[4]; 423*c4762a1bSJed Brown PetscScalar ic1 = 2.2,ic2 = -0.7; /* initial guess for TAO */ 424*c4762a1bSJed Brown Tao tao; 425*c4762a1bSJed Brown KSP ksp; 426*c4762a1bSJed Brown PC pc; 427*c4762a1bSJed Brown PetscErrorCode ierr; 428*c4762a1bSJed Brown 429*c4762a1bSJed Brown /* Initialize program */ 430*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 431*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 432*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 433*c4762a1bSJed Brown 434*c4762a1bSJed Brown /* Set runtime options */ 435*c4762a1bSJed Brown user.next_output = 0.0; 436*c4762a1bSJed Brown user.mu = 1.0e3; 437*c4762a1bSJed Brown user.steps = 0; 438*c4762a1bSJed Brown user.ftime = 0.5; 439*c4762a1bSJed Brown user.implicitform = PETSC_TRUE; 440*c4762a1bSJed Brown 441*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 442*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);CHKERRQ(ierr); 443*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-mode",&mode,NULL);CHKERRQ(ierr); 444*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-ic1",&ic1,NULL);CHKERRQ(ierr); 445*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-ic2",&ic2,NULL);CHKERRQ(ierr); 446*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-my_tao_mf",&mf,NULL);CHKERRQ(ierr); /* matrix-free hessian for optimization */ 447*c4762a1bSJed Brown ierr = PetscOptionsGetBool(NULL,NULL,"-implicitform",&user.implicitform,NULL);CHKERRQ(ierr); 448*c4762a1bSJed Brown 449*c4762a1bSJed Brown /* Create necessary matrix and vectors, solve same ODE on every process */ 450*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user.A);CHKERRQ(ierr); 451*c4762a1bSJed Brown ierr = MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 452*c4762a1bSJed Brown ierr = MatSetFromOptions(user.A);CHKERRQ(ierr); 453*c4762a1bSJed Brown ierr = MatSetUp(user.A);CHKERRQ(ierr); 454*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.U,NULL);CHKERRQ(ierr); 455*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.Dir,NULL);CHKERRQ(ierr); 456*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.Lambda[0],NULL);CHKERRQ(ierr); 457*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.Lambda2[0],NULL);CHKERRQ(ierr); 458*c4762a1bSJed Brown ierr = MatCreateVecs(user.A,&user.Ihp1[0],NULL);CHKERRQ(ierr); 459*c4762a1bSJed Brown 460*c4762a1bSJed Brown /* Create timestepping solver context */ 461*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&user.ts);CHKERRQ(ierr); 462*c4762a1bSJed Brown ierr = TSSetEquationType(user.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 463*c4762a1bSJed Brown if (user.implicitform) { 464*c4762a1bSJed Brown ierr = TSSetIFunction(user.ts,NULL,IFunction,&user);CHKERRQ(ierr); 465*c4762a1bSJed Brown ierr = TSSetIJacobian(user.ts,user.A,user.A,IJacobian,&user);CHKERRQ(ierr); 466*c4762a1bSJed Brown ierr = TSSetType(user.ts,TSCN);CHKERRQ(ierr); 467*c4762a1bSJed Brown } else { 468*c4762a1bSJed Brown ierr = TSSetRHSFunction(user.ts,NULL,RHSFunction,&user);CHKERRQ(ierr); 469*c4762a1bSJed Brown ierr = TSSetRHSJacobian(user.ts,user.A,user.A,RHSJacobian,&user);CHKERRQ(ierr); 470*c4762a1bSJed Brown ierr = TSSetType(user.ts,TSRK);CHKERRQ(ierr); 471*c4762a1bSJed Brown } 472*c4762a1bSJed Brown ierr = TSSetMaxTime(user.ts,user.ftime);CHKERRQ(ierr); 473*c4762a1bSJed Brown ierr = TSSetExactFinalTime(user.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 474*c4762a1bSJed Brown 475*c4762a1bSJed Brown if (monitor) { 476*c4762a1bSJed Brown ierr = TSMonitorSet(user.ts,Monitor,&user,NULL);CHKERRQ(ierr); 477*c4762a1bSJed Brown } 478*c4762a1bSJed Brown 479*c4762a1bSJed Brown /* Set ODE initial conditions */ 480*c4762a1bSJed Brown ierr = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr); 481*c4762a1bSJed Brown x_ptr[0] = 2.0; 482*c4762a1bSJed Brown x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu); 483*c4762a1bSJed Brown ierr = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr); 484*c4762a1bSJed Brown 485*c4762a1bSJed Brown /* Set runtime options */ 486*c4762a1bSJed Brown ierr = TSSetFromOptions(user.ts);CHKERRQ(ierr); 487*c4762a1bSJed Brown 488*c4762a1bSJed Brown /* Obtain the observation */ 489*c4762a1bSJed Brown ierr = TSSolve(user.ts,user.U);CHKERRQ(ierr); 490*c4762a1bSJed Brown ierr = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr); 491*c4762a1bSJed Brown user.ob[0] = x_ptr[0]; 492*c4762a1bSJed Brown user.ob[1] = x_ptr[1]; 493*c4762a1bSJed Brown ierr = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr); 494*c4762a1bSJed Brown 495*c4762a1bSJed Brown ierr = VecDuplicate(user.U,&x);CHKERRQ(ierr); 496*c4762a1bSJed Brown ierr = VecGetArray(x,&x_ptr);CHKERRQ(ierr); 497*c4762a1bSJed Brown x_ptr[0] = ic1; 498*c4762a1bSJed Brown x_ptr[1] = ic2; 499*c4762a1bSJed Brown ierr = VecRestoreArray(x,&x_ptr);CHKERRQ(ierr); 500*c4762a1bSJed Brown 501*c4762a1bSJed Brown /* Save trajectory of solution so that TSAdjointSolve() may be used */ 502*c4762a1bSJed Brown ierr = TSSetSaveTrajectory(user.ts);CHKERRQ(ierr); 503*c4762a1bSJed Brown 504*c4762a1bSJed Brown /* Compare finite difference and second-order adjoint. */ 505*c4762a1bSJed Brown switch (mode) { 506*c4762a1bSJed Brown case 2 : 507*c4762a1bSJed Brown ierr = FiniteDiff(x,arr,&user);CHKERRQ(ierr); 508*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Finite difference approximation of the Hessian\n");CHKERRQ(ierr); 509*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%g %g\n%g %g\n",(double)arr[0],(double)arr[1],(double)arr[2],(double)arr[3]);CHKERRQ(ierr); 510*c4762a1bSJed Brown break; 511*c4762a1bSJed Brown case 1 : /* Compute the Hessian column by column */ 512*c4762a1bSJed Brown ierr = VecCopy(x,user.U);CHKERRQ(ierr); 513*c4762a1bSJed Brown ierr = VecGetArray(user.Dir,&x_ptr);CHKERRQ(ierr); 514*c4762a1bSJed Brown x_ptr[0] = 1.; 515*c4762a1bSJed Brown x_ptr[1] = 0.; 516*c4762a1bSJed Brown ierr = VecRestoreArray(user.Dir,&x_ptr);CHKERRQ(ierr); 517*c4762a1bSJed Brown ierr = Adjoint2(user.U,arr,&user);CHKERRQ(ierr); 518*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nFirst column of the Hessian\n");CHKERRQ(ierr); 519*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%g\n%g\n",(double)arr[0],(double)arr[1]);CHKERRQ(ierr); 520*c4762a1bSJed Brown ierr = VecCopy(x,user.U);CHKERRQ(ierr); 521*c4762a1bSJed Brown ierr = VecGetArray(user.Dir,&x_ptr);CHKERRQ(ierr); 522*c4762a1bSJed Brown x_ptr[0] = 0.; 523*c4762a1bSJed Brown x_ptr[1] = 1.; 524*c4762a1bSJed Brown ierr = VecRestoreArray(user.Dir,&x_ptr);CHKERRQ(ierr); 525*c4762a1bSJed Brown ierr = Adjoint2(user.U,arr,&user);CHKERRQ(ierr); 526*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"\nSecond column of the Hessian\n");CHKERRQ(ierr); 527*c4762a1bSJed Brown ierr = PetscPrintf(PETSC_COMM_WORLD,"%g\n%g\n",(double)arr[0],(double)arr[1]);CHKERRQ(ierr); 528*c4762a1bSJed Brown break; 529*c4762a1bSJed Brown case 0 : 530*c4762a1bSJed Brown /* Create optimization context and set up */ 531*c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 532*c4762a1bSJed Brown ierr = TaoSetType(tao,TAOBLMVM);CHKERRQ(ierr); 533*c4762a1bSJed Brown ierr = TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);CHKERRQ(ierr); 534*c4762a1bSJed Brown 535*c4762a1bSJed Brown if (mf) { 536*c4762a1bSJed Brown ierr = MatCreateShell(PETSC_COMM_SELF,2,2,2,2,(void*)&user,&user.H);CHKERRQ(ierr); 537*c4762a1bSJed Brown ierr = MatShellSetOperation(user.H,MATOP_MULT,(void(*)(void))HessianProductMat);CHKERRQ(ierr); 538*c4762a1bSJed Brown ierr = MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 539*c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao,user.H,user.H,MatrixFreeHessian,(void *)&user);CHKERRQ(ierr); 540*c4762a1bSJed Brown } else { /* Create Hessian matrix */ 541*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&user.H);CHKERRQ(ierr); 542*c4762a1bSJed Brown ierr = MatSetSizes(user.H,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 543*c4762a1bSJed Brown ierr = MatSetUp(user.H);CHKERRQ(ierr); 544*c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao,user.H,user.H,FormHessian,(void *)&user);CHKERRQ(ierr); 545*c4762a1bSJed Brown } 546*c4762a1bSJed Brown 547*c4762a1bSJed Brown /* Not use any preconditioner */ 548*c4762a1bSJed Brown ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 549*c4762a1bSJed Brown if (ksp) { 550*c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 551*c4762a1bSJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 552*c4762a1bSJed Brown } 553*c4762a1bSJed Brown 554*c4762a1bSJed Brown /* Set initial solution guess */ 555*c4762a1bSJed Brown ierr = TaoSetInitialVector(tao,x);CHKERRQ(ierr); 556*c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 557*c4762a1bSJed Brown ierr = TaoSolve(tao); CHKERRQ(ierr); 558*c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 559*c4762a1bSJed Brown ierr = MatDestroy(&user.H);CHKERRQ(ierr); 560*c4762a1bSJed Brown break; 561*c4762a1bSJed Brown default: 562*c4762a1bSJed Brown break; 563*c4762a1bSJed Brown } 564*c4762a1bSJed Brown 565*c4762a1bSJed Brown /* Free work space. All PETSc objects should be destroyed when they are no longer needed. */ 566*c4762a1bSJed Brown ierr = MatDestroy(&user.A);CHKERRQ(ierr); 567*c4762a1bSJed Brown ierr = VecDestroy(&user.U);CHKERRQ(ierr); 568*c4762a1bSJed Brown ierr = VecDestroy(&user.Lambda[0]);CHKERRQ(ierr); 569*c4762a1bSJed Brown ierr = VecDestroy(&user.Lambda2[0]);CHKERRQ(ierr); 570*c4762a1bSJed Brown ierr = VecDestroy(&user.Ihp1[0]);CHKERRQ(ierr); 571*c4762a1bSJed Brown ierr = VecDestroy(&user.Dir);CHKERRQ(ierr); 572*c4762a1bSJed Brown ierr = TSDestroy(&user.ts);CHKERRQ(ierr); 573*c4762a1bSJed Brown ierr = VecDestroy(&x);CHKERRQ(ierr); 574*c4762a1bSJed Brown ierr = PetscFinalize(); 575*c4762a1bSJed Brown return(ierr); 576*c4762a1bSJed Brown } 577*c4762a1bSJed Brown 578*c4762a1bSJed Brown /*TEST 579*c4762a1bSJed Brown build: 580*c4762a1bSJed Brown requires: !complex !single 581*c4762a1bSJed Brown 582*c4762a1bSJed Brown test: 583*c4762a1bSJed Brown args: -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 1000 -ts_dt 0.03125 584*c4762a1bSJed Brown output_file: output/ex20opt_ic_1.out 585*c4762a1bSJed Brown 586*c4762a1bSJed Brown test: 587*c4762a1bSJed Brown suffix: 2 588*c4762a1bSJed Brown args: -ts_type beuler -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none 589*c4762a1bSJed Brown output_file: output/ex20opt_ic_2.out 590*c4762a1bSJed Brown 591*c4762a1bSJed Brown test: 592*c4762a1bSJed Brown suffix: 3 593*c4762a1bSJed Brown args: -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none 594*c4762a1bSJed Brown output_file: output/ex20opt_ic_3.out 595*c4762a1bSJed Brown 596*c4762a1bSJed Brown test: 597*c4762a1bSJed Brown suffix: 4 598*c4762a1bSJed Brown args: -implicitform 0 -ts_dt 0.01 -ts_max_steps 2 -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view -mode 1 -my_tao_mf 599*c4762a1bSJed Brown TEST*/ 600