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 PetscErrorCode ierr; 3282ad9101SHong Zhang User user = (User)ctx; 3382ad9101SHong Zhang const PetscScalar *x,*xdot; 3482ad9101SHong Zhang PetscScalar *f; 3582ad9101SHong Zhang 3682ad9101SHong Zhang PetscFunctionBeginUser; 3782ad9101SHong Zhang ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 3882ad9101SHong Zhang ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 3982ad9101SHong Zhang ierr = VecGetArrayWrite(F,&f);CHKERRQ(ierr); 4082ad9101SHong Zhang f[0] = user->c*xdot[0] - user->b*x[0]; 4182ad9101SHong Zhang ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 4282ad9101SHong Zhang ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr); 4382ad9101SHong Zhang ierr = VecRestoreArrayWrite(F,&f);CHKERRQ(ierr); 4482ad9101SHong Zhang PetscFunctionReturn(0); 4582ad9101SHong Zhang } 4682ad9101SHong Zhang 4782ad9101SHong Zhang static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx) 4882ad9101SHong Zhang { 4982ad9101SHong Zhang PetscErrorCode ierr; 5082ad9101SHong Zhang User user = (User)ctx; 5182ad9101SHong Zhang PetscInt rowcol[] = {0}; 5282ad9101SHong Zhang PetscScalar J[1][1]; 5382ad9101SHong Zhang const PetscScalar *x; 5482ad9101SHong Zhang 5582ad9101SHong Zhang PetscFunctionBeginUser; 5682ad9101SHong Zhang ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 5782ad9101SHong Zhang J[0][0] = user->c*a - user->b*1.0; 5882ad9101SHong Zhang ierr = MatSetValues(B,1,rowcol,1,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 5982ad9101SHong Zhang ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 6082ad9101SHong Zhang 6182ad9101SHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6282ad9101SHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6382ad9101SHong Zhang if (A != B) { 6482ad9101SHong Zhang ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6582ad9101SHong Zhang ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 6682ad9101SHong Zhang } 6782ad9101SHong Zhang PetscFunctionReturn(0); 6882ad9101SHong Zhang } 6982ad9101SHong Zhang 7082ad9101SHong Zhang static PetscErrorCode IJacobianP(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat A,void *ctx) 7182ad9101SHong Zhang { 7282ad9101SHong Zhang User user = (User)ctx; 7382ad9101SHong Zhang PetscInt row[] = {0},col[]={0}; 7482ad9101SHong Zhang PetscScalar J[1][1]; 7582ad9101SHong Zhang const PetscScalar *x,*xdot; 7682ad9101SHong Zhang PetscReal dt; 7782ad9101SHong Zhang PetscErrorCode ierr; 7882ad9101SHong Zhang 7982ad9101SHong Zhang PetscFunctionBeginUser; 8082ad9101SHong Zhang ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr); 8182ad9101SHong Zhang ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr); 8282ad9101SHong Zhang ierr = TSGetTimeStep(ts,&dt);CHKERRQ(ierr); 8382ad9101SHong Zhang if (user->der == 1) J[0][0] = xdot[0]; 8482ad9101SHong Zhang if (user->der == 2) J[0][0] = -x[0]; 8582ad9101SHong Zhang ierr = MatSetValues(A,1,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 8682ad9101SHong Zhang ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr); 8782ad9101SHong Zhang 8882ad9101SHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 8982ad9101SHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 9082ad9101SHong Zhang PetscFunctionReturn(0); 9182ad9101SHong Zhang } 9282ad9101SHong Zhang 9382ad9101SHong Zhang int main(int argc,char **argv) 9482ad9101SHong Zhang { 9582ad9101SHong Zhang TS ts; 9682ad9101SHong Zhang PetscScalar *x_ptr; 9782ad9101SHong Zhang PetscMPIInt size; 9882ad9101SHong Zhang struct _n_User user; 9982ad9101SHong Zhang PetscInt rows,cols; 10082ad9101SHong Zhang PetscErrorCode ierr; 10182ad9101SHong Zhang 10282ad9101SHong Zhang ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 10382ad9101SHong Zhang 10455b25c41SPierre Jolivet ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 105*3c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 10682ad9101SHong Zhang 10782ad9101SHong Zhang user.a = 2.0; 10882ad9101SHong Zhang user.b = 4.0; 10982ad9101SHong Zhang user.c = 3.0; 11082ad9101SHong Zhang user.steps = 0; 11182ad9101SHong Zhang user.ftime = 1.0; 11282ad9101SHong Zhang user.der = 1; 11382ad9101SHong Zhang ierr = PetscOptionsGetInt(NULL,NULL,"-der",&user.der,NULL);CHKERRQ(ierr); 11482ad9101SHong Zhang 11582ad9101SHong Zhang rows = 1; 11682ad9101SHong Zhang cols = 1; 11782ad9101SHong Zhang ierr = MatCreate(PETSC_COMM_WORLD,&user.Jac);CHKERRQ(ierr); 11882ad9101SHong Zhang ierr = MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,1,1);CHKERRQ(ierr); 11982ad9101SHong Zhang ierr = MatSetFromOptions(user.Jac);CHKERRQ(ierr); 12082ad9101SHong Zhang ierr = MatSetUp(user.Jac);CHKERRQ(ierr); 12182ad9101SHong Zhang ierr = MatCreateVecs(user.Jac,&user.x,NULL);CHKERRQ(ierr); 12282ad9101SHong Zhang 12382ad9101SHong Zhang ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 12482ad9101SHong Zhang ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr); 12582ad9101SHong Zhang ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr); 12682ad9101SHong Zhang ierr = TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user);CHKERRQ(ierr); 12782ad9101SHong Zhang ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 12882ad9101SHong Zhang ierr = TSSetMaxTime(ts,user.ftime);CHKERRQ(ierr); 12982ad9101SHong Zhang 13082ad9101SHong Zhang ierr = VecGetArrayWrite(user.x,&x_ptr);CHKERRQ(ierr); 13182ad9101SHong Zhang x_ptr[0] = user.a; 13282ad9101SHong Zhang ierr = VecRestoreArrayWrite(user.x,&x_ptr);CHKERRQ(ierr); 13382ad9101SHong Zhang ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr); 13482ad9101SHong Zhang 13582ad9101SHong Zhang /* Set up forward sensitivity */ 13682ad9101SHong Zhang ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr); 13782ad9101SHong Zhang ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr); 13882ad9101SHong Zhang ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr); 13982ad9101SHong Zhang ierr = MatSetUp(user.Jacp);CHKERRQ(ierr); 14082ad9101SHong Zhang ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp);CHKERRQ(ierr); 14182ad9101SHong Zhang ierr = MatZeroEntries(user.sp);CHKERRQ(ierr); 14282ad9101SHong Zhang ierr = TSForwardSetSensitivities(ts,cols,user.sp);CHKERRQ(ierr); 14382ad9101SHong Zhang ierr = TSSetIJacobianP(ts,user.Jacp,IJacobianP,&user);CHKERRQ(ierr); 14482ad9101SHong Zhang 14582ad9101SHong Zhang ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 14682ad9101SHong Zhang ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 14782ad9101SHong Zhang 14882ad9101SHong Zhang ierr = TSSolve(ts,user.x);CHKERRQ(ierr); 14982ad9101SHong Zhang ierr = TSGetSolveTime(ts,&user.ftime);CHKERRQ(ierr); 15082ad9101SHong Zhang ierr = TSGetStepNumber(ts,&user.steps);CHKERRQ(ierr); 15182ad9101SHong Zhang ierr = VecGetArray(user.x,&x_ptr);CHKERRQ(ierr); 15282ad9101SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n ode solution %g\n",(double)PetscRealPart(x_ptr[0]));CHKERRQ(ierr); 15382ad9101SHong Zhang ierr = VecRestoreArray(user.x,&x_ptr);CHKERRQ(ierr); 15482ad9101SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n analytical solution %g\n",(double)user.a*PetscExpReal(user.b/user.c*user.ftime));CHKERRQ(ierr); 15582ad9101SHong Zhang 15682ad9101SHong Zhang if (user.der == 1) { 15782ad9101SHong Zhang ierr = 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));CHKERRQ(ierr); 15882ad9101SHong Zhang } 15982ad9101SHong Zhang if (user.der == 2) { 16082ad9101SHong Zhang ierr = 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));CHKERRQ(ierr); 16182ad9101SHong Zhang } 16282ad9101SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity:\n");CHKERRQ(ierr); 16382ad9101SHong Zhang ierr = MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 16482ad9101SHong Zhang 16582ad9101SHong Zhang ierr = MatCreateVecs(user.Jac,&user.lambda[0],NULL);CHKERRQ(ierr); 16682ad9101SHong Zhang /* Set initial conditions for the adjoint integration */ 16782ad9101SHong Zhang ierr = VecGetArrayWrite(user.lambda[0],&x_ptr);CHKERRQ(ierr); 16882ad9101SHong Zhang x_ptr[0] = 1.0; 16982ad9101SHong Zhang ierr = VecRestoreArrayWrite(user.lambda[0],&x_ptr);CHKERRQ(ierr); 17082ad9101SHong Zhang ierr = MatCreateVecs(user.Jacp,&user.mup[0],NULL);CHKERRQ(ierr); 17182ad9101SHong Zhang ierr = VecGetArrayWrite(user.mup[0],&x_ptr);CHKERRQ(ierr); 17282ad9101SHong Zhang x_ptr[0] = 0.0; 17382ad9101SHong Zhang ierr = VecRestoreArrayWrite(user.mup[0],&x_ptr);CHKERRQ(ierr); 17482ad9101SHong Zhang 17582ad9101SHong Zhang ierr = TSSetCostGradients(ts,1,user.lambda,user.mup);CHKERRQ(ierr); 17682ad9101SHong Zhang ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 17782ad9101SHong Zhang 17882ad9101SHong Zhang ierr = PetscPrintf(PETSC_COMM_WORLD,"\n adjoint sensitivity:\n");CHKERRQ(ierr); 17982ad9101SHong Zhang ierr = VecView(user.mup[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 18082ad9101SHong Zhang 18182ad9101SHong Zhang ierr = MatDestroy(&user.Jac);CHKERRQ(ierr); 18282ad9101SHong Zhang ierr = MatDestroy(&user.sp);CHKERRQ(ierr); 18382ad9101SHong Zhang ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr); 18482ad9101SHong Zhang ierr = VecDestroy(&user.x);CHKERRQ(ierr); 18582ad9101SHong Zhang ierr = VecDestroy(&user.lambda[0]);CHKERRQ(ierr); 18682ad9101SHong Zhang ierr = VecDestroy(&user.mup[0]);CHKERRQ(ierr); 18782ad9101SHong Zhang ierr = TSDestroy(&ts);CHKERRQ(ierr); 18882ad9101SHong Zhang 18982ad9101SHong Zhang ierr = PetscFinalize(); 19082ad9101SHong Zhang return(ierr); 19182ad9101SHong Zhang } 19282ad9101SHong Zhang 19382ad9101SHong Zhang /*TEST 19482ad9101SHong Zhang 19582ad9101SHong Zhang test: 19682ad9101SHong Zhang args: -ts_type beuler 19782ad9101SHong Zhang 19882ad9101SHong Zhang test: 19982ad9101SHong Zhang suffix: 2 20082ad9101SHong Zhang args: -ts_type cn 20182ad9101SHong Zhang output_file: output/ex23fwdadj_1.out 20282ad9101SHong Zhang 20382ad9101SHong Zhang TEST*/ 204