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