xref: /petsc/src/ts/tutorials/ex23fwdadj.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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;
365f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Xdot,&xdot));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(F,&f));
3982ad9101SHong Zhang   f[0] = user->c*xdot[0] - user->b*x[0];
405f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(Xdot,&xdot));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(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;
545f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
5582ad9101SHong Zhang   J[0][0] = user->c*a - user->b*1.0;
565f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(B,1,rowcol,1,rowcol,&J[0][0],INSERT_VALUES));
575f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
5882ad9101SHong Zhang 
595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
605f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
6182ad9101SHong Zhang   if (A != B) {
625f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
635f80ce2aSJacob Faibussowitsch     CHKERRQ(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;
775f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
785f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(Xdot,&xdot));
795f80ce2aSJacob Faibussowitsch   CHKERRQ(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];
825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,1,row,1,col,&J[0][0],INSERT_VALUES));
835f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
8482ad9101SHong Zhang 
855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
865f80ce2aSJacob Faibussowitsch   CHKERRQ(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*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
9982ad9101SHong Zhang 
1005f80ce2aSJacob Faibussowitsch   CHKERRMPI(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;
1095f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL,"-der",&user.der,NULL));
11082ad9101SHong Zhang 
11182ad9101SHong Zhang   rows = 1;
11282ad9101SHong Zhang   cols = 1;
1135f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jac));
1145f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,1,1));
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.Jac));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.Jac));
1175f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jac,&user.x,NULL));
11882ad9101SHong Zhang 
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1205f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSBEULER));
1215f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,IFunction,&user));
1225f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user));
1235f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
1245f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,user.ftime));
12582ad9101SHong Zhang 
1265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(user.x,&x_ptr));
12782ad9101SHong Zhang   x_ptr[0] = user.a;
1285f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(user.x,&x_ptr));
1295f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,0.001));
13082ad9101SHong Zhang 
13182ad9101SHong Zhang   /* Set up forward sensitivity */
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(user.Jacp));
1355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(user.Jacp));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(user.sp));
1385f80ce2aSJacob Faibussowitsch   CHKERRQ(TSForwardSetSensitivities(ts,cols,user.sp));
1395f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobianP(ts,user.Jacp,IJacobianP,&user));
14082ad9101SHong Zhang 
1415f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSaveTrajectory(ts));
1425f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
14382ad9101SHong Zhang 
1445f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,user.x));
1455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&user.ftime));
1465f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&user.steps));
1475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(user.x,&x_ptr));
1485f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n ode solution %g\n",(double)PetscRealPart(x_ptr[0])));
1495f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(user.x,&x_ptr));
1505f80ce2aSJacob Faibussowitsch   CHKERRQ(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) {
1535f80ce2aSJacob Faibussowitsch     CHKERRQ(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) {
1565f80ce2aSJacob Faibussowitsch     CHKERRQ(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   }
1585f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity:\n"));
1595f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD));
16082ad9101SHong Zhang 
1615f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jac,&user.lambda[0],NULL));
16282ad9101SHong Zhang   /* Set initial conditions for the adjoint integration */
1635f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(user.lambda[0],&x_ptr));
16482ad9101SHong Zhang   x_ptr[0] = 1.0;
1655f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(user.lambda[0],&x_ptr));
1665f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(user.Jacp,&user.mup[0],NULL));
1675f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayWrite(user.mup[0],&x_ptr));
16882ad9101SHong Zhang   x_ptr[0] = 0.0;
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayWrite(user.mup[0],&x_ptr));
17082ad9101SHong Zhang 
1715f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetCostGradients(ts,1,user.lambda,user.mup));
1725f80ce2aSJacob Faibussowitsch   CHKERRQ(TSAdjointSolve(ts));
17382ad9101SHong Zhang 
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n adjoint sensitivity:\n"));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(user.mup[0],PETSC_VIEWER_STDOUT_WORLD));
17682ad9101SHong Zhang 
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.Jac));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.sp));
1795f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&user.Jacp));
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.x));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.lambda[0]));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&user.mup[0]));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
18482ad9101SHong Zhang 
185*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
186*b122ec5aSJacob 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