1c4762a1bSJed Brown static char help[] = "Solves a ODE-constrained optimization problem -- finding the optimal initial conditions for the van der Pol equation.\n"; 2c4762a1bSJed Brown 30e3d61c9SBarry Smith /* 4c4762a1bSJed Brown Concepts: TS^time-dependent nonlinear problems 5c4762a1bSJed Brown Concepts: TS^van der Pol equation DAE equivalent 6ee300463SSatish Balay Concepts: TS^Optimization using adjoint sensitivity analysis 7c4762a1bSJed Brown Processors: 1 8c4762a1bSJed Brown */ 90e3d61c9SBarry Smith /* 10c4762a1bSJed Brown Notes: 11c4762a1bSJed Brown This code demonstrates how to solve an ODE-constrained optimization problem with TAO, TSAdjoint and TS. 12c4762a1bSJed Brown The nonlinear problem is written in an ODE equivalent form. 13c4762a1bSJed Brown The objective is to minimize the difference between observation and model prediction by finding optimal values for initial conditions. 14c4762a1bSJed Brown The gradient is computed with the discrete adjoint of an implicit method or an explicit method, see ex20adj.c for details. 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown 17c4762a1bSJed Brown #include <petsctao.h> 18c4762a1bSJed Brown #include <petscts.h> 19c4762a1bSJed Brown 20c4762a1bSJed Brown typedef struct _n_User *User; 21c4762a1bSJed Brown struct _n_User { 22c4762a1bSJed Brown TS ts; 23c4762a1bSJed Brown PetscReal mu; 24c4762a1bSJed Brown PetscReal next_output; 25c4762a1bSJed Brown 26c4762a1bSJed Brown /* Sensitivity analysis support */ 27c4762a1bSJed Brown PetscInt steps; 28c4762a1bSJed Brown PetscReal ftime; 29c4762a1bSJed Brown Mat A; /* Jacobian matrix for ODE */ 30c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix for ODE*/ 31c4762a1bSJed Brown Mat H; /* Hessian matrix for optimization */ 32c4762a1bSJed Brown Vec U,Lambda[1],Mup[1]; /* first-order adjoint variables */ 33c4762a1bSJed Brown Vec Lambda2[2]; /* second-order adjoint variables */ 34c4762a1bSJed Brown Vec Ihp1[1]; /* working space for Hessian evaluations */ 35c4762a1bSJed Brown Vec Dir; /* direction vector */ 36c4762a1bSJed Brown PetscReal ob[2]; /* observation used by the cost function */ 37c4762a1bSJed Brown PetscBool implicitform; /* implicit ODE? */ 38c4762a1bSJed Brown }; 39c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec,PetscScalar[],User); 40c4762a1bSJed Brown 41c4762a1bSJed Brown /* ----------------------- Explicit form of the ODE -------------------- */ 42c4762a1bSJed Brown 43c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 44c4762a1bSJed Brown { 45c4762a1bSJed Brown User user = (User)ctx; 46c4762a1bSJed Brown PetscScalar *f; 47c4762a1bSJed Brown const PetscScalar *u; 48c4762a1bSJed Brown 49c4762a1bSJed Brown PetscFunctionBeginUser; 505f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 515f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 52c4762a1bSJed Brown f[0] = u[1]; 53c4762a1bSJed Brown f[1] = user->mu*((1.-u[0]*u[0])*u[1]-u[0]); 545f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 555f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 56c4762a1bSJed Brown PetscFunctionReturn(0); 57c4762a1bSJed Brown } 58c4762a1bSJed Brown 59c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 60c4762a1bSJed Brown { 61c4762a1bSJed Brown User user = (User)ctx; 62c4762a1bSJed Brown PetscReal mu = user->mu; 63c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 64c4762a1bSJed Brown PetscScalar J[2][2]; 65c4762a1bSJed Brown const PetscScalar *u; 66c4762a1bSJed Brown 67c4762a1bSJed Brown PetscFunctionBeginUser; 685f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 69c4762a1bSJed Brown J[0][0] = 0; 70c4762a1bSJed Brown J[1][0] = -mu*(2.0*u[1]*u[0]+1.); 71c4762a1bSJed Brown J[0][1] = 1.0; 72c4762a1bSJed Brown J[1][1] = mu*(1.0-u[0]*u[0]); 735f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 745f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 755f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 76c4762a1bSJed Brown if (A != B) { 775f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 785f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 79c4762a1bSJed Brown } 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 81c4762a1bSJed Brown PetscFunctionReturn(0); 82c4762a1bSJed Brown } 83c4762a1bSJed Brown 84c4762a1bSJed Brown static PetscErrorCode RHSHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 85c4762a1bSJed Brown { 86c4762a1bSJed Brown const PetscScalar *vl,*vr,*u; 87c4762a1bSJed Brown PetscScalar *vhv; 88c4762a1bSJed Brown PetscScalar dJdU[2][2][2]={{{0}}}; 89c4762a1bSJed Brown PetscInt i,j,k; 90c4762a1bSJed Brown User user = (User)ctx; 91c4762a1bSJed Brown 92c4762a1bSJed Brown PetscFunctionBeginUser; 935f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 945f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Vl[0],&vl)); 955f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Vr,&vr)); 965f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(VHV[0],&vhv)); 97c4762a1bSJed Brown 98c4762a1bSJed Brown dJdU[1][0][0] = -2.*user->mu*u[1]; 99c4762a1bSJed Brown dJdU[1][1][0] = -2.*user->mu*u[0]; 100c4762a1bSJed Brown dJdU[1][0][1] = -2.*user->mu*u[0]; 101c4762a1bSJed Brown for (j=0;j<2;j++) { 102c4762a1bSJed Brown vhv[j] = 0; 103c4762a1bSJed Brown for (k=0;k<2;k++) 104c4762a1bSJed Brown for (i=0;i<2;i++) 105c4762a1bSJed Brown vhv[j] += vl[i]*dJdU[i][j][k]*vr[k]; 106c4762a1bSJed Brown } 1075f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Vl[0],&vl)); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Vr,&vr)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(VHV[0],&vhv)); 111c4762a1bSJed Brown PetscFunctionReturn(0); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* ----------------------- Implicit form of the ODE -------------------- */ 115c4762a1bSJed Brown 116c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,void *ctx) 117c4762a1bSJed Brown { 118c4762a1bSJed Brown User user = (User)ctx; 119c4762a1bSJed Brown const PetscScalar *u,*udot; 120c4762a1bSJed Brown PetscScalar *f; 121c4762a1bSJed Brown 122c4762a1bSJed Brown PetscFunctionBeginUser; 1235f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 1245f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Udot,&udot)); 1255f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(F,&f)); 126c4762a1bSJed Brown f[0] = udot[0] - u[1]; 127c4762a1bSJed Brown f[1] = udot[1] - user->mu*((1.0-u[0]*u[0])*u[1] - u[0]) ; 1285f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 1295f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Udot,&udot)); 1305f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(F,&f)); 131c4762a1bSJed Brown PetscFunctionReturn(0); 132c4762a1bSJed Brown } 133c4762a1bSJed Brown 134c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,void *ctx) 135c4762a1bSJed Brown { 136c4762a1bSJed Brown User user = (User)ctx; 137c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 138c4762a1bSJed Brown PetscScalar J[2][2]; 139c4762a1bSJed Brown const PetscScalar *u; 140c4762a1bSJed Brown 141c4762a1bSJed Brown PetscFunctionBeginUser; 1425f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 143c4762a1bSJed Brown J[0][0] = a; J[0][1] = -1.0; 144c4762a1bSJed Brown J[1][0] = user->mu*(1.0 + 2.0*u[0]*u[1]); J[1][1] = a - user->mu*(1.0-u[0]*u[0]); 1455f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 1465f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 1475f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 1485f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 149c4762a1bSJed Brown if (A != B) { 1505f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 152c4762a1bSJed Brown } 153c4762a1bSJed Brown PetscFunctionReturn(0); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown 156c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */ 157c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx) 158c4762a1bSJed Brown { 159c4762a1bSJed Brown const PetscScalar *u; 160c4762a1bSJed Brown PetscReal tfinal, dt; 161c4762a1bSJed Brown User user = (User)ctx; 162c4762a1bSJed Brown Vec interpolatedU; 163c4762a1bSJed Brown 164c4762a1bSJed Brown PetscFunctionBeginUser; 1655f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetTimeStep(ts,&dt)); 1665f80ce2aSJacob Faibussowitsch CHKERRQ(TSGetMaxTime(ts,&tfinal)); 167c4762a1bSJed Brown 168c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) { 1695f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&interpolatedU)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(TSInterpolate(ts,user->next_output,interpolatedU)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(interpolatedU,&u)); 172*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"[%g] %D TS %g (dt = %g) X %g %g\n", 173c4762a1bSJed Brown (double)user->next_output,step,(double)t,(double)dt,(double)PetscRealPart(u[0]), 174*b122ec5aSJacob Faibussowitsch (double)PetscRealPart(u[1]))); 1755f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(interpolatedU,&u)); 1765f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&interpolatedU)); 177c4762a1bSJed Brown user->next_output += 0.1; 178c4762a1bSJed Brown } 179c4762a1bSJed Brown PetscFunctionReturn(0); 180c4762a1bSJed Brown } 181c4762a1bSJed Brown 182c4762a1bSJed Brown static PetscErrorCode IHessianProductUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx) 183c4762a1bSJed Brown { 184c4762a1bSJed Brown const PetscScalar *vl,*vr,*u; 185c4762a1bSJed Brown PetscScalar *vhv; 186c4762a1bSJed Brown PetscScalar dJdU[2][2][2]={{{0}}}; 187c4762a1bSJed Brown PetscInt i,j,k; 188c4762a1bSJed Brown User user = (User)ctx; 189c4762a1bSJed Brown 190c4762a1bSJed Brown PetscFunctionBeginUser; 1915f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 1925f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Vl[0],&vl)); 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(Vr,&vr)); 1945f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(VHV[0],&vhv)); 195c4762a1bSJed Brown dJdU[1][0][0] = 2.*user->mu*u[1]; 196c4762a1bSJed Brown dJdU[1][1][0] = 2.*user->mu*u[0]; 197c4762a1bSJed Brown dJdU[1][0][1] = 2.*user->mu*u[0]; 198c4762a1bSJed Brown for (j=0;j<2;j++) { 199c4762a1bSJed Brown vhv[j] = 0; 200c4762a1bSJed Brown for (k=0;k<2;k++) 201c4762a1bSJed Brown for (i=0;i<2;i++) 202c4762a1bSJed Brown vhv[j] += vl[i]*dJdU[i][j][k]*vr[k]; 203c4762a1bSJed Brown } 2045f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Vl[0],&vl)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(Vr,&vr)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(VHV[0],&vhv)); 208c4762a1bSJed Brown PetscFunctionReturn(0); 209c4762a1bSJed Brown } 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* ------------------ User-defined routines for TAO -------------------------- */ 212c4762a1bSJed Brown 213c4762a1bSJed Brown static PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx) 214c4762a1bSJed Brown { 215c4762a1bSJed Brown User user_ptr = (User)ctx; 216c4762a1bSJed Brown TS ts = user_ptr->ts; 217c4762a1bSJed Brown const PetscScalar *x_ptr; 218c4762a1bSJed Brown PetscScalar *y_ptr; 219c4762a1bSJed Brown 220c4762a1bSJed Brown PetscFunctionBeginUser; 2215f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(IC,user_ptr->U)); /* set up the initial condition */ 222c4762a1bSJed Brown 2235f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(ts,0.0)); 2245f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetStepNumber(ts,0)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,0.001)); /* can be overwritten by command line options */ 2265f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,user_ptr->U)); 228c4762a1bSJed Brown 2295f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(user_ptr->U,&x_ptr)); 2305f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user_ptr->Lambda[0],&y_ptr)); 231c4762a1bSJed Brown *f = (x_ptr[0]-user_ptr->ob[0])*(x_ptr[0]-user_ptr->ob[0])+(x_ptr[1]-user_ptr->ob[1])*(x_ptr[1]-user_ptr->ob[1]); 232c4762a1bSJed Brown y_ptr[0] = 2.*(x_ptr[0]-user_ptr->ob[0]); 233c4762a1bSJed Brown y_ptr[1] = 2.*(x_ptr[1]-user_ptr->ob[1]); 2345f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user_ptr->Lambda[0],&y_ptr)); 2355f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(user_ptr->U,&x_ptr)); 236c4762a1bSJed Brown 2375f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetCostGradients(ts,1,user_ptr->Lambda,NULL)); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointSolve(ts)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(user_ptr->Lambda[0],G)); 240c4762a1bSJed Brown PetscFunctionReturn(0); 241c4762a1bSJed Brown } 242c4762a1bSJed Brown 243c4762a1bSJed Brown static PetscErrorCode FormHessian(Tao tao,Vec U,Mat H,Mat Hpre,void *ctx) 244c4762a1bSJed Brown { 245c4762a1bSJed Brown User user_ptr = (User)ctx; 246c4762a1bSJed Brown PetscScalar harr[2]; 247c4762a1bSJed Brown PetscScalar *x_ptr; 248c4762a1bSJed Brown const PetscInt rows[2] = {0,1}; 249c4762a1bSJed Brown PetscInt col; 250c4762a1bSJed Brown 251c4762a1bSJed Brown PetscFunctionBeginUser; 2525f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,user_ptr->U)); 2535f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user_ptr->Dir,&x_ptr)); 254c4762a1bSJed Brown x_ptr[0] = 1.; 255c4762a1bSJed Brown x_ptr[1] = 0.; 2565f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user_ptr->Dir,&x_ptr)); 2575f80ce2aSJacob Faibussowitsch CHKERRQ(Adjoint2(user_ptr->U,harr,user_ptr)); 258c4762a1bSJed Brown col = 0; 2595f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(H,2,rows,1,&col,harr,INSERT_VALUES)); 260c4762a1bSJed Brown 2615f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,user_ptr->U)); 2625f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user_ptr->Dir,&x_ptr)); 263c4762a1bSJed Brown x_ptr[0] = 0.; 264c4762a1bSJed Brown x_ptr[1] = 1.; 2655f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user_ptr->Dir,&x_ptr)); 2665f80ce2aSJacob Faibussowitsch CHKERRQ(Adjoint2(user_ptr->U,harr,user_ptr)); 267c4762a1bSJed Brown col = 1; 2685f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(H,2,rows,1,&col,harr,INSERT_VALUES)); 269c4762a1bSJed Brown 2705f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(H,MAT_FINAL_ASSEMBLY)); 2715f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(H,MAT_FINAL_ASSEMBLY)); 272c4762a1bSJed Brown if (H != Hpre) { 2735f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Hpre,MAT_FINAL_ASSEMBLY)); 2745f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Hpre,MAT_FINAL_ASSEMBLY)); 275c4762a1bSJed Brown } 276c4762a1bSJed Brown PetscFunctionReturn(0); 277c4762a1bSJed Brown } 278c4762a1bSJed Brown 279c4762a1bSJed Brown static PetscErrorCode MatrixFreeHessian(Tao tao,Vec U,Mat H,Mat Hpre,void *ctx) 280c4762a1bSJed Brown { 281c4762a1bSJed Brown User user_ptr = (User)ctx; 282c4762a1bSJed Brown 283c4762a1bSJed Brown PetscFunctionBeginUser; 2845f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,user_ptr->U)); 285c4762a1bSJed Brown PetscFunctionReturn(0); 286c4762a1bSJed Brown } 287c4762a1bSJed Brown 288c4762a1bSJed Brown /* ------------ Routines calculating second-order derivatives -------------- */ 289c4762a1bSJed Brown 2900e3d61c9SBarry Smith /* 291c4762a1bSJed Brown Compute the Hessian-vector product for the cost function using Second-order adjoint 292c4762a1bSJed Brown */ 293c4762a1bSJed Brown PetscErrorCode Adjoint2(Vec U,PetscScalar arr[],User ctx) 294c4762a1bSJed Brown { 295c4762a1bSJed Brown TS ts = ctx->ts; 296c4762a1bSJed Brown PetscScalar *x_ptr,*y_ptr; 297c4762a1bSJed Brown Mat tlmsen; 298c4762a1bSJed Brown 299c4762a1bSJed Brown PetscFunctionBeginUser; 3005f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointReset(ts)); 301c4762a1bSJed Brown 3025f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTime(ts,0.0)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetStepNumber(ts,0)); 3045f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,0.001)); 3055f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 3065f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetCostHessianProducts(ts,1,ctx->Lambda2,NULL,ctx->Dir)); 3075f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointSetForward(ts,NULL)); 3085f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 309c4762a1bSJed Brown 310c4762a1bSJed Brown /* Set terminal conditions for first- and second-order adjonts */ 3115f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&x_ptr)); 3125f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ctx->Lambda[0],&y_ptr)); 313c4762a1bSJed Brown y_ptr[0] = 2.*(x_ptr[0]-ctx->ob[0]); 314c4762a1bSJed Brown y_ptr[1] = 2.*(x_ptr[1]-ctx->ob[1]); 3155f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ctx->Lambda[0],&y_ptr)); 3165f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&x_ptr)); 3175f80ce2aSJacob Faibussowitsch CHKERRQ(TSForwardGetSensitivities(ts,NULL,&tlmsen)); 3185f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumn(tlmsen,0,&x_ptr)); 3195f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ctx->Lambda2[0],&y_ptr)); 320c4762a1bSJed Brown y_ptr[0] = 2.*x_ptr[0]; 321c4762a1bSJed Brown y_ptr[1] = 2.*x_ptr[1]; 3225f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&y_ptr)); 3235f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumn(tlmsen,&x_ptr)); 324c4762a1bSJed Brown 3255f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetCostGradients(ts,1,ctx->Lambda,NULL)); 326c4762a1bSJed Brown if (ctx->implicitform) { 3275f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIHessianProduct(ts,ctx->Ihp1,IHessianProductUU,NULL,NULL,NULL,NULL,NULL,NULL,ctx)); 328c4762a1bSJed Brown } else { 3295f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSHessianProduct(ts,ctx->Ihp1,RHSHessianProductUU,NULL,NULL,NULL,NULL,NULL,NULL,ctx)); 330c4762a1bSJed Brown } 3315f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointSolve(ts)); 332c4762a1bSJed Brown 3335f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(ctx->Lambda2[0],&x_ptr)); 334c4762a1bSJed Brown arr[0] = x_ptr[0]; 335c4762a1bSJed Brown arr[1] = x_ptr[1]; 3365f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(ctx->Lambda2[0],&x_ptr)); 337c4762a1bSJed Brown 3385f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointReset(ts)); 3395f80ce2aSJacob Faibussowitsch CHKERRQ(TSAdjointResetForward(ts)); 340c4762a1bSJed Brown PetscFunctionReturn(0); 341c4762a1bSJed Brown } 342c4762a1bSJed Brown 343c4762a1bSJed Brown PetscErrorCode FiniteDiff(Vec U,PetscScalar arr[],User ctx) 344c4762a1bSJed Brown { 345c4762a1bSJed Brown Vec Up,G,Gp; 346c4762a1bSJed Brown const PetscScalar eps = PetscRealConstant(1e-7); 347c4762a1bSJed Brown PetscScalar *u; 348c4762a1bSJed Brown Tao tao = NULL; 349c4762a1bSJed Brown PetscReal f; 350c4762a1bSJed Brown 351c4762a1bSJed Brown PetscFunctionBeginUser; 3525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&Up)); 3535f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&G)); 3545f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(U,&Gp)); 355c4762a1bSJed Brown 3565f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunctionGradient(tao,U,&f,G,ctx)); 357c4762a1bSJed Brown 3585f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,Up)); 3595f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Up,&u)); 360c4762a1bSJed Brown u[0] += eps; 3615f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Up,&u)); 3625f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunctionGradient(tao,Up,&f,Gp,ctx)); 3635f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(Gp,-1,G)); 3645f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(Gp,1./eps)); 3655f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Gp,&u)); 366c4762a1bSJed Brown arr[0] = u[0]; 367c4762a1bSJed Brown arr[1] = u[1]; 3685f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Gp,&u)); 369c4762a1bSJed Brown 3705f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(U,Up)); 3715f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Up,&u)); 372c4762a1bSJed Brown u[1] += eps; 3735f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Up,&u)); 3745f80ce2aSJacob Faibussowitsch CHKERRQ(FormFunctionGradient(tao,Up,&f,Gp,ctx)); 3755f80ce2aSJacob Faibussowitsch CHKERRQ(VecAXPY(Gp,-1,G)); 3765f80ce2aSJacob Faibussowitsch CHKERRQ(VecScale(Gp,1./eps)); 3775f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(Gp,&u)); 378c4762a1bSJed Brown arr[2] = u[0]; 379c4762a1bSJed Brown arr[3] = u[1]; 3805f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(Gp,&u)); 381c4762a1bSJed Brown 3825f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&G)); 3835f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Gp)); 3845f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&Up)); 385c4762a1bSJed Brown PetscFunctionReturn(0); 386c4762a1bSJed Brown } 387c4762a1bSJed Brown 388c4762a1bSJed Brown static PetscErrorCode HessianProductMat(Mat mat,Vec svec,Vec y) 389c4762a1bSJed Brown { 390c4762a1bSJed Brown User user_ptr; 391c4762a1bSJed Brown PetscScalar *y_ptr; 392c4762a1bSJed Brown 393c4762a1bSJed Brown PetscFunctionBeginUser; 3945f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellGetContext(mat,&user_ptr)); 3955f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(svec,user_ptr->Dir)); 3965f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(y,&y_ptr)); 3975f80ce2aSJacob Faibussowitsch CHKERRQ(Adjoint2(user_ptr->U,y_ptr,user_ptr)); 3985f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(y,&y_ptr)); 399c4762a1bSJed Brown PetscFunctionReturn(0); 400c4762a1bSJed Brown } 401c4762a1bSJed Brown 402c4762a1bSJed Brown int main(int argc,char **argv) 403c4762a1bSJed Brown { 404c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE,mf = PETSC_TRUE; 405c4762a1bSJed Brown PetscInt mode = 0; 406c4762a1bSJed Brown PetscMPIInt size; 407c4762a1bSJed Brown struct _n_User user; 408c4762a1bSJed Brown Vec x; /* working vector for TAO */ 409c4762a1bSJed Brown PetscScalar *x_ptr,arr[4]; 410c4762a1bSJed Brown PetscScalar ic1 = 2.2,ic2 = -0.7; /* initial guess for TAO */ 411c4762a1bSJed Brown Tao tao; 412c4762a1bSJed Brown KSP ksp; 413c4762a1bSJed Brown PC pc; 414c4762a1bSJed Brown 415c4762a1bSJed Brown /* Initialize program */ 416*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,NULL,help)); 4175f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 4183c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 419c4762a1bSJed Brown 420c4762a1bSJed Brown /* Set runtime options */ 421c4762a1bSJed Brown user.next_output = 0.0; 422c4762a1bSJed Brown user.mu = 1.0e3; 423c4762a1bSJed Brown user.steps = 0; 424c4762a1bSJed Brown user.ftime = 0.5; 425c4762a1bSJed Brown user.implicitform = PETSC_TRUE; 426c4762a1bSJed Brown 4275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL)); 4285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL)); 4295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-mode",&mode,NULL)); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ic1",&ic1,NULL)); 4315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-ic2",&ic2,NULL)); 4325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-my_tao_mf",&mf,NULL)); /* matrix-free hessian for optimization */ 4335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-implicitform",&user.implicitform,NULL)); 434c4762a1bSJed Brown 435c4762a1bSJed Brown /* Create necessary matrix and vectors, solve same ODE on every process */ 4365f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.A)); 4375f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2)); 4385f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(user.A)); 4395f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user.A)); 4405f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(user.A,&user.U,NULL)); 4415f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(user.A,&user.Dir,NULL)); 4425f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(user.A,&user.Lambda[0],NULL)); 4435f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(user.A,&user.Lambda2[0],NULL)); 4445f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(user.A,&user.Ihp1[0],NULL)); 445c4762a1bSJed Brown 446c4762a1bSJed Brown /* Create timestepping solver context */ 4475f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&user.ts)); 4485f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetEquationType(user.ts,TS_EQ_ODE_EXPLICIT)); /* less Jacobian evaluations when adjoint BEuler is used, otherwise no effect */ 449c4762a1bSJed Brown if (user.implicitform) { 4505f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIFunction(user.ts,NULL,IFunction,&user)); 4515f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetIJacobian(user.ts,user.A,user.A,IJacobian,&user)); 4525f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(user.ts,TSCN)); 453c4762a1bSJed Brown } else { 4545f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(user.ts,NULL,RHSFunction,&user)); 4555f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(user.ts,user.A,user.A,RHSJacobian,&user)); 4565f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(user.ts,TSRK)); 457c4762a1bSJed Brown } 4585f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(user.ts,user.ftime)); 4595f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(user.ts,TS_EXACTFINALTIME_MATCHSTEP)); 460c4762a1bSJed Brown 461c4762a1bSJed Brown if (monitor) { 4625f80ce2aSJacob Faibussowitsch CHKERRQ(TSMonitorSet(user.ts,Monitor,&user,NULL)); 463c4762a1bSJed Brown } 464c4762a1bSJed Brown 465c4762a1bSJed Brown /* Set ODE initial conditions */ 4665f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user.U,&x_ptr)); 467c4762a1bSJed Brown x_ptr[0] = 2.0; 468c4762a1bSJed Brown x_ptr[1] = -2.0/3.0 + 10.0/(81.0*user.mu) - 292.0/(2187.0*user.mu*user.mu); 4695f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user.U,&x_ptr)); 470c4762a1bSJed Brown 471c4762a1bSJed Brown /* Set runtime options */ 4725f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(user.ts)); 473c4762a1bSJed Brown 474c4762a1bSJed Brown /* Obtain the observation */ 4755f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(user.ts,user.U)); 4765f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user.U,&x_ptr)); 477c4762a1bSJed Brown user.ob[0] = x_ptr[0]; 478c4762a1bSJed Brown user.ob[1] = x_ptr[1]; 4795f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user.U,&x_ptr)); 480c4762a1bSJed Brown 4815f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(user.U,&x)); 4825f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(x,&x_ptr)); 483c4762a1bSJed Brown x_ptr[0] = ic1; 484c4762a1bSJed Brown x_ptr[1] = ic2; 4855f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(x,&x_ptr)); 486c4762a1bSJed Brown 487c4762a1bSJed Brown /* Save trajectory of solution so that TSAdjointSolve() may be used */ 4885f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSaveTrajectory(user.ts)); 489c4762a1bSJed Brown 490c4762a1bSJed Brown /* Compare finite difference and second-order adjoint. */ 491c4762a1bSJed Brown switch (mode) { 492c4762a1bSJed Brown case 2 : 4935f80ce2aSJacob Faibussowitsch CHKERRQ(FiniteDiff(x,arr,&user)); 4945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n Finite difference approximation of the Hessian\n")); 4955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%g %g\n%g %g\n",(double)arr[0],(double)arr[1],(double)arr[2],(double)arr[3])); 496c4762a1bSJed Brown break; 497c4762a1bSJed Brown case 1 : /* Compute the Hessian column by column */ 4985f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,user.U)); 4995f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user.Dir,&x_ptr)); 500c4762a1bSJed Brown x_ptr[0] = 1.; 501c4762a1bSJed Brown x_ptr[1] = 0.; 5025f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user.Dir,&x_ptr)); 5035f80ce2aSJacob Faibussowitsch CHKERRQ(Adjoint2(user.U,arr,&user)); 5045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nFirst column of the Hessian\n")); 5055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%g\n%g\n",(double)arr[0],(double)arr[1])); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(VecCopy(x,user.U)); 5075f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(user.Dir,&x_ptr)); 508c4762a1bSJed Brown x_ptr[0] = 0.; 509c4762a1bSJed Brown x_ptr[1] = 1.; 5105f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(user.Dir,&x_ptr)); 5115f80ce2aSJacob Faibussowitsch CHKERRQ(Adjoint2(user.U,arr,&user)); 5125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\nSecond column of the Hessian\n")); 5135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"%g\n%g\n",(double)arr[0],(double)arr[1])); 514c4762a1bSJed Brown break; 515c4762a1bSJed Brown case 0 : 516c4762a1bSJed Brown /* Create optimization context and set up */ 5175f80ce2aSJacob Faibussowitsch CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao)); 5185f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetType(tao,TAOBLMVM)); 5195f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&user)); 520c4762a1bSJed Brown 521c4762a1bSJed Brown if (mf) { 5225f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateShell(PETSC_COMM_SELF,2,2,2,2,(void*)&user,&user.H)); 5235f80ce2aSJacob Faibussowitsch CHKERRQ(MatShellSetOperation(user.H,MATOP_MULT,(void(*)(void))HessianProductMat)); 5245f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(user.H,MAT_SYMMETRIC,PETSC_TRUE)); 5255f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,user.H,user.H,MatrixFreeHessian,(void *)&user)); 526c4762a1bSJed Brown } else { /* Create Hessian matrix */ 5275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.H)); 5285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(user.H,PETSC_DECIDE,PETSC_DECIDE,2,2)); 5295f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(user.H)); 5305f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetHessian(tao,user.H,user.H,FormHessian,(void *)&user)); 531c4762a1bSJed Brown } 532c4762a1bSJed Brown 533c4762a1bSJed Brown /* Not use any preconditioner */ 5345f80ce2aSJacob Faibussowitsch CHKERRQ(TaoGetKSP(tao,&ksp)); 535c4762a1bSJed Brown if (ksp) { 5365f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,&pc)); 5375f80ce2aSJacob Faibussowitsch CHKERRQ(PCSetType(pc,PCNONE)); 538c4762a1bSJed Brown } 539c4762a1bSJed Brown 540c4762a1bSJed Brown /* Set initial solution guess */ 5415f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetSolution(tao,x)); 5425f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSetFromOptions(tao)); 5435f80ce2aSJacob Faibussowitsch CHKERRQ(TaoSolve(tao)); 5445f80ce2aSJacob Faibussowitsch CHKERRQ(TaoDestroy(&tao)); 5455f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.H)); 546c4762a1bSJed Brown break; 547c4762a1bSJed Brown default: 548c4762a1bSJed Brown break; 549c4762a1bSJed Brown } 550c4762a1bSJed Brown 551c4762a1bSJed Brown /* Free work space. All PETSc objects should be destroyed when they are no longer needed. */ 5525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&user.A)); 5535f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.U)); 5545f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.Lambda[0])); 5555f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.Lambda2[0])); 5565f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.Ihp1[0])); 5575f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&user.Dir)); 5585f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&user.ts)); 5595f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 560*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 561*b122ec5aSJacob Faibussowitsch return 0; 562c4762a1bSJed Brown } 563c4762a1bSJed Brown 564c4762a1bSJed Brown /*TEST 565c4762a1bSJed Brown build: 566c4762a1bSJed Brown requires: !complex !single 567c4762a1bSJed Brown 568c4762a1bSJed Brown test: 569c4762a1bSJed Brown args: -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 1000 -ts_dt 0.03125 570c4762a1bSJed Brown output_file: output/ex20opt_ic_1.out 571c4762a1bSJed Brown 572c4762a1bSJed Brown test: 573c4762a1bSJed Brown suffix: 2 574c4762a1bSJed Brown args: -ts_type beuler -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none 575c4762a1bSJed Brown output_file: output/ex20opt_ic_2.out 576c4762a1bSJed Brown 577c4762a1bSJed Brown test: 578c4762a1bSJed Brown suffix: 3 579c4762a1bSJed Brown args: -ts_type cn -viewer_binary_skip_info -tao_monitor -tao_view -mu 100 -ts_dt 0.01 -tao_type bntr -tao_bnk_pc_type none 580c4762a1bSJed Brown output_file: output/ex20opt_ic_3.out 581c4762a1bSJed Brown 582c4762a1bSJed Brown test: 583c4762a1bSJed Brown suffix: 4 584c4762a1bSJed Brown args: -implicitform 0 -ts_dt 0.01 -ts_max_steps 2 -ts_rhs_jacobian_test_mult_transpose -mat_shell_test_mult_transpose_view -ts_rhs_jacobian_test_mult -mat_shell_test_mult_view -mode 1 -my_tao_mf 585c4762a1bSJed Brown TEST*/ 586