1*c4762a1bSJed Brown #include <petsctao.h> 2*c4762a1bSJed Brown #include <petscts.h> 3*c4762a1bSJed Brown 4*c4762a1bSJed Brown typedef struct _n_aircraft *Aircraft; 5*c4762a1bSJed Brown struct _n_aircraft { 6*c4762a1bSJed Brown TS ts,quadts; 7*c4762a1bSJed Brown Vec V,W; /* control variables V and W */ 8*c4762a1bSJed Brown PetscInt nsteps; /* number of time steps */ 9*c4762a1bSJed Brown PetscReal ftime; 10*c4762a1bSJed Brown Mat A,H; 11*c4762a1bSJed Brown Mat Jacp,DRDU,DRDP; 12*c4762a1bSJed Brown Vec U,Lambda[1],Mup[1],Lambda2[1],Mup2[1],Dir; 13*c4762a1bSJed Brown Vec rhshp1[1],rhshp2[1],rhshp3[1],rhshp4[1],inthp1[1],inthp2[1],inthp3[1],inthp4[1]; 14*c4762a1bSJed Brown PetscReal lv,lw; 15*c4762a1bSJed Brown PetscBool mf,eh; 16*c4762a1bSJed Brown }; 17*c4762a1bSJed Brown 18*c4762a1bSJed Brown PetscErrorCode FormObjFunctionGradient(Tao,Vec,PetscReal *,Vec,void *); 19*c4762a1bSJed Brown PetscErrorCode FormObjHessian(Tao,Vec,Mat,Mat,void *); 20*c4762a1bSJed Brown PetscErrorCode ComputeObjHessianWithSOA(Vec,PetscScalar[],Aircraft); 21*c4762a1bSJed Brown PetscErrorCode MatrixFreeObjHessian(Tao,Vec,Mat,Mat,void *); 22*c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat,Vec,Vec); 23*c4762a1bSJed Brown 24*c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 25*c4762a1bSJed Brown { 26*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 27*c4762a1bSJed Brown const PetscScalar *u,*v,*w; 28*c4762a1bSJed Brown PetscScalar *f; 29*c4762a1bSJed Brown PetscInt step; 30*c4762a1bSJed Brown PetscErrorCode ierr; 31*c4762a1bSJed Brown 32*c4762a1bSJed Brown PetscFunctionBeginUser; 33*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr); 34*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 35*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr); 36*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr); 37*c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 38*c4762a1bSJed Brown f[0] = v[step]*PetscCosReal(w[step]); 39*c4762a1bSJed Brown f[1] = v[step]*PetscSinReal(w[step]); 40*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 41*c4762a1bSJed Brown ierr = VecRestoreArrayRead(actx->V,&v);CHKERRQ(ierr); 42*c4762a1bSJed Brown ierr = VecRestoreArrayRead(actx->W,&w);CHKERRQ(ierr); 43*c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 44*c4762a1bSJed Brown PetscFunctionReturn(0); 45*c4762a1bSJed Brown } 46*c4762a1bSJed Brown 47*c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec U,Mat A,void *ctx) 48*c4762a1bSJed Brown { 49*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 50*c4762a1bSJed Brown const PetscScalar *u,*v,*w; 51*c4762a1bSJed Brown PetscInt step,rows[2] = {0,1},rowcol[2]; 52*c4762a1bSJed Brown PetscScalar Jp[2][2]; 53*c4762a1bSJed Brown PetscErrorCode ierr; 54*c4762a1bSJed Brown 55*c4762a1bSJed Brown PetscFunctionBeginUser; 56*c4762a1bSJed Brown ierr = MatZeroEntries(A);CHKERRQ(ierr); 57*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr); 58*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 59*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr); 60*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr); 61*c4762a1bSJed Brown 62*c4762a1bSJed Brown Jp[0][0] = PetscCosReal(w[step-1]); 63*c4762a1bSJed Brown Jp[0][1] = -v[step-1]*PetscSinReal(w[step-1]); 64*c4762a1bSJed Brown Jp[1][0] = PetscSinReal(w[step-1]); 65*c4762a1bSJed Brown Jp[1][1] = v[step-1]*PetscCosReal(w[step-1]); 66*c4762a1bSJed Brown 67*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 68*c4762a1bSJed Brown ierr = VecRestoreArrayRead(actx->V,&v);CHKERRQ(ierr); 69*c4762a1bSJed Brown ierr = VecRestoreArrayRead(actx->W,&w);CHKERRQ(ierr); 70*c4762a1bSJed Brown 71*c4762a1bSJed Brown rowcol[0] = 2*(step-1); 72*c4762a1bSJed Brown rowcol[1] = 2*(step-1)+1; 73*c4762a1bSJed Brown ierr = MatSetValues(A,2,rows,2,rowcol,&Jp[0][0],INSERT_VALUES);CHKERRQ(ierr); 74*c4762a1bSJed Brown 75*c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 76*c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 77*c4762a1bSJed Brown PetscFunctionReturn(0); 78*c4762a1bSJed Brown } 79*c4762a1bSJed Brown 80*c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 81*c4762a1bSJed Brown { 82*c4762a1bSJed Brown PetscFunctionBeginUser; 83*c4762a1bSJed Brown PetscFunctionReturn(0); 84*c4762a1bSJed Brown } 85*c4762a1bSJed Brown 86*c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 87*c4762a1bSJed Brown { 88*c4762a1bSJed Brown PetscFunctionBeginUser; 89*c4762a1bSJed Brown PetscFunctionReturn(0); 90*c4762a1bSJed Brown } 91*c4762a1bSJed Brown 92*c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 93*c4762a1bSJed Brown { 94*c4762a1bSJed Brown PetscFunctionBeginUser; 95*c4762a1bSJed Brown PetscFunctionReturn(0); 96*c4762a1bSJed Brown } 97*c4762a1bSJed Brown 98*c4762a1bSJed Brown static PetscErrorCode RHSHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 99*c4762a1bSJed Brown { 100*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 101*c4762a1bSJed Brown const PetscScalar *v,*w,*vl,*vr,*u; 102*c4762a1bSJed Brown PetscScalar *vhv; 103*c4762a1bSJed Brown PetscScalar dJpdP[2][2][2]={{{0}}}; 104*c4762a1bSJed Brown PetscInt step,i,j,k; 105*c4762a1bSJed Brown PetscErrorCode ierr; 106*c4762a1bSJed Brown 107*c4762a1bSJed Brown PetscFunctionBeginUser; 108*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr); 109*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 110*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr); 111*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr); 112*c4762a1bSJed Brown ierr = VecGetArrayRead(Vl[0],&vl);CHKERRQ(ierr); 113*c4762a1bSJed Brown ierr = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr); 114*c4762a1bSJed Brown ierr = VecSet(VHV[0],0.0);CHKERRQ(ierr); 115*c4762a1bSJed Brown ierr = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr); 116*c4762a1bSJed Brown 117*c4762a1bSJed Brown dJpdP[0][0][1] = -PetscSinReal(w[step-1]); 118*c4762a1bSJed Brown dJpdP[0][1][0] = -PetscSinReal(w[step-1]); 119*c4762a1bSJed Brown dJpdP[0][1][1] = -v[step-1]*PetscCosReal(w[step-1]); 120*c4762a1bSJed Brown dJpdP[1][0][1] = PetscCosReal(w[step-1]); 121*c4762a1bSJed Brown dJpdP[1][1][0] = PetscCosReal(w[step-1]); 122*c4762a1bSJed Brown dJpdP[1][1][1] = -v[step-1]*PetscSinReal(w[step-1]); 123*c4762a1bSJed Brown 124*c4762a1bSJed Brown for (j=0; j<2; j++) { 125*c4762a1bSJed Brown vhv[2*(step-1)+j] = 0; 126*c4762a1bSJed Brown for (k=0; k<2; k++) 127*c4762a1bSJed Brown for (i=0; i<2; i++) 128*c4762a1bSJed Brown vhv[2*(step-1)+j] += vl[i]*dJpdP[i][j][k]*vr[2*(step-1)+k]; 129*c4762a1bSJed Brown } 130*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 131*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Vl[0],&vl);CHKERRQ(ierr); 132*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr); 133*c4762a1bSJed Brown ierr = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr); 134*c4762a1bSJed Brown PetscFunctionReturn(0); 135*c4762a1bSJed Brown } 136*c4762a1bSJed Brown 137*c4762a1bSJed Brown /* Vl in NULL,updates to VHV must be added */ 138*c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 139*c4762a1bSJed Brown { 140*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 141*c4762a1bSJed Brown const PetscScalar *v,*w,*vr,*u; 142*c4762a1bSJed Brown PetscScalar *vhv; 143*c4762a1bSJed Brown PetscScalar dRudU[2][2]={{0}}; 144*c4762a1bSJed Brown PetscInt step,j,k; 145*c4762a1bSJed Brown PetscErrorCode ierr; 146*c4762a1bSJed Brown 147*c4762a1bSJed Brown PetscFunctionBeginUser; 148*c4762a1bSJed Brown ierr = TSGetStepNumber(ts,&step);CHKERRQ(ierr); 149*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 150*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->V,&v);CHKERRQ(ierr); 151*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->W,&w);CHKERRQ(ierr); 152*c4762a1bSJed Brown ierr = VecGetArrayRead(Vr,&vr);CHKERRQ(ierr); 153*c4762a1bSJed Brown ierr = VecGetArray(VHV[0],&vhv);CHKERRQ(ierr); 154*c4762a1bSJed Brown 155*c4762a1bSJed Brown dRudU[0][0] = 2.0; 156*c4762a1bSJed Brown dRudU[1][1] = 2.0; 157*c4762a1bSJed Brown 158*c4762a1bSJed Brown for (j=0; j<2; j++) { 159*c4762a1bSJed Brown vhv[j] = 0; 160*c4762a1bSJed Brown for (k=0; k<2; k++) 161*c4762a1bSJed Brown vhv[j] += dRudU[j][k]*vr[k]; 162*c4762a1bSJed Brown } 163*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 164*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Vr,&vr);CHKERRQ(ierr); 165*c4762a1bSJed Brown ierr = VecRestoreArray(VHV[0],&vhv);CHKERRQ(ierr); 166*c4762a1bSJed Brown PetscFunctionReturn(0); 167*c4762a1bSJed Brown } 168*c4762a1bSJed Brown 169*c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 170*c4762a1bSJed Brown { 171*c4762a1bSJed Brown PetscFunctionBeginUser; 172*c4762a1bSJed Brown PetscFunctionReturn(0); 173*c4762a1bSJed Brown } 174*c4762a1bSJed Brown 175*c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 176*c4762a1bSJed Brown { 177*c4762a1bSJed Brown PetscFunctionBeginUser; 178*c4762a1bSJed Brown PetscFunctionReturn(0); 179*c4762a1bSJed Brown } 180*c4762a1bSJed Brown 181*c4762a1bSJed Brown static PetscErrorCode IntegrandHessianProductPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 182*c4762a1bSJed Brown { 183*c4762a1bSJed Brown PetscFunctionBeginUser; 184*c4762a1bSJed Brown PetscFunctionReturn(0); 185*c4762a1bSJed Brown } 186*c4762a1bSJed Brown 187*c4762a1bSJed Brown 188*c4762a1bSJed Brown static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,void *ctx) 189*c4762a1bSJed Brown { 190*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 191*c4762a1bSJed Brown PetscScalar *r; 192*c4762a1bSJed Brown PetscReal dx,dy; 193*c4762a1bSJed Brown const PetscScalar *u; 194*c4762a1bSJed Brown PetscErrorCode ierr; 195*c4762a1bSJed Brown 196*c4762a1bSJed Brown PetscFunctionBegin; 197*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 198*c4762a1bSJed Brown ierr = VecGetArray(R,&r);CHKERRQ(ierr); 199*c4762a1bSJed Brown dx = u[0] - actx->lv*t*PetscCosReal(actx->lw); 200*c4762a1bSJed Brown dy = u[1] - actx->lv*t*PetscSinReal(actx->lw); 201*c4762a1bSJed Brown r[0] = dx*dx+dy*dy; 202*c4762a1bSJed Brown ierr = VecRestoreArray(R,&r);CHKERRQ(ierr); 203*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 204*c4762a1bSJed Brown PetscFunctionReturn(0); 205*c4762a1bSJed Brown } 206*c4762a1bSJed Brown 207*c4762a1bSJed Brown static PetscErrorCode DRDUJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDU,Mat B,void *ctx) 208*c4762a1bSJed Brown { 209*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 210*c4762a1bSJed Brown PetscScalar drdu[2][1]; 211*c4762a1bSJed Brown const PetscScalar *u; 212*c4762a1bSJed Brown PetscReal dx,dy; 213*c4762a1bSJed Brown PetscInt row[] = {0,1},col[] = {0}; 214*c4762a1bSJed Brown PetscErrorCode ierr; 215*c4762a1bSJed Brown 216*c4762a1bSJed Brown PetscFunctionBegin; 217*c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 218*c4762a1bSJed Brown dx = u[0] - actx->lv*t*PetscCosReal(actx->lw); 219*c4762a1bSJed Brown dy = u[1] - actx->lv*t*PetscSinReal(actx->lw); 220*c4762a1bSJed Brown drdu[0][0] = 2.*dx; 221*c4762a1bSJed Brown drdu[1][0] = 2.*dy; 222*c4762a1bSJed Brown ierr = MatSetValues(DRDU,2,row,1,col,&drdu[0][0],INSERT_VALUES);CHKERRQ(ierr); 223*c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 224*c4762a1bSJed Brown ierr = MatAssemblyBegin(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 225*c4762a1bSJed Brown ierr = MatAssemblyEnd(DRDU,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 226*c4762a1bSJed Brown PetscFunctionReturn(0); 227*c4762a1bSJed Brown } 228*c4762a1bSJed Brown 229*c4762a1bSJed Brown static PetscErrorCode DRDPJacobianTranspose(TS ts,PetscReal t,Vec U,Mat DRDP,void *ctx) 230*c4762a1bSJed Brown { 231*c4762a1bSJed Brown PetscErrorCode ierr; 232*c4762a1bSJed Brown 233*c4762a1bSJed Brown PetscFunctionBegin; 234*c4762a1bSJed Brown ierr = MatZeroEntries(DRDP);CHKERRQ(ierr); 235*c4762a1bSJed Brown ierr = MatAssemblyBegin(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 236*c4762a1bSJed Brown ierr = MatAssemblyEnd(DRDP,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 237*c4762a1bSJed Brown PetscFunctionReturn(0); 238*c4762a1bSJed Brown } 239*c4762a1bSJed Brown 240*c4762a1bSJed Brown int main(int argc,char **argv) 241*c4762a1bSJed Brown { 242*c4762a1bSJed Brown Vec P,PL,PU; 243*c4762a1bSJed Brown struct _n_aircraft aircraft; 244*c4762a1bSJed Brown PetscMPIInt size; 245*c4762a1bSJed Brown Tao tao; 246*c4762a1bSJed Brown KSP ksp; 247*c4762a1bSJed Brown PC pc; 248*c4762a1bSJed Brown PetscScalar *u,*p; 249*c4762a1bSJed Brown PetscInt i; 250*c4762a1bSJed Brown PetscErrorCode ierr; 251*c4762a1bSJed Brown 252*c4762a1bSJed Brown /* Initialize program */ 253*c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,NULL,NULL);if(ierr) return ierr; 254*c4762a1bSJed Brown ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 255*c4762a1bSJed Brown if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 256*c4762a1bSJed Brown 257*c4762a1bSJed Brown /* Parameter settings */ 258*c4762a1bSJed Brown aircraft.ftime = 1.; /* time interval in hour */ 259*c4762a1bSJed Brown aircraft.nsteps = 10; /* number of steps */ 260*c4762a1bSJed Brown aircraft.lv = 2.0; /* leader speed in kmph */ 261*c4762a1bSJed Brown aircraft.lw = PETSC_PI/4.; /* leader heading angle */ 262*c4762a1bSJed Brown 263*c4762a1bSJed Brown ierr = PetscOptionsGetReal(NULL,NULL,"-ftime",&aircraft.ftime,NULL);CHKERRQ(ierr); 264*c4762a1bSJed Brown ierr = PetscOptionsGetInt(NULL,NULL,"-nsteps",&aircraft.nsteps,NULL);CHKERRQ(ierr); 265*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-matrixfree",&aircraft.mf);CHKERRQ(ierr); 266*c4762a1bSJed Brown ierr = PetscOptionsHasName(NULL,NULL,"-exacthessian",&aircraft.eh);CHKERRQ(ierr); 267*c4762a1bSJed Brown 268*c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 269*c4762a1bSJed Brown ierr = TaoCreate(PETSC_COMM_WORLD,&tao);CHKERRQ(ierr); 270*c4762a1bSJed Brown ierr = TaoSetType(tao,TAOBQNLS);CHKERRQ(ierr); 271*c4762a1bSJed Brown 272*c4762a1bSJed Brown /* Create necessary matrix and vectors, solve same ODE on every process */ 273*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&aircraft.A);CHKERRQ(ierr); 274*c4762a1bSJed Brown ierr = MatSetSizes(aircraft.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 275*c4762a1bSJed Brown ierr = MatSetFromOptions(aircraft.A);CHKERRQ(ierr); 276*c4762a1bSJed Brown ierr = MatSetUp(aircraft.A);CHKERRQ(ierr); 277*c4762a1bSJed Brown ierr = MatAssemblyBegin(aircraft.A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 278*c4762a1bSJed Brown ierr = MatAssemblyEnd(aircraft.A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 279*c4762a1bSJed Brown ierr = MatShift(aircraft.A,1);CHKERRQ(ierr); 280*c4762a1bSJed Brown ierr = MatShift(aircraft.A,-1);CHKERRQ(ierr); 281*c4762a1bSJed Brown 282*c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&aircraft.Jacp);CHKERRQ(ierr); 283*c4762a1bSJed Brown ierr = MatSetSizes(aircraft.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2*aircraft.nsteps);CHKERRQ(ierr); 284*c4762a1bSJed Brown ierr = MatSetFromOptions(aircraft.Jacp);CHKERRQ(ierr); 285*c4762a1bSJed Brown ierr = MatSetUp(aircraft.Jacp);CHKERRQ(ierr); 286*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,1,NULL,&aircraft.DRDP);CHKERRQ(ierr); 287*c4762a1bSJed Brown ierr = MatSetUp(aircraft.DRDP);CHKERRQ(ierr); 288*c4762a1bSJed Brown ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&aircraft.DRDU);CHKERRQ(ierr); 289*c4762a1bSJed Brown ierr = MatSetUp(aircraft.DRDU);CHKERRQ(ierr); 290*c4762a1bSJed Brown 291*c4762a1bSJed Brown /* Create timestepping solver context */ 292*c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&aircraft.ts);CHKERRQ(ierr); 293*c4762a1bSJed Brown ierr = TSSetType(aircraft.ts,TSRK);CHKERRQ(ierr); 294*c4762a1bSJed Brown ierr = TSSetRHSFunction(aircraft.ts,NULL,RHSFunction,&aircraft);CHKERRQ(ierr); 295*c4762a1bSJed Brown ierr = TSSetRHSJacobian(aircraft.ts,aircraft.A,aircraft.A,TSComputeRHSJacobianConstant,&aircraft);CHKERRQ(ierr); 296*c4762a1bSJed Brown ierr = TSSetRHSJacobianP(aircraft.ts,aircraft.Jacp,RHSJacobianP,&aircraft);CHKERRQ(ierr); 297*c4762a1bSJed Brown ierr = TSSetExactFinalTime(aircraft.ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 298*c4762a1bSJed Brown ierr = TSSetEquationType(aircraft.ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 299*c4762a1bSJed Brown 300*c4762a1bSJed Brown /* Set initial conditions */ 301*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.A,&aircraft.U,NULL);CHKERRQ(ierr); 302*c4762a1bSJed Brown ierr = TSSetSolution(aircraft.ts,aircraft.U);CHKERRQ(ierr); 303*c4762a1bSJed Brown ierr = VecGetArray(aircraft.U,&u);CHKERRQ(ierr); 304*c4762a1bSJed Brown u[0] = 1.5; 305*c4762a1bSJed Brown u[1] = 0; 306*c4762a1bSJed Brown ierr = VecRestoreArray(aircraft.U,&u);CHKERRQ(ierr); 307*c4762a1bSJed Brown ierr = VecCreate(PETSC_COMM_WORLD,&aircraft.V);CHKERRQ(ierr); 308*c4762a1bSJed Brown ierr = VecSetSizes(aircraft.V,PETSC_DECIDE,aircraft.nsteps);CHKERRQ(ierr); 309*c4762a1bSJed Brown ierr = VecSetUp(aircraft.V);CHKERRQ(ierr); 310*c4762a1bSJed Brown ierr = VecDuplicate(aircraft.V,&aircraft.W);CHKERRQ(ierr); 311*c4762a1bSJed Brown ierr = VecSet(aircraft.V,1.);CHKERRQ(ierr); 312*c4762a1bSJed Brown ierr = VecSet(aircraft.W,PETSC_PI/4.);CHKERRQ(ierr); 313*c4762a1bSJed Brown 314*c4762a1bSJed Brown /* Save trajectory of solution so that TSAdjointSolve() may be used */ 315*c4762a1bSJed Brown ierr = TSSetSaveTrajectory(aircraft.ts);CHKERRQ(ierr); 316*c4762a1bSJed Brown 317*c4762a1bSJed Brown /* Set sensitivity context */ 318*c4762a1bSJed Brown ierr = TSCreateQuadratureTS(aircraft.ts,PETSC_FALSE,&aircraft.quadts);CHKERRQ(ierr); 319*c4762a1bSJed Brown ierr = TSSetRHSFunction(aircraft.quadts,NULL,(TSRHSFunction)CostIntegrand,&aircraft);CHKERRQ(ierr); 320*c4762a1bSJed Brown ierr = TSSetRHSJacobian(aircraft.quadts,aircraft.DRDU,aircraft.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&aircraft);CHKERRQ(ierr); 321*c4762a1bSJed Brown ierr = TSSetRHSJacobianP(aircraft.quadts,aircraft.DRDP,(TSRHSJacobianP)DRDPJacobianTranspose,&aircraft);CHKERRQ(ierr); 322*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.A,&aircraft.Lambda[0],NULL);CHKERRQ(ierr); 323*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.Jacp,&aircraft.Mup[0],NULL);CHKERRQ(ierr); 324*c4762a1bSJed Brown if (aircraft.eh) { 325*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.A,&aircraft.rhshp1[0],NULL);CHKERRQ(ierr); 326*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.A,&aircraft.rhshp2[0],NULL);CHKERRQ(ierr); 327*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.Jacp,&aircraft.rhshp3[0],NULL);CHKERRQ(ierr); 328*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.Jacp,&aircraft.rhshp4[0],NULL);CHKERRQ(ierr); 329*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.DRDU,&aircraft.inthp1[0],NULL);CHKERRQ(ierr); 330*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.DRDU,&aircraft.inthp2[0],NULL);CHKERRQ(ierr); 331*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.DRDP,&aircraft.inthp3[0],NULL);CHKERRQ(ierr); 332*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.DRDP,&aircraft.inthp4[0],NULL);CHKERRQ(ierr); 333*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.Jacp,&aircraft.Dir,NULL);CHKERRQ(ierr); 334*c4762a1bSJed Brown ierr = TSSetRHSHessianProduct(aircraft.ts,aircraft.rhshp1,RHSHessianProductUU,aircraft.rhshp2,RHSHessianProductUP,aircraft.rhshp3,RHSHessianProductPU,aircraft.rhshp4,RHSHessianProductPP,&aircraft);CHKERRQ(ierr); 335*c4762a1bSJed Brown ierr = TSSetRHSHessianProduct(aircraft.quadts,aircraft.inthp1,IntegrandHessianProductUU,aircraft.inthp2,IntegrandHessianProductUP,aircraft.inthp3,IntegrandHessianProductPU,aircraft.inthp4,IntegrandHessianProductPP,&aircraft);CHKERRQ(ierr); 336*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.A,&aircraft.Lambda2[0],NULL);CHKERRQ(ierr); 337*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.Jacp,&aircraft.Mup2[0],NULL);CHKERRQ(ierr); 338*c4762a1bSJed Brown } 339*c4762a1bSJed Brown ierr = TSSetFromOptions(aircraft.ts);CHKERRQ(ierr); 340*c4762a1bSJed Brown ierr = TSSetMaxTime(aircraft.ts,aircraft.ftime);CHKERRQ(ierr); 341*c4762a1bSJed Brown ierr = TSSetTimeStep(aircraft.ts,aircraft.ftime/aircraft.nsteps);CHKERRQ(ierr); 342*c4762a1bSJed Brown 343*c4762a1bSJed Brown /* Set initial solution guess */ 344*c4762a1bSJed Brown ierr = MatCreateVecs(aircraft.Jacp,&P,NULL);CHKERRQ(ierr); 345*c4762a1bSJed Brown ierr = VecGetArray(P,&p);CHKERRQ(ierr); 346*c4762a1bSJed Brown for (i=0; i<aircraft.nsteps; i++) { 347*c4762a1bSJed Brown p[2*i] = 2.0; 348*c4762a1bSJed Brown p[2*i+1] = PETSC_PI/2.0; 349*c4762a1bSJed Brown } 350*c4762a1bSJed Brown ierr = VecRestoreArray(P,&p);CHKERRQ(ierr); 351*c4762a1bSJed Brown ierr = VecDuplicate(P,&PU);CHKERRQ(ierr); 352*c4762a1bSJed Brown ierr = VecDuplicate(P,&PL);CHKERRQ(ierr); 353*c4762a1bSJed Brown ierr = VecGetArray(PU,&p);CHKERRQ(ierr); 354*c4762a1bSJed Brown for (i=0; i<aircraft.nsteps; i++) { 355*c4762a1bSJed Brown p[2*i] = 2.0; 356*c4762a1bSJed Brown p[2*i+1] = PETSC_PI; 357*c4762a1bSJed Brown } 358*c4762a1bSJed Brown ierr = VecRestoreArray(PU,&p);CHKERRQ(ierr); 359*c4762a1bSJed Brown ierr = VecGetArray(PL,&p);CHKERRQ(ierr); 360*c4762a1bSJed Brown for (i=0; i<aircraft.nsteps; i++) { 361*c4762a1bSJed Brown p[2*i] = 0.0; 362*c4762a1bSJed Brown p[2*i+1] = -PETSC_PI; 363*c4762a1bSJed Brown } 364*c4762a1bSJed Brown ierr = VecRestoreArray(PL,&p);CHKERRQ(ierr); 365*c4762a1bSJed Brown 366*c4762a1bSJed Brown ierr = TaoSetInitialVector(tao,P);CHKERRQ(ierr); 367*c4762a1bSJed Brown ierr = TaoSetVariableBounds(tao,PL,PU);CHKERRQ(ierr); 368*c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 369*c4762a1bSJed Brown ierr = TaoSetObjectiveAndGradientRoutine(tao,FormObjFunctionGradient,(void *)&aircraft);CHKERRQ(ierr); 370*c4762a1bSJed Brown 371*c4762a1bSJed Brown if (aircraft.eh) { 372*c4762a1bSJed Brown if (aircraft.mf) { 373*c4762a1bSJed Brown ierr = MatCreateShell(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2*aircraft.nsteps,2*aircraft.nsteps,(void*)&aircraft,&aircraft.H);CHKERRQ(ierr); 374*c4762a1bSJed Brown ierr = MatShellSetOperation(aircraft.H,MATOP_MULT,(void(*)(void))MyMatMult);CHKERRQ(ierr); 375*c4762a1bSJed Brown ierr = MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 376*c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao,aircraft.H,aircraft.H,MatrixFreeObjHessian,(void*)&aircraft);CHKERRQ(ierr); 377*c4762a1bSJed Brown } else { 378*c4762a1bSJed Brown ierr = MatCreateDense(MPI_COMM_WORLD,PETSC_DETERMINE,PETSC_DETERMINE,2*aircraft.nsteps,2*aircraft.nsteps,NULL,&(aircraft.H));CHKERRQ(ierr); 379*c4762a1bSJed Brown ierr = MatSetOption(aircraft.H,MAT_SYMMETRIC,PETSC_TRUE);CHKERRQ(ierr); 380*c4762a1bSJed Brown ierr = TaoSetHessianRoutine(tao,aircraft.H,aircraft.H,FormObjHessian,(void *)&aircraft);CHKERRQ(ierr); 381*c4762a1bSJed Brown } 382*c4762a1bSJed Brown } 383*c4762a1bSJed Brown 384*c4762a1bSJed Brown /* Check for any TAO command line options */ 385*c4762a1bSJed Brown ierr = TaoGetKSP(tao,&ksp);CHKERRQ(ierr); 386*c4762a1bSJed Brown if (ksp) { 387*c4762a1bSJed Brown ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 388*c4762a1bSJed Brown ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 389*c4762a1bSJed Brown } 390*c4762a1bSJed Brown ierr = TaoSetFromOptions(tao);CHKERRQ(ierr); 391*c4762a1bSJed Brown 392*c4762a1bSJed Brown ierr = TaoSolve(tao); CHKERRQ(ierr); 393*c4762a1bSJed Brown ierr = VecView(P,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 394*c4762a1bSJed Brown 395*c4762a1bSJed Brown /* Free TAO data structures */ 396*c4762a1bSJed Brown ierr = TaoDestroy(&tao);CHKERRQ(ierr); 397*c4762a1bSJed Brown 398*c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 399*c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they 400*c4762a1bSJed Brown are no longer needed. 401*c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 402*c4762a1bSJed Brown ierr = TSDestroy(&aircraft.ts);CHKERRQ(ierr); 403*c4762a1bSJed Brown ierr = MatDestroy(&aircraft.A);CHKERRQ(ierr); 404*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.U);CHKERRQ(ierr); 405*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.V);CHKERRQ(ierr); 406*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.W);CHKERRQ(ierr); 407*c4762a1bSJed Brown ierr = VecDestroy(&P);CHKERRQ(ierr); 408*c4762a1bSJed Brown ierr = VecDestroy(&PU);CHKERRQ(ierr); 409*c4762a1bSJed Brown ierr = VecDestroy(&PL);CHKERRQ(ierr); 410*c4762a1bSJed Brown ierr = MatDestroy(&aircraft.Jacp);CHKERRQ(ierr); 411*c4762a1bSJed Brown ierr = MatDestroy(&aircraft.DRDU);CHKERRQ(ierr); 412*c4762a1bSJed Brown ierr = MatDestroy(&aircraft.DRDP);CHKERRQ(ierr); 413*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.Lambda[0]);CHKERRQ(ierr); 414*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.Mup[0]);CHKERRQ(ierr); 415*c4762a1bSJed Brown ierr = VecDestroy(&P);CHKERRQ(ierr); 416*c4762a1bSJed Brown if (aircraft.eh) { 417*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.Lambda2[0]);CHKERRQ(ierr); 418*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.Mup2[0]);CHKERRQ(ierr); 419*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.Dir);CHKERRQ(ierr); 420*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.rhshp1[0]);CHKERRQ(ierr); 421*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.rhshp2[0]);CHKERRQ(ierr); 422*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.rhshp3[0]);CHKERRQ(ierr); 423*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.rhshp4[0]);CHKERRQ(ierr); 424*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.inthp1[0]);CHKERRQ(ierr); 425*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.inthp2[0]);CHKERRQ(ierr); 426*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.inthp3[0]);CHKERRQ(ierr); 427*c4762a1bSJed Brown ierr = VecDestroy(&aircraft.inthp4[0]);CHKERRQ(ierr); 428*c4762a1bSJed Brown ierr = MatDestroy(&aircraft.H);CHKERRQ(ierr); 429*c4762a1bSJed Brown } 430*c4762a1bSJed Brown ierr = PetscFinalize(); 431*c4762a1bSJed Brown return ierr; 432*c4762a1bSJed Brown } 433*c4762a1bSJed Brown 434*c4762a1bSJed Brown /* 435*c4762a1bSJed Brown FormObjFunctionGradient - Evaluates the function and corresponding gradient. 436*c4762a1bSJed Brown 437*c4762a1bSJed Brown Input Parameters: 438*c4762a1bSJed Brown tao - the Tao context 439*c4762a1bSJed Brown P - the input vector 440*c4762a1bSJed Brown ctx - optional aircraft-defined context, as set by TaoSetObjectiveAndGradientRoutine() 441*c4762a1bSJed Brown 442*c4762a1bSJed Brown Output Parameters: 443*c4762a1bSJed Brown f - the newly evaluated function 444*c4762a1bSJed Brown G - the newly evaluated gradient 445*c4762a1bSJed Brown */ 446*c4762a1bSJed Brown PetscErrorCode FormObjFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx) 447*c4762a1bSJed Brown { 448*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 449*c4762a1bSJed Brown TS ts = actx->ts; 450*c4762a1bSJed Brown Vec Q; 451*c4762a1bSJed Brown const PetscScalar *p,*q; 452*c4762a1bSJed Brown PetscScalar *u,*v,*w; 453*c4762a1bSJed Brown PetscInt i; 454*c4762a1bSJed Brown PetscErrorCode ierr; 455*c4762a1bSJed Brown 456*c4762a1bSJed Brown PetscFunctionBeginUser; 457*c4762a1bSJed Brown ierr = VecGetArrayRead(P,&p);CHKERRQ(ierr); 458*c4762a1bSJed Brown ierr = VecGetArray(actx->V,&v);CHKERRQ(ierr); 459*c4762a1bSJed Brown ierr = VecGetArray(actx->W,&w);CHKERRQ(ierr); 460*c4762a1bSJed Brown for (i=0; i<actx->nsteps; i++) { 461*c4762a1bSJed Brown v[i] = p[2*i]; 462*c4762a1bSJed Brown w[i] = p[2*i+1]; 463*c4762a1bSJed Brown } 464*c4762a1bSJed Brown ierr = VecRestoreArrayRead(P,&p);CHKERRQ(ierr); 465*c4762a1bSJed Brown ierr = VecRestoreArray(actx->V,&v);CHKERRQ(ierr); 466*c4762a1bSJed Brown ierr = VecRestoreArray(actx->W,&w);CHKERRQ(ierr); 467*c4762a1bSJed Brown 468*c4762a1bSJed Brown ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 469*c4762a1bSJed Brown ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 470*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 471*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,actx->ftime/actx->nsteps);CHKERRQ(ierr); 472*c4762a1bSJed Brown 473*c4762a1bSJed Brown /* reinitialize system state */ 474*c4762a1bSJed Brown ierr = VecGetArray(actx->U,&u);CHKERRQ(ierr); 475*c4762a1bSJed Brown u[0] = 2.0; 476*c4762a1bSJed Brown u[1] = 0; 477*c4762a1bSJed Brown ierr = VecRestoreArray(actx->U,&u);CHKERRQ(ierr); 478*c4762a1bSJed Brown 479*c4762a1bSJed Brown /* reinitialize the integral value */ 480*c4762a1bSJed Brown ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr); 481*c4762a1bSJed Brown ierr = VecSet(Q,0.0);CHKERRQ(ierr); 482*c4762a1bSJed Brown 483*c4762a1bSJed Brown ierr = TSSolve(ts,actx->U);CHKERRQ(ierr); 484*c4762a1bSJed Brown 485*c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 486*c4762a1bSJed Brown ierr = VecSet(actx->Lambda[0],0.0);CHKERRQ(ierr); 487*c4762a1bSJed Brown ierr = VecSet(actx->Mup[0],0.0);CHKERRQ(ierr); 488*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,actx->Lambda,actx->Mup);CHKERRQ(ierr); 489*c4762a1bSJed Brown 490*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 491*c4762a1bSJed Brown ierr = VecCopy(actx->Mup[0],G);CHKERRQ(ierr); 492*c4762a1bSJed Brown ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr); 493*c4762a1bSJed Brown ierr = VecGetArrayRead(Q,&q);CHKERRQ(ierr); 494*c4762a1bSJed Brown *f = q[0]; 495*c4762a1bSJed Brown ierr = VecRestoreArrayRead(Q,&q);CHKERRQ(ierr); 496*c4762a1bSJed Brown PetscFunctionReturn(0); 497*c4762a1bSJed Brown } 498*c4762a1bSJed Brown 499*c4762a1bSJed Brown PetscErrorCode FormObjHessian(Tao tao,Vec P,Mat H,Mat Hpre,void *ctx) 500*c4762a1bSJed Brown { 501*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 502*c4762a1bSJed Brown const PetscScalar *p; 503*c4762a1bSJed Brown PetscScalar *harr,*v,*w,one = 1.0; 504*c4762a1bSJed Brown PetscInt ind[1]; 505*c4762a1bSJed Brown PetscInt *cols,i; 506*c4762a1bSJed Brown Vec Dir; 507*c4762a1bSJed Brown PetscErrorCode ierr; 508*c4762a1bSJed Brown 509*c4762a1bSJed Brown PetscFunctionBeginUser; 510*c4762a1bSJed Brown /* set up control parameters */ 511*c4762a1bSJed Brown ierr = VecGetArrayRead(P,&p);CHKERRQ(ierr); 512*c4762a1bSJed Brown ierr = VecGetArray(actx->V,&v);CHKERRQ(ierr); 513*c4762a1bSJed Brown ierr = VecGetArray(actx->W,&w);CHKERRQ(ierr); 514*c4762a1bSJed Brown for (i=0; i<actx->nsteps; i++) { 515*c4762a1bSJed Brown v[i] = p[2*i]; 516*c4762a1bSJed Brown w[i] = p[2*i+1]; 517*c4762a1bSJed Brown } 518*c4762a1bSJed Brown ierr = VecRestoreArrayRead(P,&p);CHKERRQ(ierr); 519*c4762a1bSJed Brown ierr = VecRestoreArray(actx->V,&v);CHKERRQ(ierr); 520*c4762a1bSJed Brown ierr = VecRestoreArray(actx->W,&w);CHKERRQ(ierr); 521*c4762a1bSJed Brown 522*c4762a1bSJed Brown ierr = PetscMalloc1(2*actx->nsteps,&harr);CHKERRQ(ierr); 523*c4762a1bSJed Brown ierr = PetscMalloc1(2*actx->nsteps,&cols);CHKERRQ(ierr); 524*c4762a1bSJed Brown for (i=0; i<2*actx->nsteps; i++) cols[i] = i; 525*c4762a1bSJed Brown ierr = VecDuplicate(P,&Dir);CHKERRQ(ierr); 526*c4762a1bSJed Brown for (i=0; i<2*actx->nsteps; i++) { 527*c4762a1bSJed Brown ind[0] = i; 528*c4762a1bSJed Brown ierr = VecSet(Dir,0.0);CHKERRQ(ierr); 529*c4762a1bSJed Brown ierr = VecSetValues(Dir,1,ind,&one,INSERT_VALUES);CHKERRQ(ierr); 530*c4762a1bSJed Brown ierr = VecAssemblyBegin(Dir);CHKERRQ(ierr); 531*c4762a1bSJed Brown ierr = VecAssemblyEnd(Dir);CHKERRQ(ierr); 532*c4762a1bSJed Brown ierr = ComputeObjHessianWithSOA(Dir,harr,actx);CHKERRQ(ierr); 533*c4762a1bSJed Brown ierr = MatSetValues(H,1,ind,2*actx->nsteps,cols,harr,INSERT_VALUES);CHKERRQ(ierr); 534*c4762a1bSJed Brown ierr = MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 535*c4762a1bSJed Brown ierr = MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 536*c4762a1bSJed Brown if (H != Hpre) { 537*c4762a1bSJed Brown ierr = MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 538*c4762a1bSJed Brown ierr = MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 539*c4762a1bSJed Brown } 540*c4762a1bSJed Brown } 541*c4762a1bSJed Brown ierr = PetscFree(cols);CHKERRQ(ierr); 542*c4762a1bSJed Brown ierr = PetscFree(harr);CHKERRQ(ierr); 543*c4762a1bSJed Brown ierr = VecDestroy(&Dir);CHKERRQ(ierr); 544*c4762a1bSJed Brown PetscFunctionReturn(0); 545*c4762a1bSJed Brown } 546*c4762a1bSJed Brown 547*c4762a1bSJed Brown PetscErrorCode MatrixFreeObjHessian(Tao tao, Vec P, Mat H, Mat Hpre, void *ctx) 548*c4762a1bSJed Brown { 549*c4762a1bSJed Brown Aircraft actx = (Aircraft)ctx; 550*c4762a1bSJed Brown PetscScalar *v,*w; 551*c4762a1bSJed Brown const PetscScalar *p; 552*c4762a1bSJed Brown PetscInt i; 553*c4762a1bSJed Brown PetscErrorCode ierr; 554*c4762a1bSJed Brown 555*c4762a1bSJed Brown PetscFunctionBegin; 556*c4762a1bSJed Brown ierr = VecGetArrayRead(P,&p);CHKERRQ(ierr); 557*c4762a1bSJed Brown ierr = VecGetArray(actx->V,&v);CHKERRQ(ierr); 558*c4762a1bSJed Brown ierr = VecGetArray(actx->W,&w);CHKERRQ(ierr); 559*c4762a1bSJed Brown for (i=0; i<actx->nsteps; i++) { 560*c4762a1bSJed Brown v[i] = p[2*i]; 561*c4762a1bSJed Brown w[i] = p[2*i+1]; 562*c4762a1bSJed Brown } 563*c4762a1bSJed Brown ierr = VecRestoreArrayRead(P,&p);CHKERRQ(ierr); 564*c4762a1bSJed Brown ierr = VecRestoreArray(actx->V,&v);CHKERRQ(ierr); 565*c4762a1bSJed Brown ierr = VecRestoreArray(actx->W,&w);CHKERRQ(ierr); 566*c4762a1bSJed Brown PetscFunctionReturn(0); 567*c4762a1bSJed Brown } 568*c4762a1bSJed Brown 569*c4762a1bSJed Brown PetscErrorCode MyMatMult(Mat H_shell, Vec X, Vec Y) 570*c4762a1bSJed Brown { 571*c4762a1bSJed Brown PetscScalar *y; 572*c4762a1bSJed Brown void *ptr; 573*c4762a1bSJed Brown PetscErrorCode ierr; 574*c4762a1bSJed Brown 575*c4762a1bSJed Brown PetscFunctionBegin; 576*c4762a1bSJed Brown ierr = MatShellGetContext(H_shell,&ptr);CHKERRQ(ierr); 577*c4762a1bSJed Brown ierr = VecGetArray(Y,&y);CHKERRQ(ierr); 578*c4762a1bSJed Brown ierr = ComputeObjHessianWithSOA(X,y,(Aircraft)ptr);CHKERRQ(ierr); 579*c4762a1bSJed Brown ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr); 580*c4762a1bSJed Brown PetscFunctionReturn(0); 581*c4762a1bSJed Brown } 582*c4762a1bSJed Brown 583*c4762a1bSJed Brown PetscErrorCode ComputeObjHessianWithSOA(Vec Dir,PetscScalar arr[],Aircraft actx) 584*c4762a1bSJed Brown { 585*c4762a1bSJed Brown TS ts = actx->ts; 586*c4762a1bSJed Brown const PetscScalar *z_ptr; 587*c4762a1bSJed Brown PetscScalar *u; 588*c4762a1bSJed Brown Vec Q; 589*c4762a1bSJed Brown PetscInt i; 590*c4762a1bSJed Brown PetscErrorCode ierr; 591*c4762a1bSJed Brown 592*c4762a1bSJed Brown PetscFunctionBeginUser; 593*c4762a1bSJed Brown /* Reset TSAdjoint so that AdjointSetUp will be called again */ 594*c4762a1bSJed Brown ierr = TSAdjointReset(ts);CHKERRQ(ierr); 595*c4762a1bSJed Brown 596*c4762a1bSJed Brown ierr = TSSetTime(ts,0.0);CHKERRQ(ierr); 597*c4762a1bSJed Brown ierr = TSSetStepNumber(ts,0);CHKERRQ(ierr); 598*c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 599*c4762a1bSJed Brown ierr = TSSetTimeStep(ts,actx->ftime/actx->nsteps);CHKERRQ(ierr); 600*c4762a1bSJed Brown ierr = TSSetCostHessianProducts(actx->ts,1,actx->Lambda2,actx->Mup2,Dir);CHKERRQ(ierr); 601*c4762a1bSJed Brown 602*c4762a1bSJed Brown /* reinitialize system state */ 603*c4762a1bSJed Brown ierr = VecGetArray(actx->U,&u);CHKERRQ(ierr); 604*c4762a1bSJed Brown u[0] = 2.0; 605*c4762a1bSJed Brown u[1] = 0; 606*c4762a1bSJed Brown ierr = VecRestoreArray(actx->U,&u);CHKERRQ(ierr); 607*c4762a1bSJed Brown 608*c4762a1bSJed Brown /* reinitialize the integral value */ 609*c4762a1bSJed Brown ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr); 610*c4762a1bSJed Brown ierr = VecSet(Q,0.0);CHKERRQ(ierr); 611*c4762a1bSJed Brown 612*c4762a1bSJed Brown /* initialize tlm variable */ 613*c4762a1bSJed Brown ierr = MatZeroEntries(actx->Jacp);CHKERRQ(ierr); 614*c4762a1bSJed Brown ierr = MatAssemblyBegin(actx->Jacp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 615*c4762a1bSJed Brown ierr = MatAssemblyEnd(actx->Jacp,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 616*c4762a1bSJed Brown ierr = TSAdjointSetForward(ts,actx->Jacp);CHKERRQ(ierr); 617*c4762a1bSJed Brown 618*c4762a1bSJed Brown ierr = TSSolve(ts,actx->U);CHKERRQ(ierr); 619*c4762a1bSJed Brown 620*c4762a1bSJed Brown /* Set terminal conditions for first- and second-order adjonts */ 621*c4762a1bSJed Brown ierr = VecSet(actx->Lambda[0],0.0);CHKERRQ(ierr); 622*c4762a1bSJed Brown ierr = VecSet(actx->Mup[0],0.0);CHKERRQ(ierr); 623*c4762a1bSJed Brown ierr = VecSet(actx->Lambda2[0],0.0);CHKERRQ(ierr); 624*c4762a1bSJed Brown ierr = VecSet(actx->Mup2[0],0.0);CHKERRQ(ierr); 625*c4762a1bSJed Brown ierr = TSSetCostGradients(ts,1,actx->Lambda,actx->Mup);CHKERRQ(ierr); 626*c4762a1bSJed Brown 627*c4762a1bSJed Brown ierr = TSGetCostIntegral(ts,&Q);CHKERRQ(ierr); 628*c4762a1bSJed Brown 629*c4762a1bSJed Brown /* Reset initial conditions for the adjoint integration */ 630*c4762a1bSJed Brown ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 631*c4762a1bSJed Brown 632*c4762a1bSJed Brown /* initial condition does not depend on p, so that lambda is not needed to assemble G */ 633*c4762a1bSJed Brown ierr = VecGetArrayRead(actx->Mup2[0],&z_ptr);CHKERRQ(ierr); 634*c4762a1bSJed Brown for (i=0; i<2*actx->nsteps; i++) arr[i] = z_ptr[i]; 635*c4762a1bSJed Brown ierr = VecRestoreArrayRead(actx->Mup2[0],&z_ptr);CHKERRQ(ierr); 636*c4762a1bSJed Brown 637*c4762a1bSJed Brown /* Disable second-order adjoint mode */ 638*c4762a1bSJed Brown ierr = TSAdjointReset(ts);CHKERRQ(ierr); 639*c4762a1bSJed Brown ierr = TSAdjointResetForward(ts);CHKERRQ(ierr); 640*c4762a1bSJed Brown PetscFunctionReturn(0); 641*c4762a1bSJed Brown } 642*c4762a1bSJed Brown 643*c4762a1bSJed Brown /*TEST 644*c4762a1bSJed Brown build: 645*c4762a1bSJed Brown requires: !complex !single 646*c4762a1bSJed Brown 647*c4762a1bSJed Brown test: 648*c4762a1bSJed Brown args: -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_gatol 1e-7 649*c4762a1bSJed Brown 650*c4762a1bSJed Brown test: 651*c4762a1bSJed Brown suffix: 2 652*c4762a1bSJed Brown args: -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_view -tao_type bntr -tao_bnk_pc_type none -exacthessian 653*c4762a1bSJed Brown 654*c4762a1bSJed Brown test: 655*c4762a1bSJed Brown suffix: 3 656*c4762a1bSJed Brown args: -ts_adapt_type none -ts_type rk -ts_rk_type 3 -viewer_binary_skip_info -tao_monitor -tao_view -tao_type bntr -tao_bnk_pc_type none -exacthessian -matrixfree 657*c4762a1bSJed Brown TEST*/ 658*c4762a1bSJed Brown 659*c4762a1bSJed Brown 660