xref: /petsc/src/ts/tutorials/ex23fwdadj.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1baca6076SPierre Jolivet static char help[] = "A toy example for testing forward and adjoint sensitivity analysis of an implicit ODE with a parametrized 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 
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)29*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
30d71ae5a4SJacob Faibussowitsch {
3182ad9101SHong Zhang   User               user = (User)ctx;
3282ad9101SHong Zhang   const PetscScalar *x, *xdot;
3382ad9101SHong Zhang   PetscScalar       *f;
3482ad9101SHong Zhang 
3582ad9101SHong Zhang   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(F, &f));
3982ad9101SHong Zhang   f[0] = user->c * xdot[0] - user->b * x[0];
409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(F, &f));
433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4482ad9101SHong Zhang }
4582ad9101SHong Zhang 
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,PetscCtx ctx)46*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
47d71ae5a4SJacob Faibussowitsch {
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;
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
5582ad9101SHong Zhang   J[0][0] = user->c * a - user->b * 1.0;
569566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, rowcol, 1, rowcol, &J[0][0], INSERT_VALUES));
579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
5882ad9101SHong Zhang 
599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
609566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
6182ad9101SHong Zhang   if (A != B) {
629566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
639566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
6482ad9101SHong Zhang   }
653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6682ad9101SHong Zhang }
6782ad9101SHong Zhang 
IJacobianP(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat A,PetscCtx ctx)68*2a8381b2SBarry Smith static PetscErrorCode IJacobianP(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal shift, Mat A, PetscCtx ctx)
69d71ae5a4SJacob Faibussowitsch {
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;
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
799566063dSJacob 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];
829566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, row, 1, col, &J[0][0], INSERT_VALUES));
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
8482ad9101SHong Zhang 
859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
869566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8882ad9101SHong Zhang }
8982ad9101SHong Zhang 
main(int argc,char ** argv)90d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
91d71ae5a4SJacob Faibussowitsch {
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 
98327415f7SBarry Smith   PetscFunctionBeginUser;
999566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
10082ad9101SHong Zhang 
1019566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1023c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
10382ad9101SHong Zhang 
10482ad9101SHong Zhang   user.a     = 2.0;
10582ad9101SHong Zhang   user.b     = 4.0;
10682ad9101SHong Zhang   user.c     = 3.0;
10782ad9101SHong Zhang   user.steps = 0;
10882ad9101SHong Zhang   user.ftime = 1.0;
10982ad9101SHong Zhang   user.der   = 1;
1109566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-der", &user.der, NULL));
11182ad9101SHong Zhang 
11282ad9101SHong Zhang   rows = 1;
11382ad9101SHong Zhang   cols = 1;
1149566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jac));
1159566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jac, PETSC_DECIDE, PETSC_DECIDE, 1, 1));
1169566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jac));
1179566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jac));
1189566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jac, &user.x, NULL));
11982ad9101SHong Zhang 
1209566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1219566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
1229566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
1239566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, user.Jac, user.Jac, IJacobian, &user));
1249566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1259566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, user.ftime));
12682ad9101SHong Zhang 
1279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(user.x, &x_ptr));
12882ad9101SHong Zhang   x_ptr[0] = user.a;
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(user.x, &x_ptr));
1309566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.001));
13182ad9101SHong Zhang 
13282ad9101SHong Zhang   /* Set up forward sensitivity */
1339566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp));
1349566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, rows, cols));
1359566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
1369566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
1379566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, rows, cols, NULL, &user.sp));
1389566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(user.sp));
1399566063dSJacob Faibussowitsch   PetscCall(TSForwardSetSensitivities(ts, cols, user.sp));
1409566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobianP(ts, user.Jacp, IJacobianP, &user));
14182ad9101SHong Zhang 
1429566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
1439566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
14482ad9101SHong Zhang 
1459566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, user.x));
1469566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &user.ftime));
1479566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &user.steps));
1489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.x, &x_ptr));
1499566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n ode solution %g\n", (double)PetscRealPart(x_ptr[0])));
1509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.x, &x_ptr));
15163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n analytical solution %g\n", (double)(user.a * PetscExpReal(user.b / user.c * user.ftime))));
15282ad9101SHong Zhang 
15348a46eb9SPierre Jolivet   if (user.der == 1) 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))));
15448a46eb9SPierre Jolivet   if (user.der == 2) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n analytical derivative w.r.t. b %g\n", (double)(user.a * user.ftime / user.c * PetscExpReal(user.b / user.c * user.ftime))));
1559566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n forward sensitivity:\n"));
1569566063dSJacob Faibussowitsch   PetscCall(MatView(user.sp, PETSC_VIEWER_STDOUT_WORLD));
15782ad9101SHong Zhang 
1589566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jac, &user.lambda[0], NULL));
15982ad9101SHong Zhang   /* Set initial conditions for the adjoint integration */
1609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(user.lambda[0], &x_ptr));
16182ad9101SHong Zhang   x_ptr[0] = 1.0;
1629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(user.lambda[0], &x_ptr));
1639566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp, &user.mup[0], NULL));
1649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayWrite(user.mup[0], &x_ptr));
16582ad9101SHong Zhang   x_ptr[0] = 0.0;
1669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayWrite(user.mup[0], &x_ptr));
16782ad9101SHong Zhang 
1689566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts, 1, user.lambda, user.mup));
1699566063dSJacob Faibussowitsch   PetscCall(TSAdjointSolve(ts));
17082ad9101SHong Zhang 
1719566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n adjoint sensitivity:\n"));
1729566063dSJacob Faibussowitsch   PetscCall(VecView(user.mup[0], PETSC_VIEWER_STDOUT_WORLD));
17382ad9101SHong Zhang 
1749566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jac));
1759566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.sp));
1769566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
1779566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.x));
1789566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.lambda[0]));
1799566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.mup[0]));
1809566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
18182ad9101SHong Zhang 
1829566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
183b122ec5aSJacob Faibussowitsch   return 0;
18482ad9101SHong Zhang }
18582ad9101SHong Zhang 
18682ad9101SHong Zhang /*TEST
18782ad9101SHong Zhang 
18882ad9101SHong Zhang     test:
18982ad9101SHong Zhang       args: -ts_type beuler
19082ad9101SHong Zhang 
19182ad9101SHong Zhang     test:
19282ad9101SHong Zhang       suffix: 2
19382ad9101SHong Zhang       args: -ts_type cn
19482ad9101SHong Zhang       output_file: output/ex23fwdadj_1.out
19582ad9101SHong Zhang 
19682ad9101SHong Zhang TEST*/
197