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 Simple example to demonstrate TSAdjoint capabilities for time dependent params 1080f8ae99SSajid Ali without integral cost terms using either a tracking or global method. 1180f8ae99SSajid Ali 1280f8ae99SSajid Ali Modify the Van Der Pol Eq to : 1380f8ae99SSajid Ali [u1'] = [mu1(t)*u1] 1480f8ae99SSajid Ali [u2'] = [mu2(t)*((1-u1^2)*u2-u1)] 1580f8ae99SSajid Ali (with initial conditions & params independent) 1680f8ae99SSajid Ali 1780f8ae99SSajid Ali Define uref to be solution with initail conditions (2,-2/3), mu=(1,1e3) 1880f8ae99SSajid Ali - u_ref : (1.5967,-1.02969) 1980f8ae99SSajid Ali 2080f8ae99SSajid Ali Define const function as cost = 2-norm(u - u_ref); 2180f8ae99SSajid Ali 2280f8ae99SSajid Ali Initialization for the adjoint TS : 2380f8ae99SSajid Ali - dcost/dy|final_time = 2*(u-u_ref)|final_time 2480f8ae99SSajid Ali - dcost/dp|final_time = 0 2580f8ae99SSajid Ali 2680f8ae99SSajid Ali The tracking method only tracks local sensitivity at each time step 2780f8ae99SSajid Ali and accumulates these sensitivities in a global array. Since the structure 2880f8ae99SSajid Ali of the equations being solved at each time step does not change, the jacobian 2980f8ae99SSajid Ali wrt parameters is defined analogous to constant RHSJacobian for a liner 3080f8ae99SSajid Ali TSSolve and the size of the jacP is independent of the number of time 3180f8ae99SSajid Ali steps. Enable this mode of adjoint analysis by -sa_method track. 3280f8ae99SSajid Ali 3380f8ae99SSajid Ali The global method combines the parameters at all timesteps and tracks them 3480f8ae99SSajid Ali together. Thus, the columns of the jacP matrix are filled dependent upon the 3580f8ae99SSajid Ali time step. Also, the dimensions of the jacP matrix now depend upon the number 3680f8ae99SSajid Ali of time steps. Enable this mode of adjoint analysis by -sa_method global. 3780f8ae99SSajid Ali 3880f8ae99SSajid Ali Since the equations here have parameters at predefined time steps, this 3980f8ae99SSajid Ali example should be run with non adaptive time stepping solvers only. This 4080f8ae99SSajid Ali can be ensured by -ts_adapt_type none (which is the default behavior only 4180f8ae99SSajid Ali for certain TS solvers like TSCN. If using an explicit TS like TSRK, 4280f8ae99SSajid Ali please be sure to add the aforementioned option to disable adaptive 4380f8ae99SSajid Ali timestepping.) 4480f8ae99SSajid Ali */ 4580f8ae99SSajid Ali 4680f8ae99SSajid Ali /* 4780f8ae99SSajid Ali Include "petscts.h" so that we can use TS solvers. Note that this file 4880f8ae99SSajid Ali automatically includes: 4980f8ae99SSajid Ali petscsys.h - base PETSc routines petscvec.h - vectors 5080f8ae99SSajid Ali petscmat.h - matrices 5180f8ae99SSajid Ali petscis.h - index sets petscksp.h - Krylov subspace methods 5280f8ae99SSajid Ali petscviewer.h - viewers petscpc.h - preconditioners 5380f8ae99SSajid Ali petscksp.h - linear solvers petscsnes.h - nonlinear solvers 5480f8ae99SSajid Ali */ 5580f8ae99SSajid Ali #include <petscts.h> 5680f8ae99SSajid Ali 5780f8ae99SSajid Ali extern PetscErrorCode RHSFunction(TS, PetscReal, Vec, Vec, void *); 5880f8ae99SSajid Ali extern PetscErrorCode RHSJacobian(TS, PetscReal, Vec, Mat, Mat, void *); 5980f8ae99SSajid Ali extern PetscErrorCode RHSJacobianP_track(TS, PetscReal, Vec, Mat, void *); 6080f8ae99SSajid Ali extern PetscErrorCode RHSJacobianP_global(TS, PetscReal, Vec, Mat, void *); 6180f8ae99SSajid Ali extern PetscErrorCode Monitor(TS, PetscInt, PetscReal, Vec, void *); 6280f8ae99SSajid Ali extern PetscErrorCode AdjointMonitor(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *); 6380f8ae99SSajid Ali 6480f8ae99SSajid Ali /* 6580f8ae99SSajid Ali User-defined application context - contains data needed by the 6680f8ae99SSajid Ali application-provided call-back routines. 6780f8ae99SSajid Ali */ 6880f8ae99SSajid Ali 6980f8ae99SSajid Ali typedef struct { 7080f8ae99SSajid Ali /*------------- Forward solve data structures --------------*/ 7180f8ae99SSajid Ali PetscInt max_steps; /* number of steps to run ts for */ 7280f8ae99SSajid Ali PetscReal final_time; /* final time to integrate to*/ 7380f8ae99SSajid Ali PetscReal time_step; /* ts integration time step */ 7480f8ae99SSajid Ali Vec mu1; /* time dependent params */ 7580f8ae99SSajid Ali Vec mu2; /* time dependent params */ 7680f8ae99SSajid Ali Vec U; /* solution vector */ 7780f8ae99SSajid Ali Mat A; /* Jacobian matrix */ 7880f8ae99SSajid Ali 7980f8ae99SSajid Ali /*------------- Adjoint solve data structures --------------*/ 8080f8ae99SSajid Ali Mat Jacp; /* JacobianP matrix */ 8180f8ae99SSajid Ali Vec lambda; /* adjoint variable */ 8280f8ae99SSajid Ali Vec mup; /* adjoint variable */ 8380f8ae99SSajid Ali 8480f8ae99SSajid Ali /*------------- Global accumation vecs for monitor based tracking --------------*/ 8580f8ae99SSajid Ali Vec sens_mu1; /* global sensitivity accumulation */ 8680f8ae99SSajid Ali Vec sens_mu2; /* global sensitivity accumulation */ 8780f8ae99SSajid Ali PetscInt adj_idx; /* to keep track of adjoint solve index */ 8880f8ae99SSajid Ali } AppCtx; 8980f8ae99SSajid Ali 909371c9d4SSatish Balay typedef enum { 919371c9d4SSatish Balay SA_TRACK, 929371c9d4SSatish Balay SA_GLOBAL 939371c9d4SSatish Balay } SAMethod; 9480f8ae99SSajid Ali static const char *const SAMethods[] = {"TRACK", "GLOBAL", "SAMethod", "SA_", 0}; 9580f8ae99SSajid Ali 9680f8ae99SSajid Ali /* ----------------------- Explicit form of the ODE -------------------- */ 9780f8ae99SSajid Ali 989371c9d4SSatish Balay PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, void *ctx) { 9980f8ae99SSajid Ali AppCtx *user = (AppCtx *)ctx; 10080f8ae99SSajid Ali PetscScalar *f; 10180f8ae99SSajid Ali PetscInt curr_step; 10280f8ae99SSajid Ali const PetscScalar *u; 10380f8ae99SSajid Ali const PetscScalar *mu1; 10480f8ae99SSajid Ali const PetscScalar *mu2; 10580f8ae99SSajid Ali 10680f8ae99SSajid Ali PetscFunctionBeginUser; 1079566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &curr_step)); 1089566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 1099566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user->mu1, &mu1)); 1109566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user->mu2, &mu2)); 1119566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 11280f8ae99SSajid Ali f[0] = mu1[curr_step] * u[1]; 11380f8ae99SSajid Ali f[1] = mu2[curr_step] * ((1. - u[0] * u[0]) * u[1] - u[0]); 1149566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1159566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user->mu1, &mu1)); 1169566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user->mu2, &mu2)); 1179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 11880f8ae99SSajid Ali PetscFunctionReturn(0); 11980f8ae99SSajid Ali } 12080f8ae99SSajid Ali 1219371c9d4SSatish Balay PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, void *ctx) { 12280f8ae99SSajid Ali AppCtx *user = (AppCtx *)ctx; 12380f8ae99SSajid Ali PetscInt rowcol[] = {0, 1}; 12480f8ae99SSajid Ali PetscScalar J[2][2]; 12580f8ae99SSajid Ali PetscInt curr_step; 12680f8ae99SSajid Ali const PetscScalar *u; 12780f8ae99SSajid Ali const PetscScalar *mu1; 12880f8ae99SSajid Ali const PetscScalar *mu2; 12980f8ae99SSajid Ali 13080f8ae99SSajid Ali PetscFunctionBeginUser; 1319566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &curr_step)); 1329566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user->mu1, &mu1)); 1339566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(user->mu2, &mu2)); 1349566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 13580f8ae99SSajid Ali J[0][0] = 0; 13680f8ae99SSajid Ali J[1][0] = -mu2[curr_step] * (2.0 * u[1] * u[0] + 1.); 13780f8ae99SSajid Ali J[0][1] = mu1[curr_step]; 13880f8ae99SSajid Ali J[1][1] = mu2[curr_step] * (1.0 - u[0] * u[0]); 1399566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 1409566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1419566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1429566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 1439566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user->mu1, &mu1)); 1449566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(user->mu2, &mu2)); 14580f8ae99SSajid Ali PetscFunctionReturn(0); 14680f8ae99SSajid Ali } 14780f8ae99SSajid Ali 14880f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for tracking method ------------------ */ 14980f8ae99SSajid Ali 1509371c9d4SSatish Balay PetscErrorCode RHSJacobianP_track(TS ts, PetscReal t, Vec U, Mat A, void *ctx) { 15180f8ae99SSajid Ali PetscInt row[] = {0, 1}, col[] = {0, 1}; 15280f8ae99SSajid Ali PetscScalar J[2][2]; 15380f8ae99SSajid Ali const PetscScalar *u; 15480f8ae99SSajid Ali 15580f8ae99SSajid Ali PetscFunctionBeginUser; 1569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 15780f8ae99SSajid Ali J[0][0] = u[1]; 15880f8ae99SSajid Ali J[1][0] = 0; 15980f8ae99SSajid Ali J[0][1] = 0; 16080f8ae99SSajid Ali J[1][1] = (1. - u[0] * u[0]) * u[1] - u[0]; 1619566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 2, col, &J[0][0], INSERT_VALUES)); 1629566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1639566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1649566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 16580f8ae99SSajid Ali PetscFunctionReturn(0); 16680f8ae99SSajid Ali } 16780f8ae99SSajid Ali 16880f8ae99SSajid Ali /* ------------------ Jacobian wrt parameters for global method ------------------ */ 16980f8ae99SSajid Ali 1709371c9d4SSatish Balay PetscErrorCode RHSJacobianP_global(TS ts, PetscReal t, Vec U, Mat A, void *ctx) { 17180f8ae99SSajid Ali PetscInt row[] = {0, 1}, col[] = {0, 1}; 17280f8ae99SSajid Ali PetscScalar J[2][2]; 17380f8ae99SSajid Ali const PetscScalar *u; 17480f8ae99SSajid Ali PetscInt curr_step; 17580f8ae99SSajid Ali 17680f8ae99SSajid Ali PetscFunctionBeginUser; 1779566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &curr_step)); 1789566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 17980f8ae99SSajid Ali J[0][0] = u[1]; 18080f8ae99SSajid Ali J[1][0] = 0; 18180f8ae99SSajid Ali J[0][1] = 0; 18280f8ae99SSajid Ali J[1][1] = (1. - u[0] * u[0]) * u[1] - u[0]; 18380f8ae99SSajid Ali col[0] = (curr_step)*2; 18480f8ae99SSajid Ali col[1] = (curr_step)*2 + 1; 1859566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 2, col, &J[0][0], INSERT_VALUES)); 1869566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 1879566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 1889566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 18980f8ae99SSajid Ali PetscFunctionReturn(0); 19080f8ae99SSajid Ali } 19180f8ae99SSajid Ali 19280f8ae99SSajid Ali /* Dump solution to console if called */ 1939371c9d4SSatish Balay PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec U, void *ctx) { 19480f8ae99SSajid Ali PetscFunctionBeginUser; 19563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Solution at time %e is \n", (double)t)); 1969566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD)); 19780f8ae99SSajid Ali PetscFunctionReturn(0); 19880f8ae99SSajid Ali } 19980f8ae99SSajid Ali 20080f8ae99SSajid Ali /* Customized adjoint monitor to keep track of local 20180f8ae99SSajid Ali sensitivities by storing them in a global sensitivity array. 20280f8ae99SSajid Ali Note : This routine is only used for the tracking method. */ 2039371c9d4SSatish Balay PetscErrorCode AdjointMonitor(TS ts, PetscInt steps, PetscReal time, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *ctx) { 20480f8ae99SSajid Ali AppCtx *user = (AppCtx *)ctx; 20580f8ae99SSajid Ali PetscInt curr_step; 20680f8ae99SSajid Ali PetscScalar *sensmu1_glob; 20780f8ae99SSajid Ali PetscScalar *sensmu2_glob; 20880f8ae99SSajid Ali const PetscScalar *sensmu_loc; 20980f8ae99SSajid Ali 21080f8ae99SSajid Ali PetscFunctionBeginUser; 2119566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &curr_step)); 21280f8ae99SSajid Ali /* Note that we skip the first call to the monitor in the adjoint 21380f8ae99SSajid Ali solve since the sensitivities are already set (during 21480f8ae99SSajid Ali initialization of adjoint vectors). 21580f8ae99SSajid Ali We also note that each indvidial TSAdjointSolve calls the monitor 21680f8ae99SSajid Ali twice, once at the step it is integrating from and once at the step 21780f8ae99SSajid Ali it integrates to. Only the second call is useful for transferring 21880f8ae99SSajid Ali local sensitivities to the global array. */ 21980f8ae99SSajid Ali if (curr_step == user->adj_idx) { 22080f8ae99SSajid Ali PetscFunctionReturn(0); 22180f8ae99SSajid Ali } else { 2229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(*mu, &sensmu_loc)); 2239566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->sens_mu1, &sensmu1_glob)); 2249566063dSJacob Faibussowitsch PetscCall(VecGetArray(user->sens_mu2, &sensmu2_glob)); 22580f8ae99SSajid Ali sensmu1_glob[curr_step] = sensmu_loc[0]; 22680f8ae99SSajid Ali sensmu2_glob[curr_step] = sensmu_loc[1]; 2279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->sens_mu1, &sensmu1_glob)); 2289566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user->sens_mu2, &sensmu2_glob)); 2299566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(*mu, &sensmu_loc)); 23080f8ae99SSajid Ali PetscFunctionReturn(0); 23180f8ae99SSajid Ali } 23280f8ae99SSajid Ali } 23380f8ae99SSajid Ali 2349371c9d4SSatish Balay int main(int argc, char **argv) { 23580f8ae99SSajid Ali TS ts; 23680f8ae99SSajid Ali AppCtx user; 23780f8ae99SSajid Ali PetscScalar *x_ptr, *y_ptr, *u_ptr; 23880f8ae99SSajid Ali PetscMPIInt size; 23980f8ae99SSajid Ali PetscBool monitor = PETSC_FALSE; 24080f8ae99SSajid Ali SAMethod sa = SA_GLOBAL; 24180f8ae99SSajid Ali 24280f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 24380f8ae99SSajid Ali Initialize program 24480f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 245327415f7SBarry Smith PetscFunctionBeginUser; 2469566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 2479566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 2483c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 24980f8ae99SSajid Ali 25080f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 25180f8ae99SSajid Ali Set runtime options 25280f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2539371c9d4SSatish Balay PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "SA analysis options.", ""); 2549371c9d4SSatish Balay { 2559566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL)); 2569566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-sa_method", "Sensitivity analysis method (track or global)", "", SAMethods, (PetscEnum)sa, (PetscEnum *)&sa, NULL)); 25780f8ae99SSajid Ali } 258d0609cedSBarry Smith PetscOptionsEnd(); 25980f8ae99SSajid Ali 26080f8ae99SSajid Ali user.final_time = 0.1; 26180f8ae99SSajid Ali user.max_steps = 5; 26280f8ae99SSajid Ali user.time_step = user.final_time / user.max_steps; 26380f8ae99SSajid Ali 26480f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 26580f8ae99SSajid Ali Create necessary matrix and vectors for forward solve. 26680f8ae99SSajid Ali Create Jacp matrix for adjoint solve. 26780f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2689566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, user.max_steps, &user.mu1)); 2699566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, user.max_steps, &user.mu2)); 2709566063dSJacob Faibussowitsch PetscCall(VecSet(user.mu1, 1.25)); 2719566063dSJacob Faibussowitsch PetscCall(VecSet(user.mu2, 1.0e2)); 27280f8ae99SSajid Ali 27380f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 27480f8ae99SSajid Ali For tracking method : create the global sensitivity array to 27580f8ae99SSajid Ali accumulate sensitivity with respect to parameters at each step. 27680f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 27780f8ae99SSajid Ali if (sa == SA_TRACK) { 2789566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, user.max_steps, &user.sens_mu1)); 2799566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, user.max_steps, &user.sens_mu2)); 28080f8ae99SSajid Ali } 28180f8ae99SSajid Ali 2829566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.A)); 2839566063dSJacob Faibussowitsch PetscCall(MatSetSizes(user.A, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 2849566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.A)); 2859566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.A)); 2869566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.U, NULL)); 28780f8ae99SSajid Ali 28880f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 28980f8ae99SSajid Ali Note that the dimensions of the Jacp matrix depend upon the 29080f8ae99SSajid Ali sensitivity analysis method being used ! 29180f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2929566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp)); 293*48a46eb9SPierre Jolivet if (sa == SA_TRACK) PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 2)); 294*48a46eb9SPierre Jolivet if (sa == SA_GLOBAL) PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, user.max_steps * 2)); 2959566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(user.Jacp)); 2969566063dSJacob Faibussowitsch PetscCall(MatSetUp(user.Jacp)); 29780f8ae99SSajid Ali 29880f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 29980f8ae99SSajid Ali Create timestepping solver context 30080f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3019566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 3029566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ts, TS_EQ_ODE_EXPLICIT)); 3039566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 30480f8ae99SSajid Ali 3059566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user)); 3069566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, user.A, user.A, RHSJacobian, &user)); 307*48a46eb9SPierre Jolivet if (sa == SA_TRACK) PetscCall(TSSetRHSJacobianP(ts, user.Jacp, RHSJacobianP_track, &user)); 308*48a46eb9SPierre Jolivet if (sa == SA_GLOBAL) PetscCall(TSSetRHSJacobianP(ts, user.Jacp, RHSJacobianP_global, &user)); 30980f8ae99SSajid Ali 3109566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 3119566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, user.final_time)); 3129566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, user.final_time / user.max_steps)); 31380f8ae99SSajid Ali 314*48a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); 315*48a46eb9SPierre Jolivet if (sa == SA_TRACK) PetscCall(TSAdjointMonitorSet(ts, AdjointMonitor, &user, NULL)); 31680f8ae99SSajid Ali 31780f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 31880f8ae99SSajid Ali Set initial conditions 31980f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3209566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.U, &x_ptr)); 32180f8ae99SSajid Ali x_ptr[0] = 2.0; 32280f8ae99SSajid Ali x_ptr[1] = -2.0 / 3.0; 3239566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.U, &x_ptr)); 32480f8ae99SSajid Ali 32580f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 32680f8ae99SSajid Ali Save trajectory of solution so that TSAdjointSolve() may be used 32780f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3289566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 32980f8ae99SSajid Ali 33080f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 33180f8ae99SSajid Ali Set runtime options 33280f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3339566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 33480f8ae99SSajid Ali 33580f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 33680f8ae99SSajid Ali Execute forward model and print solution. 33780f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3389566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, user.U)); 3399566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Solution of forward TS :\n")); 3409566063dSJacob Faibussowitsch PetscCall(VecView(user.U, PETSC_VIEWER_STDOUT_WORLD)); 3416aad120cSJose E. Roman PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n Forward TS solve successful! Adjoint run begins!\n")); 34280f8ae99SSajid Ali 34380f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 34480f8ae99SSajid Ali Adjoint model starts here! Create adjoint vectors. 34580f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3469566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.A, &user.lambda, NULL)); 3479566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(user.Jacp, &user.mup, NULL)); 34880f8ae99SSajid Ali 34980f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 35080f8ae99SSajid Ali Set initial conditions for the adjoint vector 35180f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3529566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.U, &u_ptr)); 3539566063dSJacob Faibussowitsch PetscCall(VecGetArray(user.lambda, &y_ptr)); 35480f8ae99SSajid Ali y_ptr[0] = 2 * (u_ptr[0] - 1.5967); 35580f8ae99SSajid Ali y_ptr[1] = 2 * (u_ptr[1] - -(1.02969)); 3569566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.lambda, &y_ptr)); 3579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(user.U, &y_ptr)); 3589566063dSJacob Faibussowitsch PetscCall(VecSet(user.mup, 0)); 35980f8ae99SSajid Ali 36080f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 36180f8ae99SSajid Ali Set number of cost functions. 36280f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 3639566063dSJacob Faibussowitsch PetscCall(TSSetCostGradients(ts, 1, &user.lambda, &user.mup)); 36480f8ae99SSajid Ali 36580f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 36680f8ae99SSajid Ali The adjoint vector mup has to be reset for each adjoint step when 36780f8ae99SSajid Ali using the tracking method as we want to treat the parameters at each 36880f8ae99SSajid Ali time step one at a time and prevent accumulation of the sensitivities 36980f8ae99SSajid Ali from parameters at previous time steps. 37080f8ae99SSajid Ali This is not necessary for the global method as each time dependent 37180f8ae99SSajid Ali parameter is treated as an independent parameter. 37280f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 37380f8ae99SSajid Ali if (sa == SA_TRACK) { 37480f8ae99SSajid Ali for (user.adj_idx = user.max_steps; user.adj_idx > 0; user.adj_idx--) { 3759566063dSJacob Faibussowitsch PetscCall(VecSet(user.mup, 0)); 3769566063dSJacob Faibussowitsch PetscCall(TSAdjointSetSteps(ts, 1)); 3779566063dSJacob Faibussowitsch PetscCall(TSAdjointSolve(ts)); 37880f8ae99SSajid Ali } 37980f8ae99SSajid Ali } 380*48a46eb9SPierre Jolivet if (sa == SA_GLOBAL) PetscCall(TSAdjointSolve(ts)); 38180f8ae99SSajid Ali 38280f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 38380f8ae99SSajid Ali Dispaly adjoint sensitivities wrt parameters and initial conditions 38480f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 38580f8ae99SSajid Ali if (sa == SA_TRACK) { 3869566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt mu1: d[cost]/d[mu1]\n")); 3879566063dSJacob Faibussowitsch PetscCall(VecView(user.sens_mu1, PETSC_VIEWER_STDOUT_WORLD)); 3889566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt mu2: d[cost]/d[mu2]\n")); 3899566063dSJacob Faibussowitsch PetscCall(VecView(user.sens_mu2, PETSC_VIEWER_STDOUT_WORLD)); 39080f8ae99SSajid Ali } 39180f8ae99SSajid Ali 39280f8ae99SSajid Ali if (sa == SA_GLOBAL) { 393d0609cedSBarry Smith PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt params: d[cost]/d[p], where p refers to \nthe interlaced vector made by combining mu1,mu2\n")); 3949566063dSJacob Faibussowitsch PetscCall(VecView(user.mup, PETSC_VIEWER_STDOUT_WORLD)); 39580f8ae99SSajid Ali } 39680f8ae99SSajid Ali 3979566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: d[cost]/d[u(t=0)]\n")); 3989566063dSJacob Faibussowitsch PetscCall(VecView(user.lambda, PETSC_VIEWER_STDOUT_WORLD)); 39980f8ae99SSajid Ali 40080f8ae99SSajid Ali /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 40180f8ae99SSajid Ali Free work space! 40280f8ae99SSajid Ali All PETSc objects should be destroyed when they are no longer needed. 40380f8ae99SSajid Ali - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 4049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.A)); 4059566063dSJacob Faibussowitsch PetscCall(MatDestroy(&user.Jacp)); 4069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.U)); 4079566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.lambda)); 4089566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.mup)); 4099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.mu1)); 4109566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.mu2)); 41180f8ae99SSajid Ali if (sa == SA_TRACK) { 4129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.sens_mu1)); 4139566063dSJacob Faibussowitsch PetscCall(VecDestroy(&user.sens_mu2)); 41480f8ae99SSajid Ali } 4159566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 4169566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 417d0609cedSBarry Smith return (0); 41880f8ae99SSajid Ali } 41980f8ae99SSajid Ali 42080f8ae99SSajid Ali /*TEST 42180f8ae99SSajid Ali 42280f8ae99SSajid Ali test: 42380f8ae99SSajid Ali requires: !complex 42480f8ae99SSajid Ali suffix : track 42580f8ae99SSajid Ali args : -sa_method track 42680f8ae99SSajid Ali 42780f8ae99SSajid Ali test: 42880f8ae99SSajid Ali requires: !complex 42980f8ae99SSajid Ali suffix : global 43080f8ae99SSajid Ali args : -sa_method global 43180f8ae99SSajid Ali 43280f8ae99SSajid Ali TEST*/ 433