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