182ad9101SHong Zhang static char help[] = "A toy example for testing forward and adjoint sensitivity analysis of an implicit ODE with a paramerized mass matrice.\n"; 282ad9101SHong Zhang 382ad9101SHong Zhang /* 482ad9101SHong Zhang This example solves the simple ODE 582ad9101SHong Zhang c x' = b x, x(0) = a, 682ad9101SHong Zhang whose analytical solution is x(T)=a*exp(b/c*T), and calculates the derivative of x(T) w.r.t. c (by default) or w.r.t. b (can be enabled with command line option -der 2). 782ad9101SHong Zhang 882ad9101SHong Zhang */ 982ad9101SHong Zhang 1082ad9101SHong Zhang #include <petscts.h> 1182ad9101SHong Zhang 1282ad9101SHong Zhang typedef struct _n_User *User; 1382ad9101SHong Zhang struct _n_User { 1482ad9101SHong Zhang PetscReal a; 1582ad9101SHong Zhang PetscReal b; 1682ad9101SHong Zhang PetscReal c; 1782ad9101SHong Zhang /* Sensitivity analysis support */ 1882ad9101SHong Zhang PetscInt steps; 1982ad9101SHong Zhang PetscReal ftime; 2082ad9101SHong Zhang Mat Jac; /* Jacobian matrix */ 2182ad9101SHong Zhang Mat Jacp; /* JacobianP matrix */ 2282ad9101SHong Zhang Vec x; 2382ad9101SHong Zhang Mat sp; /* forward sensitivity variables */ 2482ad9101SHong Zhang Vec lambda[1]; /* adjoint sensitivity variables */ 2582ad9101SHong Zhang Vec mup[1]; /* adjoint sensitivity variables */ 2682ad9101SHong Zhang PetscInt der; 2782ad9101SHong Zhang }; 2882ad9101SHong Zhang 2982ad9101SHong Zhang static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx) 3082ad9101SHong Zhang { 3182ad9101SHong Zhang User user = (User)ctx; 3282ad9101SHong Zhang const PetscScalar *x,*xdot; 3382ad9101SHong Zhang PetscScalar *f; 3482ad9101SHong Zhang 3582ad9101SHong Zhang PetscFunctionBeginUser; 36*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 37*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot,&xdot)); 38*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(F,&f)); 3982ad9101SHong Zhang f[0] = user->c*xdot[0] - user->b*x[0]; 40*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 41*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot,&xdot)); 42*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(F,&f)); 4382ad9101SHong Zhang PetscFunctionReturn(0); 4482ad9101SHong Zhang } 4582ad9101SHong Zhang 4682ad9101SHong Zhang static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 4782ad9101SHong Zhang { 4882ad9101SHong Zhang User user = (User)ctx; 4982ad9101SHong Zhang PetscInt rowcol[] = {0}; 5082ad9101SHong Zhang PetscScalar J[1][1]; 5182ad9101SHong Zhang const PetscScalar *x; 5282ad9101SHong Zhang 5382ad9101SHong Zhang PetscFunctionBeginUser; 54*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 5582ad9101SHong Zhang J[0][0] = user->c*a - user->b*1.0; 56*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(B,1,rowcol,1,rowcol,&J[0][0],INSERT_VALUES)); 57*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 5882ad9101SHong Zhang 59*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 60*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 6182ad9101SHong Zhang if (A != B) { 62*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 63*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 6482ad9101SHong Zhang } 6582ad9101SHong Zhang PetscFunctionReturn(0); 6682ad9101SHong Zhang } 6782ad9101SHong Zhang 6882ad9101SHong Zhang static PetscErrorCode IJacobianP(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat A,void *ctx) 6982ad9101SHong Zhang { 7082ad9101SHong Zhang User user = (User)ctx; 7182ad9101SHong Zhang PetscInt row[] = {0},col[]={0}; 7282ad9101SHong Zhang PetscScalar J[1][1]; 7382ad9101SHong Zhang const PetscScalar *x,*xdot; 7482ad9101SHong Zhang PetscReal dt; 7582ad9101SHong Zhang 7682ad9101SHong Zhang PetscFunctionBeginUser; 77*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X,&x)); 78*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot,&xdot)); 79*9566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts,&dt)); 8082ad9101SHong Zhang if (user->der == 1) J[0][0] = xdot[0]; 8182ad9101SHong Zhang if (user->der == 2) J[0][0] = -x[0]; 82*9566063dSJacob Faibussowitsch PetscCall(MatSetValues(A,1,row,1,col,&J[0][0],INSERT_VALUES)); 83*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X,&x)); 8482ad9101SHong Zhang 85*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 86*9566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 8782ad9101SHong Zhang PetscFunctionReturn(0); 8882ad9101SHong Zhang } 8982ad9101SHong Zhang 9082ad9101SHong Zhang int main(int argc,char **argv) 9182ad9101SHong Zhang { 9282ad9101SHong Zhang TS ts; 9382ad9101SHong Zhang PetscScalar *x_ptr; 9482ad9101SHong Zhang PetscMPIInt size; 9582ad9101SHong Zhang struct _n_User user; 9682ad9101SHong Zhang PetscInt rows,cols; 9782ad9101SHong Zhang 98*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,NULL,help)); 9982ad9101SHong Zhang 100*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1013c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 10282ad9101SHong Zhang 10382ad9101SHong Zhang user.a = 2.0; 10482ad9101SHong Zhang user.b = 4.0; 10582ad9101SHong Zhang user.c = 3.0; 10682ad9101SHong Zhang user.steps = 0; 10782ad9101SHong Zhang user.ftime = 1.0; 10882ad9101SHong Zhang user.der = 1; 109*9566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetInt(NULL,NULL,"-der",&user.der,NULL)); 11082ad9101SHong Zhang 11182ad9101SHong Zhang rows = 1; 11282ad9101SHong Zhang cols = 1; 113*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jac)); 114*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,1,1)); 115*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.Jac)); 116*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.Jac)); 117*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jac,&user.x,NULL)); 11882ad9101SHong Zhang 119*9566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 120*9566063dSJacob Faibussowitsch PetscCall(TSSetType(ts,TSBEULER)); 121*9566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts,NULL,IFunction,&user)); 122*9566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user)); 123*9566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 124*9566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts,user.ftime)); 12582ad9101SHong Zhang 126*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(user.x,&x_ptr)); 12782ad9101SHong Zhang x_ptr[0] = user.a; 128*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(user.x,&x_ptr)); 129*9566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts,0.001)); 13082ad9101SHong Zhang 13182ad9101SHong Zhang /* Set up forward sensitivity */ 132*9566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jacp)); 133*9566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols)); 134*9566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.Jacp)); 135*9566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.Jacp)); 136*9566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp)); 137*9566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(user.sp)); 138*9566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ts,cols,user.sp)); 139*9566063dSJacob Faibussowitsch PetscCall(TSSetIJacobianP(ts,user.Jacp,IJacobianP,&user)); 14082ad9101SHong Zhang 141*9566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 142*9566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 14382ad9101SHong Zhang 144*9566063dSJacob Faibussowitsch PetscCall(TSSolve(ts,user.x)); 145*9566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts,&user.ftime)); 146*9566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts,&user.steps)); 147*9566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.x,&x_ptr)); 148*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n ode solution %g\n",(double)PetscRealPart(x_ptr[0]))); 149*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.x,&x_ptr)); 150*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n analytical solution %g\n",(double)user.a*PetscExpReal(user.b/user.c*user.ftime))); 15182ad9101SHong Zhang 15282ad9101SHong Zhang if (user.der == 1) { 153*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n analytical derivative w.r.t. c %g\n",(double)-user.a*user.ftime*user.b/(user.c*user.c)*PetscExpReal(user.b/user.c*user.ftime))); 15482ad9101SHong Zhang } 15582ad9101SHong Zhang if (user.der == 2) { 156*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n analytical derivative w.r.t. b %g\n",user.a*user.ftime/user.c*PetscExpReal(user.b/user.c*user.ftime))); 15782ad9101SHong Zhang } 158*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity:\n")); 159*9566063dSJacob Faibussowitsch PetscCall(MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD)); 16082ad9101SHong Zhang 161*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jac,&user.lambda[0],NULL)); 16282ad9101SHong Zhang /* Set initial conditions for the adjoint integration */ 163*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(user.lambda[0],&x_ptr)); 16482ad9101SHong Zhang x_ptr[0] = 1.0; 165*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(user.lambda[0],&x_ptr)); 166*9566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp,&user.mup[0],NULL)); 167*9566063dSJacob Faibussowitsch PetscCall(VecGetArrayWrite(user.mup[0],&x_ptr)); 16882ad9101SHong Zhang x_ptr[0] = 0.0; 169*9566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(user.mup[0],&x_ptr)); 17082ad9101SHong Zhang 171*9566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts,1,user.lambda,user.mup)); 172*9566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 17382ad9101SHong Zhang 174*9566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n adjoint sensitivity:\n")); 175*9566063dSJacob Faibussowitsch PetscCall(VecView(user.mup[0],PETSC_VIEWER_STDOUT_WORLD)); 17682ad9101SHong Zhang 177*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Jac)); 178*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.sp)); 179*9566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Jacp)); 180*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.x)); 181*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.lambda[0])); 182*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.mup[0])); 183*9566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 18482ad9101SHong Zhang 185*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 186b122ec5aSJacob Faibussowitsch return 0; 18782ad9101SHong Zhang } 18882ad9101SHong Zhang 18982ad9101SHong Zhang /*TEST 19082ad9101SHong Zhang 19182ad9101SHong Zhang test: 19282ad9101SHong Zhang args: -ts_type beuler 19382ad9101SHong Zhang 19482ad9101SHong Zhang test: 19582ad9101SHong Zhang suffix: 2 19682ad9101SHong Zhang args: -ts_type cn 19782ad9101SHong Zhang output_file: output/ex23fwdadj_1.out 19882ad9101SHong Zhang 19982ad9101SHong Zhang TEST*/ 200