xref: /petsc/src/ts/tutorials/ex20td.c (revision 9566063d113dddea24716c546802770db7481bc0)
180f8ae99SSajid Ali static char help[] = "Performs adjoint sensitivity analysis for a van der Pol like \n\
280f8ae99SSajid Ali equation with time dependent parameters using two approaches :  \n\
380f8ae99SSajid Ali track  : track only local sensitivities at each adjoint step \n\
480f8ae99SSajid Ali          and accumulate them in a global array \n\
580f8ae99SSajid Ali global : track parameters at all timesteps together \n\
680f8ae99SSajid Ali Choose one of the two at runtime by -sa_method {track,global}. \n";
780f8ae99SSajid Ali 
880f8ae99SSajid Ali /*
980f8ae99SSajid Ali   Concepts: TS^adjoint for time dependent parameters
1099c90e12SSatish Balay   Concepts: TS^Customized adjoint monitor based sensitivity tracking
1199c90e12SSatish Balay   Concepts: TS^All at once approach to sensitivity tracking
1280f8ae99SSajid Ali   Processors: 1
1380f8ae99SSajid Ali */
1480f8ae99SSajid Ali 
1580f8ae99SSajid Ali /*
1680f8ae99SSajid Ali    Simple example to demonstrate TSAdjoint capabilities for time dependent params
1780f8ae99SSajid Ali    without integral cost terms using either a tracking or global method.
1880f8ae99SSajid Ali 
1980f8ae99SSajid Ali    Modify the Van Der Pol Eq to :
2080f8ae99SSajid Ali    [u1'] = [mu1(t)*u1]
2180f8ae99SSajid Ali    [u2'] = [mu2(t)*((1-u1^2)*u2-u1)]
2280f8ae99SSajid Ali    (with initial conditions & params independent)
2380f8ae99SSajid Ali 
2480f8ae99SSajid Ali    Define uref to be solution with initail conditions (2,-2/3), mu=(1,1e3)
2580f8ae99SSajid Ali    - u_ref : (1.5967,-1.02969)
2680f8ae99SSajid Ali 
2780f8ae99SSajid Ali    Define const function as cost = 2-norm(u - u_ref);
2880f8ae99SSajid Ali 
2980f8ae99SSajid Ali    Initialization for the adjoint TS :
3080f8ae99SSajid Ali    - dcost/dy|final_time = 2*(u-u_ref)|final_time
3180f8ae99SSajid Ali    - dcost/dp|final_time = 0
3280f8ae99SSajid Ali 
3380f8ae99SSajid Ali    The tracking method only tracks local sensitivity at each time step
3480f8ae99SSajid Ali    and accumulates these sensitivities in a global array. Since the structure
3580f8ae99SSajid Ali    of the equations being solved at each time step does not change, the jacobian
3680f8ae99SSajid Ali    wrt parameters is defined analogous to constant RHSJacobian for a liner
3780f8ae99SSajid Ali    TSSolve and the size of the jacP is independent of the number of time
3880f8ae99SSajid Ali    steps. Enable this mode of adjoint analysis by -sa_method track.
3980f8ae99SSajid Ali 
4080f8ae99SSajid Ali    The global method combines the parameters at all timesteps and tracks them
4180f8ae99SSajid Ali    together. Thus, the columns of the jacP matrix are filled dependent upon the
4280f8ae99SSajid Ali    time step. Also, the dimensions of the jacP matrix now depend upon the number
4380f8ae99SSajid Ali    of time steps. Enable this mode of adjoint analysis by -sa_method global.
4480f8ae99SSajid Ali 
4580f8ae99SSajid Ali    Since the equations here have parameters at predefined time steps, this
4680f8ae99SSajid Ali    example should be run with non adaptive time stepping solvers only. This
4780f8ae99SSajid Ali    can be ensured by -ts_adapt_type none (which is the default behavior only
4880f8ae99SSajid Ali    for certain TS solvers like TSCN. If using an explicit TS like TSRK,
4980f8ae99SSajid Ali    please be sure to add the aforementioned option to disable adaptive
5080f8ae99SSajid Ali    timestepping.)
5180f8ae99SSajid Ali */
5280f8ae99SSajid Ali 
5380f8ae99SSajid Ali /*
5480f8ae99SSajid Ali    Include "petscts.h" so that we can use TS solvers.  Note that this file
5580f8ae99SSajid Ali    automatically includes:
5680f8ae99SSajid Ali      petscsys.h    - base PETSc routines   petscvec.h  - vectors
5780f8ae99SSajid Ali      petscmat.h    - matrices
5880f8ae99SSajid Ali      petscis.h     - index sets            petscksp.h  - Krylov subspace methods
5980f8ae99SSajid Ali      petscviewer.h - viewers               petscpc.h   - preconditioners
6080f8ae99SSajid Ali      petscksp.h    - linear solvers        petscsnes.h - nonlinear solvers
6180f8ae99SSajid Ali */
6280f8ae99SSajid Ali #include <petscts.h>
6380f8ae99SSajid Ali 
6480f8ae99SSajid Ali extern PetscErrorCode RHSFunction(TS ,PetscReal ,Vec ,Vec ,void *);
6580f8ae99SSajid Ali extern PetscErrorCode RHSJacobian(TS ,PetscReal ,Vec ,Mat ,Mat ,void *);
6680f8ae99SSajid Ali extern PetscErrorCode RHSJacobianP_track(TS ,PetscReal ,Vec ,Mat ,void *);
6780f8ae99SSajid Ali extern PetscErrorCode RHSJacobianP_global(TS ,PetscReal ,Vec ,Mat ,void *);
6880f8ae99SSajid Ali extern PetscErrorCode Monitor(TS ,PetscInt ,PetscReal ,Vec ,void *);
6980f8ae99SSajid Ali extern PetscErrorCode AdjointMonitor(TS ,PetscInt ,PetscReal ,Vec ,PetscInt ,Vec *, Vec *,void *);
7080f8ae99SSajid Ali 
7180f8ae99SSajid Ali /*
7280f8ae99SSajid Ali    User-defined application context - contains data needed by the
7380f8ae99SSajid Ali    application-provided call-back routines.
7480f8ae99SSajid Ali */
7580f8ae99SSajid Ali 
7680f8ae99SSajid Ali typedef struct {
7780f8ae99SSajid Ali  /*------------- Forward solve data structures --------------*/
7880f8ae99SSajid Ali   PetscInt  max_steps;     /* number of steps to run ts for */
7980f8ae99SSajid Ali   PetscReal final_time;    /* final time to integrate to*/
8080f8ae99SSajid Ali   PetscReal time_step;     /* ts integration time step */
8180f8ae99SSajid Ali   Vec       mu1;           /* time dependent params */
8280f8ae99SSajid Ali   Vec       mu2;           /* time dependent params */
8380f8ae99SSajid Ali   Vec       U;             /* solution vector */
8480f8ae99SSajid Ali   Mat       A;             /* Jacobian matrix */
8580f8ae99SSajid Ali 
8680f8ae99SSajid Ali   /*------------- Adjoint solve data structures --------------*/
8780f8ae99SSajid Ali   Mat       Jacp;          /* JacobianP matrix */
8880f8ae99SSajid Ali   Vec       lambda;        /* adjoint variable */
8980f8ae99SSajid Ali   Vec       mup;           /* adjoint variable */
9080f8ae99SSajid Ali 
9180f8ae99SSajid Ali   /*------------- Global accumation vecs for monitor based tracking --------------*/
9280f8ae99SSajid Ali   Vec       sens_mu1;      /* global sensitivity accumulation */
9380f8ae99SSajid Ali   Vec       sens_mu2;      /* global sensitivity accumulation */
9480f8ae99SSajid Ali   PetscInt  adj_idx;       /* to keep track of adjoint solve index */
9580f8ae99SSajid Ali } AppCtx;
9680f8ae99SSajid Ali 
9780f8ae99SSajid Ali typedef enum {SA_TRACK, SA_GLOBAL} SAMethod;
9880f8ae99SSajid Ali static const char *const SAMethods[] = {"TRACK","GLOBAL","SAMethod","SA_",0};
9980f8ae99SSajid Ali 
10080f8ae99SSajid Ali /* ----------------------- Explicit form of the ODE  -------------------- */
10180f8ae99SSajid Ali 
10280f8ae99SSajid Ali PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx)
10380f8ae99SSajid Ali {
10480f8ae99SSajid Ali   AppCtx            *user = (AppCtx*) ctx;
10580f8ae99SSajid Ali   PetscScalar       *f;
10680f8ae99SSajid Ali   PetscInt          curr_step;
10780f8ae99SSajid Ali   const PetscScalar *u;
10880f8ae99SSajid Ali   const PetscScalar *mu1;
10980f8ae99SSajid Ali   const PetscScalar *mu2;
11080f8ae99SSajid Ali 
11180f8ae99SSajid Ali   PetscFunctionBeginUser;
112*9566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&curr_step));
113*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
114*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user->mu1,&mu1));
115*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user->mu2,&mu2));
116*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
11780f8ae99SSajid Ali   f[0] = mu1[curr_step]*u[1];
11880f8ae99SSajid Ali   f[1] = mu2[curr_step]*((1.-u[0]*u[0])*u[1]-u[0]);
119*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
120*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user->mu1,&mu1));
121*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user->mu2,&mu2));
122*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
12380f8ae99SSajid Ali   PetscFunctionReturn(0);
12480f8ae99SSajid Ali }
12580f8ae99SSajid Ali 
12680f8ae99SSajid Ali PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
12780f8ae99SSajid Ali {
12880f8ae99SSajid Ali   AppCtx            *user = (AppCtx*) ctx;
12980f8ae99SSajid Ali   PetscInt          rowcol[] = {0,1};
13080f8ae99SSajid Ali   PetscScalar       J[2][2];
13180f8ae99SSajid Ali   PetscInt          curr_step;
13280f8ae99SSajid Ali   const PetscScalar *u;
13380f8ae99SSajid Ali   const PetscScalar *mu1;
13480f8ae99SSajid Ali   const PetscScalar *mu2;
13580f8ae99SSajid Ali 
13680f8ae99SSajid Ali   PetscFunctionBeginUser;
137*9566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&curr_step));
138*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user->mu1,&mu1));
139*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(user->mu2,&mu2));
140*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
14180f8ae99SSajid Ali   J[0][0] = 0;
14280f8ae99SSajid Ali   J[1][0] = -mu2[curr_step]*(2.0*u[1]*u[0]+1.);
14380f8ae99SSajid Ali   J[0][1] = mu1[curr_step];
14480f8ae99SSajid Ali   J[1][1] = mu2[curr_step]*(1.0-u[0]*u[0]);
145*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
146*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
147*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
148*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
149*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user->mu1,&mu1));
150*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(user->mu2,&mu2));
15180f8ae99SSajid Ali   PetscFunctionReturn(0);
15280f8ae99SSajid Ali }
15380f8ae99SSajid Ali 
15480f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for tracking method ------------------ */
15580f8ae99SSajid Ali 
15680f8ae99SSajid Ali PetscErrorCode RHSJacobianP_track(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
15780f8ae99SSajid Ali {
15880f8ae99SSajid Ali   PetscInt          row[] = {0,1},col[] = {0,1};
15980f8ae99SSajid Ali   PetscScalar       J[2][2];
16080f8ae99SSajid Ali   const PetscScalar *u;
16180f8ae99SSajid Ali 
16280f8ae99SSajid Ali   PetscFunctionBeginUser;
163*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
16480f8ae99SSajid Ali   J[0][0] = u[1];
16580f8ae99SSajid Ali   J[1][0] = 0;
16680f8ae99SSajid Ali   J[0][1] = 0;
16780f8ae99SSajid Ali   J[1][1] = (1.-u[0]*u[0])*u[1]-u[0];
168*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,row,2,col,&J[0][0],INSERT_VALUES));
169*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
170*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
171*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
17280f8ae99SSajid Ali   PetscFunctionReturn(0);
17380f8ae99SSajid Ali }
17480f8ae99SSajid Ali 
17580f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for global method ------------------ */
17680f8ae99SSajid Ali 
17780f8ae99SSajid Ali PetscErrorCode RHSJacobianP_global(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
17880f8ae99SSajid Ali {
17980f8ae99SSajid Ali   PetscInt          row[] = {0,1},col[] = {0,1};
18080f8ae99SSajid Ali   PetscScalar       J[2][2];
18180f8ae99SSajid Ali   const PetscScalar *u;
18280f8ae99SSajid Ali   PetscInt          curr_step;
18380f8ae99SSajid Ali 
18480f8ae99SSajid Ali   PetscFunctionBeginUser;
185*9566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&curr_step));
186*9566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(U,&u));
18780f8ae99SSajid Ali   J[0][0] = u[1];
18880f8ae99SSajid Ali   J[1][0] = 0;
18980f8ae99SSajid Ali   J[0][1] = 0;
19080f8ae99SSajid Ali   J[1][1] = (1.-u[0]*u[0])*u[1]-u[0];
19180f8ae99SSajid Ali   col[0] = (curr_step)*2;
19280f8ae99SSajid Ali   col[1] = (curr_step)*2+1;
193*9566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,row,2,col,&J[0][0],INSERT_VALUES));
194*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
195*9566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
196*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(U,&u));
19780f8ae99SSajid Ali   PetscFunctionReturn(0);
19880f8ae99SSajid Ali }
19980f8ae99SSajid Ali 
20080f8ae99SSajid Ali /* Dump solution to console if called */
20180f8ae99SSajid Ali PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx)
20280f8ae99SSajid Ali {
20380f8ae99SSajid Ali   PetscFunctionBeginUser;
204*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n Solution at time %e is \n", t));
205*9566063dSJacob Faibussowitsch   PetscCall(VecView(U,PETSC_VIEWER_STDOUT_WORLD));
20680f8ae99SSajid Ali   PetscFunctionReturn(0);
20780f8ae99SSajid Ali }
20880f8ae99SSajid Ali 
20980f8ae99SSajid Ali /* Customized adjoint monitor to keep track of local
21080f8ae99SSajid Ali    sensitivities by storing them in a global sensitivity array.
21180f8ae99SSajid Ali    Note : This routine is only used for the tracking method. */
21280f8ae99SSajid Ali PetscErrorCode AdjointMonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *ctx)
21380f8ae99SSajid Ali {
21480f8ae99SSajid Ali   AppCtx            *user = (AppCtx*) ctx;
21580f8ae99SSajid Ali   PetscInt          curr_step;
21680f8ae99SSajid Ali   PetscScalar       *sensmu1_glob;
21780f8ae99SSajid Ali   PetscScalar       *sensmu2_glob;
21880f8ae99SSajid Ali   const PetscScalar *sensmu_loc;
21980f8ae99SSajid Ali 
22080f8ae99SSajid Ali   PetscFunctionBeginUser;
221*9566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&curr_step));
22280f8ae99SSajid Ali   /* Note that we skip the first call to the monitor in the adjoint
22380f8ae99SSajid Ali      solve since the sensitivities are already set (during
22480f8ae99SSajid Ali      initialization of adjoint vectors).
22580f8ae99SSajid Ali      We also note that each indvidial TSAdjointSolve calls the monitor
22680f8ae99SSajid Ali      twice, once at the step it is integrating from and once at the step
22780f8ae99SSajid Ali      it integrates to. Only the second call is useful for transferring
22880f8ae99SSajid Ali      local sensitivities to the global array. */
22980f8ae99SSajid Ali   if (curr_step == user->adj_idx) {
23080f8ae99SSajid Ali     PetscFunctionReturn(0);
23180f8ae99SSajid Ali   } else {
232*9566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(*mu,&sensmu_loc));
233*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->sens_mu1,&sensmu1_glob));
234*9566063dSJacob Faibussowitsch     PetscCall(VecGetArray(user->sens_mu2,&sensmu2_glob));
23580f8ae99SSajid Ali     sensmu1_glob[curr_step] = sensmu_loc[0];
23680f8ae99SSajid Ali     sensmu2_glob[curr_step] = sensmu_loc[1];
237*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->sens_mu1,&sensmu1_glob));
238*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(user->sens_mu2,&sensmu2_glob));
239*9566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(*mu,&sensmu_loc));
24080f8ae99SSajid Ali     PetscFunctionReturn(0);
24180f8ae99SSajid Ali   }
24280f8ae99SSajid Ali }
24380f8ae99SSajid Ali 
24480f8ae99SSajid Ali int main(int argc,char **argv)
24580f8ae99SSajid Ali {
24680f8ae99SSajid Ali   TS             ts;
24780f8ae99SSajid Ali   AppCtx         user;
24880f8ae99SSajid Ali   PetscScalar    *x_ptr,*y_ptr,*u_ptr;
24980f8ae99SSajid Ali   PetscMPIInt    size;
25080f8ae99SSajid Ali   PetscBool      monitor = PETSC_FALSE;
25180f8ae99SSajid Ali   SAMethod       sa = SA_GLOBAL;
25280f8ae99SSajid Ali   PetscErrorCode ierr;
25380f8ae99SSajid Ali 
25480f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
25580f8ae99SSajid Ali      Initialize program
25680f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
257*9566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
258*9566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
2593c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
26080f8ae99SSajid Ali 
26180f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26280f8ae99SSajid Ali      Set runtime options
26380f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
264*9566063dSJacob Faibussowitsch   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"SA analysis options.","");PetscCall(ierr);{
265*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
266*9566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-sa_method","Sensitivity analysis method (track or global)","",SAMethods,(PetscEnum)sa,(PetscEnum*)&sa,NULL));
26780f8ae99SSajid Ali   }
268*9566063dSJacob Faibussowitsch   ierr = PetscOptionsEnd();PetscCall(ierr);
26980f8ae99SSajid Ali 
27080f8ae99SSajid Ali   user.final_time = 0.1;
27180f8ae99SSajid Ali   user.max_steps  = 5;
27280f8ae99SSajid Ali   user.time_step  = user.final_time/user.max_steps;
27380f8ae99SSajid Ali 
27480f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
27580f8ae99SSajid Ali      Create necessary matrix and vectors for forward solve.
27680f8ae99SSajid Ali      Create Jacp matrix for adjoint solve.
27780f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.mu1));
279*9566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.mu2));
280*9566063dSJacob Faibussowitsch   PetscCall(VecSet(user.mu1,1.25));
281*9566063dSJacob Faibussowitsch   PetscCall(VecSet(user.mu2,1.0e2));
28280f8ae99SSajid Ali 
28380f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28480f8ae99SSajid Ali       For tracking method : create the global sensitivity array to
28580f8ae99SSajid Ali       accumulate sensitivity with respect to parameters at each step.
28680f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
28780f8ae99SSajid Ali   if (sa == SA_TRACK) {
288*9566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.sens_mu1));
289*9566063dSJacob Faibussowitsch     PetscCall(VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.sens_mu2));
29080f8ae99SSajid Ali   }
29180f8ae99SSajid Ali 
292*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.A));
293*9566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2));
294*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.A));
295*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.A));
296*9566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.U,NULL));
29780f8ae99SSajid Ali 
29880f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29980f8ae99SSajid Ali       Note that the dimensions of the Jacp matrix depend upon the
30080f8ae99SSajid Ali       sensitivity analysis method being used !
30180f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
302*9566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
30380f8ae99SSajid Ali   if (sa == SA_TRACK) {
304*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2));
30580f8ae99SSajid Ali   }
30680f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
307*9566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,user.max_steps*2));
30880f8ae99SSajid Ali   }
309*9566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
310*9566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
31180f8ae99SSajid Ali 
31280f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
31380f8ae99SSajid Ali      Create timestepping solver context
31480f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
315*9566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
316*9566063dSJacob Faibussowitsch   PetscCall(TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT));
317*9566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSCN));
31880f8ae99SSajid Ali 
319*9566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user));
320*9566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user));
32180f8ae99SSajid Ali   if (sa == SA_TRACK) {
322*9566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP_track,&user));
32380f8ae99SSajid Ali   }
32480f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
325*9566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP_global,&user));
32680f8ae99SSajid Ali   }
32780f8ae99SSajid Ali 
328*9566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
329*9566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,user.final_time));
330*9566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,user.final_time/user.max_steps));
33180f8ae99SSajid Ali 
33280f8ae99SSajid Ali   if (monitor) {
333*9566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
33480f8ae99SSajid Ali   }
33580f8ae99SSajid Ali   if (sa == SA_TRACK) {
336*9566063dSJacob Faibussowitsch     PetscCall(TSAdjointMonitorSet(ts,AdjointMonitor,&user,NULL));
33780f8ae99SSajid Ali   }
33880f8ae99SSajid Ali 
33980f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34080f8ae99SSajid Ali      Set initial conditions
34180f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
342*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U,&x_ptr));
34380f8ae99SSajid Ali   x_ptr[0] = 2.0;
34480f8ae99SSajid Ali   x_ptr[1] = -2.0/3.0;
345*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U,&x_ptr));
34680f8ae99SSajid Ali 
34780f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34880f8ae99SSajid Ali      Save trajectory of solution so that TSAdjointSolve() may be used
34980f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
350*9566063dSJacob Faibussowitsch   PetscCall(TSSetSaveTrajectory(ts));
35180f8ae99SSajid Ali 
35280f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35380f8ae99SSajid Ali      Set runtime options
35480f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
355*9566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
35680f8ae99SSajid Ali 
35780f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35880f8ae99SSajid Ali      Execute forward model and print solution.
35980f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
360*9566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,user.U));
361*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n Solution of forward TS :\n"));
362*9566063dSJacob Faibussowitsch   PetscCall(VecView(user.U,PETSC_VIEWER_STDOUT_WORLD));
363*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n Forward TS solve successfull! Adjoint run begins!\n"));
36480f8ae99SSajid Ali 
36580f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36680f8ae99SSajid Ali      Adjoint model starts here! Create adjoint vectors.
36780f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
368*9566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.A,&user.lambda,NULL));
369*9566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jacp,&user.mup,NULL));
37080f8ae99SSajid Ali 
37180f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37280f8ae99SSajid Ali      Set initial conditions for the adjoint vector
37380f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
374*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.U,&u_ptr));
375*9566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.lambda,&y_ptr));
37680f8ae99SSajid Ali   y_ptr[0] = 2*(u_ptr[0] - 1.5967);
37780f8ae99SSajid Ali   y_ptr[1] = 2*(u_ptr[1] - -(1.02969));
378*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.lambda,&y_ptr));
379*9566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.U,&y_ptr));
380*9566063dSJacob Faibussowitsch   PetscCall(VecSet(user.mup,0));
38180f8ae99SSajid Ali 
38280f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38380f8ae99SSajid Ali      Set number of cost functions.
38480f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
385*9566063dSJacob Faibussowitsch   PetscCall(TSSetCostGradients(ts,1,&user.lambda,&user.mup));
38680f8ae99SSajid Ali 
38780f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38880f8ae99SSajid Ali      The adjoint vector mup has to be reset for each adjoint step when
38980f8ae99SSajid Ali      using the tracking method as we want to treat the parameters at each
39080f8ae99SSajid Ali      time step one at a time and prevent accumulation of the sensitivities
39180f8ae99SSajid Ali      from parameters at previous time steps.
39280f8ae99SSajid Ali      This is not necessary for the global method as each time dependent
39380f8ae99SSajid Ali      parameter is treated as an independent parameter.
39480f8ae99SSajid Ali    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39580f8ae99SSajid Ali   if (sa == SA_TRACK) {
39680f8ae99SSajid Ali     for (user.adj_idx=user.max_steps; user.adj_idx>0; user.adj_idx--) {
397*9566063dSJacob Faibussowitsch       PetscCall(VecSet(user.mup,0));
398*9566063dSJacob Faibussowitsch       PetscCall(TSAdjointSetSteps(ts, 1));
399*9566063dSJacob Faibussowitsch       PetscCall(TSAdjointSolve(ts));
40080f8ae99SSajid Ali     }
40180f8ae99SSajid Ali   }
40280f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
403*9566063dSJacob Faibussowitsch     PetscCall(TSAdjointSolve(ts));
40480f8ae99SSajid Ali   }
40580f8ae99SSajid Ali 
40680f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
40780f8ae99SSajid Ali      Dispaly adjoint sensitivities wrt parameters and initial conditions
40880f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40980f8ae99SSajid Ali   if (sa == SA_TRACK) {
410*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt  mu1: d[cost]/d[mu1]\n"));
411*9566063dSJacob Faibussowitsch     PetscCall(VecView(user.sens_mu1,PETSC_VIEWER_STDOUT_WORLD));
412*9566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt  mu2: d[cost]/d[mu2]\n"));
413*9566063dSJacob Faibussowitsch     PetscCall(VecView(user.sens_mu2,PETSC_VIEWER_STDOUT_WORLD));
41480f8ae99SSajid Ali   }
41580f8ae99SSajid Ali 
41680f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
41780f8ae99SSajid Ali     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt  params: d[cost]/d[p], where p refers to \n\
418*9566063dSJacob Faibussowitsch                     the interlaced vector made by combining mu1,mu2\n");PetscCall(ierr);
419*9566063dSJacob Faibussowitsch     PetscCall(VecView(user.mup,PETSC_VIEWER_STDOUT_WORLD));
42080f8ae99SSajid Ali   }
42180f8ae99SSajid Ali 
422*9566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[cost]/d[u(t=0)]\n"));
423*9566063dSJacob Faibussowitsch   PetscCall(VecView(user.lambda,PETSC_VIEWER_STDOUT_WORLD));
42480f8ae99SSajid Ali 
42580f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
42680f8ae99SSajid Ali      Free work space!
42780f8ae99SSajid Ali      All PETSc objects should be destroyed when they are no longer needed.
42880f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
429*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.A));
430*9566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
431*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.U));
432*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.lambda));
433*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.mup));
434*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.mu1));
435*9566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.mu2));
43680f8ae99SSajid Ali   if (sa == SA_TRACK) {
437*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&user.sens_mu1));
438*9566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&user.sens_mu2));
43980f8ae99SSajid Ali   }
440*9566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
441*9566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
44280f8ae99SSajid Ali   return(ierr);
44380f8ae99SSajid Ali }
44480f8ae99SSajid Ali 
44580f8ae99SSajid Ali /*TEST
44680f8ae99SSajid Ali 
44780f8ae99SSajid Ali   test:
44880f8ae99SSajid Ali     requires: !complex
44980f8ae99SSajid Ali     suffix : track
45080f8ae99SSajid Ali     args : -sa_method track
45180f8ae99SSajid Ali 
45280f8ae99SSajid Ali   test:
45380f8ae99SSajid Ali     requires: !complex
45480f8ae99SSajid Ali     suffix : global
45580f8ae99SSajid Ali     args : -sa_method global
45680f8ae99SSajid Ali 
45780f8ae99SSajid Ali TEST*/
458