1814e85d6SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2814e85d6SHong Zhang #include <petscdraw.h> 3814e85d6SHong Zhang 433f72c5dSHong Zhang PetscLogEvent TS_AdjointStep, TS_ForwardStep, TS_JacobianPEval; 5814e85d6SHong Zhang 67f59fb53SHong Zhang /* #define TSADJOINT_STAGE */ 77f59fb53SHong Zhang 8814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/ 9814e85d6SHong Zhang 10814e85d6SHong Zhang /*@C 11b5b4017aSHong Zhang TSSetRHSJacobianP - Sets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix. 12814e85d6SHong Zhang 13814e85d6SHong Zhang Logically Collective on TS 14814e85d6SHong Zhang 15814e85d6SHong Zhang Input Parameters: 16b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 17b5b4017aSHong Zhang . Amat - JacobianP matrix 18b5b4017aSHong Zhang . func - function 19b5b4017aSHong Zhang - ctx - [optional] user-defined function context 20814e85d6SHong Zhang 21814e85d6SHong Zhang Calling sequence of func: 22814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 23814e85d6SHong Zhang + t - current timestep 24c9ad9de0SHong Zhang . U - input vector (current ODE solution) 25814e85d6SHong Zhang . A - output matrix 26814e85d6SHong Zhang - ctx - [optional] user-defined function context 27814e85d6SHong Zhang 28814e85d6SHong Zhang Level: intermediate 29814e85d6SHong Zhang 30814e85d6SHong Zhang Notes: 31814e85d6SHong Zhang Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 32814e85d6SHong Zhang 33db781477SPatrick Sanan .seealso: `TSGetRHSJacobianP()` 34814e85d6SHong Zhang @*/ 359371c9d4SSatish Balay PetscErrorCode TSSetRHSJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) { 36814e85d6SHong Zhang PetscFunctionBegin; 37814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 38814e85d6SHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 39814e85d6SHong Zhang 40814e85d6SHong Zhang ts->rhsjacobianp = func; 41814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 42814e85d6SHong Zhang if (Amat) { 439566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacprhs)); 4533f72c5dSHong Zhang ts->Jacprhs = Amat; 46814e85d6SHong Zhang } 47814e85d6SHong Zhang PetscFunctionReturn(0); 48814e85d6SHong Zhang } 49814e85d6SHong Zhang 50814e85d6SHong Zhang /*@C 51cd4cee2dSHong Zhang TSGetRHSJacobianP - Gets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix. 52cd4cee2dSHong Zhang 53cd4cee2dSHong Zhang Logically Collective on TS 54cd4cee2dSHong Zhang 55f899ff85SJose E. Roman Input Parameter: 56cd4cee2dSHong Zhang . ts - TS context obtained from TSCreate() 57cd4cee2dSHong Zhang 58cd4cee2dSHong Zhang Output Parameters: 59cd4cee2dSHong Zhang + Amat - JacobianP matrix 60cd4cee2dSHong Zhang . func - function 61cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 62cd4cee2dSHong Zhang 63cd4cee2dSHong Zhang Calling sequence of func: 64cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 65cd4cee2dSHong Zhang + t - current timestep 66cd4cee2dSHong Zhang . U - input vector (current ODE solution) 67cd4cee2dSHong Zhang . A - output matrix 68cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 69cd4cee2dSHong Zhang 70cd4cee2dSHong Zhang Level: intermediate 71cd4cee2dSHong Zhang 72cd4cee2dSHong Zhang Notes: 73cd4cee2dSHong Zhang Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 74cd4cee2dSHong Zhang 75db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()` 76cd4cee2dSHong Zhang @*/ 779371c9d4SSatish Balay PetscErrorCode TSGetRHSJacobianP(TS ts, Mat *Amat, PetscErrorCode (**func)(TS, PetscReal, Vec, Mat, void *), void **ctx) { 78cd4cee2dSHong Zhang PetscFunctionBegin; 79cd4cee2dSHong Zhang if (func) *func = ts->rhsjacobianp; 80cd4cee2dSHong Zhang if (ctx) *ctx = ts->rhsjacobianpctx; 81cd4cee2dSHong Zhang if (Amat) *Amat = ts->Jacprhs; 82cd4cee2dSHong Zhang PetscFunctionReturn(0); 83cd4cee2dSHong Zhang } 84cd4cee2dSHong Zhang 85cd4cee2dSHong Zhang /*@C 86814e85d6SHong Zhang TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 87814e85d6SHong Zhang 88814e85d6SHong Zhang Collective on TS 89814e85d6SHong Zhang 90814e85d6SHong Zhang Input Parameters: 91814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 92814e85d6SHong Zhang 93814e85d6SHong Zhang Level: developer 94814e85d6SHong Zhang 95db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()` 96814e85d6SHong Zhang @*/ 979371c9d4SSatish Balay PetscErrorCode TSComputeRHSJacobianP(TS ts, PetscReal t, Vec U, Mat Amat) { 98814e85d6SHong Zhang PetscFunctionBegin; 9933f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 100814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 101c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 102814e85d6SHong Zhang 103792fecdfSBarry Smith PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 104814e85d6SHong Zhang PetscFunctionReturn(0); 105814e85d6SHong Zhang } 106814e85d6SHong Zhang 107814e85d6SHong Zhang /*@C 10833f72c5dSHong Zhang TSSetIJacobianP - Sets the function that computes the Jacobian of F w.r.t. the parameters P where F(Udot,U,t) = G(U,P,t), as well as the location to store the matrix. 10933f72c5dSHong Zhang 11033f72c5dSHong Zhang Logically Collective on TS 11133f72c5dSHong Zhang 11233f72c5dSHong Zhang Input Parameters: 11333f72c5dSHong Zhang + ts - TS context obtained from TSCreate() 11433f72c5dSHong Zhang . Amat - JacobianP matrix 11533f72c5dSHong Zhang . func - function 11633f72c5dSHong Zhang - ctx - [optional] user-defined function context 11733f72c5dSHong Zhang 11833f72c5dSHong Zhang Calling sequence of func: 11933f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 12033f72c5dSHong Zhang + t - current timestep 12133f72c5dSHong Zhang . U - input vector (current ODE solution) 12233f72c5dSHong Zhang . Udot - time derivative of state vector 12333f72c5dSHong Zhang . shift - shift to apply, see note below 12433f72c5dSHong Zhang . A - output matrix 12533f72c5dSHong Zhang - ctx - [optional] user-defined function context 12633f72c5dSHong Zhang 12733f72c5dSHong Zhang Level: intermediate 12833f72c5dSHong Zhang 12933f72c5dSHong Zhang Notes: 13033f72c5dSHong Zhang Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 13133f72c5dSHong Zhang 13233f72c5dSHong Zhang .seealso: 13333f72c5dSHong Zhang @*/ 1349371c9d4SSatish Balay PetscErrorCode TSSetIJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *ctx) { 13533f72c5dSHong Zhang PetscFunctionBegin; 13633f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 13733f72c5dSHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 13833f72c5dSHong Zhang 13933f72c5dSHong Zhang ts->ijacobianp = func; 14033f72c5dSHong Zhang ts->ijacobianpctx = ctx; 14133f72c5dSHong Zhang if (Amat) { 1429566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 1439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 14433f72c5dSHong Zhang ts->Jacp = Amat; 14533f72c5dSHong Zhang } 14633f72c5dSHong Zhang PetscFunctionReturn(0); 14733f72c5dSHong Zhang } 14833f72c5dSHong Zhang 14933f72c5dSHong Zhang /*@C 15033f72c5dSHong Zhang TSComputeIJacobianP - Runs the user-defined IJacobianP function. 15133f72c5dSHong Zhang 15233f72c5dSHong Zhang Collective on TS 15333f72c5dSHong Zhang 15433f72c5dSHong Zhang Input Parameters: 15533f72c5dSHong Zhang + ts - the TS context 15633f72c5dSHong Zhang . t - current timestep 15733f72c5dSHong Zhang . U - state vector 15833f72c5dSHong Zhang . Udot - time derivative of state vector 15933f72c5dSHong Zhang . shift - shift to apply, see note below 16033f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 16133f72c5dSHong Zhang 16233f72c5dSHong Zhang Output Parameters: 16333f72c5dSHong Zhang . A - Jacobian matrix 16433f72c5dSHong Zhang 16533f72c5dSHong Zhang Level: developer 16633f72c5dSHong Zhang 167db781477SPatrick Sanan .seealso: `TSSetIJacobianP()` 16833f72c5dSHong Zhang @*/ 1699371c9d4SSatish Balay PetscErrorCode TSComputeIJacobianP(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat Amat, PetscBool imex) { 17033f72c5dSHong Zhang PetscFunctionBegin; 17133f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 17233f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 17333f72c5dSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 17433f72c5dSHong Zhang PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4); 17533f72c5dSHong Zhang 1769566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_JacobianPEval, ts, U, Amat, 0)); 177*48a46eb9SPierre Jolivet if (ts->ijacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->ijacobianp)(ts, t, U, Udot, shift, Amat, ts->ijacobianpctx)); 17833f72c5dSHong Zhang if (imex) { 17933f72c5dSHong Zhang if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */ 18033f72c5dSHong Zhang PetscBool assembled; 1819566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 1829566063dSJacob Faibussowitsch PetscCall(MatAssembled(Amat, &assembled)); 18333f72c5dSHong Zhang if (!assembled) { 1849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY)); 1859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY)); 18633f72c5dSHong Zhang } 18733f72c5dSHong Zhang } 18833f72c5dSHong Zhang } else { 1891baa6e33SBarry Smith if (ts->rhsjacobianp) PetscCall(TSComputeRHSJacobianP(ts, t, U, ts->Jacprhs)); 19033f72c5dSHong Zhang if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */ 1919566063dSJacob Faibussowitsch PetscCall(MatScale(Amat, -1)); 19233f72c5dSHong Zhang } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */ 19333f72c5dSHong Zhang MatStructure axpy = DIFFERENT_NONZERO_PATTERN; 19433f72c5dSHong Zhang if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */ 1959566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 19633f72c5dSHong Zhang } 1979566063dSJacob Faibussowitsch PetscCall(MatAXPY(Amat, -1, ts->Jacprhs, axpy)); 19833f72c5dSHong Zhang } 19933f72c5dSHong Zhang } 2009566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_JacobianPEval, ts, U, Amat, 0)); 20133f72c5dSHong Zhang PetscFunctionReturn(0); 20233f72c5dSHong Zhang } 20333f72c5dSHong Zhang 20433f72c5dSHong Zhang /*@C 205814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 206814e85d6SHong Zhang 207814e85d6SHong Zhang Logically Collective on TS 208814e85d6SHong Zhang 209814e85d6SHong Zhang Input Parameters: 210814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 211814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 212814e85d6SHong Zhang . costintegral - vector that stores the integral values 213814e85d6SHong Zhang . rf - routine for evaluating the integrand function 214c9ad9de0SHong Zhang . drduf - function that computes the gradients of the r's with respect to u 215814e85d6SHong Zhang . drdpf - function that computes the gradients of the r's with respect to p, can be NULL if parametric sensitivity is not desired (mu=NULL) 216814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 217814e85d6SHong Zhang - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 218814e85d6SHong Zhang 219814e85d6SHong Zhang Calling sequence of rf: 220c9ad9de0SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx); 221814e85d6SHong Zhang 222c9ad9de0SHong Zhang Calling sequence of drduf: 223c9ad9de0SHong Zhang $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx); 224814e85d6SHong Zhang 225814e85d6SHong Zhang Calling sequence of drdpf: 226c9ad9de0SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx); 227814e85d6SHong Zhang 228cd4cee2dSHong Zhang Level: deprecated 229814e85d6SHong Zhang 230814e85d6SHong Zhang Notes: 231814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 232814e85d6SHong Zhang 233db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()`, `TSGetCostGradients()`, `TSSetCostGradients()` 234814e85d6SHong Zhang @*/ 2359371c9d4SSatish Balay PetscErrorCode TSSetCostIntegrand(TS ts, PetscInt numcost, Vec costintegral, PetscErrorCode (*rf)(TS, PetscReal, Vec, Vec, void *), PetscErrorCode (*drduf)(TS, PetscReal, Vec, Vec *, void *), PetscErrorCode (*drdpf)(TS, PetscReal, Vec, Vec *, void *), PetscBool fwd, void *ctx) { 236814e85d6SHong Zhang PetscFunctionBegin; 237814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 238814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral, VEC_CLASSID, 3); 2393c633725SBarry Smith PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostGradients() or TSForwardSetIntegralGradients()"); 240814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numcost; 241814e85d6SHong Zhang 242814e85d6SHong Zhang if (costintegral) { 2439566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)costintegral)); 2449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_costintegral)); 245814e85d6SHong Zhang ts->vec_costintegral = costintegral; 246814e85d6SHong Zhang } else { 247814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 2489566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, numcost, &ts->vec_costintegral)); 249814e85d6SHong Zhang } else { 2509566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegral, 0.0)); 251814e85d6SHong Zhang } 252814e85d6SHong Zhang } 253814e85d6SHong Zhang if (!ts->vec_costintegrand) { 2549566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_costintegral, &ts->vec_costintegrand)); 255814e85d6SHong Zhang } else { 2569566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegrand, 0.0)); 257814e85d6SHong Zhang } 258814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 259814e85d6SHong Zhang ts->costintegrand = rf; 260814e85d6SHong Zhang ts->costintegrandctx = ctx; 261c9ad9de0SHong Zhang ts->drdufunction = drduf; 262814e85d6SHong Zhang ts->drdpfunction = drdpf; 263814e85d6SHong Zhang PetscFunctionReturn(0); 264814e85d6SHong Zhang } 265814e85d6SHong Zhang 266b5b4017aSHong Zhang /*@C 267814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 268814e85d6SHong Zhang It is valid to call the routine after a backward run. 269814e85d6SHong Zhang 270814e85d6SHong Zhang Not Collective 271814e85d6SHong Zhang 272814e85d6SHong Zhang Input Parameter: 273814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 274814e85d6SHong Zhang 275814e85d6SHong Zhang Output Parameter: 276814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 277814e85d6SHong Zhang 278814e85d6SHong Zhang Level: intermediate 279814e85d6SHong Zhang 280db781477SPatrick Sanan .seealso: `TSSetCostIntegrand()` 281814e85d6SHong Zhang 282814e85d6SHong Zhang @*/ 2839371c9d4SSatish Balay PetscErrorCode TSGetCostIntegral(TS ts, Vec *v) { 284cd4cee2dSHong Zhang TS quadts; 285cd4cee2dSHong Zhang 286814e85d6SHong Zhang PetscFunctionBegin; 287814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 288814e85d6SHong Zhang PetscValidPointer(v, 2); 2899566063dSJacob Faibussowitsch PetscCall(TSGetQuadratureTS(ts, NULL, &quadts)); 290cd4cee2dSHong Zhang *v = quadts->vec_sol; 291814e85d6SHong Zhang PetscFunctionReturn(0); 292814e85d6SHong Zhang } 293814e85d6SHong Zhang 294b5b4017aSHong Zhang /*@C 295814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 296814e85d6SHong Zhang 297814e85d6SHong Zhang Input Parameters: 298814e85d6SHong Zhang + ts - the TS context 299814e85d6SHong Zhang . t - current time 300c9ad9de0SHong Zhang - U - state vector, i.e. current solution 301814e85d6SHong Zhang 302814e85d6SHong Zhang Output Parameter: 303c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs 304814e85d6SHong Zhang 305369cf35fSHong Zhang Notes: 306814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 307814e85d6SHong Zhang is used internally within the sensitivity analysis context. 308814e85d6SHong Zhang 309cd4cee2dSHong Zhang Level: deprecated 310814e85d6SHong Zhang 311db781477SPatrick Sanan .seealso: `TSSetCostIntegrand()` 312814e85d6SHong Zhang @*/ 3139371c9d4SSatish Balay PetscErrorCode TSComputeCostIntegrand(TS ts, PetscReal t, Vec U, Vec Q) { 314814e85d6SHong Zhang PetscFunctionBegin; 315814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 316c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 317c9ad9de0SHong Zhang PetscValidHeaderSpecific(Q, VEC_CLASSID, 4); 318814e85d6SHong Zhang 3199566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_FunctionEval, ts, U, Q, 0)); 320792fecdfSBarry Smith if (ts->costintegrand) PetscCallBack("TS callback integrand in the cost function", (*ts->costintegrand)(ts, t, U, Q, ts->costintegrandctx)); 321ef1023bdSBarry Smith else PetscCall(VecZeroEntries(Q)); 3229566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_FunctionEval, ts, U, Q, 0)); 323814e85d6SHong Zhang PetscFunctionReturn(0); 324814e85d6SHong Zhang } 325814e85d6SHong Zhang 326b5b4017aSHong Zhang /*@C 327d76be551SHong Zhang TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 328814e85d6SHong Zhang 329d76be551SHong Zhang Level: deprecated 330814e85d6SHong Zhang 331814e85d6SHong Zhang @*/ 3329371c9d4SSatish Balay PetscErrorCode TSComputeDRDUFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) { 333814e85d6SHong Zhang PetscFunctionBegin; 33433f72c5dSHong Zhang if (!DRDU) PetscFunctionReturn(0); 335814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 336c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 337814e85d6SHong Zhang 338792fecdfSBarry Smith PetscCallBack("TS callback DRDU for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 339814e85d6SHong Zhang PetscFunctionReturn(0); 340814e85d6SHong Zhang } 341814e85d6SHong Zhang 342b5b4017aSHong Zhang /*@C 343d76be551SHong Zhang TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 344814e85d6SHong Zhang 345d76be551SHong Zhang Level: deprecated 346814e85d6SHong Zhang 347814e85d6SHong Zhang @*/ 3489371c9d4SSatish Balay PetscErrorCode TSComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) { 349814e85d6SHong Zhang PetscFunctionBegin; 35033f72c5dSHong Zhang if (!DRDP) PetscFunctionReturn(0); 351814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 352c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 353814e85d6SHong Zhang 354792fecdfSBarry Smith PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 355814e85d6SHong Zhang PetscFunctionReturn(0); 356814e85d6SHong Zhang } 357814e85d6SHong Zhang 358b5b4017aSHong Zhang /*@C 359b5b4017aSHong Zhang TSSetIHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of F (IFunction) w.r.t. the state variable. 360b5b4017aSHong Zhang 361b5b4017aSHong Zhang Logically Collective on TS 362b5b4017aSHong Zhang 363b5b4017aSHong Zhang Input Parameters: 364b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 365b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 366b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 367b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 368b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 369b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 370b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 371b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 372f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for F_PP 373b5b4017aSHong Zhang 374b5b4017aSHong Zhang Calling sequence of ihessianproductfunc: 375369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 376b5b4017aSHong Zhang + t - current timestep 377b5b4017aSHong Zhang . U - input vector (current ODE solution) 378369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 379b5b4017aSHong Zhang . Vr - input vector to be right-multiplied with the Hessian 380369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 381b5b4017aSHong Zhang - ctx - [optional] user-defined function context 382b5b4017aSHong Zhang 383b5b4017aSHong Zhang Level: intermediate 384b5b4017aSHong Zhang 385369cf35fSHong Zhang Notes: 386369cf35fSHong Zhang The first Hessian function and the working array are required. 387369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 388369cf35fSHong Zhang $ Vl_n^T*F_UP*Vr 389369cf35fSHong Zhang where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian F_UP is of size N x N x M. 390369cf35fSHong Zhang Each entry of F_UP corresponds to the derivative 391369cf35fSHong Zhang $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}. 392369cf35fSHong Zhang The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with the j-th entry being 393369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]} 394369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 395b5b4017aSHong Zhang 396b5b4017aSHong Zhang .seealso: 397b5b4017aSHong Zhang @*/ 3989371c9d4SSatish Balay PetscErrorCode TSSetIHessianProduct(TS ts, Vec *ihp1, PetscErrorCode (*ihessianproductfunc1)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp2, PetscErrorCode (*ihessianproductfunc2)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp3, PetscErrorCode (*ihessianproductfunc3)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *ihp4, PetscErrorCode (*ihessianproductfunc4)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *ctx) { 399b5b4017aSHong Zhang PetscFunctionBegin; 400b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 401b5b4017aSHong Zhang PetscValidPointer(ihp1, 2); 402b5b4017aSHong Zhang 403b5b4017aSHong Zhang ts->ihessianproductctx = ctx; 404b5b4017aSHong Zhang if (ihp1) ts->vecs_fuu = ihp1; 405b5b4017aSHong Zhang if (ihp2) ts->vecs_fup = ihp2; 406b5b4017aSHong Zhang if (ihp3) ts->vecs_fpu = ihp3; 407b5b4017aSHong Zhang if (ihp4) ts->vecs_fpp = ihp4; 408b5b4017aSHong Zhang ts->ihessianproduct_fuu = ihessianproductfunc1; 409b5b4017aSHong Zhang ts->ihessianproduct_fup = ihessianproductfunc2; 410b5b4017aSHong Zhang ts->ihessianproduct_fpu = ihessianproductfunc3; 411b5b4017aSHong Zhang ts->ihessianproduct_fpp = ihessianproductfunc4; 412b5b4017aSHong Zhang PetscFunctionReturn(0); 413b5b4017aSHong Zhang } 414b5b4017aSHong Zhang 415b5b4017aSHong Zhang /*@C 416b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu. 417b5b4017aSHong Zhang 418b5b4017aSHong Zhang Collective on TS 419b5b4017aSHong Zhang 420b5b4017aSHong Zhang Input Parameters: 421b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 422b5b4017aSHong Zhang 423b5b4017aSHong Zhang Notes: 424b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation, 425b5b4017aSHong Zhang so most users would not generally call this routine themselves. 426b5b4017aSHong Zhang 427b5b4017aSHong Zhang Level: developer 428b5b4017aSHong Zhang 429db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 430b5b4017aSHong Zhang @*/ 4319371c9d4SSatish Balay PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 432b5b4017aSHong Zhang PetscFunctionBegin; 43333f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 434b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 435b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 436b5b4017aSHong Zhang 437792fecdfSBarry Smith if (ts->ihessianproduct_fuu) PetscCallBack("TS callback IHessianProduct 1 for sensitivity analysis", (*ts->ihessianproduct_fuu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 438ef1023bdSBarry Smith 43967633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 44067633408SHong Zhang if (ts->rhshessianproduct_guu) { 44167633408SHong Zhang PetscInt nadj; 4429566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(ts, t, U, Vl, Vr, VHV)); 443*48a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 44467633408SHong Zhang } 445b5b4017aSHong Zhang PetscFunctionReturn(0); 446b5b4017aSHong Zhang } 447b5b4017aSHong Zhang 448b5b4017aSHong Zhang /*@C 449b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup. 450b5b4017aSHong Zhang 451b5b4017aSHong Zhang Collective on TS 452b5b4017aSHong Zhang 453b5b4017aSHong Zhang Input Parameters: 454b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 455b5b4017aSHong Zhang 456b5b4017aSHong Zhang Notes: 457b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation, 458b5b4017aSHong Zhang so most users would not generally call this routine themselves. 459b5b4017aSHong Zhang 460b5b4017aSHong Zhang Level: developer 461b5b4017aSHong Zhang 462db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 463b5b4017aSHong Zhang @*/ 4649371c9d4SSatish Balay PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 465b5b4017aSHong Zhang PetscFunctionBegin; 46633f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 467b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 468b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 469b5b4017aSHong Zhang 470792fecdfSBarry Smith if (ts->ihessianproduct_fup) PetscCallBack("TS callback IHessianProduct 2 for sensitivity analysis", (*ts->ihessianproduct_fup)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 471ef1023bdSBarry Smith 47267633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 47367633408SHong Zhang if (ts->rhshessianproduct_gup) { 47467633408SHong Zhang PetscInt nadj; 4759566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(ts, t, U, Vl, Vr, VHV)); 476*48a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 47767633408SHong Zhang } 478b5b4017aSHong Zhang PetscFunctionReturn(0); 479b5b4017aSHong Zhang } 480b5b4017aSHong Zhang 481b5b4017aSHong Zhang /*@C 482b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu. 483b5b4017aSHong Zhang 484b5b4017aSHong Zhang Collective on TS 485b5b4017aSHong Zhang 486b5b4017aSHong Zhang Input Parameters: 487b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 488b5b4017aSHong Zhang 489b5b4017aSHong Zhang Notes: 490b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation, 491b5b4017aSHong Zhang so most users would not generally call this routine themselves. 492b5b4017aSHong Zhang 493b5b4017aSHong Zhang Level: developer 494b5b4017aSHong Zhang 495db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 496b5b4017aSHong Zhang @*/ 4979371c9d4SSatish Balay PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 498b5b4017aSHong Zhang PetscFunctionBegin; 49933f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 500b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 501b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 502b5b4017aSHong Zhang 503792fecdfSBarry Smith if (ts->ihessianproduct_fpu) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 504ef1023bdSBarry Smith 50567633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 50667633408SHong Zhang if (ts->rhshessianproduct_gpu) { 50767633408SHong Zhang PetscInt nadj; 5089566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(ts, t, U, Vl, Vr, VHV)); 509*48a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 51067633408SHong Zhang } 511b5b4017aSHong Zhang PetscFunctionReturn(0); 512b5b4017aSHong Zhang } 513b5b4017aSHong Zhang 514b5b4017aSHong Zhang /*@C 515b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp. 516b5b4017aSHong Zhang 517b5b4017aSHong Zhang Collective on TS 518b5b4017aSHong Zhang 519b5b4017aSHong Zhang Input Parameters: 520b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 521b5b4017aSHong Zhang 522b5b4017aSHong Zhang Notes: 523b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation, 524b5b4017aSHong Zhang so most users would not generally call this routine themselves. 525b5b4017aSHong Zhang 526b5b4017aSHong Zhang Level: developer 527b5b4017aSHong Zhang 528db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 529b5b4017aSHong Zhang @*/ 5309371c9d4SSatish Balay PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 531b5b4017aSHong Zhang PetscFunctionBegin; 53233f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 533b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 534b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 535b5b4017aSHong Zhang 536792fecdfSBarry Smith if (ts->ihessianproduct_fpp) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpp)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 537ef1023bdSBarry Smith 53867633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 53967633408SHong Zhang if (ts->rhshessianproduct_gpp) { 54067633408SHong Zhang PetscInt nadj; 5419566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(ts, t, U, Vl, Vr, VHV)); 542*48a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 54367633408SHong Zhang } 544b5b4017aSHong Zhang PetscFunctionReturn(0); 545b5b4017aSHong Zhang } 546b5b4017aSHong Zhang 54713af1a74SHong Zhang /*@C 5484c500f23SPierre Jolivet TSSetRHSHessianProduct - Sets the function that computes the vector-Hessian-vector product. The Hessian is the second-order derivative of G (RHSFunction) w.r.t. the state variable. 54913af1a74SHong Zhang 55013af1a74SHong Zhang Logically Collective on TS 55113af1a74SHong Zhang 55213af1a74SHong Zhang Input Parameters: 55313af1a74SHong Zhang + ts - TS context obtained from TSCreate() 55413af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 55513af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 55613af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 55713af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 55813af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 55913af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 56013af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 561f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for G_PP 56213af1a74SHong Zhang 56313af1a74SHong Zhang Calling sequence of ihessianproductfunc: 564369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 56513af1a74SHong Zhang + t - current timestep 56613af1a74SHong Zhang . U - input vector (current ODE solution) 567369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 56813af1a74SHong Zhang . Vr - input vector to be right-multiplied with the Hessian 569369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 57013af1a74SHong Zhang - ctx - [optional] user-defined function context 57113af1a74SHong Zhang 57213af1a74SHong Zhang Level: intermediate 57313af1a74SHong Zhang 574369cf35fSHong Zhang Notes: 575369cf35fSHong Zhang The first Hessian function and the working array are required. 576369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 577369cf35fSHong Zhang $ Vl_n^T*G_UP*Vr 578369cf35fSHong Zhang where the vector Vl_n (n-th element in the array Vl) and Vr are of size N and M respectively, and the Hessian G_UP is of size N x N x M. 579369cf35fSHong Zhang Each entry of G_UP corresponds to the derivative 580369cf35fSHong Zhang $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}. 581369cf35fSHong Zhang The result of the vector-Hessian-vector product for Vl_n needs to be stored in vector VHV_n with j-th entry being 582369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]} 583369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 58413af1a74SHong Zhang 58513af1a74SHong Zhang .seealso: 58613af1a74SHong Zhang @*/ 5879371c9d4SSatish Balay PetscErrorCode TSSetRHSHessianProduct(TS ts, Vec *rhshp1, PetscErrorCode (*rhshessianproductfunc1)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp2, PetscErrorCode (*rhshessianproductfunc2)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp3, PetscErrorCode (*rhshessianproductfunc3)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), Vec *rhshp4, PetscErrorCode (*rhshessianproductfunc4)(TS, PetscReal, Vec, Vec *, Vec, Vec *, void *), void *ctx) { 58813af1a74SHong Zhang PetscFunctionBegin; 58913af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 59013af1a74SHong Zhang PetscValidPointer(rhshp1, 2); 59113af1a74SHong Zhang 59213af1a74SHong Zhang ts->rhshessianproductctx = ctx; 59313af1a74SHong Zhang if (rhshp1) ts->vecs_guu = rhshp1; 59413af1a74SHong Zhang if (rhshp2) ts->vecs_gup = rhshp2; 59513af1a74SHong Zhang if (rhshp3) ts->vecs_gpu = rhshp3; 59613af1a74SHong Zhang if (rhshp4) ts->vecs_gpp = rhshp4; 59713af1a74SHong Zhang ts->rhshessianproduct_guu = rhshessianproductfunc1; 59813af1a74SHong Zhang ts->rhshessianproduct_gup = rhshessianproductfunc2; 59913af1a74SHong Zhang ts->rhshessianproduct_gpu = rhshessianproductfunc3; 60013af1a74SHong Zhang ts->rhshessianproduct_gpp = rhshessianproductfunc4; 60113af1a74SHong Zhang PetscFunctionReturn(0); 60213af1a74SHong Zhang } 60313af1a74SHong Zhang 60413af1a74SHong Zhang /*@C 605b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu. 60613af1a74SHong Zhang 60713af1a74SHong Zhang Collective on TS 60813af1a74SHong Zhang 60913af1a74SHong Zhang Input Parameters: 61013af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 61113af1a74SHong Zhang 61213af1a74SHong Zhang Notes: 613b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation, 61413af1a74SHong Zhang so most users would not generally call this routine themselves. 61513af1a74SHong Zhang 61613af1a74SHong Zhang Level: developer 61713af1a74SHong Zhang 618db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 61913af1a74SHong Zhang @*/ 6209371c9d4SSatish Balay PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 62113af1a74SHong Zhang PetscFunctionBegin; 62213af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 62313af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 62413af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 62513af1a74SHong Zhang 626792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 1 for sensitivity analysis", (*ts->rhshessianproduct_guu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 62713af1a74SHong Zhang PetscFunctionReturn(0); 62813af1a74SHong Zhang } 62913af1a74SHong Zhang 63013af1a74SHong Zhang /*@C 631b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup. 63213af1a74SHong Zhang 63313af1a74SHong Zhang Collective on TS 63413af1a74SHong Zhang 63513af1a74SHong Zhang Input Parameters: 63613af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 63713af1a74SHong Zhang 63813af1a74SHong Zhang Notes: 639b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation, 64013af1a74SHong Zhang so most users would not generally call this routine themselves. 64113af1a74SHong Zhang 64213af1a74SHong Zhang Level: developer 64313af1a74SHong Zhang 644db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 64513af1a74SHong Zhang @*/ 6469371c9d4SSatish Balay PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 64713af1a74SHong Zhang PetscFunctionBegin; 64813af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 64913af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 65013af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 65113af1a74SHong Zhang 652792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 2 for sensitivity analysis", (*ts->rhshessianproduct_gup)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 65313af1a74SHong Zhang PetscFunctionReturn(0); 65413af1a74SHong Zhang } 65513af1a74SHong Zhang 65613af1a74SHong Zhang /*@C 657b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu. 65813af1a74SHong Zhang 65913af1a74SHong Zhang Collective on TS 66013af1a74SHong Zhang 66113af1a74SHong Zhang Input Parameters: 66213af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 66313af1a74SHong Zhang 66413af1a74SHong Zhang Notes: 665b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation, 66613af1a74SHong Zhang so most users would not generally call this routine themselves. 66713af1a74SHong Zhang 66813af1a74SHong Zhang Level: developer 66913af1a74SHong Zhang 670db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 67113af1a74SHong Zhang @*/ 6729371c9d4SSatish Balay PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 67313af1a74SHong Zhang PetscFunctionBegin; 67413af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 67513af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 67613af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 67713af1a74SHong Zhang 678792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 67913af1a74SHong Zhang PetscFunctionReturn(0); 68013af1a74SHong Zhang } 68113af1a74SHong Zhang 68213af1a74SHong Zhang /*@C 683b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp. 68413af1a74SHong Zhang 68513af1a74SHong Zhang Collective on TS 68613af1a74SHong Zhang 68713af1a74SHong Zhang Input Parameters: 68813af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 68913af1a74SHong Zhang 69013af1a74SHong Zhang Notes: 691b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation, 69213af1a74SHong Zhang so most users would not generally call this routine themselves. 69313af1a74SHong Zhang 69413af1a74SHong Zhang Level: developer 69513af1a74SHong Zhang 696db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 69713af1a74SHong Zhang @*/ 6989371c9d4SSatish Balay PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) { 69913af1a74SHong Zhang PetscFunctionBegin; 70013af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 70113af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 70213af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 70313af1a74SHong Zhang 704792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpp)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 70513af1a74SHong Zhang PetscFunctionReturn(0); 70613af1a74SHong Zhang } 70713af1a74SHong Zhang 708814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 709814e85d6SHong Zhang 710814e85d6SHong Zhang /*@ 711814e85d6SHong Zhang TSSetCostGradients - Sets the initial value of the gradients of the cost function w.r.t. initial values and w.r.t. the problem parameters 712814e85d6SHong Zhang for use by the TSAdjoint routines. 713814e85d6SHong Zhang 714d083f849SBarry Smith Logically Collective on TS 715814e85d6SHong Zhang 716814e85d6SHong Zhang Input Parameters: 717814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 718d2fd7bfcSHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 719814e85d6SHong Zhang . lambda - gradients with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector 720814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 721814e85d6SHong Zhang 722814e85d6SHong Zhang Level: beginner 723814e85d6SHong Zhang 724814e85d6SHong Zhang Notes: 725814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 726814e85d6SHong Zhang 727814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 728814e85d6SHong Zhang 729db781477SPatrick Sanan .seealso `TSGetCostGradients()` 730814e85d6SHong Zhang @*/ 7319371c9d4SSatish Balay PetscErrorCode TSSetCostGradients(TS ts, PetscInt numcost, Vec *lambda, Vec *mu) { 732814e85d6SHong Zhang PetscFunctionBegin; 733814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 734064a246eSJacob Faibussowitsch PetscValidPointer(lambda, 3); 735814e85d6SHong Zhang ts->vecs_sensi = lambda; 736814e85d6SHong Zhang ts->vecs_sensip = mu; 7373c633725SBarry Smith PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand"); 738814e85d6SHong Zhang ts->numcost = numcost; 739814e85d6SHong Zhang PetscFunctionReturn(0); 740814e85d6SHong Zhang } 741814e85d6SHong Zhang 742814e85d6SHong Zhang /*@ 743814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 744814e85d6SHong Zhang 745814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 746814e85d6SHong Zhang 747814e85d6SHong Zhang Input Parameter: 748814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 749814e85d6SHong Zhang 750d8d19677SJose E. Roman Output Parameters: 7516b867d5aSJose E. Roman + numcost - size of returned arrays 7526b867d5aSJose E. Roman . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 753814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 754814e85d6SHong Zhang 755814e85d6SHong Zhang Level: intermediate 756814e85d6SHong Zhang 757db781477SPatrick Sanan .seealso: `TSSetCostGradients()` 758814e85d6SHong Zhang @*/ 7599371c9d4SSatish Balay PetscErrorCode TSGetCostGradients(TS ts, PetscInt *numcost, Vec **lambda, Vec **mu) { 760814e85d6SHong Zhang PetscFunctionBegin; 761814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 762814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 763814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 764814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 765814e85d6SHong Zhang PetscFunctionReturn(0); 766814e85d6SHong Zhang } 767814e85d6SHong Zhang 768814e85d6SHong Zhang /*@ 769b5b4017aSHong Zhang TSSetCostHessianProducts - Sets the initial value of the Hessian-vector products of the cost function w.r.t. initial values and w.r.t. the problem parameters 770b5b4017aSHong Zhang for use by the TSAdjoint routines. 771b5b4017aSHong Zhang 772d083f849SBarry Smith Logically Collective on TS 773b5b4017aSHong Zhang 774b5b4017aSHong Zhang Input Parameters: 775b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 776b5b4017aSHong Zhang . numcost - number of cost functions 777b5b4017aSHong Zhang . lambda2 - Hessian-vector product with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector 778b5b4017aSHong Zhang . mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 779b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 780b5b4017aSHong Zhang 781b5b4017aSHong Zhang Level: beginner 782b5b4017aSHong Zhang 783b5b4017aSHong Zhang Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system 784b5b4017aSHong Zhang 785ecf68647SHong Zhang For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve(). 786b5b4017aSHong Zhang 787b5b4017aSHong Zhang After TSAdjointSolve() is called, the lamba2 and the mu2 will contain the computed second-order adjoint sensitivities, and can be used to produce Hessian-vector product (not the full Hessian matrix). Users must provide a direction vector; it is usually generated by an optimization solver. 788b5b4017aSHong Zhang 7893fca9621SHong Zhang Passing NULL for lambda2 disables the second-order calculation. 790db781477SPatrick Sanan .seealso: `TSAdjointSetForward()` 791b5b4017aSHong Zhang @*/ 7929371c9d4SSatish Balay PetscErrorCode TSSetCostHessianProducts(TS ts, PetscInt numcost, Vec *lambda2, Vec *mu2, Vec dir) { 793b5b4017aSHong Zhang PetscFunctionBegin; 794b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 7953c633725SBarry Smith PetscCheck(!ts->numcost || ts->numcost == numcost, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand"); 796b5b4017aSHong Zhang ts->numcost = numcost; 797b5b4017aSHong Zhang ts->vecs_sensi2 = lambda2; 79833f72c5dSHong Zhang ts->vecs_sensi2p = mu2; 799b5b4017aSHong Zhang ts->vec_dir = dir; 800b5b4017aSHong Zhang PetscFunctionReturn(0); 801b5b4017aSHong Zhang } 802b5b4017aSHong Zhang 803b5b4017aSHong Zhang /*@ 804b5b4017aSHong Zhang TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve() 805b5b4017aSHong Zhang 806b5b4017aSHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 807b5b4017aSHong Zhang 808b5b4017aSHong Zhang Input Parameter: 809b5b4017aSHong Zhang . ts - the TS context obtained from TSCreate() 810b5b4017aSHong Zhang 811d8d19677SJose E. Roman Output Parameters: 812b5b4017aSHong Zhang + numcost - number of cost functions 813b5b4017aSHong Zhang . lambda2 - Hessian-vector product with respect to the initial condition variables, the dimension and parallel layout of these vectors is the same as the ODE solution vector 814b5b4017aSHong Zhang . mu2 - Hessian-vector product with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 815b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 816b5b4017aSHong Zhang 817b5b4017aSHong Zhang Level: intermediate 818b5b4017aSHong Zhang 819db781477SPatrick Sanan .seealso: `TSSetCostHessianProducts()` 820b5b4017aSHong Zhang @*/ 8219371c9d4SSatish Balay PetscErrorCode TSGetCostHessianProducts(TS ts, PetscInt *numcost, Vec **lambda2, Vec **mu2, Vec *dir) { 822b5b4017aSHong Zhang PetscFunctionBegin; 823b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 824b5b4017aSHong Zhang if (numcost) *numcost = ts->numcost; 825b5b4017aSHong Zhang if (lambda2) *lambda2 = ts->vecs_sensi2; 82633f72c5dSHong Zhang if (mu2) *mu2 = ts->vecs_sensi2p; 827b5b4017aSHong Zhang if (dir) *dir = ts->vec_dir; 828b5b4017aSHong Zhang PetscFunctionReturn(0); 829b5b4017aSHong Zhang } 830b5b4017aSHong Zhang 831b5b4017aSHong Zhang /*@ 832ecf68647SHong Zhang TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities 833b5b4017aSHong Zhang 834d083f849SBarry Smith Logically Collective on TS 835b5b4017aSHong Zhang 836b5b4017aSHong Zhang Input Parameters: 837b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 838b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters 839b5b4017aSHong Zhang 840b5b4017aSHong Zhang Level: intermediate 841b5b4017aSHong Zhang 842ecf68647SHong Zhang Notes: When computing sensitivies w.r.t. initial condition, set didp to NULL so that the solver will take it as an identity matrix mathematically. TSAdjoint does not reset the tangent linear solver automatically, TSAdjointResetForward() should be called to reset the tangent linear solver. 843b5b4017aSHong Zhang 844db781477SPatrick Sanan .seealso: `TSSetCostHessianProducts()`, `TSAdjointResetForward()` 845b5b4017aSHong Zhang @*/ 8469371c9d4SSatish Balay PetscErrorCode TSAdjointSetForward(TS ts, Mat didp) { 84733f72c5dSHong Zhang Mat A; 84833f72c5dSHong Zhang Vec sp; 84933f72c5dSHong Zhang PetscScalar *xarr; 85033f72c5dSHong Zhang PetscInt lsize; 851b5b4017aSHong Zhang 852b5b4017aSHong Zhang PetscFunctionBegin; 853b5b4017aSHong Zhang ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 8543c633725SBarry Smith PetscCheck(ts->vecs_sensi2, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetCostHessianProducts() first"); 8553c633725SBarry Smith PetscCheck(ts->vec_dir, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 85633f72c5dSHong Zhang /* create a single-column dense matrix */ 8579566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol, &lsize)); 8589566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ts), lsize, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, &A)); 85933f72c5dSHong Zhang 8609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &sp)); 8619566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(A, 0, &xarr)); 8629566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(sp, xarr)); 863ecf68647SHong Zhang if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */ 86433f72c5dSHong Zhang if (didp) { 8659566063dSJacob Faibussowitsch PetscCall(MatMult(didp, ts->vec_dir, sp)); 8669566063dSJacob Faibussowitsch PetscCall(VecScale(sp, 2.)); 86733f72c5dSHong Zhang } else { 8689566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(sp)); 86933f72c5dSHong Zhang } 870ecf68647SHong Zhang } else { /* tangent linear variable initialized as dir */ 8719566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dir, sp)); 87233f72c5dSHong Zhang } 8739566063dSJacob Faibussowitsch PetscCall(VecResetArray(sp)); 8749566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(A, &xarr)); 8759566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sp)); 87633f72c5dSHong Zhang 8779566063dSJacob Faibussowitsch PetscCall(TSForwardSetInitialSensitivities(ts, A)); /* if didp is NULL, identity matrix is assumed */ 87833f72c5dSHong Zhang 8799566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 880b5b4017aSHong Zhang PetscFunctionReturn(0); 881b5b4017aSHong Zhang } 882b5b4017aSHong Zhang 883b5b4017aSHong Zhang /*@ 884ecf68647SHong Zhang TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context 885ecf68647SHong Zhang 886d083f849SBarry Smith Logically Collective on TS 887ecf68647SHong Zhang 888ecf68647SHong Zhang Input Parameters: 889ecf68647SHong Zhang . ts - the TS context obtained from TSCreate() 890ecf68647SHong Zhang 891ecf68647SHong Zhang Level: intermediate 892ecf68647SHong Zhang 893db781477SPatrick Sanan .seealso: `TSAdjointSetForward()` 894ecf68647SHong Zhang @*/ 8959371c9d4SSatish Balay PetscErrorCode TSAdjointResetForward(TS ts) { 896ecf68647SHong Zhang PetscFunctionBegin; 897ecf68647SHong Zhang ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */ 8989566063dSJacob Faibussowitsch PetscCall(TSForwardReset(ts)); 899ecf68647SHong Zhang PetscFunctionReturn(0); 900ecf68647SHong Zhang } 901ecf68647SHong Zhang 902ecf68647SHong Zhang /*@ 903814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 904814e85d6SHong Zhang of an adjoint solver 905814e85d6SHong Zhang 906814e85d6SHong Zhang Collective on TS 907814e85d6SHong Zhang 908814e85d6SHong Zhang Input Parameter: 909814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 910814e85d6SHong Zhang 911814e85d6SHong Zhang Level: advanced 912814e85d6SHong Zhang 913db781477SPatrick Sanan .seealso: `TSCreate()`, `TSAdjointStep()`, `TSSetCostGradients()` 914814e85d6SHong Zhang @*/ 9159371c9d4SSatish Balay PetscErrorCode TSAdjointSetUp(TS ts) { 916881c1a9bSHong Zhang TSTrajectory tj; 917881c1a9bSHong Zhang PetscBool match; 918814e85d6SHong Zhang 919814e85d6SHong Zhang PetscFunctionBegin; 920814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 921814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 9223c633725SBarry Smith PetscCheck(ts->vecs_sensi, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetCostGradients() first"); 9233c633725SBarry Smith PetscCheck(!ts->vecs_sensip || ts->Jacp || ts->Jacprhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 9249566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj)); 9259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tj, TSTRAJECTORYBASIC, &match)); 926881c1a9bSHong Zhang if (match) { 927881c1a9bSHong Zhang PetscBool solution_only; 9289566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetSolutionOnly(tj, &solution_only)); 9293c633725SBarry Smith PetscCheck(!solution_only, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "TSAdjoint cannot use the solution-only mode when choosing the Basic TSTrajectory type. Turn it off with -ts_trajectory_solution_only 0"); 930881c1a9bSHong Zhang } 9319566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUseHistory(tj, PETSC_FALSE)); /* not use TSHistory */ 93233f72c5dSHong Zhang 933cd4cee2dSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 9349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vecs_sensi[0], &ts->vec_drdu_col)); 935*48a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &ts->vec_drdp_col)); 936814e85d6SHong Zhang } 937814e85d6SHong Zhang 938dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, adjointsetup); 939814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 940814e85d6SHong Zhang PetscFunctionReturn(0); 941814e85d6SHong Zhang } 942814e85d6SHong Zhang 943814e85d6SHong Zhang /*@ 944ece46509SHong Zhang TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 945ece46509SHong Zhang 946ece46509SHong Zhang Collective on TS 947ece46509SHong Zhang 948ece46509SHong Zhang Input Parameter: 949ece46509SHong Zhang . ts - the TS context obtained from TSCreate() 950ece46509SHong Zhang 951ece46509SHong Zhang Level: beginner 952ece46509SHong Zhang 953db781477SPatrick Sanan .seealso: `TSCreate()`, `TSAdjointSetUp()`, `TSADestroy()` 954ece46509SHong Zhang @*/ 9559371c9d4SSatish Balay PetscErrorCode TSAdjointReset(TS ts) { 956ece46509SHong Zhang PetscFunctionBegin; 957ece46509SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 958dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, adjointreset); 9597207e0fdSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 9609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_drdu_col)); 961*48a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDestroy(&ts->vec_drdp_col)); 9627207e0fdSHong Zhang } 963ece46509SHong Zhang ts->vecs_sensi = NULL; 964ece46509SHong Zhang ts->vecs_sensip = NULL; 965ece46509SHong Zhang ts->vecs_sensi2 = NULL; 96633f72c5dSHong Zhang ts->vecs_sensi2p = NULL; 967ece46509SHong Zhang ts->vec_dir = NULL; 968ece46509SHong Zhang ts->adjointsetupcalled = PETSC_FALSE; 969ece46509SHong Zhang PetscFunctionReturn(0); 970ece46509SHong Zhang } 971ece46509SHong Zhang 972ece46509SHong Zhang /*@ 973814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 974814e85d6SHong Zhang 975814e85d6SHong Zhang Logically Collective on TS 976814e85d6SHong Zhang 977814e85d6SHong Zhang Input Parameters: 978814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 979a2b725a8SWilliam Gropp - steps - number of steps to use 980814e85d6SHong Zhang 981814e85d6SHong Zhang Level: intermediate 982814e85d6SHong Zhang 983814e85d6SHong Zhang Notes: 984814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 985814e85d6SHong Zhang so as to integrate back to less than the original timestep 986814e85d6SHong Zhang 987db781477SPatrick Sanan .seealso: `TSSetExactFinalTime()` 988814e85d6SHong Zhang @*/ 9899371c9d4SSatish Balay PetscErrorCode TSAdjointSetSteps(TS ts, PetscInt steps) { 990814e85d6SHong Zhang PetscFunctionBegin; 991814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 992814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts, steps, 2); 9933c633725SBarry Smith PetscCheck(steps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back a negative number of steps"); 9943c633725SBarry Smith PetscCheck(steps <= ts->steps, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back more than the total number of forward steps"); 995814e85d6SHong Zhang ts->adjoint_max_steps = steps; 996814e85d6SHong Zhang PetscFunctionReturn(0); 997814e85d6SHong Zhang } 998814e85d6SHong Zhang 999814e85d6SHong Zhang /*@C 1000814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1001814e85d6SHong Zhang 1002814e85d6SHong Zhang Level: deprecated 1003814e85d6SHong Zhang 1004814e85d6SHong Zhang @*/ 10059371c9d4SSatish Balay PetscErrorCode TSAdjointSetRHSJacobian(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) { 1006814e85d6SHong Zhang PetscFunctionBegin; 1007814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1008814e85d6SHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 1009814e85d6SHong Zhang 1010814e85d6SHong Zhang ts->rhsjacobianp = func; 1011814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 1012814e85d6SHong Zhang if (Amat) { 10139566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 10149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 1015814e85d6SHong Zhang ts->Jacp = Amat; 1016814e85d6SHong Zhang } 1017814e85d6SHong Zhang PetscFunctionReturn(0); 1018814e85d6SHong Zhang } 1019814e85d6SHong Zhang 1020814e85d6SHong Zhang /*@C 1021814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1022814e85d6SHong Zhang 1023814e85d6SHong Zhang Level: deprecated 1024814e85d6SHong Zhang 1025814e85d6SHong Zhang @*/ 10269371c9d4SSatish Balay PetscErrorCode TSAdjointComputeRHSJacobian(TS ts, PetscReal t, Vec U, Mat Amat) { 1027814e85d6SHong Zhang PetscFunctionBegin; 1028814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1029c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1030814e85d6SHong Zhang PetscValidPointer(Amat, 4); 1031814e85d6SHong Zhang 1032792fecdfSBarry Smith PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 1033814e85d6SHong Zhang PetscFunctionReturn(0); 1034814e85d6SHong Zhang } 1035814e85d6SHong Zhang 1036814e85d6SHong Zhang /*@ 1037d76be551SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 1038814e85d6SHong Zhang 1039814e85d6SHong Zhang Level: deprecated 1040814e85d6SHong Zhang 1041814e85d6SHong Zhang @*/ 10429371c9d4SSatish Balay PetscErrorCode TSAdjointComputeDRDYFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) { 1043814e85d6SHong Zhang PetscFunctionBegin; 1044814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1045c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1046814e85d6SHong Zhang 1047792fecdfSBarry Smith PetscCallBack("TS callback DRDY for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 1048814e85d6SHong Zhang PetscFunctionReturn(0); 1049814e85d6SHong Zhang } 1050814e85d6SHong Zhang 1051814e85d6SHong Zhang /*@ 1052d76be551SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 1053814e85d6SHong Zhang 1054814e85d6SHong Zhang Level: deprecated 1055814e85d6SHong Zhang 1056814e85d6SHong Zhang @*/ 10579371c9d4SSatish Balay PetscErrorCode TSAdjointComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) { 1058814e85d6SHong Zhang PetscFunctionBegin; 1059814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1060c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1061814e85d6SHong Zhang 1062792fecdfSBarry Smith PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 1063814e85d6SHong Zhang PetscFunctionReturn(0); 1064814e85d6SHong Zhang } 1065814e85d6SHong Zhang 1066814e85d6SHong Zhang /*@C 1067814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 1068814e85d6SHong Zhang 1069814e85d6SHong Zhang Level: intermediate 1070814e85d6SHong Zhang 1071db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1072814e85d6SHong Zhang @*/ 10739371c9d4SSatish Balay PetscErrorCode TSAdjointMonitorSensi(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) { 1074814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1075814e85d6SHong Zhang 1076814e85d6SHong Zhang PetscFunctionBegin; 1077064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 10789566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 10799566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], viewer)); 10809566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1081814e85d6SHong Zhang PetscFunctionReturn(0); 1082814e85d6SHong Zhang } 1083814e85d6SHong Zhang 1084814e85d6SHong Zhang /*@C 1085814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1086814e85d6SHong Zhang 1087814e85d6SHong Zhang Collective on TS 1088814e85d6SHong Zhang 1089814e85d6SHong Zhang Input Parameters: 1090814e85d6SHong Zhang + ts - TS object you wish to monitor 1091814e85d6SHong Zhang . name - the monitor type one is seeking 1092814e85d6SHong Zhang . help - message indicating what monitoring is done 1093814e85d6SHong Zhang . manual - manual page for the monitor 1094814e85d6SHong Zhang . monitor - the monitor function 1095814e85d6SHong Zhang - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects 1096814e85d6SHong Zhang 1097814e85d6SHong Zhang Level: developer 1098814e85d6SHong Zhang 1099db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 1100db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 1101db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1102db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1103c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1104db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1105db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 1106814e85d6SHong Zhang @*/ 11079371c9d4SSatish Balay PetscErrorCode TSAdjointMonitorSetFromOptions(TS ts, const char name[], const char help[], const char manual[], PetscErrorCode (*monitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, PetscViewerAndFormat *), PetscErrorCode (*monitorsetup)(TS, PetscViewerAndFormat *)) { 1108814e85d6SHong Zhang PetscViewer viewer; 1109814e85d6SHong Zhang PetscViewerFormat format; 1110814e85d6SHong Zhang PetscBool flg; 1111814e85d6SHong Zhang 1112814e85d6SHong Zhang PetscFunctionBegin; 11139566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg)); 1114814e85d6SHong Zhang if (flg) { 1115814e85d6SHong Zhang PetscViewerAndFormat *vf; 11169566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 11179566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 11181baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 11199566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 1120814e85d6SHong Zhang } 1121814e85d6SHong Zhang PetscFunctionReturn(0); 1122814e85d6SHong Zhang } 1123814e85d6SHong Zhang 1124814e85d6SHong Zhang /*@C 1125814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1126814e85d6SHong Zhang timestep to display the iteration's progress. 1127814e85d6SHong Zhang 1128814e85d6SHong Zhang Logically Collective on TS 1129814e85d6SHong Zhang 1130814e85d6SHong Zhang Input Parameters: 1131814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1132814e85d6SHong Zhang . adjointmonitor - monitoring routine 1133814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 1134814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 1135814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 1136814e85d6SHong Zhang (may be NULL) 1137814e85d6SHong Zhang 1138814e85d6SHong Zhang Calling sequence of monitor: 1139814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1140814e85d6SHong Zhang 1141814e85d6SHong Zhang + ts - the TS context 1142814e85d6SHong Zhang . steps - iteration number (after the final time step the monitor routine is called with a step of -1, this is at the final time which may have 1143814e85d6SHong Zhang been interpolated to) 1144814e85d6SHong Zhang . time - current time 1145814e85d6SHong Zhang . u - current iterate 1146814e85d6SHong Zhang . numcost - number of cost functionos 1147814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1148814e85d6SHong Zhang . mu - sensitivities to parameters 1149814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 1150814e85d6SHong Zhang 1151814e85d6SHong Zhang Notes: 1152814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 1153814e85d6SHong Zhang already has been loaded. 1154814e85d6SHong Zhang 1155814e85d6SHong Zhang Fortran Notes: 1156814e85d6SHong Zhang Only a single monitor function can be set for each TS object 1157814e85d6SHong Zhang 1158814e85d6SHong Zhang Level: intermediate 1159814e85d6SHong Zhang 1160db781477SPatrick Sanan .seealso: `TSAdjointMonitorCancel()` 1161814e85d6SHong Zhang @*/ 11629371c9d4SSatish Balay PetscErrorCode TSAdjointMonitorSet(TS ts, PetscErrorCode (*adjointmonitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *adjointmctx, PetscErrorCode (*adjointmdestroy)(void **)) { 1163814e85d6SHong Zhang PetscInt i; 1164814e85d6SHong Zhang PetscBool identical; 1165814e85d6SHong Zhang 1166814e85d6SHong Zhang PetscFunctionBegin; 1167814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1168814e85d6SHong Zhang for (i = 0; i < ts->numbermonitors; i++) { 11699566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))adjointmonitor, adjointmctx, adjointmdestroy, (PetscErrorCode(*)(void))ts->adjointmonitor[i], ts->adjointmonitorcontext[i], ts->adjointmonitordestroy[i], &identical)); 1170814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 1171814e85d6SHong Zhang } 11723c633725SBarry Smith PetscCheck(ts->numberadjointmonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many adjoint monitors set"); 1173814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1174814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1175814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void *)adjointmctx; 1176814e85d6SHong Zhang PetscFunctionReturn(0); 1177814e85d6SHong Zhang } 1178814e85d6SHong Zhang 1179814e85d6SHong Zhang /*@C 1180814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1181814e85d6SHong Zhang 1182814e85d6SHong Zhang Logically Collective on TS 1183814e85d6SHong Zhang 1184814e85d6SHong Zhang Input Parameters: 1185814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1186814e85d6SHong Zhang 1187814e85d6SHong Zhang Notes: 1188814e85d6SHong Zhang There is no way to remove a single, specific monitor. 1189814e85d6SHong Zhang 1190814e85d6SHong Zhang Level: intermediate 1191814e85d6SHong Zhang 1192db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1193814e85d6SHong Zhang @*/ 11949371c9d4SSatish Balay PetscErrorCode TSAdjointMonitorCancel(TS ts) { 1195814e85d6SHong Zhang PetscInt i; 1196814e85d6SHong Zhang 1197814e85d6SHong Zhang PetscFunctionBegin; 1198814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1199814e85d6SHong Zhang for (i = 0; i < ts->numberadjointmonitors; i++) { 1200*48a46eb9SPierre Jolivet if (ts->adjointmonitordestroy[i]) PetscCall((*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i])); 1201814e85d6SHong Zhang } 1202814e85d6SHong Zhang ts->numberadjointmonitors = 0; 1203814e85d6SHong Zhang PetscFunctionReturn(0); 1204814e85d6SHong Zhang } 1205814e85d6SHong Zhang 1206814e85d6SHong Zhang /*@C 1207814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 1208814e85d6SHong Zhang 1209814e85d6SHong Zhang Level: intermediate 1210814e85d6SHong Zhang 1211db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1212814e85d6SHong Zhang @*/ 12139371c9d4SSatish Balay PetscErrorCode TSAdjointMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) { 1214814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1215814e85d6SHong Zhang 1216814e85d6SHong Zhang PetscFunctionBegin; 1217064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 12189566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 12199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 122063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 12219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 12229566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1223814e85d6SHong Zhang PetscFunctionReturn(0); 1224814e85d6SHong Zhang } 1225814e85d6SHong Zhang 1226814e85d6SHong Zhang /*@C 1227814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1228814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 1229814e85d6SHong Zhang 1230814e85d6SHong Zhang Collective on TS 1231814e85d6SHong Zhang 1232814e85d6SHong Zhang Input Parameters: 1233814e85d6SHong Zhang + ts - the TS context 1234814e85d6SHong Zhang . step - current time-step 1235814e85d6SHong Zhang . ptime - current time 1236814e85d6SHong Zhang . u - current state 1237814e85d6SHong Zhang . numcost - number of cost functions 1238814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1239814e85d6SHong Zhang . mu - sensitivities to parameters 1240814e85d6SHong Zhang - dummy - either a viewer or NULL 1241814e85d6SHong Zhang 1242814e85d6SHong Zhang Level: intermediate 1243814e85d6SHong Zhang 1244db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()`, `TSAdjointMonitorDefault()`, `VecView()` 1245814e85d6SHong Zhang @*/ 12469371c9d4SSatish Balay PetscErrorCode TSAdjointMonitorDrawSensi(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *dummy) { 1247814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1248814e85d6SHong Zhang PetscDraw draw; 1249814e85d6SHong Zhang PetscReal xl, yl, xr, yr, h; 1250814e85d6SHong Zhang char time[32]; 1251814e85d6SHong Zhang 1252814e85d6SHong Zhang PetscFunctionBegin; 1253814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1254814e85d6SHong Zhang 12559566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], ictx->viewer)); 12569566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 12579566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 12589566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 1259814e85d6SHong Zhang h = yl + .95 * (yr - yl); 12609566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 12619566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1262814e85d6SHong Zhang PetscFunctionReturn(0); 1263814e85d6SHong Zhang } 1264814e85d6SHong Zhang 1265814e85d6SHong Zhang /* 1266814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1267814e85d6SHong Zhang 1268814e85d6SHong Zhang Collective on TSAdjoint 1269814e85d6SHong Zhang 1270814e85d6SHong Zhang Input Parameter: 1271814e85d6SHong Zhang . ts - the TS context 1272814e85d6SHong Zhang 1273814e85d6SHong Zhang Options Database Keys: 1274814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1275814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 1276814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1277814e85d6SHong Zhang 1278814e85d6SHong Zhang Level: developer 1279814e85d6SHong Zhang 1280814e85d6SHong Zhang Notes: 1281814e85d6SHong Zhang This is not normally called directly by users 1282814e85d6SHong Zhang 1283db781477SPatrick Sanan .seealso: `TSSetSaveTrajectory()`, `TSTrajectorySetUp()` 1284814e85d6SHong Zhang */ 12859371c9d4SSatish Balay PetscErrorCode TSAdjointSetFromOptions(TS ts, PetscOptionItems *PetscOptionsObject) { 1286814e85d6SHong Zhang PetscBool tflg, opt; 1287814e85d6SHong Zhang 1288814e85d6SHong Zhang PetscFunctionBegin; 1289dbbe0bcdSBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1290d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS Adjoint options"); 1291814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 12929566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adjoint_solve", "Solve the adjoint problem immediately after solving the forward problem", "", tflg, &tflg, &opt)); 1293814e85d6SHong Zhang if (opt) { 12949566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 1295814e85d6SHong Zhang ts->adjoint_solve = tflg; 1296814e85d6SHong Zhang } 12979566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor", "Monitor adjoint timestep size", "TSAdjointMonitorDefault", TSAdjointMonitorDefault, NULL)); 12989566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor_sensi", "Monitor sensitivity in the adjoint computation", "TSAdjointMonitorSensi", TSAdjointMonitorSensi, NULL)); 1299814e85d6SHong Zhang opt = PETSC_FALSE; 13009566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", &opt)); 1301814e85d6SHong Zhang if (opt) { 1302814e85d6SHong Zhang TSMonitorDrawCtx ctx; 1303814e85d6SHong Zhang PetscInt howoften = 1; 1304814e85d6SHong Zhang 13059566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", howoften, &howoften, NULL)); 13069566063dSJacob Faibussowitsch PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx)); 13079566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts, TSAdjointMonitorDrawSensi, ctx, (PetscErrorCode(*)(void **))TSMonitorDrawCtxDestroy)); 1308814e85d6SHong Zhang } 1309814e85d6SHong Zhang PetscFunctionReturn(0); 1310814e85d6SHong Zhang } 1311814e85d6SHong Zhang 1312814e85d6SHong Zhang /*@ 1313814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 1314814e85d6SHong Zhang 1315814e85d6SHong Zhang Collective on TS 1316814e85d6SHong Zhang 1317814e85d6SHong Zhang Input Parameter: 1318814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1319814e85d6SHong Zhang 1320814e85d6SHong Zhang Level: intermediate 1321814e85d6SHong Zhang 1322db781477SPatrick Sanan .seealso: `TSAdjointSetUp()`, `TSAdjointSolve()` 1323814e85d6SHong Zhang @*/ 13249371c9d4SSatish Balay PetscErrorCode TSAdjointStep(TS ts) { 1325814e85d6SHong Zhang DM dm; 1326814e85d6SHong Zhang 1327814e85d6SHong Zhang PetscFunctionBegin; 1328814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 13299566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 13309566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 13317b0e2f17SHong Zhang ts->steps--; /* must decrease the step index before the adjoint step is taken. */ 1332814e85d6SHong Zhang 1333814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1334814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 13359566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_AdjointStep, ts, 0, 0, 0)); 1336dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, adjointstep); 13379566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_AdjointStep, ts, 0, 0, 0)); 13387b0e2f17SHong Zhang ts->adjoint_steps++; 1339814e85d6SHong Zhang 1340814e85d6SHong Zhang if (ts->reason < 0) { 13413c633725SBarry Smith PetscCheck(!ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSAdjointStep has failed due to %s", TSConvergedReasons[ts->reason]); 1342814e85d6SHong Zhang } else if (!ts->reason) { 1343814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1344814e85d6SHong Zhang } 1345814e85d6SHong Zhang PetscFunctionReturn(0); 1346814e85d6SHong Zhang } 1347814e85d6SHong Zhang 1348814e85d6SHong Zhang /*@ 1349814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1350814e85d6SHong Zhang 1351814e85d6SHong Zhang Collective on TS 1352814e85d6SHong Zhang 1353814e85d6SHong Zhang Input Parameter: 1354814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1355814e85d6SHong Zhang 1356814e85d6SHong Zhang Options Database: 1357814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1358814e85d6SHong Zhang 1359814e85d6SHong Zhang Level: intermediate 1360814e85d6SHong Zhang 1361814e85d6SHong Zhang Notes: 1362814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 1363814e85d6SHong Zhang 1364814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1365814e85d6SHong Zhang 1366db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetCostGradients()`, `TSSetSolution()`, `TSAdjointStep()` 1367814e85d6SHong Zhang @*/ 13689371c9d4SSatish Balay PetscErrorCode TSAdjointSolve(TS ts) { 1369f1d62c27SHong Zhang static PetscBool cite = PETSC_FALSE; 13707f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 13717f59fb53SHong Zhang PetscLogStage adjoint_stage; 13727f59fb53SHong Zhang #endif 1373814e85d6SHong Zhang 1374814e85d6SHong Zhang PetscFunctionBegin; 1375814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1376421b783aSHong Zhang PetscCall(PetscCitationsRegister("@article{Zhang2022tsadjoint,\n" 1377421b783aSHong Zhang " title = {{PETSc TSAdjoint: A Discrete Adjoint ODE Solver for First-Order and Second-Order Sensitivity Analysis}},\n" 1378f1d62c27SHong Zhang " author = {Zhang, Hong and Constantinescu, Emil M. and Smith, Barry F.},\n" 1379421b783aSHong Zhang " journal = {SIAM Journal on Scientific Computing},\n" 1380421b783aSHong Zhang " volume = {44},\n" 1381421b783aSHong Zhang " number = {1},\n" 1382421b783aSHong Zhang " pages = {C1-C24},\n" 1383421b783aSHong Zhang " doi = {10.1137/21M140078X},\n" 13849371c9d4SSatish Balay " year = {2022}\n}\n", 13859371c9d4SSatish Balay &cite)); 13867f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 13879566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("TSAdjoint", &adjoint_stage)); 13889566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(adjoint_stage)); 13897f59fb53SHong Zhang #endif 13909566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 1391814e85d6SHong Zhang 1392814e85d6SHong Zhang /* reset time step and iteration counters */ 1393814e85d6SHong Zhang ts->adjoint_steps = 0; 1394814e85d6SHong Zhang ts->ksp_its = 0; 1395814e85d6SHong Zhang ts->snes_its = 0; 1396814e85d6SHong Zhang ts->num_snes_failures = 0; 1397814e85d6SHong Zhang ts->reject = 0; 1398814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1399814e85d6SHong Zhang 1400814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1401814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1402814e85d6SHong Zhang 1403814e85d6SHong Zhang while (!ts->reason) { 14049566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 14059566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 14069566063dSJacob Faibussowitsch PetscCall(TSAdjointEventHandler(ts)); 14079566063dSJacob Faibussowitsch PetscCall(TSAdjointStep(ts)); 1408*48a46eb9SPierre Jolivet if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) PetscCall(TSAdjointCostIntegral(ts)); 1409814e85d6SHong Zhang } 1410bdbff044SHong Zhang if (!ts->steps) { 14119566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 14129566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 1413bdbff044SHong Zhang } 1414814e85d6SHong Zhang ts->solvetime = ts->ptime; 14159566063dSJacob Faibussowitsch PetscCall(TSTrajectoryViewFromOptions(ts->trajectory, NULL, "-ts_trajectory_view")); 14169566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ts->vecs_sensi[0], (PetscObject)ts, "-ts_adjoint_view_solution")); 1417814e85d6SHong Zhang ts->adjoint_max_steps = 0; 14187f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 14199566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 14207f59fb53SHong Zhang #endif 1421814e85d6SHong Zhang PetscFunctionReturn(0); 1422814e85d6SHong Zhang } 1423814e85d6SHong Zhang 1424814e85d6SHong Zhang /*@C 1425814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1426814e85d6SHong Zhang 1427814e85d6SHong Zhang Collective on TS 1428814e85d6SHong Zhang 1429814e85d6SHong Zhang Input Parameters: 1430814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 1431814e85d6SHong Zhang . step - step number that has just completed 1432814e85d6SHong Zhang . ptime - model time of the state 1433814e85d6SHong Zhang . u - state at the current model time 1434814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 1435814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1436814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1437814e85d6SHong Zhang 1438814e85d6SHong Zhang Notes: 1439814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1440814e85d6SHong Zhang Users would almost never call this routine directly. 1441814e85d6SHong Zhang 1442814e85d6SHong Zhang Level: developer 1443814e85d6SHong Zhang 1444814e85d6SHong Zhang @*/ 14459371c9d4SSatish Balay PetscErrorCode TSAdjointMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu) { 1446814e85d6SHong Zhang PetscInt i, n = ts->numberadjointmonitors; 1447814e85d6SHong Zhang 1448814e85d6SHong Zhang PetscFunctionBegin; 1449814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1450814e85d6SHong Zhang PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 14519566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(u)); 1452*48a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall((*ts->adjointmonitor[i])(ts, step, ptime, u, numcost, lambda, mu, ts->adjointmonitorcontext[i])); 14539566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(u)); 1454814e85d6SHong Zhang PetscFunctionReturn(0); 1455814e85d6SHong Zhang } 1456814e85d6SHong Zhang 1457814e85d6SHong Zhang /*@ 1458814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1459814e85d6SHong Zhang 1460814e85d6SHong Zhang Collective on TS 1461814e85d6SHong Zhang 14624165533cSJose E. Roman Input Parameter: 1463814e85d6SHong Zhang . ts - time stepping context 1464814e85d6SHong Zhang 1465814e85d6SHong Zhang Level: advanced 1466814e85d6SHong Zhang 1467814e85d6SHong Zhang Notes: 1468814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 1469814e85d6SHong Zhang 1470db781477SPatrick Sanan .seealso: `TSAdjointSolve()`, `TSAdjointStep` 1471814e85d6SHong Zhang @*/ 14729371c9d4SSatish Balay PetscErrorCode TSAdjointCostIntegral(TS ts) { 1473362febeeSStefano Zampini PetscFunctionBegin; 1474814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1475dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, adjointintegral); 1476814e85d6SHong Zhang PetscFunctionReturn(0); 1477814e85d6SHong Zhang } 1478814e85d6SHong Zhang 1479814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1480814e85d6SHong Zhang 1481814e85d6SHong Zhang /*@ 1482814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 1483814e85d6SHong Zhang of forward sensitivity analysis 1484814e85d6SHong Zhang 1485814e85d6SHong Zhang Collective on TS 1486814e85d6SHong Zhang 1487814e85d6SHong Zhang Input Parameter: 1488814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1489814e85d6SHong Zhang 1490814e85d6SHong Zhang Level: advanced 1491814e85d6SHong Zhang 1492db781477SPatrick Sanan .seealso: `TSCreate()`, `TSDestroy()`, `TSSetUp()` 1493814e85d6SHong Zhang @*/ 14949371c9d4SSatish Balay PetscErrorCode TSForwardSetUp(TS ts) { 1495814e85d6SHong Zhang PetscFunctionBegin; 1496814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1497814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1498dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, forwardsetup); 14999566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ts->vec_sensip_col)); 1500814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 1501814e85d6SHong Zhang PetscFunctionReturn(0); 1502814e85d6SHong Zhang } 1503814e85d6SHong Zhang 1504814e85d6SHong Zhang /*@ 15057adebddeSHong Zhang TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 15067adebddeSHong Zhang 15077adebddeSHong Zhang Collective on TS 15087adebddeSHong Zhang 15097adebddeSHong Zhang Input Parameter: 15107adebddeSHong Zhang . ts - the TS context obtained from TSCreate() 15117adebddeSHong Zhang 15127adebddeSHong Zhang Level: advanced 15137adebddeSHong Zhang 1514db781477SPatrick Sanan .seealso: `TSCreate()`, `TSDestroy()`, `TSForwardSetUp()` 15157adebddeSHong Zhang @*/ 15169371c9d4SSatish Balay PetscErrorCode TSForwardReset(TS ts) { 15177207e0fdSHong Zhang TS quadts = ts->quadraturets; 15187adebddeSHong Zhang 15197adebddeSHong Zhang PetscFunctionBegin; 15207adebddeSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1521dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, forwardreset); 15229566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 1523*48a46eb9SPierre Jolivet if (quadts) PetscCall(MatDestroy(&quadts->mat_sensip)); 15249566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_sensip_col)); 15257adebddeSHong Zhang ts->forward_solve = PETSC_FALSE; 15267adebddeSHong Zhang ts->forwardsetupcalled = PETSC_FALSE; 15277adebddeSHong Zhang PetscFunctionReturn(0); 15287adebddeSHong Zhang } 15297adebddeSHong Zhang 15307adebddeSHong Zhang /*@ 1531814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1532814e85d6SHong Zhang 1533d8d19677SJose E. Roman Input Parameters: 1534a2b725a8SWilliam Gropp + ts - the TS context obtained from TSCreate() 1535814e85d6SHong Zhang . numfwdint - number of integrals 153667b8a455SSatish Balay - vp - the vectors containing the gradients for each integral w.r.t. parameters 1537814e85d6SHong Zhang 15387207e0fdSHong Zhang Level: deprecated 1539814e85d6SHong Zhang 1540db781477SPatrick Sanan .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1541814e85d6SHong Zhang @*/ 15429371c9d4SSatish Balay PetscErrorCode TSForwardSetIntegralGradients(TS ts, PetscInt numfwdint, Vec *vp) { 1543814e85d6SHong Zhang PetscFunctionBegin; 1544814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 15453c633725SBarry Smith PetscCheck(!ts->numcost || ts->numcost == numfwdint, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "The number of cost functions (2nd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand()"); 1546814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1547814e85d6SHong Zhang 1548814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 1549814e85d6SHong Zhang PetscFunctionReturn(0); 1550814e85d6SHong Zhang } 1551814e85d6SHong Zhang 1552814e85d6SHong Zhang /*@ 1553814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1554814e85d6SHong Zhang 1555814e85d6SHong Zhang Input Parameter: 1556814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1557814e85d6SHong Zhang 1558814e85d6SHong Zhang Output Parameter: 155967b8a455SSatish Balay . vp - the vectors containing the gradients for each integral w.r.t. parameters 1560814e85d6SHong Zhang 15617207e0fdSHong Zhang Level: deprecated 1562814e85d6SHong Zhang 1563db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1564814e85d6SHong Zhang @*/ 15659371c9d4SSatish Balay PetscErrorCode TSForwardGetIntegralGradients(TS ts, PetscInt *numfwdint, Vec **vp) { 1566814e85d6SHong Zhang PetscFunctionBegin; 1567814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1568814e85d6SHong Zhang PetscValidPointer(vp, 3); 1569814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1570814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1571814e85d6SHong Zhang PetscFunctionReturn(0); 1572814e85d6SHong Zhang } 1573814e85d6SHong Zhang 1574814e85d6SHong Zhang /*@ 1575814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1576814e85d6SHong Zhang 1577814e85d6SHong Zhang Collective on TS 1578814e85d6SHong Zhang 15794165533cSJose E. Roman Input Parameter: 1580814e85d6SHong Zhang . ts - time stepping context 1581814e85d6SHong Zhang 1582814e85d6SHong Zhang Level: advanced 1583814e85d6SHong Zhang 1584814e85d6SHong Zhang Notes: 1585814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1586814e85d6SHong Zhang 1587db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardSetUp()` 1588814e85d6SHong Zhang @*/ 15899371c9d4SSatish Balay PetscErrorCode TSForwardStep(TS ts) { 1590362febeeSStefano Zampini PetscFunctionBegin; 1591814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 15929566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_ForwardStep, ts, 0, 0, 0)); 1593dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, forwardstep); 15949566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_ForwardStep, ts, 0, 0, 0)); 15953c633725SBarry Smith PetscCheck(ts->reason >= 0 || !ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSFowardStep has failed due to %s", TSConvergedReasons[ts->reason]); 1596814e85d6SHong Zhang PetscFunctionReturn(0); 1597814e85d6SHong Zhang } 1598814e85d6SHong Zhang 1599814e85d6SHong Zhang /*@ 1600814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1601814e85d6SHong Zhang 1602d083f849SBarry Smith Logically Collective on TS 1603814e85d6SHong Zhang 1604814e85d6SHong Zhang Input Parameters: 1605814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1606814e85d6SHong Zhang . nump - number of parameters 1607814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1608814e85d6SHong Zhang 1609814e85d6SHong Zhang Level: beginner 1610814e85d6SHong Zhang 1611814e85d6SHong Zhang Notes: 1612814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1613814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1614814e85d6SHong Zhang You must call this function before TSSolve(). 1615814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1616814e85d6SHong Zhang 1617db781477SPatrick Sanan .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1618814e85d6SHong Zhang @*/ 16199371c9d4SSatish Balay PetscErrorCode TSForwardSetSensitivities(TS ts, PetscInt nump, Mat Smat) { 1620814e85d6SHong Zhang PetscFunctionBegin; 1621814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1622814e85d6SHong Zhang PetscValidHeaderSpecific(Smat, MAT_CLASSID, 3); 1623814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1624b5b4017aSHong Zhang if (nump == PETSC_DEFAULT) { 16259566063dSJacob Faibussowitsch PetscCall(MatGetSize(Smat, NULL, &ts->num_parameters)); 1626b5b4017aSHong Zhang } else ts->num_parameters = nump; 16279566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Smat)); 16289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 1629814e85d6SHong Zhang ts->mat_sensip = Smat; 1630814e85d6SHong Zhang PetscFunctionReturn(0); 1631814e85d6SHong Zhang } 1632814e85d6SHong Zhang 1633814e85d6SHong Zhang /*@ 1634814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1635814e85d6SHong Zhang 1636814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1637814e85d6SHong Zhang 1638d8d19677SJose E. Roman Output Parameters: 1639814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1640814e85d6SHong Zhang . nump - number of parameters 1641814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1642814e85d6SHong Zhang 1643814e85d6SHong Zhang Level: intermediate 1644814e85d6SHong Zhang 1645db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1646814e85d6SHong Zhang @*/ 16479371c9d4SSatish Balay PetscErrorCode TSForwardGetSensitivities(TS ts, PetscInt *nump, Mat *Smat) { 1648814e85d6SHong Zhang PetscFunctionBegin; 1649814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1650814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1651814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1652814e85d6SHong Zhang PetscFunctionReturn(0); 1653814e85d6SHong Zhang } 1654814e85d6SHong Zhang 1655814e85d6SHong Zhang /*@ 1656814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1657814e85d6SHong Zhang 1658814e85d6SHong Zhang Collective on TS 1659814e85d6SHong Zhang 16604165533cSJose E. Roman Input Parameter: 1661814e85d6SHong Zhang . ts - time stepping context 1662814e85d6SHong Zhang 1663814e85d6SHong Zhang Level: advanced 1664814e85d6SHong Zhang 1665814e85d6SHong Zhang Notes: 1666814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1667814e85d6SHong Zhang 1668db781477SPatrick Sanan .seealso: `TSSolve()`, `TSAdjointCostIntegral()` 1669814e85d6SHong Zhang @*/ 16709371c9d4SSatish Balay PetscErrorCode TSForwardCostIntegral(TS ts) { 1671362febeeSStefano Zampini PetscFunctionBegin; 1672814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1673dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, forwardintegral); 1674814e85d6SHong Zhang PetscFunctionReturn(0); 1675814e85d6SHong Zhang } 1676b5b4017aSHong Zhang 1677b5b4017aSHong Zhang /*@ 1678b5b4017aSHong Zhang TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1679b5b4017aSHong Zhang 1680d083f849SBarry Smith Collective on TS 1681b5b4017aSHong Zhang 1682d8d19677SJose E. Roman Input Parameters: 1683b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 1684b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition 1685b5b4017aSHong Zhang 1686b5b4017aSHong Zhang Level: intermediate 1687b5b4017aSHong Zhang 1688b5b4017aSHong Zhang Notes: TSSolve() allows users to pass the initial solution directly to TS. But the tangent linear variables cannot be initialized in this way. This function is used to set initial values for tangent linear variables. 1689b5b4017aSHong Zhang 1690db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()` 1691b5b4017aSHong Zhang @*/ 16929371c9d4SSatish Balay PetscErrorCode TSForwardSetInitialSensitivities(TS ts, Mat didp) { 1693362febeeSStefano Zampini PetscFunctionBegin; 1694b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1695b5b4017aSHong Zhang PetscValidHeaderSpecific(didp, MAT_CLASSID, 2); 1696*48a46eb9SPierre Jolivet if (!ts->mat_sensip) PetscCall(TSForwardSetSensitivities(ts, PETSC_DEFAULT, didp)); 1697b5b4017aSHong Zhang PetscFunctionReturn(0); 1698b5b4017aSHong Zhang } 16996affb6f8SHong Zhang 17006affb6f8SHong Zhang /*@ 17016affb6f8SHong Zhang TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 17026affb6f8SHong Zhang 17036affb6f8SHong Zhang Input Parameter: 17046affb6f8SHong Zhang . ts - the TS context obtained from TSCreate() 17056affb6f8SHong Zhang 17066affb6f8SHong Zhang Output Parameters: 1707cd4cee2dSHong Zhang + ns - number of stages 17086affb6f8SHong Zhang - S - tangent linear sensitivities at the intermediate stages 17096affb6f8SHong Zhang 17106affb6f8SHong Zhang Level: advanced 17116affb6f8SHong Zhang 17126affb6f8SHong Zhang @*/ 17139371c9d4SSatish Balay PetscErrorCode TSForwardGetStages(TS ts, PetscInt *ns, Mat **S) { 17146affb6f8SHong Zhang PetscFunctionBegin; 17156affb6f8SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 17166affb6f8SHong Zhang 17176affb6f8SHong Zhang if (!ts->ops->getstages) *S = NULL; 1718dbbe0bcdSBarry Smith else PetscUseTypeMethod(ts, forwardgetstages, ns, S); 17196affb6f8SHong Zhang PetscFunctionReturn(0); 17206affb6f8SHong Zhang } 1721cd4cee2dSHong Zhang 1722cd4cee2dSHong Zhang /*@ 1723cd4cee2dSHong Zhang TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 1724cd4cee2dSHong Zhang 1725d8d19677SJose E. Roman Input Parameters: 1726cd4cee2dSHong Zhang + ts - the TS context obtained from TSCreate() 1727cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1728cd4cee2dSHong Zhang 1729cd4cee2dSHong Zhang Output Parameters: 1730cd4cee2dSHong Zhang . quadts - the child TS context 1731cd4cee2dSHong Zhang 1732cd4cee2dSHong Zhang Level: intermediate 1733cd4cee2dSHong Zhang 1734db781477SPatrick Sanan .seealso: `TSGetQuadratureTS()` 1735cd4cee2dSHong Zhang @*/ 17369371c9d4SSatish Balay PetscErrorCode TSCreateQuadratureTS(TS ts, PetscBool fwd, TS *quadts) { 1737cd4cee2dSHong Zhang char prefix[128]; 1738cd4cee2dSHong Zhang 1739cd4cee2dSHong Zhang PetscFunctionBegin; 1740cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1741064a246eSJacob Faibussowitsch PetscValidPointer(quadts, 3); 17429566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts->quadraturets)); 17439566063dSJacob Faibussowitsch PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &ts->quadraturets)); 17449566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets, (PetscObject)ts, 1)); 17459566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts, (PetscObject)ts->quadraturets)); 17469566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%squad_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "")); 17479566063dSJacob Faibussowitsch PetscCall(TSSetOptionsPrefix(ts->quadraturets, prefix)); 1748cd4cee2dSHong Zhang *quadts = ts->quadraturets; 1749cd4cee2dSHong Zhang 1750cd4cee2dSHong Zhang if (ts->numcost) { 17519566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, ts->numcost, &(*quadts)->vec_sol)); 1752cd4cee2dSHong Zhang } else { 17539566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &(*quadts)->vec_sol)); 1754cd4cee2dSHong Zhang } 1755cd4cee2dSHong Zhang ts->costintegralfwd = fwd; 1756cd4cee2dSHong Zhang PetscFunctionReturn(0); 1757cd4cee2dSHong Zhang } 1758cd4cee2dSHong Zhang 1759cd4cee2dSHong Zhang /*@ 1760cd4cee2dSHong Zhang TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 1761cd4cee2dSHong Zhang 1762cd4cee2dSHong Zhang Input Parameter: 1763cd4cee2dSHong Zhang . ts - the TS context obtained from TSCreate() 1764cd4cee2dSHong Zhang 1765cd4cee2dSHong Zhang Output Parameters: 1766cd4cee2dSHong Zhang + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1767cd4cee2dSHong Zhang - quadts - the child TS context 1768cd4cee2dSHong Zhang 1769cd4cee2dSHong Zhang Level: intermediate 1770cd4cee2dSHong Zhang 1771db781477SPatrick Sanan .seealso: `TSCreateQuadratureTS()` 1772cd4cee2dSHong Zhang @*/ 17739371c9d4SSatish Balay PetscErrorCode TSGetQuadratureTS(TS ts, PetscBool *fwd, TS *quadts) { 1774cd4cee2dSHong Zhang PetscFunctionBegin; 1775cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1776cd4cee2dSHong Zhang if (fwd) *fwd = ts->costintegralfwd; 1777cd4cee2dSHong Zhang if (quadts) *quadts = ts->quadraturets; 1778cd4cee2dSHong Zhang PetscFunctionReturn(0); 1779cd4cee2dSHong Zhang } 1780ba0675f6SHong Zhang 1781ba0675f6SHong Zhang /*@ 1782ba0675f6SHong Zhang TSComputeSNESJacobian - Compute the SNESJacobian 1783ba0675f6SHong Zhang 1784ba0675f6SHong Zhang Input Parameters: 1785ba0675f6SHong Zhang + ts - the TS context obtained from TSCreate() 1786ba0675f6SHong Zhang - x - state vector 1787ba0675f6SHong Zhang 1788ba0675f6SHong Zhang Output Parameters: 1789ba0675f6SHong Zhang + J - Jacobian matrix 1790ba0675f6SHong Zhang - Jpre - preconditioning matrix for J (may be same as J) 1791ba0675f6SHong Zhang 1792ba0675f6SHong Zhang Level: developer 1793ba0675f6SHong Zhang 1794ba0675f6SHong Zhang Notes: 1795ba0675f6SHong Zhang Using SNES to compute the Jacobian enables finite differencing when TS Jacobian is not available. 1796ba0675f6SHong Zhang @*/ 17979371c9d4SSatish Balay PetscErrorCode TSComputeSNESJacobian(TS ts, Vec x, Mat J, Mat Jpre) { 1798ba0675f6SHong Zhang SNES snes = ts->snes; 1799ba0675f6SHong Zhang PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *) = NULL; 1800ba0675f6SHong Zhang 1801ba0675f6SHong Zhang PetscFunctionBegin; 1802ba0675f6SHong Zhang /* 1803ba0675f6SHong Zhang Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object 1804ba0675f6SHong Zhang because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for 1805ba0675f6SHong Zhang explicit methods. Instead, we check the Jacobian compute function directly to determin if FD 1806ba0675f6SHong Zhang coloring is used. 1807ba0675f6SHong Zhang */ 18089566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, NULL, NULL, &jac, NULL)); 1809ba0675f6SHong Zhang if (jac == SNESComputeJacobianDefaultColor) { 1810ba0675f6SHong Zhang Vec f; 18119566063dSJacob Faibussowitsch PetscCall(SNESSetSolution(snes, x)); 18129566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &f, NULL, NULL)); 1813ba0675f6SHong Zhang /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */ 18149566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, x, f)); 1815ba0675f6SHong Zhang } 18169566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, x, J, Jpre)); 1817ba0675f6SHong Zhang PetscFunctionReturn(0); 1818ba0675f6SHong Zhang } 1819