xref: /petsc/src/ts/tutorials/ex20td.c (revision 3c633725528e547aaaa9b672a746f5d686a276e1)
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   PetscErrorCode    ierr;
10580f8ae99SSajid Ali   AppCtx            *user = (AppCtx*) ctx;
10680f8ae99SSajid Ali   PetscScalar       *f;
10780f8ae99SSajid Ali   PetscInt          curr_step;
10880f8ae99SSajid Ali   const PetscScalar *u;
10980f8ae99SSajid Ali   const PetscScalar *mu1;
11080f8ae99SSajid Ali   const PetscScalar *mu2;
11180f8ae99SSajid Ali 
11280f8ae99SSajid Ali   PetscFunctionBeginUser;
11380f8ae99SSajid Ali   ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr);
11480f8ae99SSajid Ali   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
11580f8ae99SSajid Ali   ierr = VecGetArrayRead(user->mu1,&mu1);CHKERRQ(ierr);
11680f8ae99SSajid Ali   ierr = VecGetArrayRead(user->mu2,&mu2);CHKERRQ(ierr);
11780f8ae99SSajid Ali   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
11880f8ae99SSajid Ali   f[0] = mu1[curr_step]*u[1];
11980f8ae99SSajid Ali   f[1] = mu2[curr_step]*((1.-u[0]*u[0])*u[1]-u[0]);
12080f8ae99SSajid Ali   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
12180f8ae99SSajid Ali   ierr = VecRestoreArrayRead(user->mu1,&mu1);CHKERRQ(ierr);
12280f8ae99SSajid Ali   ierr = VecRestoreArrayRead(user->mu2,&mu2);CHKERRQ(ierr);
12380f8ae99SSajid Ali   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
12480f8ae99SSajid Ali   PetscFunctionReturn(0);
12580f8ae99SSajid Ali }
12680f8ae99SSajid Ali 
12780f8ae99SSajid Ali PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx)
12880f8ae99SSajid Ali {
12980f8ae99SSajid Ali   PetscErrorCode    ierr;
13080f8ae99SSajid Ali   AppCtx            *user = (AppCtx*) ctx;
13180f8ae99SSajid Ali   PetscInt          rowcol[] = {0,1};
13280f8ae99SSajid Ali   PetscScalar       J[2][2];
13380f8ae99SSajid Ali   PetscInt          curr_step;
13480f8ae99SSajid Ali   const PetscScalar *u;
13580f8ae99SSajid Ali   const PetscScalar *mu1;
13680f8ae99SSajid Ali   const PetscScalar *mu2;
13780f8ae99SSajid Ali 
13880f8ae99SSajid Ali   PetscFunctionBeginUser;
13980f8ae99SSajid Ali   ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr);
14080f8ae99SSajid Ali   ierr = VecGetArrayRead(user->mu1,&mu1);CHKERRQ(ierr);
14180f8ae99SSajid Ali   ierr = VecGetArrayRead(user->mu2,&mu2);CHKERRQ(ierr);
14280f8ae99SSajid Ali   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
14380f8ae99SSajid Ali   J[0][0] = 0;
14480f8ae99SSajid Ali   J[1][0] = -mu2[curr_step]*(2.0*u[1]*u[0]+1.);
14580f8ae99SSajid Ali   J[0][1] = mu1[curr_step];
14680f8ae99SSajid Ali   J[1][1] = mu2[curr_step]*(1.0-u[0]*u[0]);
14780f8ae99SSajid Ali   ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
14880f8ae99SSajid Ali   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
14980f8ae99SSajid Ali   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15080f8ae99SSajid Ali   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
15180f8ae99SSajid Ali   ierr = VecRestoreArrayRead(user->mu1,&mu1);CHKERRQ(ierr);
15280f8ae99SSajid Ali   ierr = VecRestoreArrayRead(user->mu2,&mu2);CHKERRQ(ierr);
15380f8ae99SSajid Ali   PetscFunctionReturn(0);
15480f8ae99SSajid Ali }
15580f8ae99SSajid Ali 
15680f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for tracking method ------------------ */
15780f8ae99SSajid Ali 
15880f8ae99SSajid Ali PetscErrorCode RHSJacobianP_track(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
15980f8ae99SSajid Ali {
16080f8ae99SSajid Ali   PetscErrorCode    ierr;
16180f8ae99SSajid Ali   PetscInt          row[] = {0,1},col[] = {0,1};
16280f8ae99SSajid Ali   PetscScalar       J[2][2];
16380f8ae99SSajid Ali   const PetscScalar *u;
16480f8ae99SSajid Ali 
16580f8ae99SSajid Ali   PetscFunctionBeginUser;
16680f8ae99SSajid Ali   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
16780f8ae99SSajid Ali   J[0][0] = u[1];
16880f8ae99SSajid Ali   J[1][0] = 0;
16980f8ae99SSajid Ali   J[0][1] = 0;
17080f8ae99SSajid Ali   J[1][1] = (1.-u[0]*u[0])*u[1]-u[0];
17180f8ae99SSajid Ali   ierr = MatSetValues(A,2,row,2,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
17280f8ae99SSajid Ali   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17380f8ae99SSajid Ali   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
17480f8ae99SSajid Ali   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
17580f8ae99SSajid Ali   PetscFunctionReturn(0);
17680f8ae99SSajid Ali }
17780f8ae99SSajid Ali 
17880f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for global method ------------------ */
17980f8ae99SSajid Ali 
18080f8ae99SSajid Ali PetscErrorCode RHSJacobianP_global(TS ts,PetscReal t,Vec U,Mat A,void *ctx)
18180f8ae99SSajid Ali {
18280f8ae99SSajid Ali   PetscErrorCode    ierr;
18380f8ae99SSajid Ali   PetscInt          row[] = {0,1},col[] = {0,1};
18480f8ae99SSajid Ali   PetscScalar       J[2][2];
18580f8ae99SSajid Ali   const PetscScalar *u;
18680f8ae99SSajid Ali   PetscInt          curr_step;
18780f8ae99SSajid Ali 
18880f8ae99SSajid Ali   PetscFunctionBeginUser;
18980f8ae99SSajid Ali   ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr);
19080f8ae99SSajid Ali   ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr);
19180f8ae99SSajid Ali   J[0][0] = u[1];
19280f8ae99SSajid Ali   J[1][0] = 0;
19380f8ae99SSajid Ali   J[0][1] = 0;
19480f8ae99SSajid Ali   J[1][1] = (1.-u[0]*u[0])*u[1]-u[0];
19580f8ae99SSajid Ali   col[0] = (curr_step)*2;
19680f8ae99SSajid Ali   col[1] = (curr_step)*2+1;
19780f8ae99SSajid Ali   ierr = MatSetValues(A,2,row,2,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
19880f8ae99SSajid Ali   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
19980f8ae99SSajid Ali   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
20080f8ae99SSajid Ali   ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr);
20180f8ae99SSajid Ali   PetscFunctionReturn(0);
20280f8ae99SSajid Ali }
20380f8ae99SSajid Ali 
20480f8ae99SSajid Ali /* Dump solution to console if called */
20580f8ae99SSajid Ali PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx)
20680f8ae99SSajid Ali {
20780f8ae99SSajid Ali   PetscErrorCode    ierr;
20880f8ae99SSajid Ali 
20980f8ae99SSajid Ali   PetscFunctionBeginUser;
21080f8ae99SSajid Ali   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Solution at time %e is \n", t);CHKERRQ(ierr);
21180f8ae99SSajid Ali   ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
21280f8ae99SSajid Ali   PetscFunctionReturn(0);
21380f8ae99SSajid Ali }
21480f8ae99SSajid Ali 
21580f8ae99SSajid Ali /* Customized adjoint monitor to keep track of local
21680f8ae99SSajid Ali    sensitivities by storing them in a global sensitivity array.
21780f8ae99SSajid Ali    Note : This routine is only used for the tracking method. */
21880f8ae99SSajid Ali PetscErrorCode AdjointMonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *ctx)
21980f8ae99SSajid Ali {
22080f8ae99SSajid Ali   PetscErrorCode    ierr;
22180f8ae99SSajid Ali   AppCtx            *user = (AppCtx*) ctx;
22280f8ae99SSajid Ali   PetscInt          curr_step;
22380f8ae99SSajid Ali   PetscScalar       *sensmu1_glob;
22480f8ae99SSajid Ali   PetscScalar       *sensmu2_glob;
22580f8ae99SSajid Ali   const PetscScalar *sensmu_loc;
22680f8ae99SSajid Ali 
22780f8ae99SSajid Ali   PetscFunctionBeginUser;
22880f8ae99SSajid Ali   ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr);
22980f8ae99SSajid Ali   /* Note that we skip the first call to the monitor in the adjoint
23080f8ae99SSajid Ali      solve since the sensitivities are already set (during
23180f8ae99SSajid Ali      initialization of adjoint vectors).
23280f8ae99SSajid Ali      We also note that each indvidial TSAdjointSolve calls the monitor
23380f8ae99SSajid Ali      twice, once at the step it is integrating from and once at the step
23480f8ae99SSajid Ali      it integrates to. Only the second call is useful for transferring
23580f8ae99SSajid Ali      local sensitivities to the global array. */
23680f8ae99SSajid Ali   if (curr_step == user->adj_idx) {
23780f8ae99SSajid Ali     PetscFunctionReturn(0);
23880f8ae99SSajid Ali   } else {
23980f8ae99SSajid Ali     ierr = VecGetArrayRead(*mu,&sensmu_loc);CHKERRQ(ierr);
24080f8ae99SSajid Ali     ierr = VecGetArray(user->sens_mu1,&sensmu1_glob);CHKERRQ(ierr);
24180f8ae99SSajid Ali     ierr = VecGetArray(user->sens_mu2,&sensmu2_glob);CHKERRQ(ierr);
24280f8ae99SSajid Ali     sensmu1_glob[curr_step] = sensmu_loc[0];
24380f8ae99SSajid Ali     sensmu2_glob[curr_step] = sensmu_loc[1];
24480f8ae99SSajid Ali     ierr = VecRestoreArray(user->sens_mu1,&sensmu1_glob);CHKERRQ(ierr);
24580f8ae99SSajid Ali     ierr = VecRestoreArray(user->sens_mu2,&sensmu2_glob);CHKERRQ(ierr);
24680f8ae99SSajid Ali     ierr = VecRestoreArrayRead(*mu,&sensmu_loc);CHKERRQ(ierr);
24780f8ae99SSajid Ali     PetscFunctionReturn(0);
24880f8ae99SSajid Ali   }
24980f8ae99SSajid Ali }
25080f8ae99SSajid Ali 
25180f8ae99SSajid Ali int main(int argc,char **argv)
25280f8ae99SSajid Ali {
25380f8ae99SSajid Ali   TS             ts;
25480f8ae99SSajid Ali   AppCtx         user;
25580f8ae99SSajid Ali   PetscScalar    *x_ptr,*y_ptr,*u_ptr;
25680f8ae99SSajid Ali   PetscMPIInt    size;
25780f8ae99SSajid Ali   PetscBool      monitor = PETSC_FALSE;
25880f8ae99SSajid Ali   SAMethod       sa = SA_GLOBAL;
25980f8ae99SSajid Ali   PetscErrorCode ierr;
26080f8ae99SSajid Ali 
26180f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26280f8ae99SSajid Ali      Initialize program
26380f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
26480f8ae99SSajid Ali   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
265ffc4695bSBarry Smith   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr);
266*3c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
26780f8ae99SSajid Ali 
26880f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
26980f8ae99SSajid Ali      Set runtime options
27080f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
27180f8ae99SSajid Ali   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"SA analysis options.","");CHKERRQ(ierr);{
27280f8ae99SSajid Ali   ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr);
27380f8ae99SSajid Ali   ierr = PetscOptionsEnum("-sa_method","Sensitivity analysis method (track or global)","",SAMethods,(PetscEnum)sa,(PetscEnum*)&sa,NULL);CHKERRQ(ierr);
27480f8ae99SSajid Ali   }
27580f8ae99SSajid Ali   ierr = PetscOptionsEnd();CHKERRQ(ierr);
27680f8ae99SSajid Ali 
27780f8ae99SSajid Ali   user.final_time = 0.1;
27880f8ae99SSajid Ali   user.max_steps  = 5;
27980f8ae99SSajid Ali   user.time_step  = user.final_time/user.max_steps;
28080f8ae99SSajid Ali 
28180f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
28280f8ae99SSajid Ali      Create necessary matrix and vectors for forward solve.
28380f8ae99SSajid Ali      Create Jacp matrix for adjoint solve.
28480f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
28580f8ae99SSajid Ali   ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.mu1);CHKERRQ(ierr);
28680f8ae99SSajid Ali   ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.mu2);CHKERRQ(ierr);
28780f8ae99SSajid Ali   ierr = VecSet(user.mu1,1.25);CHKERRQ(ierr);
28880f8ae99SSajid Ali   ierr = VecSet(user.mu2,1.0e2);CHKERRQ(ierr);
28980f8ae99SSajid Ali 
29080f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
29180f8ae99SSajid Ali       For tracking method : create the global sensitivity array to
29280f8ae99SSajid Ali       accumulate sensitivity with respect to parameters at each step.
29380f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
29480f8ae99SSajid Ali   if (sa == SA_TRACK) {
29580f8ae99SSajid Ali     ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.sens_mu1);CHKERRQ(ierr);
29680f8ae99SSajid Ali     ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.sens_mu2);CHKERRQ(ierr);
29780f8ae99SSajid Ali   }
29880f8ae99SSajid Ali 
29980f8ae99SSajid Ali   ierr = MatCreate(PETSC_COMM_WORLD,&user.A);CHKERRQ(ierr);
30080f8ae99SSajid Ali   ierr = MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
30180f8ae99SSajid Ali   ierr = MatSetFromOptions(user.A);CHKERRQ(ierr);
30280f8ae99SSajid Ali   ierr = MatSetUp(user.A);CHKERRQ(ierr);
30380f8ae99SSajid Ali   ierr = MatCreateVecs(user.A,&user.U,NULL);CHKERRQ(ierr);
30480f8ae99SSajid Ali 
30580f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
30680f8ae99SSajid Ali       Note that the dimensions of the Jacp matrix depend upon the
30780f8ae99SSajid Ali       sensitivity analysis method being used !
30880f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
30980f8ae99SSajid Ali   ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr);
31080f8ae99SSajid Ali   if (sa == SA_TRACK) {
31180f8ae99SSajid Ali     ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr);
31280f8ae99SSajid Ali   }
31380f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
31480f8ae99SSajid Ali     ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,user.max_steps*2);CHKERRQ(ierr);
31580f8ae99SSajid Ali   }
31680f8ae99SSajid Ali   ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr);
31780f8ae99SSajid Ali   ierr = MatSetUp(user.Jacp);CHKERRQ(ierr);
31880f8ae99SSajid Ali 
31980f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
32080f8ae99SSajid Ali      Create timestepping solver context
32180f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
32280f8ae99SSajid Ali   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
32380f8ae99SSajid Ali   ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr);
32480f8ae99SSajid Ali   ierr = TSSetType(ts,TSCN);CHKERRQ(ierr);
32580f8ae99SSajid Ali 
32680f8ae99SSajid Ali   ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr);
32780f8ae99SSajid Ali   ierr = TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user);CHKERRQ(ierr);
32880f8ae99SSajid Ali   if (sa == SA_TRACK) {
32980f8ae99SSajid Ali     ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP_track,&user);CHKERRQ(ierr);
33080f8ae99SSajid Ali   }
33180f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
33280f8ae99SSajid Ali     ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP_global,&user);CHKERRQ(ierr);
33380f8ae99SSajid Ali   }
33480f8ae99SSajid Ali 
33580f8ae99SSajid Ali   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
33680f8ae99SSajid Ali   ierr = TSSetMaxTime(ts,user.final_time);CHKERRQ(ierr);
33780f8ae99SSajid Ali   ierr = TSSetTimeStep(ts,user.final_time/user.max_steps);CHKERRQ(ierr);
33880f8ae99SSajid Ali 
33980f8ae99SSajid Ali   if (monitor) {
34080f8ae99SSajid Ali     ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr);
34180f8ae99SSajid Ali   }
34280f8ae99SSajid Ali   if (sa == SA_TRACK) {
34380f8ae99SSajid Ali     ierr = TSAdjointMonitorSet(ts,AdjointMonitor,&user,NULL);CHKERRQ(ierr);
34480f8ae99SSajid Ali   }
34580f8ae99SSajid Ali 
34680f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
34780f8ae99SSajid Ali      Set initial conditions
34880f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
34980f8ae99SSajid Ali   ierr = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr);
35080f8ae99SSajid Ali   x_ptr[0] = 2.0;
35180f8ae99SSajid Ali   x_ptr[1] = -2.0/3.0;
35280f8ae99SSajid Ali   ierr = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr);
35380f8ae99SSajid Ali 
35480f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35580f8ae99SSajid Ali      Save trajectory of solution so that TSAdjointSolve() may be used
35680f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
35780f8ae99SSajid Ali   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
35880f8ae99SSajid Ali 
35980f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36080f8ae99SSajid Ali      Set runtime options
36180f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36280f8ae99SSajid Ali   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
36380f8ae99SSajid Ali 
36480f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36580f8ae99SSajid Ali      Execute forward model and print solution.
36680f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
36780f8ae99SSajid Ali   ierr = TSSolve(ts,user.U);CHKERRQ(ierr);
36880f8ae99SSajid Ali   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Solution of forward TS :\n");CHKERRQ(ierr);
36980f8ae99SSajid Ali   ierr = VecView(user.U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
37080f8ae99SSajid Ali   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Forward TS solve successfull! Adjoint run begins!\n");CHKERRQ(ierr);
37180f8ae99SSajid Ali 
37280f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37380f8ae99SSajid Ali      Adjoint model starts here! Create adjoint vectors.
37480f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
37580f8ae99SSajid Ali   ierr = MatCreateVecs(user.A,&user.lambda,NULL);CHKERRQ(ierr);
37680f8ae99SSajid Ali   ierr = MatCreateVecs(user.Jacp,&user.mup,NULL);CHKERRQ(ierr);
37780f8ae99SSajid Ali 
37880f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
37980f8ae99SSajid Ali      Set initial conditions for the adjoint vector
38080f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
38180f8ae99SSajid Ali   ierr = VecGetArray(user.U,&u_ptr);CHKERRQ(ierr);
38280f8ae99SSajid Ali   ierr = VecGetArray(user.lambda,&y_ptr);CHKERRQ(ierr);
38380f8ae99SSajid Ali   y_ptr[0] = 2*(u_ptr[0] - 1.5967);
38480f8ae99SSajid Ali   y_ptr[1] = 2*(u_ptr[1] - -(1.02969));
38580f8ae99SSajid Ali   ierr = VecRestoreArray(user.lambda,&y_ptr);CHKERRQ(ierr);
38680f8ae99SSajid Ali   ierr = VecRestoreArray(user.U,&y_ptr);CHKERRQ(ierr);
38780f8ae99SSajid Ali   ierr = VecSet(user.mup,0);CHKERRQ(ierr);
38880f8ae99SSajid Ali 
38980f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39080f8ae99SSajid Ali      Set number of cost functions.
39180f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
39280f8ae99SSajid Ali   ierr = TSSetCostGradients(ts,1,&user.lambda,&user.mup);CHKERRQ(ierr);
39380f8ae99SSajid Ali 
39480f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
39580f8ae99SSajid Ali      The adjoint vector mup has to be reset for each adjoint step when
39680f8ae99SSajid Ali      using the tracking method as we want to treat the parameters at each
39780f8ae99SSajid Ali      time step one at a time and prevent accumulation of the sensitivities
39880f8ae99SSajid Ali      from parameters at previous time steps.
39980f8ae99SSajid Ali      This is not necessary for the global method as each time dependent
40080f8ae99SSajid Ali      parameter is treated as an independent parameter.
40180f8ae99SSajid Ali    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
40280f8ae99SSajid Ali   if (sa == SA_TRACK) {
40380f8ae99SSajid Ali     for (user.adj_idx=user.max_steps; user.adj_idx>0; user.adj_idx--) {
40480f8ae99SSajid Ali       ierr = VecSet(user.mup,0);CHKERRQ(ierr);
40580f8ae99SSajid Ali       ierr = TSAdjointSetSteps(ts, 1);CHKERRQ(ierr);
40680f8ae99SSajid Ali       ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
40780f8ae99SSajid Ali     }
40880f8ae99SSajid Ali   }
40980f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
41080f8ae99SSajid Ali     ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
41180f8ae99SSajid Ali   }
41280f8ae99SSajid Ali 
41380f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41480f8ae99SSajid Ali      Dispaly adjoint sensitivities wrt parameters and initial conditions
41580f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
41680f8ae99SSajid Ali   if (sa == SA_TRACK) {
41780f8ae99SSajid Ali     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt  mu1: d[cost]/d[mu1]\n");CHKERRQ(ierr);
41880f8ae99SSajid Ali     ierr = VecView(user.sens_mu1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
41980f8ae99SSajid Ali     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt  mu2: d[cost]/d[mu2]\n");CHKERRQ(ierr);
42080f8ae99SSajid Ali     ierr = VecView(user.sens_mu2,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
42180f8ae99SSajid Ali   }
42280f8ae99SSajid Ali 
42380f8ae99SSajid Ali   if (sa == SA_GLOBAL) {
42480f8ae99SSajid Ali     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt  params: d[cost]/d[p], where p refers to \n\
42580f8ae99SSajid Ali                     the interlaced vector made by combining mu1,mu2\n");CHKERRQ(ierr);
42680f8ae99SSajid Ali     ierr = VecView(user.mup,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
42780f8ae99SSajid Ali   }
42880f8ae99SSajid Ali 
42980f8ae99SSajid Ali   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[cost]/d[u(t=0)]\n");CHKERRQ(ierr);
43080f8ae99SSajid Ali   ierr = VecView(user.lambda,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
43180f8ae99SSajid Ali 
43280f8ae99SSajid Ali   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43380f8ae99SSajid Ali      Free work space!
43480f8ae99SSajid Ali      All PETSc objects should be destroyed when they are no longer needed.
43580f8ae99SSajid Ali      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
43680f8ae99SSajid Ali   ierr = MatDestroy(&user.A);CHKERRQ(ierr);
43780f8ae99SSajid Ali   ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr);
43880f8ae99SSajid Ali   ierr = VecDestroy(&user.U);CHKERRQ(ierr);
43980f8ae99SSajid Ali   ierr = VecDestroy(&user.lambda);CHKERRQ(ierr);
44080f8ae99SSajid Ali   ierr = VecDestroy(&user.mup);CHKERRQ(ierr);
44180f8ae99SSajid Ali   ierr = VecDestroy(&user.mu1);CHKERRQ(ierr);
44280f8ae99SSajid Ali   ierr = VecDestroy(&user.mu2);CHKERRQ(ierr);
44380f8ae99SSajid Ali   if (sa == SA_TRACK) {
44480f8ae99SSajid Ali     ierr = VecDestroy(&user.sens_mu1);CHKERRQ(ierr);
44580f8ae99SSajid Ali     ierr = VecDestroy(&user.sens_mu2);CHKERRQ(ierr);
44680f8ae99SSajid Ali   }
44780f8ae99SSajid Ali   ierr = TSDestroy(&ts);CHKERRQ(ierr);
44880f8ae99SSajid Ali   ierr = PetscFinalize();
44980f8ae99SSajid Ali   return(ierr);
45080f8ae99SSajid Ali }
45180f8ae99SSajid Ali 
45280f8ae99SSajid Ali /*TEST
45380f8ae99SSajid Ali 
45480f8ae99SSajid Ali   test:
45580f8ae99SSajid Ali     requires: !complex
45680f8ae99SSajid Ali     suffix : track
45780f8ae99SSajid Ali     args : -sa_method track
45880f8ae99SSajid Ali 
45980f8ae99SSajid Ali   test:
46080f8ae99SSajid Ali     requires: !complex
46180f8ae99SSajid Ali     suffix : global
46280f8ae99SSajid Ali     args : -sa_method global
46380f8ae99SSajid Ali 
46480f8ae99SSajid Ali TEST*/
46580f8ae99SSajid Ali 
466