1*80f8ae99SSajid Ali static char help[] = "Performs adjoint sensitivity analysis for a van der Pol like \n\ 2*80f8ae99SSajid Ali equation with time dependent parameters using two approaches : \n\ 3*80f8ae99SSajid Ali track : track only local sensitivities at each adjoint step \n\ 4*80f8ae99SSajid Ali and accumulate them in a global array \n\ 5*80f8ae99SSajid Ali global : track parameters at all timesteps together \n\ 6*80f8ae99SSajid Ali Choose one of the two at runtime by -sa_method {track,global}. \n"; 7*80f8ae99SSajid Ali 8*80f8ae99SSajid Ali /* 9*80f8ae99SSajid Ali Concepts: TS^adjoint for time dependent parameters 10*80f8ae99SSajid Ali Concepts: Customized adjoint monitor based sensitivity tracking 11*80f8ae99SSajid Ali Concepts: All at once approach to sensitivity tracking 12*80f8ae99SSajid Ali Processors: 1 13*80f8ae99SSajid Ali */ 14*80f8ae99SSajid Ali 15*80f8ae99SSajid Ali 16*80f8ae99SSajid Ali /* 17*80f8ae99SSajid Ali Simple example to demonstrate TSAdjoint capabilities for time dependent params 18*80f8ae99SSajid Ali without integral cost terms using either a tracking or global method. 19*80f8ae99SSajid Ali 20*80f8ae99SSajid Ali Modify the Van Der Pol Eq to : 21*80f8ae99SSajid Ali [u1'] = [mu1(t)*u1] 22*80f8ae99SSajid Ali [u2'] = [mu2(t)*((1-u1^2)*u2-u1)] 23*80f8ae99SSajid Ali (with initial conditions & params independent) 24*80f8ae99SSajid Ali 25*80f8ae99SSajid Ali Define uref to be solution with initail conditions (2,-2/3), mu=(1,1e3) 26*80f8ae99SSajid Ali - u_ref : (1.5967,-1.02969) 27*80f8ae99SSajid Ali 28*80f8ae99SSajid Ali Define const function as cost = 2-norm(u - u_ref); 29*80f8ae99SSajid Ali 30*80f8ae99SSajid Ali Initialization for the adjoint TS : 31*80f8ae99SSajid Ali - dcost/dy|final_time = 2*(u-u_ref)|final_time 32*80f8ae99SSajid Ali - dcost/dp|final_time = 0 33*80f8ae99SSajid Ali 34*80f8ae99SSajid Ali The tracking method only tracks local sensitivity at each time step 35*80f8ae99SSajid Ali and accumulates these sensitivities in a global array. Since the structure 36*80f8ae99SSajid Ali of the equations being solved at each time step does not change, the jacobian 37*80f8ae99SSajid Ali wrt parameters is defined analogous to constant RHSJacobian for a liner 38*80f8ae99SSajid Ali TSSolve and the size of the jacP is independent of the number of time 39*80f8ae99SSajid Ali steps. Enable this mode of adjoint analysis by -sa_method track. 40*80f8ae99SSajid Ali 41*80f8ae99SSajid Ali The global method combines the parameters at all timesteps and tracks them 42*80f8ae99SSajid Ali together. Thus, the columns of the jacP matrix are filled dependent upon the 43*80f8ae99SSajid Ali time step. Also, the dimensions of the jacP matrix now depend upon the number 44*80f8ae99SSajid Ali of time steps. Enable this mode of adjoint analysis by -sa_method global. 45*80f8ae99SSajid Ali 46*80f8ae99SSajid Ali Since the equations here have parameters at predefined time steps, this 47*80f8ae99SSajid Ali example should be run with non adaptive time stepping solvers only. This 48*80f8ae99SSajid Ali can be ensured by -ts_adapt_type none (which is the default behavior only 49*80f8ae99SSajid Ali for certain TS solvers like TSCN. If using an explicit TS like TSRK, 50*80f8ae99SSajid Ali please be sure to add the aforementioned option to disable adaptive 51*80f8ae99SSajid Ali timestepping.) 52*80f8ae99SSajid Ali */ 53*80f8ae99SSajid Ali 54*80f8ae99SSajid Ali /* 55*80f8ae99SSajid Ali Include "petscts.h" so that we can use TS solvers. Note that this file 56*80f8ae99SSajid Ali automatically includes: 57*80f8ae99SSajid Ali petscsys.h - base PETSc routines petscvec.h - vectors 58*80f8ae99SSajid Ali petscmat.h - matrices 59*80f8ae99SSajid Ali petscis.h - index sets petscksp.h - Krylov subspace methods 60*80f8ae99SSajid Ali petscviewer.h - viewers petscpc.h - preconditioners 61*80f8ae99SSajid Ali petscksp.h - linear solvers petscsnes.h - nonlinear solvers 62*80f8ae99SSajid Ali */ 63*80f8ae99SSajid Ali #include <petscts.h> 64*80f8ae99SSajid Ali 65*80f8ae99SSajid Ali extern PetscErrorCode RHSFunction(TS ,PetscReal ,Vec ,Vec ,void *); 66*80f8ae99SSajid Ali extern PetscErrorCode RHSJacobian(TS ,PetscReal ,Vec ,Mat ,Mat ,void *); 67*80f8ae99SSajid Ali extern PetscErrorCode RHSJacobianP_track(TS ,PetscReal ,Vec ,Mat ,void *); 68*80f8ae99SSajid Ali extern PetscErrorCode RHSJacobianP_global(TS ,PetscReal ,Vec ,Mat ,void *); 69*80f8ae99SSajid Ali extern PetscErrorCode Monitor(TS ,PetscInt ,PetscReal ,Vec ,void *); 70*80f8ae99SSajid Ali extern PetscErrorCode AdjointMonitor(TS ,PetscInt ,PetscReal ,Vec ,PetscInt ,Vec *, Vec *,void *); 71*80f8ae99SSajid Ali 72*80f8ae99SSajid Ali /* 73*80f8ae99SSajid Ali User-defined application context - contains data needed by the 74*80f8ae99SSajid Ali application-provided call-back routines. 75*80f8ae99SSajid Ali */ 76*80f8ae99SSajid Ali 77*80f8ae99SSajid Ali typedef struct { 78*80f8ae99SSajid Ali /*------------- Forward solve data structures --------------*/ 79*80f8ae99SSajid Ali PetscInt max_steps; /* number of steps to run ts for */ 80*80f8ae99SSajid Ali PetscReal final_time; /* final time to integrate to*/ 81*80f8ae99SSajid Ali PetscReal time_step; /* ts integration time step */ 82*80f8ae99SSajid Ali Vec mu1; /* time dependent params */ 83*80f8ae99SSajid Ali Vec mu2; /* time dependent params */ 84*80f8ae99SSajid Ali Vec U; /* solution vector */ 85*80f8ae99SSajid Ali Mat A; /* Jacobian matrix */ 86*80f8ae99SSajid Ali 87*80f8ae99SSajid Ali /*------------- Adjoint solve data structures --------------*/ 88*80f8ae99SSajid Ali Mat Jacp; /* JacobianP matrix */ 89*80f8ae99SSajid Ali Vec lambda; /* adjoint variable */ 90*80f8ae99SSajid Ali Vec mup; /* adjoint variable */ 91*80f8ae99SSajid Ali 92*80f8ae99SSajid Ali /*------------- Global accumation vecs for monitor based tracking --------------*/ 93*80f8ae99SSajid Ali Vec sens_mu1; /* global sensitivity accumulation */ 94*80f8ae99SSajid Ali Vec sens_mu2; /* global sensitivity accumulation */ 95*80f8ae99SSajid Ali PetscInt adj_idx; /* to keep track of adjoint solve index */ 96*80f8ae99SSajid Ali } AppCtx; 97*80f8ae99SSajid Ali 98*80f8ae99SSajid Ali typedef enum {SA_TRACK, SA_GLOBAL} SAMethod; 99*80f8ae99SSajid Ali static const char *const SAMethods[] = {"TRACK","GLOBAL","SAMethod","SA_",0}; 100*80f8ae99SSajid Ali 101*80f8ae99SSajid Ali /* ----------------------- Explicit form of the ODE -------------------- */ 102*80f8ae99SSajid Ali 103*80f8ae99SSajid Ali PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,void *ctx) 104*80f8ae99SSajid Ali { 105*80f8ae99SSajid Ali PetscErrorCode ierr; 106*80f8ae99SSajid Ali AppCtx *user = (AppCtx*) ctx; 107*80f8ae99SSajid Ali PetscScalar *f; 108*80f8ae99SSajid Ali PetscInt curr_step; 109*80f8ae99SSajid Ali const PetscScalar *u; 110*80f8ae99SSajid Ali const PetscScalar *mu1; 111*80f8ae99SSajid Ali const PetscScalar *mu2; 112*80f8ae99SSajid Ali 113*80f8ae99SSajid Ali PetscFunctionBeginUser; 114*80f8ae99SSajid Ali ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr); 115*80f8ae99SSajid Ali ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 116*80f8ae99SSajid Ali ierr = VecGetArrayRead(user->mu1,&mu1);CHKERRQ(ierr); 117*80f8ae99SSajid Ali ierr = VecGetArrayRead(user->mu2,&mu2);CHKERRQ(ierr); 118*80f8ae99SSajid Ali ierr = VecGetArray(F,&f);CHKERRQ(ierr); 119*80f8ae99SSajid Ali f[0] = mu1[curr_step]*u[1]; 120*80f8ae99SSajid Ali f[1] = mu2[curr_step]*((1.-u[0]*u[0])*u[1]-u[0]); 121*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 122*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(user->mu1,&mu1);CHKERRQ(ierr); 123*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(user->mu2,&mu2);CHKERRQ(ierr); 124*80f8ae99SSajid Ali ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 125*80f8ae99SSajid Ali PetscFunctionReturn(0); 126*80f8ae99SSajid Ali } 127*80f8ae99SSajid Ali 128*80f8ae99SSajid Ali PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,void *ctx) 129*80f8ae99SSajid Ali { 130*80f8ae99SSajid Ali PetscErrorCode ierr; 131*80f8ae99SSajid Ali AppCtx *user = (AppCtx*) ctx; 132*80f8ae99SSajid Ali PetscInt rowcol[] = {0,1}; 133*80f8ae99SSajid Ali PetscScalar J[2][2]; 134*80f8ae99SSajid Ali PetscInt curr_step; 135*80f8ae99SSajid Ali const PetscScalar *u; 136*80f8ae99SSajid Ali const PetscScalar *mu1; 137*80f8ae99SSajid Ali const PetscScalar *mu2; 138*80f8ae99SSajid Ali 139*80f8ae99SSajid Ali PetscFunctionBeginUser; 140*80f8ae99SSajid Ali ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr); 141*80f8ae99SSajid Ali ierr = VecGetArrayRead(user->mu1,&mu1);CHKERRQ(ierr); 142*80f8ae99SSajid Ali ierr = VecGetArrayRead(user->mu2,&mu2);CHKERRQ(ierr); 143*80f8ae99SSajid Ali ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 144*80f8ae99SSajid Ali J[0][0] = 0; 145*80f8ae99SSajid Ali J[1][0] = -mu2[curr_step]*(2.0*u[1]*u[0]+1.); 146*80f8ae99SSajid Ali J[0][1] = mu1[curr_step]; 147*80f8ae99SSajid Ali J[1][1] = mu2[curr_step]*(1.0-u[0]*u[0]); 148*80f8ae99SSajid Ali ierr = MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 149*80f8ae99SSajid Ali ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 150*80f8ae99SSajid Ali ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 151*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 152*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(user->mu1,&mu1);CHKERRQ(ierr); 153*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(user->mu2,&mu2);CHKERRQ(ierr); 154*80f8ae99SSajid Ali PetscFunctionReturn(0); 155*80f8ae99SSajid Ali } 156*80f8ae99SSajid Ali 157*80f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for tracking method ------------------ */ 158*80f8ae99SSajid Ali 159*80f8ae99SSajid Ali PetscErrorCode RHSJacobianP_track(TS ts,PetscReal t,Vec U,Mat A,void *ctx) 160*80f8ae99SSajid Ali { 161*80f8ae99SSajid Ali PetscErrorCode ierr; 162*80f8ae99SSajid Ali PetscInt row[] = {0,1},col[] = {0,1}; 163*80f8ae99SSajid Ali PetscScalar J[2][2]; 164*80f8ae99SSajid Ali const PetscScalar *u; 165*80f8ae99SSajid Ali 166*80f8ae99SSajid Ali PetscFunctionBeginUser; 167*80f8ae99SSajid Ali ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 168*80f8ae99SSajid Ali J[0][0] = u[1]; 169*80f8ae99SSajid Ali J[1][0] = 0; 170*80f8ae99SSajid Ali J[0][1] = 0; 171*80f8ae99SSajid Ali J[1][1] = (1.-u[0]*u[0])*u[1]-u[0]; 172*80f8ae99SSajid Ali ierr = MatSetValues(A,2,row,2,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 173*80f8ae99SSajid Ali ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 174*80f8ae99SSajid Ali ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 175*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 176*80f8ae99SSajid Ali PetscFunctionReturn(0); 177*80f8ae99SSajid Ali } 178*80f8ae99SSajid Ali 179*80f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for global method ------------------ */ 180*80f8ae99SSajid Ali 181*80f8ae99SSajid Ali PetscErrorCode RHSJacobianP_global(TS ts,PetscReal t,Vec U,Mat A,void *ctx) 182*80f8ae99SSajid Ali { 183*80f8ae99SSajid Ali PetscErrorCode ierr; 184*80f8ae99SSajid Ali PetscInt row[] = {0,1},col[] = {0,1}; 185*80f8ae99SSajid Ali PetscScalar J[2][2]; 186*80f8ae99SSajid Ali const PetscScalar *u; 187*80f8ae99SSajid Ali PetscInt curr_step; 188*80f8ae99SSajid Ali 189*80f8ae99SSajid Ali PetscFunctionBeginUser; 190*80f8ae99SSajid Ali ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr); 191*80f8ae99SSajid Ali ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 192*80f8ae99SSajid Ali J[0][0] = u[1]; 193*80f8ae99SSajid Ali J[1][0] = 0; 194*80f8ae99SSajid Ali J[0][1] = 0; 195*80f8ae99SSajid Ali J[1][1] = (1.-u[0]*u[0])*u[1]-u[0]; 196*80f8ae99SSajid Ali col[0] = (curr_step)*2; 197*80f8ae99SSajid Ali col[1] = (curr_step)*2+1; 198*80f8ae99SSajid Ali ierr = MatSetValues(A,2,row,2,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 199*80f8ae99SSajid Ali ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 200*80f8ae99SSajid Ali ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 201*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 202*80f8ae99SSajid Ali PetscFunctionReturn(0); 203*80f8ae99SSajid Ali } 204*80f8ae99SSajid Ali 205*80f8ae99SSajid Ali /* Dump solution to console if called */ 206*80f8ae99SSajid Ali PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec U,void *ctx) 207*80f8ae99SSajid Ali { 208*80f8ae99SSajid Ali PetscErrorCode ierr; 209*80f8ae99SSajid Ali 210*80f8ae99SSajid Ali PetscFunctionBeginUser; 211*80f8ae99SSajid Ali ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Solution at time %e is \n", t);CHKERRQ(ierr); 212*80f8ae99SSajid Ali ierr = VecView(U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 213*80f8ae99SSajid Ali PetscFunctionReturn(0); 214*80f8ae99SSajid Ali } 215*80f8ae99SSajid Ali 216*80f8ae99SSajid Ali /* Customized adjoint monitor to keep track of local 217*80f8ae99SSajid Ali sensitivities by storing them in a global sensitivity array. 218*80f8ae99SSajid Ali Note : This routine is only used for the tracking method. */ 219*80f8ae99SSajid Ali PetscErrorCode AdjointMonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *ctx) 220*80f8ae99SSajid Ali { 221*80f8ae99SSajid Ali PetscErrorCode ierr; 222*80f8ae99SSajid Ali AppCtx *user = (AppCtx*) ctx; 223*80f8ae99SSajid Ali PetscInt curr_step; 224*80f8ae99SSajid Ali PetscScalar *sensmu1_glob; 225*80f8ae99SSajid Ali PetscScalar *sensmu2_glob; 226*80f8ae99SSajid Ali const PetscScalar *sensmu_loc; 227*80f8ae99SSajid Ali 228*80f8ae99SSajid Ali PetscFunctionBeginUser; 229*80f8ae99SSajid Ali ierr = TSGetStepNumber(ts,&curr_step);CHKERRQ(ierr); 230*80f8ae99SSajid Ali /* Note that we skip the first call to the monitor in the adjoint 231*80f8ae99SSajid Ali solve since the sensitivities are already set (during 232*80f8ae99SSajid Ali initialization of adjoint vectors). 233*80f8ae99SSajid Ali We also note that each indvidial TSAdjointSolve calls the monitor 234*80f8ae99SSajid Ali twice, once at the step it is integrating from and once at the step 235*80f8ae99SSajid Ali it integrates to. Only the second call is useful for transferring 236*80f8ae99SSajid Ali local sensitivities to the global array. */ 237*80f8ae99SSajid Ali if (curr_step == user->adj_idx) { 238*80f8ae99SSajid Ali PetscFunctionReturn(0); 239*80f8ae99SSajid Ali } else { 240*80f8ae99SSajid Ali ierr = VecGetArrayRead(*mu,&sensmu_loc);CHKERRQ(ierr); 241*80f8ae99SSajid Ali ierr = VecGetArray(user->sens_mu1,&sensmu1_glob);CHKERRQ(ierr); 242*80f8ae99SSajid Ali ierr = VecGetArray(user->sens_mu2,&sensmu2_glob);CHKERRQ(ierr); 243*80f8ae99SSajid Ali sensmu1_glob[curr_step] = sensmu_loc[0]; 244*80f8ae99SSajid Ali sensmu2_glob[curr_step] = sensmu_loc[1]; 245*80f8ae99SSajid Ali ierr = VecRestoreArray(user->sens_mu1,&sensmu1_glob);CHKERRQ(ierr); 246*80f8ae99SSajid Ali ierr = VecRestoreArray(user->sens_mu2,&sensmu2_glob);CHKERRQ(ierr); 247*80f8ae99SSajid Ali ierr = VecRestoreArrayRead(*mu,&sensmu_loc);CHKERRQ(ierr); 248*80f8ae99SSajid Ali PetscFunctionReturn(0); 249*80f8ae99SSajid Ali } 250*80f8ae99SSajid Ali } 251*80f8ae99SSajid Ali 252*80f8ae99SSajid Ali 253*80f8ae99SSajid Ali int main(int argc,char **argv) 254*80f8ae99SSajid Ali { 255*80f8ae99SSajid Ali TS ts; 256*80f8ae99SSajid Ali AppCtx user; 257*80f8ae99SSajid Ali PetscScalar *x_ptr,*y_ptr,*u_ptr; 258*80f8ae99SSajid Ali PetscMPIInt size; 259*80f8ae99SSajid Ali PetscBool monitor = PETSC_FALSE; 260*80f8ae99SSajid Ali SAMethod sa = SA_GLOBAL; 261*80f8ae99SSajid Ali PetscErrorCode ierr; 262*80f8ae99SSajid Ali 263*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 264*80f8ae99SSajid Ali Initialize program 265*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 266*80f8ae99SSajid Ali ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr; 267*80f8ae99SSajid Ali ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr); 268*80f8ae99SSajid Ali if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!"); 269*80f8ae99SSajid Ali 270*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 271*80f8ae99SSajid Ali Set runtime options 272*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 273*80f8ae99SSajid Ali ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"SA analysis options.","");CHKERRQ(ierr);{ 274*80f8ae99SSajid Ali ierr = PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);CHKERRQ(ierr); 275*80f8ae99SSajid Ali ierr = PetscOptionsEnum("-sa_method","Sensitivity analysis method (track or global)","",SAMethods,(PetscEnum)sa,(PetscEnum*)&sa,NULL);CHKERRQ(ierr); 276*80f8ae99SSajid Ali } 277*80f8ae99SSajid Ali ierr = PetscOptionsEnd();CHKERRQ(ierr); 278*80f8ae99SSajid Ali 279*80f8ae99SSajid Ali user.final_time = 0.1; 280*80f8ae99SSajid Ali user.max_steps = 5; 281*80f8ae99SSajid Ali user.time_step = user.final_time/user.max_steps; 282*80f8ae99SSajid Ali 283*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 284*80f8ae99SSajid Ali Create necessary matrix and vectors for forward solve. 285*80f8ae99SSajid Ali Create Jacp matrix for adjoint solve. 286*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 287*80f8ae99SSajid Ali ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.mu1);CHKERRQ(ierr); 288*80f8ae99SSajid Ali ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.mu2);CHKERRQ(ierr); 289*80f8ae99SSajid Ali ierr = VecSet(user.mu1,1.25);CHKERRQ(ierr); 290*80f8ae99SSajid Ali ierr = VecSet(user.mu2,1.0e2);CHKERRQ(ierr); 291*80f8ae99SSajid Ali 292*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 293*80f8ae99SSajid Ali For tracking method : create the global sensitivity array to 294*80f8ae99SSajid Ali accumulate sensitivity with respect to parameters at each step. 295*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 296*80f8ae99SSajid Ali if (sa == SA_TRACK) { 297*80f8ae99SSajid Ali ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.sens_mu1);CHKERRQ(ierr); 298*80f8ae99SSajid Ali ierr = VecCreateSeq(PETSC_COMM_WORLD,user.max_steps,&user.sens_mu2);CHKERRQ(ierr); 299*80f8ae99SSajid Ali } 300*80f8ae99SSajid Ali 301*80f8ae99SSajid Ali ierr = MatCreate(PETSC_COMM_WORLD,&user.A);CHKERRQ(ierr); 302*80f8ae99SSajid Ali ierr = MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 303*80f8ae99SSajid Ali ierr = MatSetFromOptions(user.A);CHKERRQ(ierr); 304*80f8ae99SSajid Ali ierr = MatSetUp(user.A);CHKERRQ(ierr); 305*80f8ae99SSajid Ali ierr = MatCreateVecs(user.A,&user.U,NULL);CHKERRQ(ierr); 306*80f8ae99SSajid Ali 307*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 308*80f8ae99SSajid Ali Note that the dimensions of the Jacp matrix depend upon the 309*80f8ae99SSajid Ali sensitivity analysis method being used ! 310*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 311*80f8ae99SSajid Ali ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr); 312*80f8ae99SSajid Ali if (sa == SA_TRACK) { 313*80f8ae99SSajid Ali ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,2);CHKERRQ(ierr); 314*80f8ae99SSajid Ali } 315*80f8ae99SSajid Ali if (sa == SA_GLOBAL) { 316*80f8ae99SSajid Ali ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,user.max_steps*2);CHKERRQ(ierr); 317*80f8ae99SSajid Ali } 318*80f8ae99SSajid Ali ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr); 319*80f8ae99SSajid Ali ierr = MatSetUp(user.Jacp);CHKERRQ(ierr); 320*80f8ae99SSajid Ali 321*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 322*80f8ae99SSajid Ali Create timestepping solver context 323*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 324*80f8ae99SSajid Ali ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 325*80f8ae99SSajid Ali ierr = TSSetEquationType(ts,TS_EQ_ODE_EXPLICIT);CHKERRQ(ierr); 326*80f8ae99SSajid Ali ierr = TSSetType(ts,TSCN);CHKERRQ(ierr); 327*80f8ae99SSajid Ali 328*80f8ae99SSajid Ali ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&user);CHKERRQ(ierr); 329*80f8ae99SSajid Ali ierr = TSSetRHSJacobian(ts,user.A,user.A,RHSJacobian,&user);CHKERRQ(ierr); 330*80f8ae99SSajid Ali if (sa == SA_TRACK) { 331*80f8ae99SSajid Ali ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP_track,&user);CHKERRQ(ierr); 332*80f8ae99SSajid Ali } 333*80f8ae99SSajid Ali if (sa == SA_GLOBAL) { 334*80f8ae99SSajid Ali ierr = TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP_global,&user);CHKERRQ(ierr); 335*80f8ae99SSajid Ali } 336*80f8ae99SSajid Ali 337*80f8ae99SSajid Ali ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 338*80f8ae99SSajid Ali ierr = TSSetMaxTime(ts,user.final_time);CHKERRQ(ierr); 339*80f8ae99SSajid Ali ierr = TSSetTimeStep(ts,user.final_time/user.max_steps);CHKERRQ(ierr); 340*80f8ae99SSajid Ali 341*80f8ae99SSajid Ali if (monitor) { 342*80f8ae99SSajid Ali ierr = TSMonitorSet(ts,Monitor,&user,NULL);CHKERRQ(ierr); 343*80f8ae99SSajid Ali } 344*80f8ae99SSajid Ali if (sa == SA_TRACK) { 345*80f8ae99SSajid Ali ierr = TSAdjointMonitorSet(ts,AdjointMonitor,&user,NULL);CHKERRQ(ierr); 346*80f8ae99SSajid Ali } 347*80f8ae99SSajid Ali 348*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 349*80f8ae99SSajid Ali Set initial conditions 350*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 351*80f8ae99SSajid Ali ierr = VecGetArray(user.U,&x_ptr);CHKERRQ(ierr); 352*80f8ae99SSajid Ali x_ptr[0] = 2.0; 353*80f8ae99SSajid Ali x_ptr[1] = -2.0/3.0; 354*80f8ae99SSajid Ali ierr = VecRestoreArray(user.U,&x_ptr);CHKERRQ(ierr); 355*80f8ae99SSajid Ali 356*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 357*80f8ae99SSajid Ali Save trajectory of solution so that TSAdjointSolve() may be used 358*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 359*80f8ae99SSajid Ali ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 360*80f8ae99SSajid Ali 361*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 362*80f8ae99SSajid Ali Set runtime options 363*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 364*80f8ae99SSajid Ali ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 365*80f8ae99SSajid Ali 366*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 367*80f8ae99SSajid Ali Execute forward model and print solution. 368*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 369*80f8ae99SSajid Ali ierr = TSSolve(ts,user.U);CHKERRQ(ierr); 370*80f8ae99SSajid Ali ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Solution of forward TS :\n");CHKERRQ(ierr); 371*80f8ae99SSajid Ali ierr = VecView(user.U,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 372*80f8ae99SSajid Ali ierr = PetscPrintf(PETSC_COMM_WORLD,"\n Forward TS solve successfull! Adjoint run begins!\n");CHKERRQ(ierr); 373*80f8ae99SSajid Ali 374*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 375*80f8ae99SSajid Ali Adjoint model starts here! Create adjoint vectors. 376*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 377*80f8ae99SSajid Ali ierr = MatCreateVecs(user.A,&user.lambda,NULL);CHKERRQ(ierr); 378*80f8ae99SSajid Ali ierr = MatCreateVecs(user.Jacp,&user.mup,NULL);CHKERRQ(ierr); 379*80f8ae99SSajid Ali 380*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 381*80f8ae99SSajid Ali Set initial conditions for the adjoint vector 382*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 383*80f8ae99SSajid Ali ierr = VecGetArray(user.U,&u_ptr);CHKERRQ(ierr); 384*80f8ae99SSajid Ali ierr = VecGetArray(user.lambda,&y_ptr);CHKERRQ(ierr); 385*80f8ae99SSajid Ali y_ptr[0] = 2*(u_ptr[0] - 1.5967); 386*80f8ae99SSajid Ali y_ptr[1] = 2*(u_ptr[1] - -(1.02969)); 387*80f8ae99SSajid Ali ierr = VecRestoreArray(user.lambda,&y_ptr);CHKERRQ(ierr); 388*80f8ae99SSajid Ali ierr = VecRestoreArray(user.U,&y_ptr);CHKERRQ(ierr); 389*80f8ae99SSajid Ali ierr = VecSet(user.mup,0);CHKERRQ(ierr); 390*80f8ae99SSajid Ali 391*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 392*80f8ae99SSajid Ali Set number of cost functions. 393*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 394*80f8ae99SSajid Ali ierr = TSSetCostGradients(ts,1,&user.lambda,&user.mup);CHKERRQ(ierr); 395*80f8ae99SSajid Ali 396*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 397*80f8ae99SSajid Ali The adjoint vector mup has to be reset for each adjoint step when 398*80f8ae99SSajid Ali using the tracking method as we want to treat the parameters at each 399*80f8ae99SSajid Ali time step one at a time and prevent accumulation of the sensitivities 400*80f8ae99SSajid Ali from parameters at previous time steps. 401*80f8ae99SSajid Ali This is not necessary for the global method as each time dependent 402*80f8ae99SSajid Ali parameter is treated as an independent parameter. 403*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 404*80f8ae99SSajid Ali if (sa == SA_TRACK) { 405*80f8ae99SSajid Ali for (user.adj_idx=user.max_steps; user.adj_idx>0; user.adj_idx--) { 406*80f8ae99SSajid Ali ierr = VecSet(user.mup,0);CHKERRQ(ierr); 407*80f8ae99SSajid Ali ierr = TSAdjointSetSteps(ts, 1);CHKERRQ(ierr); 408*80f8ae99SSajid Ali ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 409*80f8ae99SSajid Ali } 410*80f8ae99SSajid Ali } 411*80f8ae99SSajid Ali if (sa == SA_GLOBAL) { 412*80f8ae99SSajid Ali ierr = TSAdjointSolve(ts);CHKERRQ(ierr); 413*80f8ae99SSajid Ali } 414*80f8ae99SSajid Ali 415*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 416*80f8ae99SSajid Ali Dispaly adjoint sensitivities wrt parameters and initial conditions 417*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 418*80f8ae99SSajid Ali if (sa == SA_TRACK) { 419*80f8ae99SSajid Ali ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt mu1: d[cost]/d[mu1]\n");CHKERRQ(ierr); 420*80f8ae99SSajid Ali ierr = VecView(user.sens_mu1,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 421*80f8ae99SSajid Ali ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt mu2: d[cost]/d[mu2]\n");CHKERRQ(ierr); 422*80f8ae99SSajid Ali ierr = VecView(user.sens_mu2,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 423*80f8ae99SSajid Ali } 424*80f8ae99SSajid Ali 425*80f8ae99SSajid Ali if (sa == SA_GLOBAL) { 426*80f8ae99SSajid Ali ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt params: d[cost]/d[p], where p refers to \n\ 427*80f8ae99SSajid Ali the interlaced vector made by combining mu1,mu2\n");CHKERRQ(ierr); 428*80f8ae99SSajid Ali ierr = VecView(user.mup,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 429*80f8ae99SSajid Ali } 430*80f8ae99SSajid Ali 431*80f8ae99SSajid Ali ierr = PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[cost]/d[u(t=0)]\n");CHKERRQ(ierr); 432*80f8ae99SSajid Ali ierr = VecView(user.lambda,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 433*80f8ae99SSajid Ali 434*80f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 435*80f8ae99SSajid Ali Free work space! 436*80f8ae99SSajid Ali All PETSc objects should be destroyed when they are no longer needed. 437*80f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 438*80f8ae99SSajid Ali ierr = MatDestroy(&user.A);CHKERRQ(ierr); 439*80f8ae99SSajid Ali ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr); 440*80f8ae99SSajid Ali ierr = VecDestroy(&user.U);CHKERRQ(ierr); 441*80f8ae99SSajid Ali ierr = VecDestroy(&user.lambda);CHKERRQ(ierr); 442*80f8ae99SSajid Ali ierr = VecDestroy(&user.mup);CHKERRQ(ierr); 443*80f8ae99SSajid Ali ierr = VecDestroy(&user.mu1);CHKERRQ(ierr); 444*80f8ae99SSajid Ali ierr = VecDestroy(&user.mu2);CHKERRQ(ierr); 445*80f8ae99SSajid Ali if (sa == SA_TRACK) { 446*80f8ae99SSajid Ali ierr = VecDestroy(&user.sens_mu1);CHKERRQ(ierr); 447*80f8ae99SSajid Ali ierr = VecDestroy(&user.sens_mu2);CHKERRQ(ierr); 448*80f8ae99SSajid Ali } 449*80f8ae99SSajid Ali ierr = TSDestroy(&ts);CHKERRQ(ierr); 450*80f8ae99SSajid Ali ierr = PetscFinalize(); 451*80f8ae99SSajid Ali return(ierr); 452*80f8ae99SSajid Ali } 453*80f8ae99SSajid Ali 454*80f8ae99SSajid Ali 455*80f8ae99SSajid Ali /*TEST 456*80f8ae99SSajid Ali 457*80f8ae99SSajid Ali test: 458*80f8ae99SSajid Ali requires: !complex 459*80f8ae99SSajid Ali suffix : track 460*80f8ae99SSajid Ali args : -sa_method track 461*80f8ae99SSajid Ali 462*80f8ae99SSajid Ali test: 463*80f8ae99SSajid Ali requires: !complex 464*80f8ae99SSajid Ali suffix : global 465*80f8ae99SSajid Ali args : -sa_method global 466*80f8ae99SSajid Ali 467*80f8ae99SSajid Ali TEST*/ 468*80f8ae99SSajid Ali 469