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 @*/ 35*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetRHSJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) 36*d71ae5a4SJacob Faibussowitsch { 37814e85d6SHong Zhang PetscFunctionBegin; 38814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 39814e85d6SHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 40814e85d6SHong Zhang 41814e85d6SHong Zhang ts->rhsjacobianp = func; 42814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 43814e85d6SHong Zhang if (Amat) { 449566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacprhs)); 4633f72c5dSHong Zhang ts->Jacprhs = Amat; 47814e85d6SHong Zhang } 48814e85d6SHong Zhang PetscFunctionReturn(0); 49814e85d6SHong Zhang } 50814e85d6SHong Zhang 51814e85d6SHong Zhang /*@C 52cd4cee2dSHong 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. 53cd4cee2dSHong Zhang 54cd4cee2dSHong Zhang Logically Collective on TS 55cd4cee2dSHong Zhang 56f899ff85SJose E. Roman Input Parameter: 57cd4cee2dSHong Zhang . ts - TS context obtained from TSCreate() 58cd4cee2dSHong Zhang 59cd4cee2dSHong Zhang Output Parameters: 60cd4cee2dSHong Zhang + Amat - JacobianP matrix 61cd4cee2dSHong Zhang . func - function 62cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 63cd4cee2dSHong Zhang 64cd4cee2dSHong Zhang Calling sequence of func: 65cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 66cd4cee2dSHong Zhang + t - current timestep 67cd4cee2dSHong Zhang . U - input vector (current ODE solution) 68cd4cee2dSHong Zhang . A - output matrix 69cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 70cd4cee2dSHong Zhang 71cd4cee2dSHong Zhang Level: intermediate 72cd4cee2dSHong Zhang 73cd4cee2dSHong Zhang Notes: 74cd4cee2dSHong 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 75cd4cee2dSHong Zhang 76db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()` 77cd4cee2dSHong Zhang @*/ 78*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetRHSJacobianP(TS ts, Mat *Amat, PetscErrorCode (**func)(TS, PetscReal, Vec, Mat, void *), void **ctx) 79*d71ae5a4SJacob Faibussowitsch { 80cd4cee2dSHong Zhang PetscFunctionBegin; 81cd4cee2dSHong Zhang if (func) *func = ts->rhsjacobianp; 82cd4cee2dSHong Zhang if (ctx) *ctx = ts->rhsjacobianpctx; 83cd4cee2dSHong Zhang if (Amat) *Amat = ts->Jacprhs; 84cd4cee2dSHong Zhang PetscFunctionReturn(0); 85cd4cee2dSHong Zhang } 86cd4cee2dSHong Zhang 87cd4cee2dSHong Zhang /*@C 88814e85d6SHong Zhang TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 89814e85d6SHong Zhang 90814e85d6SHong Zhang Collective on TS 91814e85d6SHong Zhang 92814e85d6SHong Zhang Input Parameters: 93814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 94814e85d6SHong Zhang 95814e85d6SHong Zhang Level: developer 96814e85d6SHong Zhang 97db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()` 98814e85d6SHong Zhang @*/ 99*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSJacobianP(TS ts, PetscReal t, Vec U, Mat Amat) 100*d71ae5a4SJacob Faibussowitsch { 101814e85d6SHong Zhang PetscFunctionBegin; 10233f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 103814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 104c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 105814e85d6SHong Zhang 106792fecdfSBarry Smith PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 107814e85d6SHong Zhang PetscFunctionReturn(0); 108814e85d6SHong Zhang } 109814e85d6SHong Zhang 110814e85d6SHong Zhang /*@C 11133f72c5dSHong 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. 11233f72c5dSHong Zhang 11333f72c5dSHong Zhang Logically Collective on TS 11433f72c5dSHong Zhang 11533f72c5dSHong Zhang Input Parameters: 11633f72c5dSHong Zhang + ts - TS context obtained from TSCreate() 11733f72c5dSHong Zhang . Amat - JacobianP matrix 11833f72c5dSHong Zhang . func - function 11933f72c5dSHong Zhang - ctx - [optional] user-defined function context 12033f72c5dSHong Zhang 12133f72c5dSHong Zhang Calling sequence of func: 12233f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 12333f72c5dSHong Zhang + t - current timestep 12433f72c5dSHong Zhang . U - input vector (current ODE solution) 12533f72c5dSHong Zhang . Udot - time derivative of state vector 12633f72c5dSHong Zhang . shift - shift to apply, see note below 12733f72c5dSHong Zhang . A - output matrix 12833f72c5dSHong Zhang - ctx - [optional] user-defined function context 12933f72c5dSHong Zhang 13033f72c5dSHong Zhang Level: intermediate 13133f72c5dSHong Zhang 13233f72c5dSHong Zhang Notes: 13333f72c5dSHong 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 13433f72c5dSHong Zhang 13533f72c5dSHong Zhang .seealso: 13633f72c5dSHong Zhang @*/ 137*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetIJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *ctx) 138*d71ae5a4SJacob Faibussowitsch { 13933f72c5dSHong Zhang PetscFunctionBegin; 14033f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 14133f72c5dSHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 14233f72c5dSHong Zhang 14333f72c5dSHong Zhang ts->ijacobianp = func; 14433f72c5dSHong Zhang ts->ijacobianpctx = ctx; 14533f72c5dSHong Zhang if (Amat) { 1469566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 1479566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 14833f72c5dSHong Zhang ts->Jacp = Amat; 14933f72c5dSHong Zhang } 15033f72c5dSHong Zhang PetscFunctionReturn(0); 15133f72c5dSHong Zhang } 15233f72c5dSHong Zhang 15333f72c5dSHong Zhang /*@C 15433f72c5dSHong Zhang TSComputeIJacobianP - Runs the user-defined IJacobianP function. 15533f72c5dSHong Zhang 15633f72c5dSHong Zhang Collective on TS 15733f72c5dSHong Zhang 15833f72c5dSHong Zhang Input Parameters: 15933f72c5dSHong Zhang + ts - the TS context 16033f72c5dSHong Zhang . t - current timestep 16133f72c5dSHong Zhang . U - state vector 16233f72c5dSHong Zhang . Udot - time derivative of state vector 16333f72c5dSHong Zhang . shift - shift to apply, see note below 16433f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 16533f72c5dSHong Zhang 16633f72c5dSHong Zhang Output Parameters: 16733f72c5dSHong Zhang . A - Jacobian matrix 16833f72c5dSHong Zhang 16933f72c5dSHong Zhang Level: developer 17033f72c5dSHong Zhang 171db781477SPatrick Sanan .seealso: `TSSetIJacobianP()` 17233f72c5dSHong Zhang @*/ 173*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIJacobianP(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat Amat, PetscBool imex) 174*d71ae5a4SJacob Faibussowitsch { 17533f72c5dSHong Zhang PetscFunctionBegin; 17633f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 17733f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 17833f72c5dSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 17933f72c5dSHong Zhang PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4); 18033f72c5dSHong Zhang 1819566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_JacobianPEval, ts, U, Amat, 0)); 18248a46eb9SPierre Jolivet if (ts->ijacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->ijacobianp)(ts, t, U, Udot, shift, Amat, ts->ijacobianpctx)); 18333f72c5dSHong Zhang if (imex) { 18433f72c5dSHong Zhang if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */ 18533f72c5dSHong Zhang PetscBool assembled; 1869566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 1879566063dSJacob Faibussowitsch PetscCall(MatAssembled(Amat, &assembled)); 18833f72c5dSHong Zhang if (!assembled) { 1899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY)); 1909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY)); 19133f72c5dSHong Zhang } 19233f72c5dSHong Zhang } 19333f72c5dSHong Zhang } else { 1941baa6e33SBarry Smith if (ts->rhsjacobianp) PetscCall(TSComputeRHSJacobianP(ts, t, U, ts->Jacprhs)); 19533f72c5dSHong Zhang if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */ 1969566063dSJacob Faibussowitsch PetscCall(MatScale(Amat, -1)); 19733f72c5dSHong Zhang } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */ 19833f72c5dSHong Zhang MatStructure axpy = DIFFERENT_NONZERO_PATTERN; 19933f72c5dSHong Zhang if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */ 2009566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 20133f72c5dSHong Zhang } 2029566063dSJacob Faibussowitsch PetscCall(MatAXPY(Amat, -1, ts->Jacprhs, axpy)); 20333f72c5dSHong Zhang } 20433f72c5dSHong Zhang } 2059566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_JacobianPEval, ts, U, Amat, 0)); 20633f72c5dSHong Zhang PetscFunctionReturn(0); 20733f72c5dSHong Zhang } 20833f72c5dSHong Zhang 20933f72c5dSHong Zhang /*@C 210814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 211814e85d6SHong Zhang 212814e85d6SHong Zhang Logically Collective on TS 213814e85d6SHong Zhang 214814e85d6SHong Zhang Input Parameters: 215814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 216814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 217814e85d6SHong Zhang . costintegral - vector that stores the integral values 218814e85d6SHong Zhang . rf - routine for evaluating the integrand function 219c9ad9de0SHong Zhang . drduf - function that computes the gradients of the r's with respect to u 220814e85d6SHong 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) 221814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 222814e85d6SHong Zhang - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 223814e85d6SHong Zhang 224814e85d6SHong Zhang Calling sequence of rf: 225c9ad9de0SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx); 226814e85d6SHong Zhang 227c9ad9de0SHong Zhang Calling sequence of drduf: 228c9ad9de0SHong Zhang $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx); 229814e85d6SHong Zhang 230814e85d6SHong Zhang Calling sequence of drdpf: 231c9ad9de0SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx); 232814e85d6SHong Zhang 233cd4cee2dSHong Zhang Level: deprecated 234814e85d6SHong Zhang 235814e85d6SHong Zhang Notes: 236814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 237814e85d6SHong Zhang 238db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()`, `TSGetCostGradients()`, `TSSetCostGradients()` 239814e85d6SHong Zhang @*/ 240*d71ae5a4SJacob Faibussowitsch 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) 241*d71ae5a4SJacob Faibussowitsch { 242814e85d6SHong Zhang PetscFunctionBegin; 243814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 244814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral, VEC_CLASSID, 3); 2453c633725SBarry 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()"); 246814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numcost; 247814e85d6SHong Zhang 248814e85d6SHong Zhang if (costintegral) { 2499566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)costintegral)); 2509566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_costintegral)); 251814e85d6SHong Zhang ts->vec_costintegral = costintegral; 252814e85d6SHong Zhang } else { 253814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 2549566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, numcost, &ts->vec_costintegral)); 255814e85d6SHong Zhang } else { 2569566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegral, 0.0)); 257814e85d6SHong Zhang } 258814e85d6SHong Zhang } 259814e85d6SHong Zhang if (!ts->vec_costintegrand) { 2609566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_costintegral, &ts->vec_costintegrand)); 261814e85d6SHong Zhang } else { 2629566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegrand, 0.0)); 263814e85d6SHong Zhang } 264814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 265814e85d6SHong Zhang ts->costintegrand = rf; 266814e85d6SHong Zhang ts->costintegrandctx = ctx; 267c9ad9de0SHong Zhang ts->drdufunction = drduf; 268814e85d6SHong Zhang ts->drdpfunction = drdpf; 269814e85d6SHong Zhang PetscFunctionReturn(0); 270814e85d6SHong Zhang } 271814e85d6SHong Zhang 272b5b4017aSHong Zhang /*@C 273814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 274814e85d6SHong Zhang It is valid to call the routine after a backward run. 275814e85d6SHong Zhang 276814e85d6SHong Zhang Not Collective 277814e85d6SHong Zhang 278814e85d6SHong Zhang Input Parameter: 279814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 280814e85d6SHong Zhang 281814e85d6SHong Zhang Output Parameter: 282814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 283814e85d6SHong Zhang 284814e85d6SHong Zhang Level: intermediate 285814e85d6SHong Zhang 286db781477SPatrick Sanan .seealso: `TSSetCostIntegrand()` 287814e85d6SHong Zhang 288814e85d6SHong Zhang @*/ 289*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostIntegral(TS ts, Vec *v) 290*d71ae5a4SJacob Faibussowitsch { 291cd4cee2dSHong Zhang TS quadts; 292cd4cee2dSHong Zhang 293814e85d6SHong Zhang PetscFunctionBegin; 294814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 295814e85d6SHong Zhang PetscValidPointer(v, 2); 2969566063dSJacob Faibussowitsch PetscCall(TSGetQuadratureTS(ts, NULL, &quadts)); 297cd4cee2dSHong Zhang *v = quadts->vec_sol; 298814e85d6SHong Zhang PetscFunctionReturn(0); 299814e85d6SHong Zhang } 300814e85d6SHong Zhang 301b5b4017aSHong Zhang /*@C 302814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 303814e85d6SHong Zhang 304814e85d6SHong Zhang Input Parameters: 305814e85d6SHong Zhang + ts - the TS context 306814e85d6SHong Zhang . t - current time 307c9ad9de0SHong Zhang - U - state vector, i.e. current solution 308814e85d6SHong Zhang 309814e85d6SHong Zhang Output Parameter: 310c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs 311814e85d6SHong Zhang 312369cf35fSHong Zhang Notes: 313814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 314814e85d6SHong Zhang is used internally within the sensitivity analysis context. 315814e85d6SHong Zhang 316cd4cee2dSHong Zhang Level: deprecated 317814e85d6SHong Zhang 318db781477SPatrick Sanan .seealso: `TSSetCostIntegrand()` 319814e85d6SHong Zhang @*/ 320*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeCostIntegrand(TS ts, PetscReal t, Vec U, Vec Q) 321*d71ae5a4SJacob Faibussowitsch { 322814e85d6SHong Zhang PetscFunctionBegin; 323814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 324c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 325c9ad9de0SHong Zhang PetscValidHeaderSpecific(Q, VEC_CLASSID, 4); 326814e85d6SHong Zhang 3279566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_FunctionEval, ts, U, Q, 0)); 328792fecdfSBarry Smith if (ts->costintegrand) PetscCallBack("TS callback integrand in the cost function", (*ts->costintegrand)(ts, t, U, Q, ts->costintegrandctx)); 329ef1023bdSBarry Smith else PetscCall(VecZeroEntries(Q)); 3309566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_FunctionEval, ts, U, Q, 0)); 331814e85d6SHong Zhang PetscFunctionReturn(0); 332814e85d6SHong Zhang } 333814e85d6SHong Zhang 334b5b4017aSHong Zhang /*@C 335d76be551SHong Zhang TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 336814e85d6SHong Zhang 337d76be551SHong Zhang Level: deprecated 338814e85d6SHong Zhang 339814e85d6SHong Zhang @*/ 340*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDUFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) 341*d71ae5a4SJacob Faibussowitsch { 342814e85d6SHong Zhang PetscFunctionBegin; 34333f72c5dSHong Zhang if (!DRDU) PetscFunctionReturn(0); 344814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 345c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 346814e85d6SHong Zhang 347792fecdfSBarry Smith PetscCallBack("TS callback DRDU for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 348814e85d6SHong Zhang PetscFunctionReturn(0); 349814e85d6SHong Zhang } 350814e85d6SHong Zhang 351b5b4017aSHong Zhang /*@C 352d76be551SHong Zhang TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 353814e85d6SHong Zhang 354d76be551SHong Zhang Level: deprecated 355814e85d6SHong Zhang 356814e85d6SHong Zhang @*/ 357*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) 358*d71ae5a4SJacob Faibussowitsch { 359814e85d6SHong Zhang PetscFunctionBegin; 36033f72c5dSHong Zhang if (!DRDP) PetscFunctionReturn(0); 361814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 362c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 363814e85d6SHong Zhang 364792fecdfSBarry Smith PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 365814e85d6SHong Zhang PetscFunctionReturn(0); 366814e85d6SHong Zhang } 367814e85d6SHong Zhang 368b5b4017aSHong Zhang /*@C 369b5b4017aSHong 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. 370b5b4017aSHong Zhang 371b5b4017aSHong Zhang Logically Collective on TS 372b5b4017aSHong Zhang 373b5b4017aSHong Zhang Input Parameters: 374b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 375b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 376b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 377b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 378b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 379b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 380b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 381b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 382f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for F_PP 383b5b4017aSHong Zhang 384b5b4017aSHong Zhang Calling sequence of ihessianproductfunc: 385369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 386b5b4017aSHong Zhang + t - current timestep 387b5b4017aSHong Zhang . U - input vector (current ODE solution) 388369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 389b5b4017aSHong Zhang . Vr - input vector to be right-multiplied with the Hessian 390369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 391b5b4017aSHong Zhang - ctx - [optional] user-defined function context 392b5b4017aSHong Zhang 393b5b4017aSHong Zhang Level: intermediate 394b5b4017aSHong Zhang 395369cf35fSHong Zhang Notes: 396369cf35fSHong Zhang The first Hessian function and the working array are required. 397369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 398369cf35fSHong Zhang $ Vl_n^T*F_UP*Vr 399369cf35fSHong 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. 400369cf35fSHong Zhang Each entry of F_UP corresponds to the derivative 401369cf35fSHong Zhang $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}. 402369cf35fSHong 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 403369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]} 404369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 405b5b4017aSHong Zhang 406b5b4017aSHong Zhang .seealso: 407b5b4017aSHong Zhang @*/ 408*d71ae5a4SJacob Faibussowitsch 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) 409*d71ae5a4SJacob Faibussowitsch { 410b5b4017aSHong Zhang PetscFunctionBegin; 411b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 412b5b4017aSHong Zhang PetscValidPointer(ihp1, 2); 413b5b4017aSHong Zhang 414b5b4017aSHong Zhang ts->ihessianproductctx = ctx; 415b5b4017aSHong Zhang if (ihp1) ts->vecs_fuu = ihp1; 416b5b4017aSHong Zhang if (ihp2) ts->vecs_fup = ihp2; 417b5b4017aSHong Zhang if (ihp3) ts->vecs_fpu = ihp3; 418b5b4017aSHong Zhang if (ihp4) ts->vecs_fpp = ihp4; 419b5b4017aSHong Zhang ts->ihessianproduct_fuu = ihessianproductfunc1; 420b5b4017aSHong Zhang ts->ihessianproduct_fup = ihessianproductfunc2; 421b5b4017aSHong Zhang ts->ihessianproduct_fpu = ihessianproductfunc3; 422b5b4017aSHong Zhang ts->ihessianproduct_fpp = ihessianproductfunc4; 423b5b4017aSHong Zhang PetscFunctionReturn(0); 424b5b4017aSHong Zhang } 425b5b4017aSHong Zhang 426b5b4017aSHong Zhang /*@C 427b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu. 428b5b4017aSHong Zhang 429b5b4017aSHong Zhang Collective on TS 430b5b4017aSHong Zhang 431b5b4017aSHong Zhang Input Parameters: 432b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 433b5b4017aSHong Zhang 434b5b4017aSHong Zhang Notes: 435b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation, 436b5b4017aSHong Zhang so most users would not generally call this routine themselves. 437b5b4017aSHong Zhang 438b5b4017aSHong Zhang Level: developer 439b5b4017aSHong Zhang 440db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 441b5b4017aSHong Zhang @*/ 442*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 443*d71ae5a4SJacob Faibussowitsch { 444b5b4017aSHong Zhang PetscFunctionBegin; 44533f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 446b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 447b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 448b5b4017aSHong Zhang 449792fecdfSBarry Smith if (ts->ihessianproduct_fuu) PetscCallBack("TS callback IHessianProduct 1 for sensitivity analysis", (*ts->ihessianproduct_fuu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 450ef1023bdSBarry Smith 45167633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 45267633408SHong Zhang if (ts->rhshessianproduct_guu) { 45367633408SHong Zhang PetscInt nadj; 4549566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(ts, t, U, Vl, Vr, VHV)); 45548a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 45667633408SHong Zhang } 457b5b4017aSHong Zhang PetscFunctionReturn(0); 458b5b4017aSHong Zhang } 459b5b4017aSHong Zhang 460b5b4017aSHong Zhang /*@C 461b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup. 462b5b4017aSHong Zhang 463b5b4017aSHong Zhang Collective on TS 464b5b4017aSHong Zhang 465b5b4017aSHong Zhang Input Parameters: 466b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 467b5b4017aSHong Zhang 468b5b4017aSHong Zhang Notes: 469b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation, 470b5b4017aSHong Zhang so most users would not generally call this routine themselves. 471b5b4017aSHong Zhang 472b5b4017aSHong Zhang Level: developer 473b5b4017aSHong Zhang 474db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 475b5b4017aSHong Zhang @*/ 476*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 477*d71ae5a4SJacob Faibussowitsch { 478b5b4017aSHong Zhang PetscFunctionBegin; 47933f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 480b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 481b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 482b5b4017aSHong Zhang 483792fecdfSBarry Smith if (ts->ihessianproduct_fup) PetscCallBack("TS callback IHessianProduct 2 for sensitivity analysis", (*ts->ihessianproduct_fup)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 484ef1023bdSBarry Smith 48567633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 48667633408SHong Zhang if (ts->rhshessianproduct_gup) { 48767633408SHong Zhang PetscInt nadj; 4889566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(ts, t, U, Vl, Vr, VHV)); 48948a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 49067633408SHong Zhang } 491b5b4017aSHong Zhang PetscFunctionReturn(0); 492b5b4017aSHong Zhang } 493b5b4017aSHong Zhang 494b5b4017aSHong Zhang /*@C 495b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu. 496b5b4017aSHong Zhang 497b5b4017aSHong Zhang Collective on TS 498b5b4017aSHong Zhang 499b5b4017aSHong Zhang Input Parameters: 500b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 501b5b4017aSHong Zhang 502b5b4017aSHong Zhang Notes: 503b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation, 504b5b4017aSHong Zhang so most users would not generally call this routine themselves. 505b5b4017aSHong Zhang 506b5b4017aSHong Zhang Level: developer 507b5b4017aSHong Zhang 508db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 509b5b4017aSHong Zhang @*/ 510*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 511*d71ae5a4SJacob Faibussowitsch { 512b5b4017aSHong Zhang PetscFunctionBegin; 51333f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 514b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 515b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 516b5b4017aSHong Zhang 517792fecdfSBarry Smith if (ts->ihessianproduct_fpu) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 518ef1023bdSBarry Smith 51967633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 52067633408SHong Zhang if (ts->rhshessianproduct_gpu) { 52167633408SHong Zhang PetscInt nadj; 5229566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(ts, t, U, Vl, Vr, VHV)); 52348a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 52467633408SHong Zhang } 525b5b4017aSHong Zhang PetscFunctionReturn(0); 526b5b4017aSHong Zhang } 527b5b4017aSHong Zhang 528b5b4017aSHong Zhang /*@C 529b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp. 530b5b4017aSHong Zhang 531b5b4017aSHong Zhang Collective on TS 532b5b4017aSHong Zhang 533b5b4017aSHong Zhang Input Parameters: 534b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 535b5b4017aSHong Zhang 536b5b4017aSHong Zhang Notes: 537b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation, 538b5b4017aSHong Zhang so most users would not generally call this routine themselves. 539b5b4017aSHong Zhang 540b5b4017aSHong Zhang Level: developer 541b5b4017aSHong Zhang 542db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 543b5b4017aSHong Zhang @*/ 544*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 545*d71ae5a4SJacob Faibussowitsch { 546b5b4017aSHong Zhang PetscFunctionBegin; 54733f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 548b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 549b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 550b5b4017aSHong Zhang 551792fecdfSBarry Smith if (ts->ihessianproduct_fpp) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpp)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 552ef1023bdSBarry Smith 55367633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 55467633408SHong Zhang if (ts->rhshessianproduct_gpp) { 55567633408SHong Zhang PetscInt nadj; 5569566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(ts, t, U, Vl, Vr, VHV)); 55748a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 55867633408SHong Zhang } 559b5b4017aSHong Zhang PetscFunctionReturn(0); 560b5b4017aSHong Zhang } 561b5b4017aSHong Zhang 56213af1a74SHong Zhang /*@C 5634c500f23SPierre 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. 56413af1a74SHong Zhang 56513af1a74SHong Zhang Logically Collective on TS 56613af1a74SHong Zhang 56713af1a74SHong Zhang Input Parameters: 56813af1a74SHong Zhang + ts - TS context obtained from TSCreate() 56913af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 57013af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 57113af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 57213af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 57313af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 57413af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 57513af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 576f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for G_PP 57713af1a74SHong Zhang 57813af1a74SHong Zhang Calling sequence of ihessianproductfunc: 579369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 58013af1a74SHong Zhang + t - current timestep 58113af1a74SHong Zhang . U - input vector (current ODE solution) 582369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 58313af1a74SHong Zhang . Vr - input vector to be right-multiplied with the Hessian 584369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 58513af1a74SHong Zhang - ctx - [optional] user-defined function context 58613af1a74SHong Zhang 58713af1a74SHong Zhang Level: intermediate 58813af1a74SHong Zhang 589369cf35fSHong Zhang Notes: 590369cf35fSHong Zhang The first Hessian function and the working array are required. 591369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 592369cf35fSHong Zhang $ Vl_n^T*G_UP*Vr 593369cf35fSHong 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. 594369cf35fSHong Zhang Each entry of G_UP corresponds to the derivative 595369cf35fSHong Zhang $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}. 596369cf35fSHong 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 597369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]} 598369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 59913af1a74SHong Zhang 60013af1a74SHong Zhang .seealso: 60113af1a74SHong Zhang @*/ 602*d71ae5a4SJacob Faibussowitsch 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) 603*d71ae5a4SJacob Faibussowitsch { 60413af1a74SHong Zhang PetscFunctionBegin; 60513af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 60613af1a74SHong Zhang PetscValidPointer(rhshp1, 2); 60713af1a74SHong Zhang 60813af1a74SHong Zhang ts->rhshessianproductctx = ctx; 60913af1a74SHong Zhang if (rhshp1) ts->vecs_guu = rhshp1; 61013af1a74SHong Zhang if (rhshp2) ts->vecs_gup = rhshp2; 61113af1a74SHong Zhang if (rhshp3) ts->vecs_gpu = rhshp3; 61213af1a74SHong Zhang if (rhshp4) ts->vecs_gpp = rhshp4; 61313af1a74SHong Zhang ts->rhshessianproduct_guu = rhshessianproductfunc1; 61413af1a74SHong Zhang ts->rhshessianproduct_gup = rhshessianproductfunc2; 61513af1a74SHong Zhang ts->rhshessianproduct_gpu = rhshessianproductfunc3; 61613af1a74SHong Zhang ts->rhshessianproduct_gpp = rhshessianproductfunc4; 61713af1a74SHong Zhang PetscFunctionReturn(0); 61813af1a74SHong Zhang } 61913af1a74SHong Zhang 62013af1a74SHong Zhang /*@C 621b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu. 62213af1a74SHong Zhang 62313af1a74SHong Zhang Collective on TS 62413af1a74SHong Zhang 62513af1a74SHong Zhang Input Parameters: 62613af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 62713af1a74SHong Zhang 62813af1a74SHong Zhang Notes: 629b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation, 63013af1a74SHong Zhang so most users would not generally call this routine themselves. 63113af1a74SHong Zhang 63213af1a74SHong Zhang Level: developer 63313af1a74SHong Zhang 634db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 63513af1a74SHong Zhang @*/ 636*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 637*d71ae5a4SJacob Faibussowitsch { 63813af1a74SHong Zhang PetscFunctionBegin; 63913af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 64013af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 64113af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 64213af1a74SHong Zhang 643792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 1 for sensitivity analysis", (*ts->rhshessianproduct_guu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 64413af1a74SHong Zhang PetscFunctionReturn(0); 64513af1a74SHong Zhang } 64613af1a74SHong Zhang 64713af1a74SHong Zhang /*@C 648b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup. 64913af1a74SHong Zhang 65013af1a74SHong Zhang Collective on TS 65113af1a74SHong Zhang 65213af1a74SHong Zhang Input Parameters: 65313af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 65413af1a74SHong Zhang 65513af1a74SHong Zhang Notes: 656b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation, 65713af1a74SHong Zhang so most users would not generally call this routine themselves. 65813af1a74SHong Zhang 65913af1a74SHong Zhang Level: developer 66013af1a74SHong Zhang 661db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 66213af1a74SHong Zhang @*/ 663*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 664*d71ae5a4SJacob Faibussowitsch { 66513af1a74SHong Zhang PetscFunctionBegin; 66613af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 66713af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 66813af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 66913af1a74SHong Zhang 670792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 2 for sensitivity analysis", (*ts->rhshessianproduct_gup)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 67113af1a74SHong Zhang PetscFunctionReturn(0); 67213af1a74SHong Zhang } 67313af1a74SHong Zhang 67413af1a74SHong Zhang /*@C 675b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu. 67613af1a74SHong Zhang 67713af1a74SHong Zhang Collective on TS 67813af1a74SHong Zhang 67913af1a74SHong Zhang Input Parameters: 68013af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 68113af1a74SHong Zhang 68213af1a74SHong Zhang Notes: 683b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation, 68413af1a74SHong Zhang so most users would not generally call this routine themselves. 68513af1a74SHong Zhang 68613af1a74SHong Zhang Level: developer 68713af1a74SHong Zhang 688db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 68913af1a74SHong Zhang @*/ 690*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 691*d71ae5a4SJacob Faibussowitsch { 69213af1a74SHong Zhang PetscFunctionBegin; 69313af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 69413af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 69513af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 69613af1a74SHong Zhang 697792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 69813af1a74SHong Zhang PetscFunctionReturn(0); 69913af1a74SHong Zhang } 70013af1a74SHong Zhang 70113af1a74SHong Zhang /*@C 702b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp. 70313af1a74SHong Zhang 70413af1a74SHong Zhang Collective on TS 70513af1a74SHong Zhang 70613af1a74SHong Zhang Input Parameters: 70713af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 70813af1a74SHong Zhang 70913af1a74SHong Zhang Notes: 710b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation, 71113af1a74SHong Zhang so most users would not generally call this routine themselves. 71213af1a74SHong Zhang 71313af1a74SHong Zhang Level: developer 71413af1a74SHong Zhang 715db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 71613af1a74SHong Zhang @*/ 717*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 718*d71ae5a4SJacob Faibussowitsch { 71913af1a74SHong Zhang PetscFunctionBegin; 72013af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 72113af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 72213af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 72313af1a74SHong Zhang 724792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpp)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 72513af1a74SHong Zhang PetscFunctionReturn(0); 72613af1a74SHong Zhang } 72713af1a74SHong Zhang 728814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 729814e85d6SHong Zhang 730814e85d6SHong Zhang /*@ 731814e85d6SHong 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 732814e85d6SHong Zhang for use by the TSAdjoint routines. 733814e85d6SHong Zhang 734d083f849SBarry Smith Logically Collective on TS 735814e85d6SHong Zhang 736814e85d6SHong Zhang Input Parameters: 737814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 738d2fd7bfcSHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 739814e85d6SHong 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 740814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 741814e85d6SHong Zhang 742814e85d6SHong Zhang Level: beginner 743814e85d6SHong Zhang 744814e85d6SHong Zhang Notes: 745814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 746814e85d6SHong Zhang 747814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 748814e85d6SHong Zhang 749db781477SPatrick Sanan .seealso `TSGetCostGradients()` 750814e85d6SHong Zhang @*/ 751*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostGradients(TS ts, PetscInt numcost, Vec *lambda, Vec *mu) 752*d71ae5a4SJacob Faibussowitsch { 753814e85d6SHong Zhang PetscFunctionBegin; 754814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 755064a246eSJacob Faibussowitsch PetscValidPointer(lambda, 3); 756814e85d6SHong Zhang ts->vecs_sensi = lambda; 757814e85d6SHong Zhang ts->vecs_sensip = mu; 7583c633725SBarry 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"); 759814e85d6SHong Zhang ts->numcost = numcost; 760814e85d6SHong Zhang PetscFunctionReturn(0); 761814e85d6SHong Zhang } 762814e85d6SHong Zhang 763814e85d6SHong Zhang /*@ 764814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 765814e85d6SHong Zhang 766814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 767814e85d6SHong Zhang 768814e85d6SHong Zhang Input Parameter: 769814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 770814e85d6SHong Zhang 771d8d19677SJose E. Roman Output Parameters: 7726b867d5aSJose E. Roman + numcost - size of returned arrays 7736b867d5aSJose E. Roman . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 774814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 775814e85d6SHong Zhang 776814e85d6SHong Zhang Level: intermediate 777814e85d6SHong Zhang 778db781477SPatrick Sanan .seealso: `TSSetCostGradients()` 779814e85d6SHong Zhang @*/ 780*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostGradients(TS ts, PetscInt *numcost, Vec **lambda, Vec **mu) 781*d71ae5a4SJacob Faibussowitsch { 782814e85d6SHong Zhang PetscFunctionBegin; 783814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 784814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 785814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 786814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 787814e85d6SHong Zhang PetscFunctionReturn(0); 788814e85d6SHong Zhang } 789814e85d6SHong Zhang 790814e85d6SHong Zhang /*@ 791b5b4017aSHong 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 792b5b4017aSHong Zhang for use by the TSAdjoint routines. 793b5b4017aSHong Zhang 794d083f849SBarry Smith Logically Collective on TS 795b5b4017aSHong Zhang 796b5b4017aSHong Zhang Input Parameters: 797b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 798b5b4017aSHong Zhang . numcost - number of cost functions 799b5b4017aSHong 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 800b5b4017aSHong 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 801b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 802b5b4017aSHong Zhang 803b5b4017aSHong Zhang Level: beginner 804b5b4017aSHong Zhang 805b5b4017aSHong Zhang Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system 806b5b4017aSHong Zhang 807ecf68647SHong Zhang For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve(). 808b5b4017aSHong Zhang 809b5b4017aSHong 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. 810b5b4017aSHong Zhang 8113fca9621SHong Zhang Passing NULL for lambda2 disables the second-order calculation. 812db781477SPatrick Sanan .seealso: `TSAdjointSetForward()` 813b5b4017aSHong Zhang @*/ 814*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostHessianProducts(TS ts, PetscInt numcost, Vec *lambda2, Vec *mu2, Vec dir) 815*d71ae5a4SJacob Faibussowitsch { 816b5b4017aSHong Zhang PetscFunctionBegin; 817b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 8183c633725SBarry 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"); 819b5b4017aSHong Zhang ts->numcost = numcost; 820b5b4017aSHong Zhang ts->vecs_sensi2 = lambda2; 82133f72c5dSHong Zhang ts->vecs_sensi2p = mu2; 822b5b4017aSHong Zhang ts->vec_dir = dir; 823b5b4017aSHong Zhang PetscFunctionReturn(0); 824b5b4017aSHong Zhang } 825b5b4017aSHong Zhang 826b5b4017aSHong Zhang /*@ 827b5b4017aSHong Zhang TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve() 828b5b4017aSHong Zhang 829b5b4017aSHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 830b5b4017aSHong Zhang 831b5b4017aSHong Zhang Input Parameter: 832b5b4017aSHong Zhang . ts - the TS context obtained from TSCreate() 833b5b4017aSHong Zhang 834d8d19677SJose E. Roman Output Parameters: 835b5b4017aSHong Zhang + numcost - number of cost functions 836b5b4017aSHong 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 837b5b4017aSHong 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 838b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 839b5b4017aSHong Zhang 840b5b4017aSHong Zhang Level: intermediate 841b5b4017aSHong Zhang 842db781477SPatrick Sanan .seealso: `TSSetCostHessianProducts()` 843b5b4017aSHong Zhang @*/ 844*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostHessianProducts(TS ts, PetscInt *numcost, Vec **lambda2, Vec **mu2, Vec *dir) 845*d71ae5a4SJacob Faibussowitsch { 846b5b4017aSHong Zhang PetscFunctionBegin; 847b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 848b5b4017aSHong Zhang if (numcost) *numcost = ts->numcost; 849b5b4017aSHong Zhang if (lambda2) *lambda2 = ts->vecs_sensi2; 85033f72c5dSHong Zhang if (mu2) *mu2 = ts->vecs_sensi2p; 851b5b4017aSHong Zhang if (dir) *dir = ts->vec_dir; 852b5b4017aSHong Zhang PetscFunctionReturn(0); 853b5b4017aSHong Zhang } 854b5b4017aSHong Zhang 855b5b4017aSHong Zhang /*@ 856ecf68647SHong Zhang TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities 857b5b4017aSHong Zhang 858d083f849SBarry Smith Logically Collective on TS 859b5b4017aSHong Zhang 860b5b4017aSHong Zhang Input Parameters: 861b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 862b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters 863b5b4017aSHong Zhang 864b5b4017aSHong Zhang Level: intermediate 865b5b4017aSHong Zhang 866ecf68647SHong 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. 867b5b4017aSHong Zhang 868db781477SPatrick Sanan .seealso: `TSSetCostHessianProducts()`, `TSAdjointResetForward()` 869b5b4017aSHong Zhang @*/ 870*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetForward(TS ts, Mat didp) 871*d71ae5a4SJacob Faibussowitsch { 87233f72c5dSHong Zhang Mat A; 87333f72c5dSHong Zhang Vec sp; 87433f72c5dSHong Zhang PetscScalar *xarr; 87533f72c5dSHong Zhang PetscInt lsize; 876b5b4017aSHong Zhang 877b5b4017aSHong Zhang PetscFunctionBegin; 878b5b4017aSHong Zhang ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 8793c633725SBarry Smith PetscCheck(ts->vecs_sensi2, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetCostHessianProducts() first"); 8803c633725SBarry Smith PetscCheck(ts->vec_dir, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 88133f72c5dSHong Zhang /* create a single-column dense matrix */ 8829566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol, &lsize)); 8839566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ts), lsize, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, &A)); 88433f72c5dSHong Zhang 8859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &sp)); 8869566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(A, 0, &xarr)); 8879566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(sp, xarr)); 888ecf68647SHong Zhang if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */ 88933f72c5dSHong Zhang if (didp) { 8909566063dSJacob Faibussowitsch PetscCall(MatMult(didp, ts->vec_dir, sp)); 8919566063dSJacob Faibussowitsch PetscCall(VecScale(sp, 2.)); 89233f72c5dSHong Zhang } else { 8939566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(sp)); 89433f72c5dSHong Zhang } 895ecf68647SHong Zhang } else { /* tangent linear variable initialized as dir */ 8969566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dir, sp)); 89733f72c5dSHong Zhang } 8989566063dSJacob Faibussowitsch PetscCall(VecResetArray(sp)); 8999566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(A, &xarr)); 9009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sp)); 90133f72c5dSHong Zhang 9029566063dSJacob Faibussowitsch PetscCall(TSForwardSetInitialSensitivities(ts, A)); /* if didp is NULL, identity matrix is assumed */ 90333f72c5dSHong Zhang 9049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 905b5b4017aSHong Zhang PetscFunctionReturn(0); 906b5b4017aSHong Zhang } 907b5b4017aSHong Zhang 908b5b4017aSHong Zhang /*@ 909ecf68647SHong Zhang TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context 910ecf68647SHong Zhang 911d083f849SBarry Smith Logically Collective on TS 912ecf68647SHong Zhang 913ecf68647SHong Zhang Input Parameters: 914ecf68647SHong Zhang . ts - the TS context obtained from TSCreate() 915ecf68647SHong Zhang 916ecf68647SHong Zhang Level: intermediate 917ecf68647SHong Zhang 918db781477SPatrick Sanan .seealso: `TSAdjointSetForward()` 919ecf68647SHong Zhang @*/ 920*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointResetForward(TS ts) 921*d71ae5a4SJacob Faibussowitsch { 922ecf68647SHong Zhang PetscFunctionBegin; 923ecf68647SHong Zhang ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */ 9249566063dSJacob Faibussowitsch PetscCall(TSForwardReset(ts)); 925ecf68647SHong Zhang PetscFunctionReturn(0); 926ecf68647SHong Zhang } 927ecf68647SHong Zhang 928ecf68647SHong Zhang /*@ 929814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 930814e85d6SHong Zhang of an adjoint solver 931814e85d6SHong Zhang 932814e85d6SHong Zhang Collective on TS 933814e85d6SHong Zhang 934814e85d6SHong Zhang Input Parameter: 935814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 936814e85d6SHong Zhang 937814e85d6SHong Zhang Level: advanced 938814e85d6SHong Zhang 939db781477SPatrick Sanan .seealso: `TSCreate()`, `TSAdjointStep()`, `TSSetCostGradients()` 940814e85d6SHong Zhang @*/ 941*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetUp(TS ts) 942*d71ae5a4SJacob Faibussowitsch { 943881c1a9bSHong Zhang TSTrajectory tj; 944881c1a9bSHong Zhang PetscBool match; 945814e85d6SHong Zhang 946814e85d6SHong Zhang PetscFunctionBegin; 947814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 948814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 9493c633725SBarry Smith PetscCheck(ts->vecs_sensi, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetCostGradients() first"); 9503c633725SBarry Smith PetscCheck(!ts->vecs_sensip || ts->Jacp || ts->Jacprhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 9519566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj)); 9529566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tj, TSTRAJECTORYBASIC, &match)); 953881c1a9bSHong Zhang if (match) { 954881c1a9bSHong Zhang PetscBool solution_only; 9559566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetSolutionOnly(tj, &solution_only)); 9563c633725SBarry 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"); 957881c1a9bSHong Zhang } 9589566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUseHistory(tj, PETSC_FALSE)); /* not use TSHistory */ 95933f72c5dSHong Zhang 960cd4cee2dSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 9619566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vecs_sensi[0], &ts->vec_drdu_col)); 96248a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &ts->vec_drdp_col)); 963814e85d6SHong Zhang } 964814e85d6SHong Zhang 965dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, adjointsetup); 966814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 967814e85d6SHong Zhang PetscFunctionReturn(0); 968814e85d6SHong Zhang } 969814e85d6SHong Zhang 970814e85d6SHong Zhang /*@ 971ece46509SHong Zhang TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 972ece46509SHong Zhang 973ece46509SHong Zhang Collective on TS 974ece46509SHong Zhang 975ece46509SHong Zhang Input Parameter: 976ece46509SHong Zhang . ts - the TS context obtained from TSCreate() 977ece46509SHong Zhang 978ece46509SHong Zhang Level: beginner 979ece46509SHong Zhang 980db781477SPatrick Sanan .seealso: `TSCreate()`, `TSAdjointSetUp()`, `TSADestroy()` 981ece46509SHong Zhang @*/ 982*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointReset(TS ts) 983*d71ae5a4SJacob Faibussowitsch { 984ece46509SHong Zhang PetscFunctionBegin; 985ece46509SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 986dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, adjointreset); 9877207e0fdSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 9889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_drdu_col)); 98948a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDestroy(&ts->vec_drdp_col)); 9907207e0fdSHong Zhang } 991ece46509SHong Zhang ts->vecs_sensi = NULL; 992ece46509SHong Zhang ts->vecs_sensip = NULL; 993ece46509SHong Zhang ts->vecs_sensi2 = NULL; 99433f72c5dSHong Zhang ts->vecs_sensi2p = NULL; 995ece46509SHong Zhang ts->vec_dir = NULL; 996ece46509SHong Zhang ts->adjointsetupcalled = PETSC_FALSE; 997ece46509SHong Zhang PetscFunctionReturn(0); 998ece46509SHong Zhang } 999ece46509SHong Zhang 1000ece46509SHong Zhang /*@ 1001814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1002814e85d6SHong Zhang 1003814e85d6SHong Zhang Logically Collective on TS 1004814e85d6SHong Zhang 1005814e85d6SHong Zhang Input Parameters: 1006814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1007a2b725a8SWilliam Gropp - steps - number of steps to use 1008814e85d6SHong Zhang 1009814e85d6SHong Zhang Level: intermediate 1010814e85d6SHong Zhang 1011814e85d6SHong Zhang Notes: 1012814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 1013814e85d6SHong Zhang so as to integrate back to less than the original timestep 1014814e85d6SHong Zhang 1015db781477SPatrick Sanan .seealso: `TSSetExactFinalTime()` 1016814e85d6SHong Zhang @*/ 1017*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetSteps(TS ts, PetscInt steps) 1018*d71ae5a4SJacob Faibussowitsch { 1019814e85d6SHong Zhang PetscFunctionBegin; 1020814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1021814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts, steps, 2); 10223c633725SBarry Smith PetscCheck(steps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back a negative number of steps"); 10233c633725SBarry Smith PetscCheck(steps <= ts->steps, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back more than the total number of forward steps"); 1024814e85d6SHong Zhang ts->adjoint_max_steps = steps; 1025814e85d6SHong Zhang PetscFunctionReturn(0); 1026814e85d6SHong Zhang } 1027814e85d6SHong Zhang 1028814e85d6SHong Zhang /*@C 1029814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1030814e85d6SHong Zhang 1031814e85d6SHong Zhang Level: deprecated 1032814e85d6SHong Zhang 1033814e85d6SHong Zhang @*/ 1034*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetRHSJacobian(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) 1035*d71ae5a4SJacob Faibussowitsch { 1036814e85d6SHong Zhang PetscFunctionBegin; 1037814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1038814e85d6SHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 1039814e85d6SHong Zhang 1040814e85d6SHong Zhang ts->rhsjacobianp = func; 1041814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 1042814e85d6SHong Zhang if (Amat) { 10439566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 10449566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 1045814e85d6SHong Zhang ts->Jacp = Amat; 1046814e85d6SHong Zhang } 1047814e85d6SHong Zhang PetscFunctionReturn(0); 1048814e85d6SHong Zhang } 1049814e85d6SHong Zhang 1050814e85d6SHong Zhang /*@C 1051814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1052814e85d6SHong Zhang 1053814e85d6SHong Zhang Level: deprecated 1054814e85d6SHong Zhang 1055814e85d6SHong Zhang @*/ 1056*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeRHSJacobian(TS ts, PetscReal t, Vec U, Mat Amat) 1057*d71ae5a4SJacob Faibussowitsch { 1058814e85d6SHong Zhang PetscFunctionBegin; 1059814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1060c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1061814e85d6SHong Zhang PetscValidPointer(Amat, 4); 1062814e85d6SHong Zhang 1063792fecdfSBarry Smith PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 1064814e85d6SHong Zhang PetscFunctionReturn(0); 1065814e85d6SHong Zhang } 1066814e85d6SHong Zhang 1067814e85d6SHong Zhang /*@ 1068d76be551SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 1069814e85d6SHong Zhang 1070814e85d6SHong Zhang Level: deprecated 1071814e85d6SHong Zhang 1072814e85d6SHong Zhang @*/ 1073*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDYFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) 1074*d71ae5a4SJacob Faibussowitsch { 1075814e85d6SHong Zhang PetscFunctionBegin; 1076814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1077c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1078814e85d6SHong Zhang 1079792fecdfSBarry Smith PetscCallBack("TS callback DRDY for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 1080814e85d6SHong Zhang PetscFunctionReturn(0); 1081814e85d6SHong Zhang } 1082814e85d6SHong Zhang 1083814e85d6SHong Zhang /*@ 1084d76be551SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 1085814e85d6SHong Zhang 1086814e85d6SHong Zhang Level: deprecated 1087814e85d6SHong Zhang 1088814e85d6SHong Zhang @*/ 1089*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) 1090*d71ae5a4SJacob Faibussowitsch { 1091814e85d6SHong Zhang PetscFunctionBegin; 1092814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1093c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1094814e85d6SHong Zhang 1095792fecdfSBarry Smith PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 1096814e85d6SHong Zhang PetscFunctionReturn(0); 1097814e85d6SHong Zhang } 1098814e85d6SHong Zhang 1099814e85d6SHong Zhang /*@C 1100814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 1101814e85d6SHong Zhang 1102814e85d6SHong Zhang Level: intermediate 1103814e85d6SHong Zhang 1104db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1105814e85d6SHong Zhang @*/ 1106*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorSensi(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) 1107*d71ae5a4SJacob Faibussowitsch { 1108814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1109814e85d6SHong Zhang 1110814e85d6SHong Zhang PetscFunctionBegin; 1111064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 11129566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 11139566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], viewer)); 11149566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1115814e85d6SHong Zhang PetscFunctionReturn(0); 1116814e85d6SHong Zhang } 1117814e85d6SHong Zhang 1118814e85d6SHong Zhang /*@C 1119814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1120814e85d6SHong Zhang 1121814e85d6SHong Zhang Collective on TS 1122814e85d6SHong Zhang 1123814e85d6SHong Zhang Input Parameters: 1124814e85d6SHong Zhang + ts - TS object you wish to monitor 1125814e85d6SHong Zhang . name - the monitor type one is seeking 1126814e85d6SHong Zhang . help - message indicating what monitoring is done 1127814e85d6SHong Zhang . manual - manual page for the monitor 1128814e85d6SHong Zhang . monitor - the monitor function 1129814e85d6SHong 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 1130814e85d6SHong Zhang 1131814e85d6SHong Zhang Level: developer 1132814e85d6SHong Zhang 1133db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 1134db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 1135db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1136db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1137c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1138db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1139db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 1140814e85d6SHong Zhang @*/ 1141*d71ae5a4SJacob Faibussowitsch 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 *)) 1142*d71ae5a4SJacob Faibussowitsch { 1143814e85d6SHong Zhang PetscViewer viewer; 1144814e85d6SHong Zhang PetscViewerFormat format; 1145814e85d6SHong Zhang PetscBool flg; 1146814e85d6SHong Zhang 1147814e85d6SHong Zhang PetscFunctionBegin; 11489566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg)); 1149814e85d6SHong Zhang if (flg) { 1150814e85d6SHong Zhang PetscViewerAndFormat *vf; 11519566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 11529566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 11531baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 11549566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 1155814e85d6SHong Zhang } 1156814e85d6SHong Zhang PetscFunctionReturn(0); 1157814e85d6SHong Zhang } 1158814e85d6SHong Zhang 1159814e85d6SHong Zhang /*@C 1160814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1161814e85d6SHong Zhang timestep to display the iteration's progress. 1162814e85d6SHong Zhang 1163814e85d6SHong Zhang Logically Collective on TS 1164814e85d6SHong Zhang 1165814e85d6SHong Zhang Input Parameters: 1166814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1167814e85d6SHong Zhang . adjointmonitor - monitoring routine 1168814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 1169814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 1170814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 1171814e85d6SHong Zhang (may be NULL) 1172814e85d6SHong Zhang 1173814e85d6SHong Zhang Calling sequence of monitor: 1174814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1175814e85d6SHong Zhang 1176814e85d6SHong Zhang + ts - the TS context 1177814e85d6SHong 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 1178814e85d6SHong Zhang been interpolated to) 1179814e85d6SHong Zhang . time - current time 1180814e85d6SHong Zhang . u - current iterate 1181814e85d6SHong Zhang . numcost - number of cost functionos 1182814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1183814e85d6SHong Zhang . mu - sensitivities to parameters 1184814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 1185814e85d6SHong Zhang 1186814e85d6SHong Zhang Notes: 1187814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 1188814e85d6SHong Zhang already has been loaded. 1189814e85d6SHong Zhang 1190814e85d6SHong Zhang Fortran Notes: 1191814e85d6SHong Zhang Only a single monitor function can be set for each TS object 1192814e85d6SHong Zhang 1193814e85d6SHong Zhang Level: intermediate 1194814e85d6SHong Zhang 1195db781477SPatrick Sanan .seealso: `TSAdjointMonitorCancel()` 1196814e85d6SHong Zhang @*/ 1197*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorSet(TS ts, PetscErrorCode (*adjointmonitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *adjointmctx, PetscErrorCode (*adjointmdestroy)(void **)) 1198*d71ae5a4SJacob Faibussowitsch { 1199814e85d6SHong Zhang PetscInt i; 1200814e85d6SHong Zhang PetscBool identical; 1201814e85d6SHong Zhang 1202814e85d6SHong Zhang PetscFunctionBegin; 1203814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1204814e85d6SHong Zhang for (i = 0; i < ts->numbermonitors; i++) { 12059566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))adjointmonitor, adjointmctx, adjointmdestroy, (PetscErrorCode(*)(void))ts->adjointmonitor[i], ts->adjointmonitorcontext[i], ts->adjointmonitordestroy[i], &identical)); 1206814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 1207814e85d6SHong Zhang } 12083c633725SBarry Smith PetscCheck(ts->numberadjointmonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many adjoint monitors set"); 1209814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1210814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1211814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void *)adjointmctx; 1212814e85d6SHong Zhang PetscFunctionReturn(0); 1213814e85d6SHong Zhang } 1214814e85d6SHong Zhang 1215814e85d6SHong Zhang /*@C 1216814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1217814e85d6SHong Zhang 1218814e85d6SHong Zhang Logically Collective on TS 1219814e85d6SHong Zhang 1220814e85d6SHong Zhang Input Parameters: 1221814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1222814e85d6SHong Zhang 1223814e85d6SHong Zhang Notes: 1224814e85d6SHong Zhang There is no way to remove a single, specific monitor. 1225814e85d6SHong Zhang 1226814e85d6SHong Zhang Level: intermediate 1227814e85d6SHong Zhang 1228db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1229814e85d6SHong Zhang @*/ 1230*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorCancel(TS ts) 1231*d71ae5a4SJacob Faibussowitsch { 1232814e85d6SHong Zhang PetscInt i; 1233814e85d6SHong Zhang 1234814e85d6SHong Zhang PetscFunctionBegin; 1235814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1236814e85d6SHong Zhang for (i = 0; i < ts->numberadjointmonitors; i++) { 123748a46eb9SPierre Jolivet if (ts->adjointmonitordestroy[i]) PetscCall((*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i])); 1238814e85d6SHong Zhang } 1239814e85d6SHong Zhang ts->numberadjointmonitors = 0; 1240814e85d6SHong Zhang PetscFunctionReturn(0); 1241814e85d6SHong Zhang } 1242814e85d6SHong Zhang 1243814e85d6SHong Zhang /*@C 1244814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 1245814e85d6SHong Zhang 1246814e85d6SHong Zhang Level: intermediate 1247814e85d6SHong Zhang 1248db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1249814e85d6SHong Zhang @*/ 1250*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) 1251*d71ae5a4SJacob Faibussowitsch { 1252814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1253814e85d6SHong Zhang 1254814e85d6SHong Zhang PetscFunctionBegin; 1255064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 12569566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 12579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 125863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 12599566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 12609566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1261814e85d6SHong Zhang PetscFunctionReturn(0); 1262814e85d6SHong Zhang } 1263814e85d6SHong Zhang 1264814e85d6SHong Zhang /*@C 1265814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1266814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 1267814e85d6SHong Zhang 1268814e85d6SHong Zhang Collective on TS 1269814e85d6SHong Zhang 1270814e85d6SHong Zhang Input Parameters: 1271814e85d6SHong Zhang + ts - the TS context 1272814e85d6SHong Zhang . step - current time-step 1273814e85d6SHong Zhang . ptime - current time 1274814e85d6SHong Zhang . u - current state 1275814e85d6SHong Zhang . numcost - number of cost functions 1276814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1277814e85d6SHong Zhang . mu - sensitivities to parameters 1278814e85d6SHong Zhang - dummy - either a viewer or NULL 1279814e85d6SHong Zhang 1280814e85d6SHong Zhang Level: intermediate 1281814e85d6SHong Zhang 1282db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()`, `TSAdjointMonitorDefault()`, `VecView()` 1283814e85d6SHong Zhang @*/ 1284*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDrawSensi(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *dummy) 1285*d71ae5a4SJacob Faibussowitsch { 1286814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1287814e85d6SHong Zhang PetscDraw draw; 1288814e85d6SHong Zhang PetscReal xl, yl, xr, yr, h; 1289814e85d6SHong Zhang char time[32]; 1290814e85d6SHong Zhang 1291814e85d6SHong Zhang PetscFunctionBegin; 1292814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1293814e85d6SHong Zhang 12949566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], ictx->viewer)); 12959566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 12969566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 12979566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 1298814e85d6SHong Zhang h = yl + .95 * (yr - yl); 12999566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 13009566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1301814e85d6SHong Zhang PetscFunctionReturn(0); 1302814e85d6SHong Zhang } 1303814e85d6SHong Zhang 1304814e85d6SHong Zhang /* 1305814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1306814e85d6SHong Zhang 1307814e85d6SHong Zhang Collective on TSAdjoint 1308814e85d6SHong Zhang 1309814e85d6SHong Zhang Input Parameter: 1310814e85d6SHong Zhang . ts - the TS context 1311814e85d6SHong Zhang 1312814e85d6SHong Zhang Options Database Keys: 1313814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1314814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 1315814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1316814e85d6SHong Zhang 1317814e85d6SHong Zhang Level: developer 1318814e85d6SHong Zhang 1319814e85d6SHong Zhang Notes: 1320814e85d6SHong Zhang This is not normally called directly by users 1321814e85d6SHong Zhang 1322db781477SPatrick Sanan .seealso: `TSSetSaveTrajectory()`, `TSTrajectorySetUp()` 1323814e85d6SHong Zhang */ 1324*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetFromOptions(TS ts, PetscOptionItems *PetscOptionsObject) 1325*d71ae5a4SJacob Faibussowitsch { 1326814e85d6SHong Zhang PetscBool tflg, opt; 1327814e85d6SHong Zhang 1328814e85d6SHong Zhang PetscFunctionBegin; 1329dbbe0bcdSBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1330d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS Adjoint options"); 1331814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 13329566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adjoint_solve", "Solve the adjoint problem immediately after solving the forward problem", "", tflg, &tflg, &opt)); 1333814e85d6SHong Zhang if (opt) { 13349566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 1335814e85d6SHong Zhang ts->adjoint_solve = tflg; 1336814e85d6SHong Zhang } 13379566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor", "Monitor adjoint timestep size", "TSAdjointMonitorDefault", TSAdjointMonitorDefault, NULL)); 13389566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor_sensi", "Monitor sensitivity in the adjoint computation", "TSAdjointMonitorSensi", TSAdjointMonitorSensi, NULL)); 1339814e85d6SHong Zhang opt = PETSC_FALSE; 13409566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", &opt)); 1341814e85d6SHong Zhang if (opt) { 1342814e85d6SHong Zhang TSMonitorDrawCtx ctx; 1343814e85d6SHong Zhang PetscInt howoften = 1; 1344814e85d6SHong Zhang 13459566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", howoften, &howoften, NULL)); 13469566063dSJacob Faibussowitsch PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx)); 13479566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts, TSAdjointMonitorDrawSensi, ctx, (PetscErrorCode(*)(void **))TSMonitorDrawCtxDestroy)); 1348814e85d6SHong Zhang } 1349814e85d6SHong Zhang PetscFunctionReturn(0); 1350814e85d6SHong Zhang } 1351814e85d6SHong Zhang 1352814e85d6SHong Zhang /*@ 1353814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 1354814e85d6SHong Zhang 1355814e85d6SHong Zhang Collective on TS 1356814e85d6SHong Zhang 1357814e85d6SHong Zhang Input Parameter: 1358814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1359814e85d6SHong Zhang 1360814e85d6SHong Zhang Level: intermediate 1361814e85d6SHong Zhang 1362db781477SPatrick Sanan .seealso: `TSAdjointSetUp()`, `TSAdjointSolve()` 1363814e85d6SHong Zhang @*/ 1364*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointStep(TS ts) 1365*d71ae5a4SJacob Faibussowitsch { 1366814e85d6SHong Zhang DM dm; 1367814e85d6SHong Zhang 1368814e85d6SHong Zhang PetscFunctionBegin; 1369814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 13709566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 13719566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 13727b0e2f17SHong Zhang ts->steps--; /* must decrease the step index before the adjoint step is taken. */ 1373814e85d6SHong Zhang 1374814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1375814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 13769566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_AdjointStep, ts, 0, 0, 0)); 1377dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, adjointstep); 13789566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_AdjointStep, ts, 0, 0, 0)); 13797b0e2f17SHong Zhang ts->adjoint_steps++; 1380814e85d6SHong Zhang 1381814e85d6SHong Zhang if (ts->reason < 0) { 13823c633725SBarry Smith PetscCheck(!ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSAdjointStep has failed due to %s", TSConvergedReasons[ts->reason]); 1383814e85d6SHong Zhang } else if (!ts->reason) { 1384814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1385814e85d6SHong Zhang } 1386814e85d6SHong Zhang PetscFunctionReturn(0); 1387814e85d6SHong Zhang } 1388814e85d6SHong Zhang 1389814e85d6SHong Zhang /*@ 1390814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1391814e85d6SHong Zhang 1392814e85d6SHong Zhang Collective on TS 1393814e85d6SHong Zhang 1394814e85d6SHong Zhang Input Parameter: 1395814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1396814e85d6SHong Zhang 1397814e85d6SHong Zhang Options Database: 1398814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1399814e85d6SHong Zhang 1400814e85d6SHong Zhang Level: intermediate 1401814e85d6SHong Zhang 1402814e85d6SHong Zhang Notes: 1403814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 1404814e85d6SHong Zhang 1405814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1406814e85d6SHong Zhang 1407db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetCostGradients()`, `TSSetSolution()`, `TSAdjointStep()` 1408814e85d6SHong Zhang @*/ 1409*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSolve(TS ts) 1410*d71ae5a4SJacob Faibussowitsch { 1411f1d62c27SHong Zhang static PetscBool cite = PETSC_FALSE; 14127f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 14137f59fb53SHong Zhang PetscLogStage adjoint_stage; 14147f59fb53SHong Zhang #endif 1415814e85d6SHong Zhang 1416814e85d6SHong Zhang PetscFunctionBegin; 1417814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1418421b783aSHong Zhang PetscCall(PetscCitationsRegister("@article{Zhang2022tsadjoint,\n" 1419421b783aSHong Zhang " title = {{PETSc TSAdjoint: A Discrete Adjoint ODE Solver for First-Order and Second-Order Sensitivity Analysis}},\n" 1420f1d62c27SHong Zhang " author = {Zhang, Hong and Constantinescu, Emil M. and Smith, Barry F.},\n" 1421421b783aSHong Zhang " journal = {SIAM Journal on Scientific Computing},\n" 1422421b783aSHong Zhang " volume = {44},\n" 1423421b783aSHong Zhang " number = {1},\n" 1424421b783aSHong Zhang " pages = {C1-C24},\n" 1425421b783aSHong Zhang " doi = {10.1137/21M140078X},\n" 14269371c9d4SSatish Balay " year = {2022}\n}\n", 14279371c9d4SSatish Balay &cite)); 14287f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 14299566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("TSAdjoint", &adjoint_stage)); 14309566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(adjoint_stage)); 14317f59fb53SHong Zhang #endif 14329566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 1433814e85d6SHong Zhang 1434814e85d6SHong Zhang /* reset time step and iteration counters */ 1435814e85d6SHong Zhang ts->adjoint_steps = 0; 1436814e85d6SHong Zhang ts->ksp_its = 0; 1437814e85d6SHong Zhang ts->snes_its = 0; 1438814e85d6SHong Zhang ts->num_snes_failures = 0; 1439814e85d6SHong Zhang ts->reject = 0; 1440814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1441814e85d6SHong Zhang 1442814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1443814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1444814e85d6SHong Zhang 1445814e85d6SHong Zhang while (!ts->reason) { 14469566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 14479566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 14489566063dSJacob Faibussowitsch PetscCall(TSAdjointEventHandler(ts)); 14499566063dSJacob Faibussowitsch PetscCall(TSAdjointStep(ts)); 145048a46eb9SPierre Jolivet if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) PetscCall(TSAdjointCostIntegral(ts)); 1451814e85d6SHong Zhang } 1452bdbff044SHong Zhang if (!ts->steps) { 14539566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 14549566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 1455bdbff044SHong Zhang } 1456814e85d6SHong Zhang ts->solvetime = ts->ptime; 14579566063dSJacob Faibussowitsch PetscCall(TSTrajectoryViewFromOptions(ts->trajectory, NULL, "-ts_trajectory_view")); 14589566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ts->vecs_sensi[0], (PetscObject)ts, "-ts_adjoint_view_solution")); 1459814e85d6SHong Zhang ts->adjoint_max_steps = 0; 14607f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 14619566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 14627f59fb53SHong Zhang #endif 1463814e85d6SHong Zhang PetscFunctionReturn(0); 1464814e85d6SHong Zhang } 1465814e85d6SHong Zhang 1466814e85d6SHong Zhang /*@C 1467814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1468814e85d6SHong Zhang 1469814e85d6SHong Zhang Collective on TS 1470814e85d6SHong Zhang 1471814e85d6SHong Zhang Input Parameters: 1472814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 1473814e85d6SHong Zhang . step - step number that has just completed 1474814e85d6SHong Zhang . ptime - model time of the state 1475814e85d6SHong Zhang . u - state at the current model time 1476814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 1477814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1478814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1479814e85d6SHong Zhang 1480814e85d6SHong Zhang Notes: 1481814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1482814e85d6SHong Zhang Users would almost never call this routine directly. 1483814e85d6SHong Zhang 1484814e85d6SHong Zhang Level: developer 1485814e85d6SHong Zhang 1486814e85d6SHong Zhang @*/ 1487*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu) 1488*d71ae5a4SJacob Faibussowitsch { 1489814e85d6SHong Zhang PetscInt i, n = ts->numberadjointmonitors; 1490814e85d6SHong Zhang 1491814e85d6SHong Zhang PetscFunctionBegin; 1492814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1493814e85d6SHong Zhang PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 14949566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(u)); 149548a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall((*ts->adjointmonitor[i])(ts, step, ptime, u, numcost, lambda, mu, ts->adjointmonitorcontext[i])); 14969566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(u)); 1497814e85d6SHong Zhang PetscFunctionReturn(0); 1498814e85d6SHong Zhang } 1499814e85d6SHong Zhang 1500814e85d6SHong Zhang /*@ 1501814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1502814e85d6SHong Zhang 1503814e85d6SHong Zhang Collective on TS 1504814e85d6SHong Zhang 15054165533cSJose E. Roman Input Parameter: 1506814e85d6SHong Zhang . ts - time stepping context 1507814e85d6SHong Zhang 1508814e85d6SHong Zhang Level: advanced 1509814e85d6SHong Zhang 1510814e85d6SHong Zhang Notes: 1511814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 1512814e85d6SHong Zhang 1513db781477SPatrick Sanan .seealso: `TSAdjointSolve()`, `TSAdjointStep` 1514814e85d6SHong Zhang @*/ 1515*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointCostIntegral(TS ts) 1516*d71ae5a4SJacob Faibussowitsch { 1517362febeeSStefano Zampini PetscFunctionBegin; 1518814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1519dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, adjointintegral); 1520814e85d6SHong Zhang PetscFunctionReturn(0); 1521814e85d6SHong Zhang } 1522814e85d6SHong Zhang 1523814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1524814e85d6SHong Zhang 1525814e85d6SHong Zhang /*@ 1526814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 1527814e85d6SHong Zhang of forward sensitivity analysis 1528814e85d6SHong Zhang 1529814e85d6SHong Zhang Collective on TS 1530814e85d6SHong Zhang 1531814e85d6SHong Zhang Input Parameter: 1532814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1533814e85d6SHong Zhang 1534814e85d6SHong Zhang Level: advanced 1535814e85d6SHong Zhang 1536db781477SPatrick Sanan .seealso: `TSCreate()`, `TSDestroy()`, `TSSetUp()` 1537814e85d6SHong Zhang @*/ 1538*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetUp(TS ts) 1539*d71ae5a4SJacob Faibussowitsch { 1540814e85d6SHong Zhang PetscFunctionBegin; 1541814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1542814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1543dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, forwardsetup); 15449566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ts->vec_sensip_col)); 1545814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 1546814e85d6SHong Zhang PetscFunctionReturn(0); 1547814e85d6SHong Zhang } 1548814e85d6SHong Zhang 1549814e85d6SHong Zhang /*@ 15507adebddeSHong Zhang TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 15517adebddeSHong Zhang 15527adebddeSHong Zhang Collective on TS 15537adebddeSHong Zhang 15547adebddeSHong Zhang Input Parameter: 15557adebddeSHong Zhang . ts - the TS context obtained from TSCreate() 15567adebddeSHong Zhang 15577adebddeSHong Zhang Level: advanced 15587adebddeSHong Zhang 1559db781477SPatrick Sanan .seealso: `TSCreate()`, `TSDestroy()`, `TSForwardSetUp()` 15607adebddeSHong Zhang @*/ 1561*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardReset(TS ts) 1562*d71ae5a4SJacob Faibussowitsch { 15637207e0fdSHong Zhang TS quadts = ts->quadraturets; 15647adebddeSHong Zhang 15657adebddeSHong Zhang PetscFunctionBegin; 15667adebddeSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1567dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, forwardreset); 15689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 156948a46eb9SPierre Jolivet if (quadts) PetscCall(MatDestroy(&quadts->mat_sensip)); 15709566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_sensip_col)); 15717adebddeSHong Zhang ts->forward_solve = PETSC_FALSE; 15727adebddeSHong Zhang ts->forwardsetupcalled = PETSC_FALSE; 15737adebddeSHong Zhang PetscFunctionReturn(0); 15747adebddeSHong Zhang } 15757adebddeSHong Zhang 15767adebddeSHong Zhang /*@ 1577814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1578814e85d6SHong Zhang 1579d8d19677SJose E. Roman Input Parameters: 1580a2b725a8SWilliam Gropp + ts - the TS context obtained from TSCreate() 1581814e85d6SHong Zhang . numfwdint - number of integrals 158267b8a455SSatish Balay - vp - the vectors containing the gradients for each integral w.r.t. parameters 1583814e85d6SHong Zhang 15847207e0fdSHong Zhang Level: deprecated 1585814e85d6SHong Zhang 1586db781477SPatrick Sanan .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1587814e85d6SHong Zhang @*/ 1588*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetIntegralGradients(TS ts, PetscInt numfwdint, Vec *vp) 1589*d71ae5a4SJacob Faibussowitsch { 1590814e85d6SHong Zhang PetscFunctionBegin; 1591814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 15923c633725SBarry 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()"); 1593814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1594814e85d6SHong Zhang 1595814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 1596814e85d6SHong Zhang PetscFunctionReturn(0); 1597814e85d6SHong Zhang } 1598814e85d6SHong Zhang 1599814e85d6SHong Zhang /*@ 1600814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1601814e85d6SHong Zhang 1602814e85d6SHong Zhang Input Parameter: 1603814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1604814e85d6SHong Zhang 1605814e85d6SHong Zhang Output Parameter: 160667b8a455SSatish Balay . vp - the vectors containing the gradients for each integral w.r.t. parameters 1607814e85d6SHong Zhang 16087207e0fdSHong Zhang Level: deprecated 1609814e85d6SHong Zhang 1610db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1611814e85d6SHong Zhang @*/ 1612*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetIntegralGradients(TS ts, PetscInt *numfwdint, Vec **vp) 1613*d71ae5a4SJacob Faibussowitsch { 1614814e85d6SHong Zhang PetscFunctionBegin; 1615814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1616814e85d6SHong Zhang PetscValidPointer(vp, 3); 1617814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1618814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1619814e85d6SHong Zhang PetscFunctionReturn(0); 1620814e85d6SHong Zhang } 1621814e85d6SHong Zhang 1622814e85d6SHong Zhang /*@ 1623814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1624814e85d6SHong Zhang 1625814e85d6SHong Zhang Collective on TS 1626814e85d6SHong Zhang 16274165533cSJose E. Roman Input Parameter: 1628814e85d6SHong Zhang . ts - time stepping context 1629814e85d6SHong Zhang 1630814e85d6SHong Zhang Level: advanced 1631814e85d6SHong Zhang 1632814e85d6SHong Zhang Notes: 1633814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1634814e85d6SHong Zhang 1635db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardSetUp()` 1636814e85d6SHong Zhang @*/ 1637*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardStep(TS ts) 1638*d71ae5a4SJacob Faibussowitsch { 1639362febeeSStefano Zampini PetscFunctionBegin; 1640814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 16419566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_ForwardStep, ts, 0, 0, 0)); 1642dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, forwardstep); 16439566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_ForwardStep, ts, 0, 0, 0)); 16443c633725SBarry Smith PetscCheck(ts->reason >= 0 || !ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSFowardStep has failed due to %s", TSConvergedReasons[ts->reason]); 1645814e85d6SHong Zhang PetscFunctionReturn(0); 1646814e85d6SHong Zhang } 1647814e85d6SHong Zhang 1648814e85d6SHong Zhang /*@ 1649814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1650814e85d6SHong Zhang 1651d083f849SBarry Smith Logically Collective on TS 1652814e85d6SHong Zhang 1653814e85d6SHong Zhang Input Parameters: 1654814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1655814e85d6SHong Zhang . nump - number of parameters 1656814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1657814e85d6SHong Zhang 1658814e85d6SHong Zhang Level: beginner 1659814e85d6SHong Zhang 1660814e85d6SHong Zhang Notes: 1661814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1662814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1663814e85d6SHong Zhang You must call this function before TSSolve(). 1664814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1665814e85d6SHong Zhang 1666db781477SPatrick Sanan .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1667814e85d6SHong Zhang @*/ 1668*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetSensitivities(TS ts, PetscInt nump, Mat Smat) 1669*d71ae5a4SJacob Faibussowitsch { 1670814e85d6SHong Zhang PetscFunctionBegin; 1671814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1672814e85d6SHong Zhang PetscValidHeaderSpecific(Smat, MAT_CLASSID, 3); 1673814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1674b5b4017aSHong Zhang if (nump == PETSC_DEFAULT) { 16759566063dSJacob Faibussowitsch PetscCall(MatGetSize(Smat, NULL, &ts->num_parameters)); 1676b5b4017aSHong Zhang } else ts->num_parameters = nump; 16779566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Smat)); 16789566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 1679814e85d6SHong Zhang ts->mat_sensip = Smat; 1680814e85d6SHong Zhang PetscFunctionReturn(0); 1681814e85d6SHong Zhang } 1682814e85d6SHong Zhang 1683814e85d6SHong Zhang /*@ 1684814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1685814e85d6SHong Zhang 1686814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1687814e85d6SHong Zhang 1688d8d19677SJose E. Roman Output Parameters: 1689814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1690814e85d6SHong Zhang . nump - number of parameters 1691814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1692814e85d6SHong Zhang 1693814e85d6SHong Zhang Level: intermediate 1694814e85d6SHong Zhang 1695db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1696814e85d6SHong Zhang @*/ 1697*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetSensitivities(TS ts, PetscInt *nump, Mat *Smat) 1698*d71ae5a4SJacob Faibussowitsch { 1699814e85d6SHong Zhang PetscFunctionBegin; 1700814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1701814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1702814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1703814e85d6SHong Zhang PetscFunctionReturn(0); 1704814e85d6SHong Zhang } 1705814e85d6SHong Zhang 1706814e85d6SHong Zhang /*@ 1707814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1708814e85d6SHong Zhang 1709814e85d6SHong Zhang Collective on TS 1710814e85d6SHong Zhang 17114165533cSJose E. Roman Input Parameter: 1712814e85d6SHong Zhang . ts - time stepping context 1713814e85d6SHong Zhang 1714814e85d6SHong Zhang Level: advanced 1715814e85d6SHong Zhang 1716814e85d6SHong Zhang Notes: 1717814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1718814e85d6SHong Zhang 1719db781477SPatrick Sanan .seealso: `TSSolve()`, `TSAdjointCostIntegral()` 1720814e85d6SHong Zhang @*/ 1721*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardCostIntegral(TS ts) 1722*d71ae5a4SJacob Faibussowitsch { 1723362febeeSStefano Zampini PetscFunctionBegin; 1724814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1725dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, forwardintegral); 1726814e85d6SHong Zhang PetscFunctionReturn(0); 1727814e85d6SHong Zhang } 1728b5b4017aSHong Zhang 1729b5b4017aSHong Zhang /*@ 1730b5b4017aSHong Zhang TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1731b5b4017aSHong Zhang 1732d083f849SBarry Smith Collective on TS 1733b5b4017aSHong Zhang 1734d8d19677SJose E. Roman Input Parameters: 1735b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 1736b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition 1737b5b4017aSHong Zhang 1738b5b4017aSHong Zhang Level: intermediate 1739b5b4017aSHong Zhang 1740b5b4017aSHong 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. 1741b5b4017aSHong Zhang 1742db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()` 1743b5b4017aSHong Zhang @*/ 1744*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetInitialSensitivities(TS ts, Mat didp) 1745*d71ae5a4SJacob Faibussowitsch { 1746362febeeSStefano Zampini PetscFunctionBegin; 1747b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1748b5b4017aSHong Zhang PetscValidHeaderSpecific(didp, MAT_CLASSID, 2); 174948a46eb9SPierre Jolivet if (!ts->mat_sensip) PetscCall(TSForwardSetSensitivities(ts, PETSC_DEFAULT, didp)); 1750b5b4017aSHong Zhang PetscFunctionReturn(0); 1751b5b4017aSHong Zhang } 17526affb6f8SHong Zhang 17536affb6f8SHong Zhang /*@ 17546affb6f8SHong Zhang TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 17556affb6f8SHong Zhang 17566affb6f8SHong Zhang Input Parameter: 17576affb6f8SHong Zhang . ts - the TS context obtained from TSCreate() 17586affb6f8SHong Zhang 17596affb6f8SHong Zhang Output Parameters: 1760cd4cee2dSHong Zhang + ns - number of stages 17616affb6f8SHong Zhang - S - tangent linear sensitivities at the intermediate stages 17626affb6f8SHong Zhang 17636affb6f8SHong Zhang Level: advanced 17646affb6f8SHong Zhang 17656affb6f8SHong Zhang @*/ 1766*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetStages(TS ts, PetscInt *ns, Mat **S) 1767*d71ae5a4SJacob Faibussowitsch { 17686affb6f8SHong Zhang PetscFunctionBegin; 17696affb6f8SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 17706affb6f8SHong Zhang 17716affb6f8SHong Zhang if (!ts->ops->getstages) *S = NULL; 1772dbbe0bcdSBarry Smith else PetscUseTypeMethod(ts, forwardgetstages, ns, S); 17736affb6f8SHong Zhang PetscFunctionReturn(0); 17746affb6f8SHong Zhang } 1775cd4cee2dSHong Zhang 1776cd4cee2dSHong Zhang /*@ 1777cd4cee2dSHong Zhang TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 1778cd4cee2dSHong Zhang 1779d8d19677SJose E. Roman Input Parameters: 1780cd4cee2dSHong Zhang + ts - the TS context obtained from TSCreate() 1781cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1782cd4cee2dSHong Zhang 1783cd4cee2dSHong Zhang Output Parameters: 1784cd4cee2dSHong Zhang . quadts - the child TS context 1785cd4cee2dSHong Zhang 1786cd4cee2dSHong Zhang Level: intermediate 1787cd4cee2dSHong Zhang 1788db781477SPatrick Sanan .seealso: `TSGetQuadratureTS()` 1789cd4cee2dSHong Zhang @*/ 1790*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSCreateQuadratureTS(TS ts, PetscBool fwd, TS *quadts) 1791*d71ae5a4SJacob Faibussowitsch { 1792cd4cee2dSHong Zhang char prefix[128]; 1793cd4cee2dSHong Zhang 1794cd4cee2dSHong Zhang PetscFunctionBegin; 1795cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1796064a246eSJacob Faibussowitsch PetscValidPointer(quadts, 3); 17979566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts->quadraturets)); 17989566063dSJacob Faibussowitsch PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &ts->quadraturets)); 17999566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets, (PetscObject)ts, 1)); 18009566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%squad_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "")); 18019566063dSJacob Faibussowitsch PetscCall(TSSetOptionsPrefix(ts->quadraturets, prefix)); 1802cd4cee2dSHong Zhang *quadts = ts->quadraturets; 1803cd4cee2dSHong Zhang 1804cd4cee2dSHong Zhang if (ts->numcost) { 18059566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, ts->numcost, &(*quadts)->vec_sol)); 1806cd4cee2dSHong Zhang } else { 18079566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &(*quadts)->vec_sol)); 1808cd4cee2dSHong Zhang } 1809cd4cee2dSHong Zhang ts->costintegralfwd = fwd; 1810cd4cee2dSHong Zhang PetscFunctionReturn(0); 1811cd4cee2dSHong Zhang } 1812cd4cee2dSHong Zhang 1813cd4cee2dSHong Zhang /*@ 1814cd4cee2dSHong Zhang TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 1815cd4cee2dSHong Zhang 1816cd4cee2dSHong Zhang Input Parameter: 1817cd4cee2dSHong Zhang . ts - the TS context obtained from TSCreate() 1818cd4cee2dSHong Zhang 1819cd4cee2dSHong Zhang Output Parameters: 1820cd4cee2dSHong Zhang + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1821cd4cee2dSHong Zhang - quadts - the child TS context 1822cd4cee2dSHong Zhang 1823cd4cee2dSHong Zhang Level: intermediate 1824cd4cee2dSHong Zhang 1825db781477SPatrick Sanan .seealso: `TSCreateQuadratureTS()` 1826cd4cee2dSHong Zhang @*/ 1827*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetQuadratureTS(TS ts, PetscBool *fwd, TS *quadts) 1828*d71ae5a4SJacob Faibussowitsch { 1829cd4cee2dSHong Zhang PetscFunctionBegin; 1830cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1831cd4cee2dSHong Zhang if (fwd) *fwd = ts->costintegralfwd; 1832cd4cee2dSHong Zhang if (quadts) *quadts = ts->quadraturets; 1833cd4cee2dSHong Zhang PetscFunctionReturn(0); 1834cd4cee2dSHong Zhang } 1835ba0675f6SHong Zhang 1836ba0675f6SHong Zhang /*@ 1837ba0675f6SHong Zhang TSComputeSNESJacobian - Compute the SNESJacobian 1838ba0675f6SHong Zhang 1839ba0675f6SHong Zhang Input Parameters: 1840ba0675f6SHong Zhang + ts - the TS context obtained from TSCreate() 1841ba0675f6SHong Zhang - x - state vector 1842ba0675f6SHong Zhang 1843ba0675f6SHong Zhang Output Parameters: 1844ba0675f6SHong Zhang + J - Jacobian matrix 1845ba0675f6SHong Zhang - Jpre - preconditioning matrix for J (may be same as J) 1846ba0675f6SHong Zhang 1847ba0675f6SHong Zhang Level: developer 1848ba0675f6SHong Zhang 1849ba0675f6SHong Zhang Notes: 1850ba0675f6SHong Zhang Using SNES to compute the Jacobian enables finite differencing when TS Jacobian is not available. 1851ba0675f6SHong Zhang @*/ 1852*d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeSNESJacobian(TS ts, Vec x, Mat J, Mat Jpre) 1853*d71ae5a4SJacob Faibussowitsch { 1854ba0675f6SHong Zhang SNES snes = ts->snes; 1855ba0675f6SHong Zhang PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *) = NULL; 1856ba0675f6SHong Zhang 1857ba0675f6SHong Zhang PetscFunctionBegin; 1858ba0675f6SHong Zhang /* 1859ba0675f6SHong Zhang Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object 1860ba0675f6SHong Zhang because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for 1861ba0675f6SHong Zhang explicit methods. Instead, we check the Jacobian compute function directly to determin if FD 1862ba0675f6SHong Zhang coloring is used. 1863ba0675f6SHong Zhang */ 18649566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, NULL, NULL, &jac, NULL)); 1865ba0675f6SHong Zhang if (jac == SNESComputeJacobianDefaultColor) { 1866ba0675f6SHong Zhang Vec f; 18679566063dSJacob Faibussowitsch PetscCall(SNESSetSolution(snes, x)); 18689566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &f, NULL, NULL)); 1869ba0675f6SHong Zhang /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */ 18709566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, x, f)); 1871ba0675f6SHong Zhang } 18729566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, x, J, Jpre)); 1873ba0675f6SHong Zhang PetscFunctionReturn(0); 1874ba0675f6SHong Zhang } 1875