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 13c3339decSBarry Smith Logically Collective 14814e85d6SHong Zhang 15814e85d6SHong Zhang Input Parameters: 16bcf0153eSBarry Smith + 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 2120f4b53cSBarry Smith Calling sequence of `func`: 2220f4b53cSBarry Smith $ PetscErrorCode 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 30bcf0153eSBarry Smith Note: 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 33*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSGetRHSJacobianP()` 34814e85d6SHong Zhang @*/ 35d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetRHSJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) 36d71ae5a4SJacob 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 } 483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 54c3339decSBarry Smith Logically Collective 55cd4cee2dSHong Zhang 56f899ff85SJose E. Roman Input Parameter: 57bcf0153eSBarry Smith . 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 6420f4b53cSBarry Smith Calling sequence of `func`: 6520f4b53cSBarry Smith $ PetscErrorCode 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 73bcf0153eSBarry Smith Note: 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 76*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSJacobianP()`, `TS`, `TSGetRHSJacobianP()` 77cd4cee2dSHong Zhang @*/ 78d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetRHSJacobianP(TS ts, Mat *Amat, PetscErrorCode (**func)(TS, PetscReal, Vec, Mat, void *), void **ctx) 79d71ae5a4SJacob 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; 843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85cd4cee2dSHong Zhang } 86cd4cee2dSHong Zhang 87cd4cee2dSHong Zhang /*@C 88814e85d6SHong Zhang TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 89814e85d6SHong Zhang 90c3339decSBarry Smith Collective 91814e85d6SHong Zhang 92814e85d6SHong Zhang Input Parameters: 932fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 942fe279fdSBarry Smith . t - the time 952fe279fdSBarry Smith - U - the solution at which to compute the Jacobian 962fe279fdSBarry Smith 972fe279fdSBarry Smith Output Parameter: 982fe279fdSBarry Smith . Amat - the computed Jacobian 99814e85d6SHong Zhang 100814e85d6SHong Zhang Level: developer 101814e85d6SHong Zhang 102*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSJacobianP()`, `TS` 103814e85d6SHong Zhang @*/ 104d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSJacobianP(TS ts, PetscReal t, Vec U, Mat Amat) 105d71ae5a4SJacob Faibussowitsch { 106814e85d6SHong Zhang PetscFunctionBegin; 1073ba16761SJacob Faibussowitsch if (!Amat) PetscFunctionReturn(PETSC_SUCCESS); 108814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 109c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 110814e85d6SHong Zhang 11180ab5e31SHong Zhang if (ts->rhsjacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 11280ab5e31SHong Zhang else { 11380ab5e31SHong Zhang PetscBool assembled; 11480ab5e31SHong Zhang PetscCall(MatZeroEntries(Amat)); 11580ab5e31SHong Zhang PetscCall(MatAssembled(Amat, &assembled)); 11680ab5e31SHong Zhang if (!assembled) { 11780ab5e31SHong Zhang PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY)); 11880ab5e31SHong Zhang PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY)); 11980ab5e31SHong Zhang } 12080ab5e31SHong Zhang } 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122814e85d6SHong Zhang } 123814e85d6SHong Zhang 124814e85d6SHong Zhang /*@C 12533f72c5dSHong 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. 12633f72c5dSHong Zhang 127c3339decSBarry Smith Logically Collective 12833f72c5dSHong Zhang 12933f72c5dSHong Zhang Input Parameters: 130bcf0153eSBarry Smith + ts - `TS` context obtained from `TSCreate()` 13133f72c5dSHong Zhang . Amat - JacobianP matrix 13233f72c5dSHong Zhang . func - function 13333f72c5dSHong Zhang - ctx - [optional] user-defined function context 13433f72c5dSHong Zhang 13520f4b53cSBarry Smith Calling sequence of `func`: 13620f4b53cSBarry Smith $ PetscErrorCode func(TS ts, PetscReal t, Vec y, Mat A, void *ctx) 13733f72c5dSHong Zhang + t - current timestep 13833f72c5dSHong Zhang . U - input vector (current ODE solution) 13933f72c5dSHong Zhang . Udot - time derivative of state vector 14033f72c5dSHong Zhang . shift - shift to apply, see note below 14133f72c5dSHong Zhang . A - output matrix 14233f72c5dSHong Zhang - ctx - [optional] user-defined function context 14333f72c5dSHong Zhang 14433f72c5dSHong Zhang Level: intermediate 14533f72c5dSHong Zhang 146bcf0153eSBarry Smith Note: 14733f72c5dSHong 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 14833f72c5dSHong Zhang 149*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSJacobianP()`, `TS` 15033f72c5dSHong Zhang @*/ 151d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetIJacobianP(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Vec, PetscReal, Mat, void *), void *ctx) 152d71ae5a4SJacob Faibussowitsch { 15333f72c5dSHong Zhang PetscFunctionBegin; 15433f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 15533f72c5dSHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 15633f72c5dSHong Zhang 15733f72c5dSHong Zhang ts->ijacobianp = func; 15833f72c5dSHong Zhang ts->ijacobianpctx = ctx; 15933f72c5dSHong Zhang if (Amat) { 1609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 1619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 16233f72c5dSHong Zhang ts->Jacp = Amat; 16333f72c5dSHong Zhang } 1643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16533f72c5dSHong Zhang } 16633f72c5dSHong Zhang 16733f72c5dSHong Zhang /*@C 16833f72c5dSHong Zhang TSComputeIJacobianP - Runs the user-defined IJacobianP function. 16933f72c5dSHong Zhang 170c3339decSBarry Smith Collective 17133f72c5dSHong Zhang 17233f72c5dSHong Zhang Input Parameters: 173bcf0153eSBarry Smith + ts - the `TS` context 17433f72c5dSHong Zhang . t - current timestep 17533f72c5dSHong Zhang . U - state vector 17633f72c5dSHong Zhang . Udot - time derivative of state vector 17733f72c5dSHong Zhang . shift - shift to apply, see note below 17880ab5e31SHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobianP should be kept separate 17933f72c5dSHong Zhang 1802fe279fdSBarry Smith Output Parameter: 18133f72c5dSHong Zhang . A - Jacobian matrix 18233f72c5dSHong Zhang 18333f72c5dSHong Zhang Level: developer 18433f72c5dSHong Zhang 185*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetIJacobianP()` 18633f72c5dSHong Zhang @*/ 187d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIJacobianP(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal shift, Mat Amat, PetscBool imex) 188d71ae5a4SJacob Faibussowitsch { 18933f72c5dSHong Zhang PetscFunctionBegin; 1903ba16761SJacob Faibussowitsch if (!Amat) PetscFunctionReturn(PETSC_SUCCESS); 19133f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 19233f72c5dSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 19333f72c5dSHong Zhang PetscValidHeaderSpecific(Udot, VEC_CLASSID, 4); 19433f72c5dSHong Zhang 1959566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_JacobianPEval, ts, U, Amat, 0)); 19648a46eb9SPierre Jolivet if (ts->ijacobianp) PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->ijacobianp)(ts, t, U, Udot, shift, Amat, ts->ijacobianpctx)); 19733f72c5dSHong Zhang if (imex) { 19833f72c5dSHong Zhang if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */ 19933f72c5dSHong Zhang PetscBool assembled; 2009566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 2019566063dSJacob Faibussowitsch PetscCall(MatAssembled(Amat, &assembled)); 20233f72c5dSHong Zhang if (!assembled) { 2039566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Amat, MAT_FINAL_ASSEMBLY)); 2049566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Amat, MAT_FINAL_ASSEMBLY)); 20533f72c5dSHong Zhang } 20633f72c5dSHong Zhang } 20733f72c5dSHong Zhang } else { 2081baa6e33SBarry Smith if (ts->rhsjacobianp) PetscCall(TSComputeRHSJacobianP(ts, t, U, ts->Jacprhs)); 20933f72c5dSHong Zhang if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */ 2109566063dSJacob Faibussowitsch PetscCall(MatScale(Amat, -1)); 21133f72c5dSHong Zhang } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */ 21233f72c5dSHong Zhang MatStructure axpy = DIFFERENT_NONZERO_PATTERN; 21333f72c5dSHong Zhang if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */ 2149566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 21533f72c5dSHong Zhang } 2169566063dSJacob Faibussowitsch PetscCall(MatAXPY(Amat, -1, ts->Jacprhs, axpy)); 21733f72c5dSHong Zhang } 21833f72c5dSHong Zhang } 2199566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_JacobianPEval, ts, U, Amat, 0)); 2203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 22133f72c5dSHong Zhang } 22233f72c5dSHong Zhang 22333f72c5dSHong Zhang /*@C 224814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 225814e85d6SHong Zhang 226c3339decSBarry Smith Logically Collective 227814e85d6SHong Zhang 228814e85d6SHong Zhang Input Parameters: 229bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 230814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 231814e85d6SHong Zhang . costintegral - vector that stores the integral values 232814e85d6SHong Zhang . rf - routine for evaluating the integrand function 233c9ad9de0SHong Zhang . drduf - function that computes the gradients of the r's with respect to u 2342fe279fdSBarry Smith . 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`) 235814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 2362fe279fdSBarry Smith - ctx - [optional] user-defined context for private data for the function evaluation routine (may be `NULL`) 237814e85d6SHong Zhang 23820f4b53cSBarry Smith Calling sequence of `rf`: 23920f4b53cSBarry Smith $ PetscErrorCode rf(TS ts, PetscReal t, Vec U, Vec F, oid *ctx) 240814e85d6SHong Zhang 24120f4b53cSBarry Smith Calling sequence of `drduf`: 24220f4b53cSBarry Smith $ PetscErroCode drduf(TS ts, PetscReal t, Vec U, Vec *dRdU, void *ctx) 243814e85d6SHong Zhang 24420f4b53cSBarry Smith Calling sequence of `drdpf`: 24520f4b53cSBarry Smith $ PetscErroCode drdpf(TS ts, PetscReal t, Vec U, Vec *dRdP, void *ctx) 246814e85d6SHong Zhang 247cd4cee2dSHong Zhang Level: deprecated 248814e85d6SHong Zhang 249bcf0153eSBarry Smith Note: 250814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 251814e85d6SHong Zhang 252*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSJacobianP()`, `TSGetCostGradients()`, `TSSetCostGradients()` 253814e85d6SHong Zhang @*/ 254d71ae5a4SJacob 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) 255d71ae5a4SJacob Faibussowitsch { 256814e85d6SHong Zhang PetscFunctionBegin; 257814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 258814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral, VEC_CLASSID, 3); 2593c633725SBarry 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()"); 260814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numcost; 261814e85d6SHong Zhang 262814e85d6SHong Zhang if (costintegral) { 2639566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)costintegral)); 2649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_costintegral)); 265814e85d6SHong Zhang ts->vec_costintegral = costintegral; 266814e85d6SHong Zhang } else { 267814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 2689566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, numcost, &ts->vec_costintegral)); 269814e85d6SHong Zhang } else { 2709566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegral, 0.0)); 271814e85d6SHong Zhang } 272814e85d6SHong Zhang } 273814e85d6SHong Zhang if (!ts->vec_costintegrand) { 2749566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_costintegral, &ts->vec_costintegrand)); 275814e85d6SHong Zhang } else { 2769566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegrand, 0.0)); 277814e85d6SHong Zhang } 278814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 279814e85d6SHong Zhang ts->costintegrand = rf; 280814e85d6SHong Zhang ts->costintegrandctx = ctx; 281c9ad9de0SHong Zhang ts->drdufunction = drduf; 282814e85d6SHong Zhang ts->drdpfunction = drdpf; 2833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 284814e85d6SHong Zhang } 285814e85d6SHong Zhang 286b5b4017aSHong Zhang /*@C 287814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 288814e85d6SHong Zhang It is valid to call the routine after a backward run. 289814e85d6SHong Zhang 290814e85d6SHong Zhang Not Collective 291814e85d6SHong Zhang 292814e85d6SHong Zhang Input Parameter: 293bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 294814e85d6SHong Zhang 295814e85d6SHong Zhang Output Parameter: 296814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 297814e85d6SHong Zhang 298814e85d6SHong Zhang Level: intermediate 299814e85d6SHong Zhang 300*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, ``TSSetCostIntegrand()` 301814e85d6SHong Zhang @*/ 302d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostIntegral(TS ts, Vec *v) 303d71ae5a4SJacob Faibussowitsch { 304cd4cee2dSHong Zhang TS quadts; 305cd4cee2dSHong Zhang 306814e85d6SHong Zhang PetscFunctionBegin; 307814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 308814e85d6SHong Zhang PetscValidPointer(v, 2); 3099566063dSJacob Faibussowitsch PetscCall(TSGetQuadratureTS(ts, NULL, &quadts)); 310cd4cee2dSHong Zhang *v = quadts->vec_sol; 3113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 312814e85d6SHong Zhang } 313814e85d6SHong Zhang 314b5b4017aSHong Zhang /*@C 315814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 316814e85d6SHong Zhang 317814e85d6SHong Zhang Input Parameters: 318bcf0153eSBarry Smith + ts - the `TS` context 319814e85d6SHong Zhang . t - current time 320c9ad9de0SHong Zhang - U - state vector, i.e. current solution 321814e85d6SHong Zhang 322814e85d6SHong Zhang Output Parameter: 323c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs 324814e85d6SHong Zhang 325bcf0153eSBarry Smith Level: deprecated 326bcf0153eSBarry Smith 327bcf0153eSBarry Smith Note: 328814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 329814e85d6SHong Zhang is used internally within the sensitivity analysis context. 330814e85d6SHong Zhang 331*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSSetCostIntegrand()` 332814e85d6SHong Zhang @*/ 333d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeCostIntegrand(TS ts, PetscReal t, Vec U, Vec Q) 334d71ae5a4SJacob Faibussowitsch { 335814e85d6SHong Zhang PetscFunctionBegin; 336814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 337c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 338c9ad9de0SHong Zhang PetscValidHeaderSpecific(Q, VEC_CLASSID, 4); 339814e85d6SHong Zhang 3409566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_FunctionEval, ts, U, Q, 0)); 341792fecdfSBarry Smith if (ts->costintegrand) PetscCallBack("TS callback integrand in the cost function", (*ts->costintegrand)(ts, t, U, Q, ts->costintegrandctx)); 342ef1023bdSBarry Smith else PetscCall(VecZeroEntries(Q)); 3439566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_FunctionEval, ts, U, Q, 0)); 3443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 345814e85d6SHong Zhang } 346814e85d6SHong Zhang 347b5b4017aSHong Zhang /*@C 348bcf0153eSBarry Smith TSComputeDRDUFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobian()` 349814e85d6SHong Zhang 350d76be551SHong Zhang Level: deprecated 351814e85d6SHong Zhang 352814e85d6SHong Zhang @*/ 353d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDUFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) 354d71ae5a4SJacob Faibussowitsch { 355814e85d6SHong Zhang PetscFunctionBegin; 3563ba16761SJacob Faibussowitsch if (!DRDU) PetscFunctionReturn(PETSC_SUCCESS); 357814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 358c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 359814e85d6SHong Zhang 360792fecdfSBarry Smith PetscCallBack("TS callback DRDU for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 3613ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 362814e85d6SHong Zhang } 363814e85d6SHong Zhang 364b5b4017aSHong Zhang /*@C 365bcf0153eSBarry Smith TSComputeDRDPFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobianP()` 366814e85d6SHong Zhang 367d76be551SHong Zhang Level: deprecated 368814e85d6SHong Zhang 369814e85d6SHong Zhang @*/ 370d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) 371d71ae5a4SJacob Faibussowitsch { 372814e85d6SHong Zhang PetscFunctionBegin; 3733ba16761SJacob Faibussowitsch if (!DRDP) PetscFunctionReturn(PETSC_SUCCESS); 374814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 375c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 376814e85d6SHong Zhang 377792fecdfSBarry Smith PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 3783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 379814e85d6SHong Zhang } 380814e85d6SHong Zhang 381b5b4017aSHong Zhang /*@C 382b5b4017aSHong 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. 383b5b4017aSHong Zhang 384c3339decSBarry Smith Logically Collective 385b5b4017aSHong Zhang 386b5b4017aSHong Zhang Input Parameters: 387bcf0153eSBarry Smith + ts - `TS` context obtained from `TSCreate()` 388b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 389b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 390b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 391b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 392b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 393b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 394b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 395f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for F_PP 396b5b4017aSHong Zhang 39720f4b53cSBarry Smith Calling sequence of `ihessianproductfunc`: 39820f4b53cSBarry Smith $ PetscErrorCode ihessianproductfunc(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx); 399b5b4017aSHong Zhang + t - current timestep 400b5b4017aSHong Zhang . U - input vector (current ODE solution) 401369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 402b5b4017aSHong Zhang . Vr - input vector to be right-multiplied with the Hessian 403369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 404b5b4017aSHong Zhang - ctx - [optional] user-defined function context 405b5b4017aSHong Zhang 406b5b4017aSHong Zhang Level: intermediate 407b5b4017aSHong Zhang 408369cf35fSHong Zhang Notes: 409369cf35fSHong Zhang The first Hessian function and the working array are required. 410369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 411369cf35fSHong Zhang $ Vl_n^T*F_UP*Vr 412369cf35fSHong 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. 413369cf35fSHong Zhang Each entry of F_UP corresponds to the derivative 414369cf35fSHong Zhang $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}. 415369cf35fSHong 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 416369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]} 417369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 418b5b4017aSHong Zhang 419*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS` 420b5b4017aSHong Zhang @*/ 421d71ae5a4SJacob 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) 422d71ae5a4SJacob Faibussowitsch { 423b5b4017aSHong Zhang PetscFunctionBegin; 424b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 425b5b4017aSHong Zhang PetscValidPointer(ihp1, 2); 426b5b4017aSHong Zhang 427b5b4017aSHong Zhang ts->ihessianproductctx = ctx; 428b5b4017aSHong Zhang if (ihp1) ts->vecs_fuu = ihp1; 429b5b4017aSHong Zhang if (ihp2) ts->vecs_fup = ihp2; 430b5b4017aSHong Zhang if (ihp3) ts->vecs_fpu = ihp3; 431b5b4017aSHong Zhang if (ihp4) ts->vecs_fpp = ihp4; 432b5b4017aSHong Zhang ts->ihessianproduct_fuu = ihessianproductfunc1; 433b5b4017aSHong Zhang ts->ihessianproduct_fup = ihessianproductfunc2; 434b5b4017aSHong Zhang ts->ihessianproduct_fpu = ihessianproductfunc3; 435b5b4017aSHong Zhang ts->ihessianproduct_fpp = ihessianproductfunc4; 4363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 437b5b4017aSHong Zhang } 438b5b4017aSHong Zhang 439b5b4017aSHong Zhang /*@C 440b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu. 441b5b4017aSHong Zhang 442c3339decSBarry Smith Collective 443b5b4017aSHong Zhang 444b5b4017aSHong Zhang Input Parameters: 4452fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 4462fe279fdSBarry Smith . t - the time 4472fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 4482fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 4492fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 4502fe279fdSBarry Smith 4512fe279fdSBarry Smith Output Parameter: 4522fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 453b5b4017aSHong Zhang 454b5b4017aSHong Zhang Level: developer 455b5b4017aSHong Zhang 456bcf0153eSBarry Smith Note: 457bcf0153eSBarry Smith `TSComputeIHessianProductFunctionUU()` is typically used for sensitivity implementation, 458bcf0153eSBarry Smith so most users would not generally call this routine themselves. 459bcf0153eSBarry Smith 460*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()` 461b5b4017aSHong Zhang @*/ 462d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 463d71ae5a4SJacob Faibussowitsch { 464b5b4017aSHong Zhang PetscFunctionBegin; 4653ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 466b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 467b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 468b5b4017aSHong Zhang 469792fecdfSBarry Smith if (ts->ihessianproduct_fuu) PetscCallBack("TS callback IHessianProduct 1 for sensitivity analysis", (*ts->ihessianproduct_fuu)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 470ef1023bdSBarry Smith 47167633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 47267633408SHong Zhang if (ts->rhshessianproduct_guu) { 47367633408SHong Zhang PetscInt nadj; 4749566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(ts, t, U, Vl, Vr, VHV)); 47548a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 47667633408SHong Zhang } 4773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 478b5b4017aSHong Zhang } 479b5b4017aSHong Zhang 480b5b4017aSHong Zhang /*@C 481b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup. 482b5b4017aSHong Zhang 483c3339decSBarry Smith Collective 484b5b4017aSHong Zhang 485b5b4017aSHong Zhang Input Parameters: 4862fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 4872fe279fdSBarry Smith . t - the time 4882fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 4892fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 4902fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 4912fe279fdSBarry Smith 4922fe279fdSBarry Smith Output Parameter: 4932fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 494b5b4017aSHong Zhang 495b5b4017aSHong Zhang Level: developer 496b5b4017aSHong Zhang 497bcf0153eSBarry Smith Note: 498bcf0153eSBarry Smith `TSComputeIHessianProductFunctionUP()` is typically used for sensitivity implementation, 499bcf0153eSBarry Smith so most users would not generally call this routine themselves. 500bcf0153eSBarry Smith 501*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()` 502b5b4017aSHong Zhang @*/ 503d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 504d71ae5a4SJacob Faibussowitsch { 505b5b4017aSHong Zhang PetscFunctionBegin; 5063ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 507b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 508b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 509b5b4017aSHong Zhang 510792fecdfSBarry Smith if (ts->ihessianproduct_fup) PetscCallBack("TS callback IHessianProduct 2 for sensitivity analysis", (*ts->ihessianproduct_fup)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 511ef1023bdSBarry Smith 51267633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 51367633408SHong Zhang if (ts->rhshessianproduct_gup) { 51467633408SHong Zhang PetscInt nadj; 5159566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(ts, t, U, Vl, Vr, VHV)); 51648a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 51767633408SHong Zhang } 5183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 519b5b4017aSHong Zhang } 520b5b4017aSHong Zhang 521b5b4017aSHong Zhang /*@C 522b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu. 523b5b4017aSHong Zhang 524c3339decSBarry Smith Collective 525b5b4017aSHong Zhang 526b5b4017aSHong Zhang Input Parameters: 5272fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 5282fe279fdSBarry Smith . t - the time 5292fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 5302fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 5312fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 5322fe279fdSBarry Smith 5332fe279fdSBarry Smith Output Parameter: 5342fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 535b5b4017aSHong Zhang 536b5b4017aSHong Zhang Level: developer 537b5b4017aSHong Zhang 538bcf0153eSBarry Smith Note: 539bcf0153eSBarry Smith `TSComputeIHessianProductFunctionPU()` is typically used for sensitivity implementation, 540bcf0153eSBarry Smith so most users would not generally call this routine themselves. 541bcf0153eSBarry Smith 542*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()` 543b5b4017aSHong Zhang @*/ 544d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 545d71ae5a4SJacob Faibussowitsch { 546b5b4017aSHong Zhang PetscFunctionBegin; 5473ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 548b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 549b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 550b5b4017aSHong Zhang 551792fecdfSBarry Smith if (ts->ihessianproduct_fpu) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpu)(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_gpu) { 55567633408SHong Zhang PetscInt nadj; 5569566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(ts, t, U, Vl, Vr, VHV)); 55748a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 55867633408SHong Zhang } 5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 560b5b4017aSHong Zhang } 561b5b4017aSHong Zhang 562b5b4017aSHong Zhang /*@C 563b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp. 564b5b4017aSHong Zhang 565c3339decSBarry Smith Collective 566b5b4017aSHong Zhang 567b5b4017aSHong Zhang Input Parameters: 5682fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 5692fe279fdSBarry Smith . t - the time 5702fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 5712fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 5722fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 5732fe279fdSBarry Smith 5742fe279fdSBarry Smith Output Parameter: 5752fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 576b5b4017aSHong Zhang 577b5b4017aSHong Zhang Level: developer 578b5b4017aSHong Zhang 579bcf0153eSBarry Smith Note: 580bcf0153eSBarry Smith `TSComputeIHessianProductFunctionPP()` is typically used for sensitivity implementation, 581bcf0153eSBarry Smith so most users would not generally call this routine themselves. 582bcf0153eSBarry Smith 583*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetIHessianProduct()` 584b5b4017aSHong Zhang @*/ 585d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 586d71ae5a4SJacob Faibussowitsch { 587b5b4017aSHong Zhang PetscFunctionBegin; 5883ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 589b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 590b5b4017aSHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 591b5b4017aSHong Zhang 592792fecdfSBarry Smith if (ts->ihessianproduct_fpp) PetscCallBack("TS callback IHessianProduct 3 for sensitivity analysis", (*ts->ihessianproduct_fpp)(ts, t, U, Vl, Vr, VHV, ts->ihessianproductctx)); 593ef1023bdSBarry Smith 59467633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 59567633408SHong Zhang if (ts->rhshessianproduct_gpp) { 59667633408SHong Zhang PetscInt nadj; 5979566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(ts, t, U, Vl, Vr, VHV)); 59848a46eb9SPierre Jolivet for (nadj = 0; nadj < ts->numcost; nadj++) PetscCall(VecScale(VHV[nadj], -1)); 59967633408SHong Zhang } 6003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 601b5b4017aSHong Zhang } 602b5b4017aSHong Zhang 60313af1a74SHong Zhang /*@C 6044c500f23SPierre 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. 60513af1a74SHong Zhang 606c3339decSBarry Smith Logically Collective 60713af1a74SHong Zhang 60813af1a74SHong Zhang Input Parameters: 609bcf0153eSBarry Smith + ts - `TS` context obtained from `TSCreate()` 61013af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 61113af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 61213af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 61313af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 61413af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 61513af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 61613af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 617f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for G_PP 61813af1a74SHong Zhang 61920f4b53cSBarry Smith Calling sequence of `ihessianproductfunc`: 62020f4b53cSBarry Smith $ PetscErrorCode rhshessianproductfunc(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV, void *ctx); 62113af1a74SHong Zhang + t - current timestep 62213af1a74SHong Zhang . U - input vector (current ODE solution) 623369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 62413af1a74SHong Zhang . Vr - input vector to be right-multiplied with the Hessian 625369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 62613af1a74SHong Zhang - ctx - [optional] user-defined function context 62713af1a74SHong Zhang 62813af1a74SHong Zhang Level: intermediate 62913af1a74SHong Zhang 630369cf35fSHong Zhang Notes: 631369cf35fSHong Zhang The first Hessian function and the working array are required. 632369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 633369cf35fSHong Zhang $ Vl_n^T*G_UP*Vr 634369cf35fSHong 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. 635369cf35fSHong Zhang Each entry of G_UP corresponds to the derivative 636369cf35fSHong Zhang $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}. 637369cf35fSHong 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 638369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]} 639369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 64013af1a74SHong Zhang 6412fe279fdSBarry Smith .seealso: `TS`, `TSAdjoint` 64213af1a74SHong Zhang @*/ 643d71ae5a4SJacob 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) 644d71ae5a4SJacob Faibussowitsch { 64513af1a74SHong Zhang PetscFunctionBegin; 64613af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 64713af1a74SHong Zhang PetscValidPointer(rhshp1, 2); 64813af1a74SHong Zhang 64913af1a74SHong Zhang ts->rhshessianproductctx = ctx; 65013af1a74SHong Zhang if (rhshp1) ts->vecs_guu = rhshp1; 65113af1a74SHong Zhang if (rhshp2) ts->vecs_gup = rhshp2; 65213af1a74SHong Zhang if (rhshp3) ts->vecs_gpu = rhshp3; 65313af1a74SHong Zhang if (rhshp4) ts->vecs_gpp = rhshp4; 65413af1a74SHong Zhang ts->rhshessianproduct_guu = rhshessianproductfunc1; 65513af1a74SHong Zhang ts->rhshessianproduct_gup = rhshessianproductfunc2; 65613af1a74SHong Zhang ts->rhshessianproduct_gpu = rhshessianproductfunc3; 65713af1a74SHong Zhang ts->rhshessianproduct_gpp = rhshessianproductfunc4; 6583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 65913af1a74SHong Zhang } 66013af1a74SHong Zhang 66113af1a74SHong Zhang /*@C 662b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu. 66313af1a74SHong Zhang 664c3339decSBarry Smith Collective 66513af1a74SHong Zhang 66613af1a74SHong Zhang Input Parameters: 6672fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 6682fe279fdSBarry Smith . t - the time 6692fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 6702fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 6712fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 6722fe279fdSBarry Smith 6732fe279fdSBarry Smith Output Parameter: 6742fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 67513af1a74SHong Zhang 67613af1a74SHong Zhang Level: developer 67713af1a74SHong Zhang 678bcf0153eSBarry Smith Note: 679bcf0153eSBarry Smith `TSComputeRHSHessianProductFunctionUU()` is typically used for sensitivity implementation, 680bcf0153eSBarry Smith so most users would not generally call this routine themselves. 681bcf0153eSBarry Smith 682*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSHessianProduct()` 68313af1a74SHong Zhang @*/ 684d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 685d71ae5a4SJacob Faibussowitsch { 68613af1a74SHong Zhang PetscFunctionBegin; 6873ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 68813af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 68913af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 69013af1a74SHong Zhang 691792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 1 for sensitivity analysis", (*ts->rhshessianproduct_guu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 6923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 69313af1a74SHong Zhang } 69413af1a74SHong Zhang 69513af1a74SHong Zhang /*@C 696b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup. 69713af1a74SHong Zhang 698c3339decSBarry Smith Collective 69913af1a74SHong Zhang 70013af1a74SHong Zhang Input Parameters: 7012fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 7022fe279fdSBarry Smith . t - the time 7032fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 7042fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 7052fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 7062fe279fdSBarry Smith 7072fe279fdSBarry Smith Output Parameter: 7082fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 70913af1a74SHong Zhang 71013af1a74SHong Zhang Level: developer 71113af1a74SHong Zhang 712bcf0153eSBarry Smith Note: 713bcf0153eSBarry Smith `TSComputeRHSHessianProductFunctionUP()` is typically used for sensitivity implementation, 714bcf0153eSBarry Smith so most users would not generally call this routine themselves. 715bcf0153eSBarry Smith 716*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSetRHSHessianProduct()` 71713af1a74SHong Zhang @*/ 718d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 719d71ae5a4SJacob Faibussowitsch { 72013af1a74SHong Zhang PetscFunctionBegin; 7213ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 72213af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 72313af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 72413af1a74SHong Zhang 725792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 2 for sensitivity analysis", (*ts->rhshessianproduct_gup)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 7263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 72713af1a74SHong Zhang } 72813af1a74SHong Zhang 72913af1a74SHong Zhang /*@C 730b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu. 73113af1a74SHong Zhang 732c3339decSBarry Smith Collective 73313af1a74SHong Zhang 73413af1a74SHong Zhang Input Parameters: 7352fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 7362fe279fdSBarry Smith . t - the time 7372fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 7382fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 7392fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 7402fe279fdSBarry Smith 7412fe279fdSBarry Smith Output Parameter: 7422fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 74313af1a74SHong Zhang 74413af1a74SHong Zhang Level: developer 74513af1a74SHong Zhang 746bcf0153eSBarry Smith Note: 747bcf0153eSBarry Smith `TSComputeRHSHessianProductFunctionPU()` is typically used for sensitivity implementation, 748bcf0153eSBarry Smith so most users would not generally call this routine themselves. 749bcf0153eSBarry Smith 750*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSHessianProduct()` 75113af1a74SHong Zhang @*/ 752d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 753d71ae5a4SJacob Faibussowitsch { 75413af1a74SHong Zhang PetscFunctionBegin; 7553ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 75613af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 75713af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 75813af1a74SHong Zhang 759792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpu)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 7603ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 76113af1a74SHong Zhang } 76213af1a74SHong Zhang 76313af1a74SHong Zhang /*@C 764b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp. 76513af1a74SHong Zhang 766c3339decSBarry Smith Collective 76713af1a74SHong Zhang 76813af1a74SHong Zhang Input Parameters: 7692fe279fdSBarry Smith + ts - The `TS` context obtained from `TSCreate()` 7702fe279fdSBarry Smith . t - the time 7712fe279fdSBarry Smith . U - the solution at which to compute the Hessian product 7722fe279fdSBarry Smith . Vl - the array of input vectors to be multiplied with the Hessian from the left 7732fe279fdSBarry Smith - Vr - the input vector to be multiplied with the Hessian from the right 7742fe279fdSBarry Smith 7752fe279fdSBarry Smith Output Parameter: 7762fe279fdSBarry Smith . VhV - the array of output vectors that store the Hessian product 77713af1a74SHong Zhang 77813af1a74SHong Zhang Level: developer 77913af1a74SHong Zhang 780bcf0153eSBarry Smith Note: 781bcf0153eSBarry Smith `TSComputeRHSHessianProductFunctionPP()` is typically used for sensitivity implementation, 782bcf0153eSBarry Smith so most users would not generally call this routine themselves. 783bcf0153eSBarry Smith 784*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetRHSHessianProduct()` 78513af1a74SHong Zhang @*/ 786d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts, PetscReal t, Vec U, Vec *Vl, Vec Vr, Vec *VHV) 787d71ae5a4SJacob Faibussowitsch { 78813af1a74SHong Zhang PetscFunctionBegin; 7893ba16761SJacob Faibussowitsch if (!VHV) PetscFunctionReturn(PETSC_SUCCESS); 79013af1a74SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 79113af1a74SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 79213af1a74SHong Zhang 793792fecdfSBarry Smith PetscCallBack("TS callback RHSHessianProduct 3 for sensitivity analysis", (*ts->rhshessianproduct_gpp)(ts, t, U, Vl, Vr, VHV, ts->rhshessianproductctx)); 7943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 79513af1a74SHong Zhang } 79613af1a74SHong Zhang 797814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 798814e85d6SHong Zhang 799814e85d6SHong Zhang /*@ 800814e85d6SHong 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 801bcf0153eSBarry Smith for use by the `TS` adjoint routines. 802814e85d6SHong Zhang 803c3339decSBarry Smith Logically Collective 804814e85d6SHong Zhang 805814e85d6SHong Zhang Input Parameters: 806bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 807d2fd7bfcSHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 808814e85d6SHong 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 809814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 810814e85d6SHong Zhang 811814e85d6SHong Zhang Level: beginner 812814e85d6SHong Zhang 813814e85d6SHong Zhang Notes: 81435cb6cd3SPierre Jolivet the entries in these vectors must be correctly initialized with the values lambda_i = df/dy|finaltime mu_i = df/dp|finaltime 815814e85d6SHong Zhang 81635cb6cd3SPierre Jolivet After `TSAdjointSolve()` is called the lambda and the mu contain the computed sensitivities 817814e85d6SHong Zhang 818bcf0153eSBarry Smith .seealso: `TS`, `TSAdjointSolve()`, `TSGetCostGradients()` 819814e85d6SHong Zhang @*/ 820d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostGradients(TS ts, PetscInt numcost, Vec *lambda, Vec *mu) 821d71ae5a4SJacob Faibussowitsch { 822814e85d6SHong Zhang PetscFunctionBegin; 823814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 824064a246eSJacob Faibussowitsch PetscValidPointer(lambda, 3); 825814e85d6SHong Zhang ts->vecs_sensi = lambda; 826814e85d6SHong Zhang ts->vecs_sensip = mu; 8273c633725SBarry 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"); 828814e85d6SHong Zhang ts->numcost = numcost; 8293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 830814e85d6SHong Zhang } 831814e85d6SHong Zhang 832814e85d6SHong Zhang /*@ 833bcf0153eSBarry Smith TSGetCostGradients - Returns the gradients from the `TSAdjointSolve()` 834814e85d6SHong Zhang 835bcf0153eSBarry Smith Not Collective, but the vectors returned are parallel if `TS` is parallel 836814e85d6SHong Zhang 837814e85d6SHong Zhang Input Parameter: 838bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 839814e85d6SHong Zhang 840d8d19677SJose E. Roman Output Parameters: 8416b867d5aSJose E. Roman + numcost - size of returned arrays 8426b867d5aSJose E. Roman . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 843814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 844814e85d6SHong Zhang 845814e85d6SHong Zhang Level: intermediate 846814e85d6SHong Zhang 847*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSSetCostGradients()` 848814e85d6SHong Zhang @*/ 849d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostGradients(TS ts, PetscInt *numcost, Vec **lambda, Vec **mu) 850d71ae5a4SJacob Faibussowitsch { 851814e85d6SHong Zhang PetscFunctionBegin; 852814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 853814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 854814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 855814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 8563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 857814e85d6SHong Zhang } 858814e85d6SHong Zhang 859814e85d6SHong Zhang /*@ 860b5b4017aSHong 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 861bcf0153eSBarry Smith for use by the `TS` adjoint routines. 862b5b4017aSHong Zhang 863c3339decSBarry Smith Logically Collective 864b5b4017aSHong Zhang 865b5b4017aSHong Zhang Input Parameters: 866bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 867b5b4017aSHong Zhang . numcost - number of cost functions 868b5b4017aSHong 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 869b5b4017aSHong 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 870b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 871b5b4017aSHong Zhang 872b5b4017aSHong Zhang Level: beginner 873b5b4017aSHong Zhang 874bcf0153eSBarry Smith Notes: 875bcf0153eSBarry Smith Hessian of the cost function is completely different from Hessian of the ODE/DAE system 876b5b4017aSHong Zhang 877bcf0153eSBarry Smith For second-order adjoint, one needs to call this function and then `TSAdjointSetForward()` before `TSSolve()`. 878b5b4017aSHong Zhang 87935cb6cd3SPierre Jolivet After `TSAdjointSolve()` is called, the lambda2 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. 880b5b4017aSHong Zhang 8812fe279fdSBarry Smith Passing `NULL` for `lambda2` disables the second-order calculation. 882bcf0153eSBarry Smith 883*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointSetForward()` 884b5b4017aSHong Zhang @*/ 885d71ae5a4SJacob Faibussowitsch PetscErrorCode TSSetCostHessianProducts(TS ts, PetscInt numcost, Vec *lambda2, Vec *mu2, Vec dir) 886d71ae5a4SJacob Faibussowitsch { 887b5b4017aSHong Zhang PetscFunctionBegin; 888b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 8893c633725SBarry 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"); 890b5b4017aSHong Zhang ts->numcost = numcost; 891b5b4017aSHong Zhang ts->vecs_sensi2 = lambda2; 89233f72c5dSHong Zhang ts->vecs_sensi2p = mu2; 893b5b4017aSHong Zhang ts->vec_dir = dir; 8943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 895b5b4017aSHong Zhang } 896b5b4017aSHong Zhang 897b5b4017aSHong Zhang /*@ 898bcf0153eSBarry Smith TSGetCostHessianProducts - Returns the gradients from the `TSAdjointSolve()` 899b5b4017aSHong Zhang 900bcf0153eSBarry Smith Not Collective, but vectors returned are parallel if `TS` is parallel 901b5b4017aSHong Zhang 902b5b4017aSHong Zhang Input Parameter: 903bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 904b5b4017aSHong Zhang 905d8d19677SJose E. Roman Output Parameters: 906b5b4017aSHong Zhang + numcost - number of cost functions 907b5b4017aSHong 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 908b5b4017aSHong 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 909b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 910b5b4017aSHong Zhang 911b5b4017aSHong Zhang Level: intermediate 912b5b4017aSHong Zhang 913*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSSetCostHessianProducts()` 914b5b4017aSHong Zhang @*/ 915d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetCostHessianProducts(TS ts, PetscInt *numcost, Vec **lambda2, Vec **mu2, Vec *dir) 916d71ae5a4SJacob Faibussowitsch { 917b5b4017aSHong Zhang PetscFunctionBegin; 918b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 919b5b4017aSHong Zhang if (numcost) *numcost = ts->numcost; 920b5b4017aSHong Zhang if (lambda2) *lambda2 = ts->vecs_sensi2; 92133f72c5dSHong Zhang if (mu2) *mu2 = ts->vecs_sensi2p; 922b5b4017aSHong Zhang if (dir) *dir = ts->vec_dir; 9233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 924b5b4017aSHong Zhang } 925b5b4017aSHong Zhang 926b5b4017aSHong Zhang /*@ 927ecf68647SHong Zhang TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities 928b5b4017aSHong Zhang 929c3339decSBarry Smith Logically Collective 930b5b4017aSHong Zhang 931b5b4017aSHong Zhang Input Parameters: 932bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 933b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters 934b5b4017aSHong Zhang 935b5b4017aSHong Zhang Level: intermediate 936b5b4017aSHong Zhang 937bcf0153eSBarry Smith Notes: 9382fe279fdSBarry Smith When computing sensitivies w.r.t. initial condition, set didp to `NULL` so that the solver will take it as an identity matrix mathematically. 9392fe279fdSBarry Smith `TSAdjoint` does not reset the tangent linear solver automatically, `TSAdjointResetForward()` should be called to reset the tangent linear solver. 940b5b4017aSHong Zhang 941*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSSetCostHessianProducts()`, `TSAdjointResetForward()` 942b5b4017aSHong Zhang @*/ 943d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetForward(TS ts, Mat didp) 944d71ae5a4SJacob Faibussowitsch { 94533f72c5dSHong Zhang Mat A; 94633f72c5dSHong Zhang Vec sp; 94733f72c5dSHong Zhang PetscScalar *xarr; 94833f72c5dSHong Zhang PetscInt lsize; 949b5b4017aSHong Zhang 950b5b4017aSHong Zhang PetscFunctionBegin; 951b5b4017aSHong Zhang ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 9523c633725SBarry Smith PetscCheck(ts->vecs_sensi2, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Must call TSSetCostHessianProducts() first"); 9533c633725SBarry Smith PetscCheck(ts->vec_dir, PetscObjectComm((PetscObject)ts), PETSC_ERR_USER, "Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 95433f72c5dSHong Zhang /* create a single-column dense matrix */ 9559566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol, &lsize)); 9569566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ts), lsize, PETSC_DECIDE, PETSC_DECIDE, 1, NULL, &A)); 95733f72c5dSHong Zhang 9589566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &sp)); 9599566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(A, 0, &xarr)); 9609566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(sp, xarr)); 961ecf68647SHong Zhang if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */ 96233f72c5dSHong Zhang if (didp) { 9639566063dSJacob Faibussowitsch PetscCall(MatMult(didp, ts->vec_dir, sp)); 9649566063dSJacob Faibussowitsch PetscCall(VecScale(sp, 2.)); 96533f72c5dSHong Zhang } else { 9669566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(sp)); 96733f72c5dSHong Zhang } 968ecf68647SHong Zhang } else { /* tangent linear variable initialized as dir */ 9699566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dir, sp)); 97033f72c5dSHong Zhang } 9719566063dSJacob Faibussowitsch PetscCall(VecResetArray(sp)); 9729566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(A, &xarr)); 9739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sp)); 97433f72c5dSHong Zhang 9759566063dSJacob Faibussowitsch PetscCall(TSForwardSetInitialSensitivities(ts, A)); /* if didp is NULL, identity matrix is assumed */ 97633f72c5dSHong Zhang 9779566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 9783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 979b5b4017aSHong Zhang } 980b5b4017aSHong Zhang 981b5b4017aSHong Zhang /*@ 982ecf68647SHong Zhang TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context 983ecf68647SHong Zhang 984c3339decSBarry Smith Logically Collective 985ecf68647SHong Zhang 9862fe279fdSBarry Smith Input Parameter: 987bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 988ecf68647SHong Zhang 989ecf68647SHong Zhang Level: intermediate 990ecf68647SHong Zhang 991*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSetForward()` 992ecf68647SHong Zhang @*/ 993d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointResetForward(TS ts) 994d71ae5a4SJacob Faibussowitsch { 995ecf68647SHong Zhang PetscFunctionBegin; 996ecf68647SHong Zhang ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */ 9979566063dSJacob Faibussowitsch PetscCall(TSForwardReset(ts)); 9983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 999ecf68647SHong Zhang } 1000ecf68647SHong Zhang 1001ecf68647SHong Zhang /*@ 1002814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 1003814e85d6SHong Zhang of an adjoint solver 1004814e85d6SHong Zhang 1005c3339decSBarry Smith Collective 1006814e85d6SHong Zhang 1007814e85d6SHong Zhang Input Parameter: 1008bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1009814e85d6SHong Zhang 1010814e85d6SHong Zhang Level: advanced 1011814e85d6SHong Zhang 1012*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSAdjointStep()`, `TSSetCostGradients()` 1013814e85d6SHong Zhang @*/ 1014d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetUp(TS ts) 1015d71ae5a4SJacob Faibussowitsch { 1016881c1a9bSHong Zhang TSTrajectory tj; 1017881c1a9bSHong Zhang PetscBool match; 1018814e85d6SHong Zhang 1019814e85d6SHong Zhang PetscFunctionBegin; 1020814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 10213ba16761SJacob Faibussowitsch if (ts->adjointsetupcalled) PetscFunctionReturn(PETSC_SUCCESS); 10223c633725SBarry Smith PetscCheck(ts->vecs_sensi, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetCostGradients() first"); 10233c633725SBarry Smith PetscCheck(!ts->vecs_sensip || ts->Jacp || ts->Jacprhs, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_WRONGSTATE, "Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 10249566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts, &tj)); 10259566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tj, TSTRAJECTORYBASIC, &match)); 1026881c1a9bSHong Zhang if (match) { 1027881c1a9bSHong Zhang PetscBool solution_only; 10289566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetSolutionOnly(tj, &solution_only)); 10293c633725SBarry 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"); 1030881c1a9bSHong Zhang } 10319566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUseHistory(tj, PETSC_FALSE)); /* not use TSHistory */ 103233f72c5dSHong Zhang 1033cd4cee2dSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 10349566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vecs_sensi[0], &ts->vec_drdu_col)); 103548a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDuplicate(ts->vecs_sensip[0], &ts->vec_drdp_col)); 1036814e85d6SHong Zhang } 1037814e85d6SHong Zhang 1038dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, adjointsetup); 1039814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 10403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1041814e85d6SHong Zhang } 1042814e85d6SHong Zhang 1043814e85d6SHong Zhang /*@ 1044bcf0153eSBarry Smith TSAdjointReset - Resets a `TS` adjoint context and removes any allocated `Vec`s and `Mat`s. 1045ece46509SHong Zhang 1046c3339decSBarry Smith Collective 1047ece46509SHong Zhang 1048ece46509SHong Zhang Input Parameter: 1049bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1050ece46509SHong Zhang 1051ece46509SHong Zhang Level: beginner 1052ece46509SHong Zhang 1053*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSAdjointSetUp()`, `TSADestroy()` 1054ece46509SHong Zhang @*/ 1055d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointReset(TS ts) 1056d71ae5a4SJacob Faibussowitsch { 1057ece46509SHong Zhang PetscFunctionBegin; 1058ece46509SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1059dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, adjointreset); 10607207e0fdSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 10619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_drdu_col)); 106248a46eb9SPierre Jolivet if (ts->vecs_sensip) PetscCall(VecDestroy(&ts->vec_drdp_col)); 10637207e0fdSHong Zhang } 1064ece46509SHong Zhang ts->vecs_sensi = NULL; 1065ece46509SHong Zhang ts->vecs_sensip = NULL; 1066ece46509SHong Zhang ts->vecs_sensi2 = NULL; 106733f72c5dSHong Zhang ts->vecs_sensi2p = NULL; 1068ece46509SHong Zhang ts->vec_dir = NULL; 1069ece46509SHong Zhang ts->adjointsetupcalled = PETSC_FALSE; 10703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1071ece46509SHong Zhang } 1072ece46509SHong Zhang 1073ece46509SHong Zhang /*@ 1074814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1075814e85d6SHong Zhang 1076c3339decSBarry Smith Logically Collective 1077814e85d6SHong Zhang 1078814e85d6SHong Zhang Input Parameters: 1079bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1080a2b725a8SWilliam Gropp - steps - number of steps to use 1081814e85d6SHong Zhang 1082814e85d6SHong Zhang Level: intermediate 1083814e85d6SHong Zhang 1084814e85d6SHong Zhang Notes: 1085bcf0153eSBarry Smith Normally one does not call this and `TSAdjointSolve()` integrates back to the original timestep. One can call this 1086814e85d6SHong Zhang so as to integrate back to less than the original timestep 1087814e85d6SHong Zhang 1088*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TS`, `TSSetExactFinalTime()` 1089814e85d6SHong Zhang @*/ 1090d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetSteps(TS ts, PetscInt steps) 1091d71ae5a4SJacob Faibussowitsch { 1092814e85d6SHong Zhang PetscFunctionBegin; 1093814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1094814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts, steps, 2); 10953c633725SBarry Smith PetscCheck(steps >= 0, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back a negative number of steps"); 10963c633725SBarry Smith PetscCheck(steps <= ts->steps, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Cannot step back more than the total number of forward steps"); 1097814e85d6SHong Zhang ts->adjoint_max_steps = steps; 10983ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1099814e85d6SHong Zhang } 1100814e85d6SHong Zhang 1101814e85d6SHong Zhang /*@C 1102bcf0153eSBarry Smith TSAdjointSetRHSJacobian - Deprecated, use `TSSetRHSJacobianP()` 1103814e85d6SHong Zhang 1104814e85d6SHong Zhang Level: deprecated 1105814e85d6SHong Zhang 1106814e85d6SHong Zhang @*/ 1107d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetRHSJacobian(TS ts, Mat Amat, PetscErrorCode (*func)(TS, PetscReal, Vec, Mat, void *), void *ctx) 1108d71ae5a4SJacob Faibussowitsch { 1109814e85d6SHong Zhang PetscFunctionBegin; 1110814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1111814e85d6SHong Zhang PetscValidHeaderSpecific(Amat, MAT_CLASSID, 2); 1112814e85d6SHong Zhang 1113814e85d6SHong Zhang ts->rhsjacobianp = func; 1114814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 1115814e85d6SHong Zhang if (Amat) { 11169566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 11179566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 1118814e85d6SHong Zhang ts->Jacp = Amat; 1119814e85d6SHong Zhang } 11203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1121814e85d6SHong Zhang } 1122814e85d6SHong Zhang 1123814e85d6SHong Zhang /*@C 1124bcf0153eSBarry Smith TSAdjointComputeRHSJacobian - Deprecated, use `TSComputeRHSJacobianP()` 1125814e85d6SHong Zhang 1126814e85d6SHong Zhang Level: deprecated 1127814e85d6SHong Zhang 1128814e85d6SHong Zhang @*/ 1129d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeRHSJacobian(TS ts, PetscReal t, Vec U, Mat Amat) 1130d71ae5a4SJacob Faibussowitsch { 1131814e85d6SHong Zhang PetscFunctionBegin; 1132814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1133c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1134814e85d6SHong Zhang PetscValidPointer(Amat, 4); 1135814e85d6SHong Zhang 1136792fecdfSBarry Smith PetscCallBack("TS callback JacobianP for sensitivity analysis", (*ts->rhsjacobianp)(ts, t, U, Amat, ts->rhsjacobianpctx)); 11373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1138814e85d6SHong Zhang } 1139814e85d6SHong Zhang 1140814e85d6SHong Zhang /*@ 1141bcf0153eSBarry Smith TSAdjointComputeDRDYFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobian()` 1142814e85d6SHong Zhang 1143814e85d6SHong Zhang Level: deprecated 1144814e85d6SHong Zhang 1145814e85d6SHong Zhang @*/ 1146d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDYFunction(TS ts, PetscReal t, Vec U, Vec *DRDU) 1147d71ae5a4SJacob Faibussowitsch { 1148814e85d6SHong Zhang PetscFunctionBegin; 1149814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1150c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1151814e85d6SHong Zhang 1152792fecdfSBarry Smith PetscCallBack("TS callback DRDY for sensitivity analysis", (*ts->drdufunction)(ts, t, U, DRDU, ts->costintegrandctx)); 11533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1154814e85d6SHong Zhang } 1155814e85d6SHong Zhang 1156814e85d6SHong Zhang /*@ 1157bcf0153eSBarry Smith TSAdjointComputeDRDPFunction - Deprecated, use `TSGetQuadratureTS()` then `TSComputeRHSJacobianP()` 1158814e85d6SHong Zhang 1159814e85d6SHong Zhang Level: deprecated 1160814e85d6SHong Zhang 1161814e85d6SHong Zhang @*/ 1162d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointComputeDRDPFunction(TS ts, PetscReal t, Vec U, Vec *DRDP) 1163d71ae5a4SJacob Faibussowitsch { 1164814e85d6SHong Zhang PetscFunctionBegin; 1165814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1166c9ad9de0SHong Zhang PetscValidHeaderSpecific(U, VEC_CLASSID, 3); 1167814e85d6SHong Zhang 1168792fecdfSBarry Smith PetscCallBack("TS callback DRDP for sensitivity analysis", (*ts->drdpfunction)(ts, t, U, DRDP, ts->costintegrandctx)); 11693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1170814e85d6SHong Zhang } 1171814e85d6SHong Zhang 1172814e85d6SHong Zhang /*@C 1173814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 1174814e85d6SHong Zhang 1175814e85d6SHong Zhang Level: intermediate 1176814e85d6SHong Zhang 1177*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointMonitorSet()` 1178814e85d6SHong Zhang @*/ 1179d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorSensi(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) 1180d71ae5a4SJacob Faibussowitsch { 1181814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1182814e85d6SHong Zhang 1183814e85d6SHong Zhang PetscFunctionBegin; 1184064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 11859566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 11869566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], viewer)); 11879566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 11883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1189814e85d6SHong Zhang } 1190814e85d6SHong Zhang 1191814e85d6SHong Zhang /*@C 1192814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1193814e85d6SHong Zhang 1194c3339decSBarry Smith Collective 1195814e85d6SHong Zhang 1196814e85d6SHong Zhang Input Parameters: 1197bcf0153eSBarry Smith + ts - `TS` object you wish to monitor 1198814e85d6SHong Zhang . name - the monitor type one is seeking 1199814e85d6SHong Zhang . help - message indicating what monitoring is done 1200814e85d6SHong Zhang . manual - manual page for the monitor 1201814e85d6SHong Zhang . monitor - the monitor function 1202bcf0153eSBarry Smith - 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 1203814e85d6SHong Zhang 1204814e85d6SHong Zhang Level: developer 1205814e85d6SHong Zhang 1206*1cc06b55SBarry Smith .seealso: [](ch_ts), `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 1207db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 1208db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1209db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1210c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1211db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1212db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 1213814e85d6SHong Zhang @*/ 1214d71ae5a4SJacob 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 *)) 1215d71ae5a4SJacob Faibussowitsch { 1216814e85d6SHong Zhang PetscViewer viewer; 1217814e85d6SHong Zhang PetscViewerFormat format; 1218814e85d6SHong Zhang PetscBool flg; 1219814e85d6SHong Zhang 1220814e85d6SHong Zhang PetscFunctionBegin; 12219566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts), ((PetscObject)ts)->options, ((PetscObject)ts)->prefix, name, &viewer, &format, &flg)); 1222814e85d6SHong Zhang if (flg) { 1223814e85d6SHong Zhang PetscViewerAndFormat *vf; 12249566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 12259566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 12261baa6e33SBarry Smith if (monitorsetup) PetscCall((*monitorsetup)(ts, vf)); 12279566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts, (PetscErrorCode(*)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *))monitor, vf, (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy)); 1228814e85d6SHong Zhang } 12293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1230814e85d6SHong Zhang } 1231814e85d6SHong Zhang 1232814e85d6SHong Zhang /*@C 1233814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1234814e85d6SHong Zhang timestep to display the iteration's progress. 1235814e85d6SHong Zhang 1236c3339decSBarry Smith Logically Collective 1237814e85d6SHong Zhang 1238814e85d6SHong Zhang Input Parameters: 1239bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1240814e85d6SHong Zhang . adjointmonitor - monitoring routine 1241814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 12422fe279fdSBarry Smith monitor routine (use `NULL` if no context is desired) 1243814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 12442fe279fdSBarry Smith (may be `NULL`) 1245814e85d6SHong Zhang 124620f4b53cSBarry Smith Calling sequence of `adjointmonitor`: 124720f4b53cSBarry Smith $ PetscErrorCode adjointmonitor(TS ts, PetscInt steps, PetscReal time, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *adjointmctx) 1248bcf0153eSBarry Smith + ts - the `TS` context 1249814e85d6SHong 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 1250814e85d6SHong Zhang been interpolated to) 1251814e85d6SHong Zhang . time - current time 1252814e85d6SHong Zhang . u - current iterate 1253814e85d6SHong Zhang . numcost - number of cost functionos 1254814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1255814e85d6SHong Zhang . mu - sensitivities to parameters 1256814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 1257814e85d6SHong Zhang 1258bcf0153eSBarry Smith Level: intermediate 1259bcf0153eSBarry Smith 1260bcf0153eSBarry Smith Note: 1261814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 1262814e85d6SHong Zhang already has been loaded. 1263814e85d6SHong Zhang 1264bcf0153eSBarry Smith Fortran Note: 126520f4b53cSBarry Smith Only a single monitor function can be set for each `TS` object 1266814e85d6SHong Zhang 1267*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorCancel()` 1268814e85d6SHong Zhang @*/ 1269d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorSet(TS ts, PetscErrorCode (*adjointmonitor)(TS, PetscInt, PetscReal, Vec, PetscInt, Vec *, Vec *, void *), void *adjointmctx, PetscErrorCode (*adjointmdestroy)(void **)) 1270d71ae5a4SJacob Faibussowitsch { 1271814e85d6SHong Zhang PetscInt i; 1272814e85d6SHong Zhang PetscBool identical; 1273814e85d6SHong Zhang 1274814e85d6SHong Zhang PetscFunctionBegin; 1275814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1276814e85d6SHong Zhang for (i = 0; i < ts->numbermonitors; i++) { 12779566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode(*)(void))adjointmonitor, adjointmctx, adjointmdestroy, (PetscErrorCode(*)(void))ts->adjointmonitor[i], ts->adjointmonitorcontext[i], ts->adjointmonitordestroy[i], &identical)); 12783ba16761SJacob Faibussowitsch if (identical) PetscFunctionReturn(PETSC_SUCCESS); 1279814e85d6SHong Zhang } 12803c633725SBarry Smith PetscCheck(ts->numberadjointmonitors < MAXTSMONITORS, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many adjoint monitors set"); 1281814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1282814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1283814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void *)adjointmctx; 12843ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1285814e85d6SHong Zhang } 1286814e85d6SHong Zhang 1287814e85d6SHong Zhang /*@C 1288814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1289814e85d6SHong Zhang 1290c3339decSBarry Smith Logically Collective 1291814e85d6SHong Zhang 12922fe279fdSBarry Smith Input Parameter: 1293bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1294814e85d6SHong Zhang 1295814e85d6SHong Zhang Notes: 1296814e85d6SHong Zhang There is no way to remove a single, specific monitor. 1297814e85d6SHong Zhang 1298814e85d6SHong Zhang Level: intermediate 1299814e85d6SHong Zhang 1300*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorSet()` 1301814e85d6SHong Zhang @*/ 1302d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorCancel(TS ts) 1303d71ae5a4SJacob Faibussowitsch { 1304814e85d6SHong Zhang PetscInt i; 1305814e85d6SHong Zhang 1306814e85d6SHong Zhang PetscFunctionBegin; 1307814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1308814e85d6SHong Zhang for (i = 0; i < ts->numberadjointmonitors; i++) { 130948a46eb9SPierre Jolivet if (ts->adjointmonitordestroy[i]) PetscCall((*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i])); 1310814e85d6SHong Zhang } 1311814e85d6SHong Zhang ts->numberadjointmonitors = 0; 13123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1313814e85d6SHong Zhang } 1314814e85d6SHong Zhang 1315814e85d6SHong Zhang /*@C 1316814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 1317814e85d6SHong Zhang 1318814e85d6SHong Zhang Level: intermediate 1319814e85d6SHong Zhang 1320*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSAdjointSolve()`, `TSAdjointMonitorSet()` 1321814e85d6SHong Zhang @*/ 1322d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDefault(TS ts, PetscInt step, PetscReal ptime, Vec v, PetscInt numcost, Vec *lambda, Vec *mu, PetscViewerAndFormat *vf) 1323d71ae5a4SJacob Faibussowitsch { 1324814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1325814e85d6SHong Zhang 1326814e85d6SHong Zhang PetscFunctionBegin; 1327064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 8); 13289566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer, vf->format)); 13299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer, ((PetscObject)ts)->tablevel)); 133063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "%" PetscInt_FMT " TS dt %g time %g%s", step, (double)ts->time_step, (double)ptime, ts->steprollback ? " (r)\n" : "\n")); 13319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer, ((PetscObject)ts)->tablevel)); 13329566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 13333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1334814e85d6SHong Zhang } 1335814e85d6SHong Zhang 1336814e85d6SHong Zhang /*@C 1337bcf0153eSBarry Smith TSAdjointMonitorDrawSensi - Monitors progress of the adjoint `TS` solvers by calling 1338bcf0153eSBarry Smith `VecView()` for the sensitivities to initial states at each timestep 1339814e85d6SHong Zhang 1340c3339decSBarry Smith Collective 1341814e85d6SHong Zhang 1342814e85d6SHong Zhang Input Parameters: 1343bcf0153eSBarry Smith + ts - the `TS` context 1344814e85d6SHong Zhang . step - current time-step 1345814e85d6SHong Zhang . ptime - current time 1346814e85d6SHong Zhang . u - current state 1347814e85d6SHong Zhang . numcost - number of cost functions 1348814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1349814e85d6SHong Zhang . mu - sensitivities to parameters 13502fe279fdSBarry Smith - dummy - either a viewer or `NULL` 1351814e85d6SHong Zhang 1352814e85d6SHong Zhang Level: intermediate 1353814e85d6SHong Zhang 1354*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSAdjointMonitorSet()`, `TSAdjointMonitorDefault()`, `VecView()` 1355814e85d6SHong Zhang @*/ 1356d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitorDrawSensi(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu, void *dummy) 1357d71ae5a4SJacob Faibussowitsch { 1358814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1359814e85d6SHong Zhang PetscDraw draw; 1360814e85d6SHong Zhang PetscReal xl, yl, xr, yr, h; 1361814e85d6SHong Zhang char time[32]; 1362814e85d6SHong Zhang 1363814e85d6SHong Zhang PetscFunctionBegin; 13643ba16761SJacob Faibussowitsch if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(PETSC_SUCCESS); 1365814e85d6SHong Zhang 13669566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0], ictx->viewer)); 13679566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer, 0, &draw)); 13689566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(time, 32, "Timestep %d Time %g", (int)step, (double)ptime)); 13699566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw, &xl, &yl, &xr, &yr)); 1370814e85d6SHong Zhang h = yl + .95 * (yr - yl); 13719566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw, .5 * (xl + xr), h, PETSC_DRAW_BLACK, time)); 13729566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 13733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1374814e85d6SHong Zhang } 1375814e85d6SHong Zhang 1376814e85d6SHong Zhang /* 1377bcf0153eSBarry Smith TSAdjointSetFromOptions - Sets various `TS` adjoint parameters from user options. 1378814e85d6SHong Zhang 1379c3339decSBarry Smith Collective 1380814e85d6SHong Zhang 1381814e85d6SHong Zhang Input Parameter: 1382bcf0153eSBarry Smith . ts - the `TS` context 1383814e85d6SHong Zhang 1384814e85d6SHong Zhang Options Database Keys: 1385814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1386814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 1387814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1388814e85d6SHong Zhang 1389814e85d6SHong Zhang Level: developer 1390814e85d6SHong Zhang 1391bcf0153eSBarry Smith Note: 1392814e85d6SHong Zhang This is not normally called directly by users 1393814e85d6SHong Zhang 1394*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSSetSaveTrajectory()`, `TSTrajectorySetUp()` 1395814e85d6SHong Zhang */ 1396d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSetFromOptions(TS ts, PetscOptionItems *PetscOptionsObject) 1397d71ae5a4SJacob Faibussowitsch { 1398814e85d6SHong Zhang PetscBool tflg, opt; 1399814e85d6SHong Zhang 1400814e85d6SHong Zhang PetscFunctionBegin; 1401dbbe0bcdSBarry Smith PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1402d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "TS Adjoint options"); 1403814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 14049566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adjoint_solve", "Solve the adjoint problem immediately after solving the forward problem", "", tflg, &tflg, &opt)); 1405814e85d6SHong Zhang if (opt) { 14069566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 1407814e85d6SHong Zhang ts->adjoint_solve = tflg; 1408814e85d6SHong Zhang } 14099566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor", "Monitor adjoint timestep size", "TSAdjointMonitorDefault", TSAdjointMonitorDefault, NULL)); 14109566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts, "-ts_adjoint_monitor_sensi", "Monitor sensitivity in the adjoint computation", "TSAdjointMonitorSensi", TSAdjointMonitorSensi, NULL)); 1411814e85d6SHong Zhang opt = PETSC_FALSE; 14129566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", &opt)); 1413814e85d6SHong Zhang if (opt) { 1414814e85d6SHong Zhang TSMonitorDrawCtx ctx; 1415814e85d6SHong Zhang PetscInt howoften = 1; 1416814e85d6SHong Zhang 14179566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_adjoint_monitor_draw_sensi", "Monitor adjoint sensitivities (lambda only) graphically", "TSAdjointMonitorDrawSensi", howoften, &howoften, NULL)); 14189566063dSJacob Faibussowitsch PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts), NULL, NULL, PETSC_DECIDE, PETSC_DECIDE, 300, 300, howoften, &ctx)); 14199566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts, TSAdjointMonitorDrawSensi, ctx, (PetscErrorCode(*)(void **))TSMonitorDrawCtxDestroy)); 1420814e85d6SHong Zhang } 14213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1422814e85d6SHong Zhang } 1423814e85d6SHong Zhang 1424814e85d6SHong Zhang /*@ 1425814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 1426814e85d6SHong Zhang 1427c3339decSBarry Smith Collective 1428814e85d6SHong Zhang 1429814e85d6SHong Zhang Input Parameter: 1430bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1431814e85d6SHong Zhang 1432814e85d6SHong Zhang Level: intermediate 1433814e85d6SHong Zhang 1434*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSetUp()`, `TSAdjointSolve()` 1435814e85d6SHong Zhang @*/ 1436d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointStep(TS ts) 1437d71ae5a4SJacob Faibussowitsch { 1438814e85d6SHong Zhang DM dm; 1439814e85d6SHong Zhang 1440814e85d6SHong Zhang PetscFunctionBegin; 1441814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 14429566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts, &dm)); 14439566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 14447b0e2f17SHong Zhang ts->steps--; /* must decrease the step index before the adjoint step is taken. */ 1445814e85d6SHong Zhang 1446814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1447814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 14489566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_AdjointStep, ts, 0, 0, 0)); 1449dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, adjointstep); 14509566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_AdjointStep, ts, 0, 0, 0)); 14517b0e2f17SHong Zhang ts->adjoint_steps++; 1452814e85d6SHong Zhang 1453814e85d6SHong Zhang if (ts->reason < 0) { 14543c633725SBarry Smith PetscCheck(!ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSAdjointStep has failed due to %s", TSConvergedReasons[ts->reason]); 1455814e85d6SHong Zhang } else if (!ts->reason) { 1456814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1457814e85d6SHong Zhang } 14583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1459814e85d6SHong Zhang } 1460814e85d6SHong Zhang 1461814e85d6SHong Zhang /*@ 1462814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1463814e85d6SHong Zhang 1464c3339decSBarry Smith Collective 1465bcf0153eSBarry Smith ` 1466814e85d6SHong Zhang Input Parameter: 1467bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1468814e85d6SHong Zhang 1469bcf0153eSBarry Smith Options Database Key: 1470814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1471814e85d6SHong Zhang 1472814e85d6SHong Zhang Level: intermediate 1473814e85d6SHong Zhang 1474814e85d6SHong Zhang Notes: 1475bcf0153eSBarry Smith This must be called after a call to `TSSolve()` that solves the forward problem 1476814e85d6SHong Zhang 1477bcf0153eSBarry Smith By default this will integrate back to the initial time, one can use `TSAdjointSetSteps()` to step back to a later time 1478814e85d6SHong Zhang 1479*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSCreate()`, `TSSetCostGradients()`, `TSSetSolution()`, `TSAdjointStep()` 1480814e85d6SHong Zhang @*/ 1481d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointSolve(TS ts) 1482d71ae5a4SJacob Faibussowitsch { 1483f1d62c27SHong Zhang static PetscBool cite = PETSC_FALSE; 14847f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 14857f59fb53SHong Zhang PetscLogStage adjoint_stage; 14867f59fb53SHong Zhang #endif 1487814e85d6SHong Zhang 1488814e85d6SHong Zhang PetscFunctionBegin; 1489814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1490421b783aSHong Zhang PetscCall(PetscCitationsRegister("@article{Zhang2022tsadjoint,\n" 1491421b783aSHong Zhang " title = {{PETSc TSAdjoint: A Discrete Adjoint ODE Solver for First-Order and Second-Order Sensitivity Analysis}},\n" 1492f1d62c27SHong Zhang " author = {Zhang, Hong and Constantinescu, Emil M. and Smith, Barry F.},\n" 1493421b783aSHong Zhang " journal = {SIAM Journal on Scientific Computing},\n" 1494421b783aSHong Zhang " volume = {44},\n" 1495421b783aSHong Zhang " number = {1},\n" 1496421b783aSHong Zhang " pages = {C1-C24},\n" 1497421b783aSHong Zhang " doi = {10.1137/21M140078X},\n" 14989371c9d4SSatish Balay " year = {2022}\n}\n", 14999371c9d4SSatish Balay &cite)); 15007f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 15019566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("TSAdjoint", &adjoint_stage)); 15029566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(adjoint_stage)); 15037f59fb53SHong Zhang #endif 15049566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 1505814e85d6SHong Zhang 1506814e85d6SHong Zhang /* reset time step and iteration counters */ 1507814e85d6SHong Zhang ts->adjoint_steps = 0; 1508814e85d6SHong Zhang ts->ksp_its = 0; 1509814e85d6SHong Zhang ts->snes_its = 0; 1510814e85d6SHong Zhang ts->num_snes_failures = 0; 1511814e85d6SHong Zhang ts->reject = 0; 1512814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1513814e85d6SHong Zhang 1514814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1515814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1516814e85d6SHong Zhang 1517814e85d6SHong Zhang while (!ts->reason) { 15189566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 15199566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 15209566063dSJacob Faibussowitsch PetscCall(TSAdjointEventHandler(ts)); 15219566063dSJacob Faibussowitsch PetscCall(TSAdjointStep(ts)); 152248a46eb9SPierre Jolivet if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) PetscCall(TSAdjointCostIntegral(ts)); 1523814e85d6SHong Zhang } 1524bdbff044SHong Zhang if (!ts->steps) { 15259566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory, ts, ts->steps, &ts->ptime)); 15269566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts, ts->steps, ts->ptime, ts->vec_sol, ts->numcost, ts->vecs_sensi, ts->vecs_sensip)); 1527bdbff044SHong Zhang } 1528814e85d6SHong Zhang ts->solvetime = ts->ptime; 15299566063dSJacob Faibussowitsch PetscCall(TSTrajectoryViewFromOptions(ts->trajectory, NULL, "-ts_trajectory_view")); 15309566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ts->vecs_sensi[0], (PetscObject)ts, "-ts_adjoint_view_solution")); 1531814e85d6SHong Zhang ts->adjoint_max_steps = 0; 15327f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 15339566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 15347f59fb53SHong Zhang #endif 15353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1536814e85d6SHong Zhang } 1537814e85d6SHong Zhang 1538814e85d6SHong Zhang /*@C 1539bcf0153eSBarry Smith TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using `TSAdjointMonitorSet()` 1540814e85d6SHong Zhang 1541c3339decSBarry Smith Collective 1542814e85d6SHong Zhang 1543814e85d6SHong Zhang Input Parameters: 1544bcf0153eSBarry Smith + ts - time stepping context obtained from `TSCreate()` 1545814e85d6SHong Zhang . step - step number that has just completed 1546814e85d6SHong Zhang . ptime - model time of the state 1547814e85d6SHong Zhang . u - state at the current model time 1548814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 1549814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1550814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1551814e85d6SHong Zhang 1552814e85d6SHong Zhang Level: developer 1553814e85d6SHong Zhang 1554bcf0153eSBarry Smith Note: 1555bcf0153eSBarry Smith `TSAdjointMonitor()` is typically used automatically within the time stepping implementations. 1556bcf0153eSBarry Smith Users would almost never call this routine directly. 1557bcf0153eSBarry Smith 1558bcf0153eSBarry Smith .seealso: `TSAdjointMonitorSet()`, `TSAdjointSolve()` 1559814e85d6SHong Zhang @*/ 1560d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointMonitor(TS ts, PetscInt step, PetscReal ptime, Vec u, PetscInt numcost, Vec *lambda, Vec *mu) 1561d71ae5a4SJacob Faibussowitsch { 1562814e85d6SHong Zhang PetscInt i, n = ts->numberadjointmonitors; 1563814e85d6SHong Zhang 1564814e85d6SHong Zhang PetscFunctionBegin; 1565814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1566814e85d6SHong Zhang PetscValidHeaderSpecific(u, VEC_CLASSID, 4); 15679566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(u)); 156848a46eb9SPierre Jolivet for (i = 0; i < n; i++) PetscCall((*ts->adjointmonitor[i])(ts, step, ptime, u, numcost, lambda, mu, ts->adjointmonitorcontext[i])); 15699566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(u)); 15703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1571814e85d6SHong Zhang } 1572814e85d6SHong Zhang 1573814e85d6SHong Zhang /*@ 1574814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1575814e85d6SHong Zhang 1576c3339decSBarry Smith Collective 1577814e85d6SHong Zhang 15784165533cSJose E. Roman Input Parameter: 1579814e85d6SHong Zhang . ts - time stepping context 1580814e85d6SHong Zhang 1581814e85d6SHong Zhang Level: advanced 1582814e85d6SHong Zhang 1583814e85d6SHong Zhang Notes: 1584bcf0153eSBarry Smith This function cannot be called until `TSAdjointStep()` has been completed. 1585814e85d6SHong Zhang 1586*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSAdjointSolve()`, `TSAdjointStep()` 1587814e85d6SHong Zhang @*/ 1588d71ae5a4SJacob Faibussowitsch PetscErrorCode TSAdjointCostIntegral(TS ts) 1589d71ae5a4SJacob Faibussowitsch { 1590362febeeSStefano Zampini PetscFunctionBegin; 1591814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1592dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, adjointintegral); 15933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1594814e85d6SHong Zhang } 1595814e85d6SHong Zhang 1596814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1597814e85d6SHong Zhang 1598814e85d6SHong Zhang /*@ 1599814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 1600814e85d6SHong Zhang of forward sensitivity analysis 1601814e85d6SHong Zhang 1602c3339decSBarry Smith Collective 1603814e85d6SHong Zhang 1604814e85d6SHong Zhang Input Parameter: 1605bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1606814e85d6SHong Zhang 1607814e85d6SHong Zhang Level: advanced 1608814e85d6SHong Zhang 1609*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSCreate()`, `TSDestroy()`, `TSSetUp()` 1610814e85d6SHong Zhang @*/ 1611d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetUp(TS ts) 1612d71ae5a4SJacob Faibussowitsch { 1613814e85d6SHong Zhang PetscFunctionBegin; 1614814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 16153ba16761SJacob Faibussowitsch if (ts->forwardsetupcalled) PetscFunctionReturn(PETSC_SUCCESS); 1616dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, forwardsetup); 16179566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol, &ts->vec_sensip_col)); 1618814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 16193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1620814e85d6SHong Zhang } 1621814e85d6SHong Zhang 1622814e85d6SHong Zhang /*@ 16237adebddeSHong Zhang TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 16247adebddeSHong Zhang 1625c3339decSBarry Smith Collective 16267adebddeSHong Zhang 16277adebddeSHong Zhang Input Parameter: 1628bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 16297adebddeSHong Zhang 16307adebddeSHong Zhang Level: advanced 16317adebddeSHong Zhang 1632*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreate()`, `TSDestroy()`, `TSForwardSetUp()` 16337adebddeSHong Zhang @*/ 1634d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardReset(TS ts) 1635d71ae5a4SJacob Faibussowitsch { 16367207e0fdSHong Zhang TS quadts = ts->quadraturets; 16377adebddeSHong Zhang 16387adebddeSHong Zhang PetscFunctionBegin; 16397adebddeSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1640dbbe0bcdSBarry Smith PetscTryTypeMethod(ts, forwardreset); 16419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 164248a46eb9SPierre Jolivet if (quadts) PetscCall(MatDestroy(&quadts->mat_sensip)); 16439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_sensip_col)); 16447adebddeSHong Zhang ts->forward_solve = PETSC_FALSE; 16457adebddeSHong Zhang ts->forwardsetupcalled = PETSC_FALSE; 16463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16477adebddeSHong Zhang } 16487adebddeSHong Zhang 16497adebddeSHong Zhang /*@ 1650814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1651814e85d6SHong Zhang 1652d8d19677SJose E. Roman Input Parameters: 1653bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1654814e85d6SHong Zhang . numfwdint - number of integrals 165567b8a455SSatish Balay - vp - the vectors containing the gradients for each integral w.r.t. parameters 1656814e85d6SHong Zhang 16577207e0fdSHong Zhang Level: deprecated 1658814e85d6SHong Zhang 1659*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1660814e85d6SHong Zhang @*/ 1661d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetIntegralGradients(TS ts, PetscInt numfwdint, Vec *vp) 1662d71ae5a4SJacob Faibussowitsch { 1663814e85d6SHong Zhang PetscFunctionBegin; 1664814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 16653c633725SBarry 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()"); 1666814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1667814e85d6SHong Zhang 1668814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 16693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1670814e85d6SHong Zhang } 1671814e85d6SHong Zhang 1672814e85d6SHong Zhang /*@ 1673814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities of the integral term. 1674814e85d6SHong Zhang 1675814e85d6SHong Zhang Input Parameter: 1676bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1677814e85d6SHong Zhang 1678814e85d6SHong Zhang Output Parameter: 167967b8a455SSatish Balay . vp - the vectors containing the gradients for each integral w.r.t. parameters 1680814e85d6SHong Zhang 16817207e0fdSHong Zhang Level: deprecated 1682814e85d6SHong Zhang 1683*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1684814e85d6SHong Zhang @*/ 1685d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetIntegralGradients(TS ts, PetscInt *numfwdint, Vec **vp) 1686d71ae5a4SJacob Faibussowitsch { 1687814e85d6SHong Zhang PetscFunctionBegin; 1688814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1689814e85d6SHong Zhang PetscValidPointer(vp, 3); 1690814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1691814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 16923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1693814e85d6SHong Zhang } 1694814e85d6SHong Zhang 1695814e85d6SHong Zhang /*@ 1696814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1697814e85d6SHong Zhang 1698c3339decSBarry Smith Collective 1699814e85d6SHong Zhang 17004165533cSJose E. Roman Input Parameter: 1701814e85d6SHong Zhang . ts - time stepping context 1702814e85d6SHong Zhang 1703814e85d6SHong Zhang Level: advanced 1704814e85d6SHong Zhang 1705814e85d6SHong Zhang Notes: 1706bcf0153eSBarry Smith This function cannot be called until `TSStep()` has been completed. 1707814e85d6SHong Zhang 1708*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardSetSensitivities()`, `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardSetUp()` 1709814e85d6SHong Zhang @*/ 1710d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardStep(TS ts) 1711d71ae5a4SJacob Faibussowitsch { 1712362febeeSStefano Zampini PetscFunctionBegin; 1713814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 17149566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_ForwardStep, ts, 0, 0, 0)); 1715dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, forwardstep); 17169566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_ForwardStep, ts, 0, 0, 0)); 17173c633725SBarry Smith PetscCheck(ts->reason >= 0 || !ts->errorifstepfailed, PetscObjectComm((PetscObject)ts), PETSC_ERR_NOT_CONVERGED, "TSFowardStep has failed due to %s", TSConvergedReasons[ts->reason]); 17183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1719814e85d6SHong Zhang } 1720814e85d6SHong Zhang 1721814e85d6SHong Zhang /*@ 1722814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1723814e85d6SHong Zhang 1724c3339decSBarry Smith Logically Collective 1725814e85d6SHong Zhang 1726814e85d6SHong Zhang Input Parameters: 1727bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1728814e85d6SHong Zhang . nump - number of parameters 1729814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1730814e85d6SHong Zhang 1731814e85d6SHong Zhang Level: beginner 1732814e85d6SHong Zhang 1733814e85d6SHong Zhang Notes: 1734814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1735bcf0153eSBarry Smith This function turns on a flag to trigger `TSSolve()` to compute forward sensitivities automatically. 1736bcf0153eSBarry Smith You must call this function before `TSSolve()`. 1737814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1738814e85d6SHong Zhang 1739*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1740814e85d6SHong Zhang @*/ 1741d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetSensitivities(TS ts, PetscInt nump, Mat Smat) 1742d71ae5a4SJacob Faibussowitsch { 1743814e85d6SHong Zhang PetscFunctionBegin; 1744814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1745814e85d6SHong Zhang PetscValidHeaderSpecific(Smat, MAT_CLASSID, 3); 1746814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1747b5b4017aSHong Zhang if (nump == PETSC_DEFAULT) { 17489566063dSJacob Faibussowitsch PetscCall(MatGetSize(Smat, NULL, &ts->num_parameters)); 1749b5b4017aSHong Zhang } else ts->num_parameters = nump; 17509566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Smat)); 17519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 1752814e85d6SHong Zhang ts->mat_sensip = Smat; 17533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1754814e85d6SHong Zhang } 1755814e85d6SHong Zhang 1756814e85d6SHong Zhang /*@ 1757814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1758814e85d6SHong Zhang 1759bcf0153eSBarry Smith Not Collective, but Smat returned is parallel if ts is parallel 1760814e85d6SHong Zhang 1761d8d19677SJose E. Roman Output Parameters: 1762bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1763814e85d6SHong Zhang . nump - number of parameters 1764814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1765814e85d6SHong Zhang 1766814e85d6SHong Zhang Level: intermediate 1767814e85d6SHong Zhang 1768*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1769814e85d6SHong Zhang @*/ 1770d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetSensitivities(TS ts, PetscInt *nump, Mat *Smat) 1771d71ae5a4SJacob Faibussowitsch { 1772814e85d6SHong Zhang PetscFunctionBegin; 1773814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1774814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1775814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 17763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1777814e85d6SHong Zhang } 1778814e85d6SHong Zhang 1779814e85d6SHong Zhang /*@ 1780814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1781814e85d6SHong Zhang 1782c3339decSBarry Smith Collective 1783814e85d6SHong Zhang 17844165533cSJose E. Roman Input Parameter: 1785814e85d6SHong Zhang . ts - time stepping context 1786814e85d6SHong Zhang 1787814e85d6SHong Zhang Level: advanced 1788814e85d6SHong Zhang 1789bcf0153eSBarry Smith Note: 1790bcf0153eSBarry Smith This function cannot be called until `TSStep()` has been completed. 1791814e85d6SHong Zhang 1792*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSSolve()`, `TSAdjointCostIntegral()` 1793814e85d6SHong Zhang @*/ 1794d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardCostIntegral(TS ts) 1795d71ae5a4SJacob Faibussowitsch { 1796362febeeSStefano Zampini PetscFunctionBegin; 1797814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1798dbbe0bcdSBarry Smith PetscUseTypeMethod(ts, forwardintegral); 17993ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1800814e85d6SHong Zhang } 1801b5b4017aSHong Zhang 1802b5b4017aSHong Zhang /*@ 1803b5b4017aSHong Zhang TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1804b5b4017aSHong Zhang 1805c3339decSBarry Smith Collective 1806b5b4017aSHong Zhang 1807d8d19677SJose E. Roman Input Parameters: 1808bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1809b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition 1810b5b4017aSHong Zhang 1811b5b4017aSHong Zhang Level: intermediate 1812b5b4017aSHong Zhang 1813bcf0153eSBarry Smith Notes: 1814bcf0153eSBarry Smith `TSSolve()` allows users to pass the initial solution directly to `TS`. But the tangent linear variables cannot be initialized in this way. 1815bcf0153eSBarry Smith This function is used to set initial values for tangent linear variables. 1816b5b4017aSHong Zhang 1817*1cc06b55SBarry Smith .seealso: [](ch_ts), `TS`, `TSForwardSetSensitivities()` 1818b5b4017aSHong Zhang @*/ 1819d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardSetInitialSensitivities(TS ts, Mat didp) 1820d71ae5a4SJacob Faibussowitsch { 1821362febeeSStefano Zampini PetscFunctionBegin; 1822b5b4017aSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1823b5b4017aSHong Zhang PetscValidHeaderSpecific(didp, MAT_CLASSID, 2); 182448a46eb9SPierre Jolivet if (!ts->mat_sensip) PetscCall(TSForwardSetSensitivities(ts, PETSC_DEFAULT, didp)); 18253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1826b5b4017aSHong Zhang } 18276affb6f8SHong Zhang 18286affb6f8SHong Zhang /*@ 18296affb6f8SHong Zhang TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 18306affb6f8SHong Zhang 18316affb6f8SHong Zhang Input Parameter: 1832bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 18336affb6f8SHong Zhang 18346affb6f8SHong Zhang Output Parameters: 1835cd4cee2dSHong Zhang + ns - number of stages 18366affb6f8SHong Zhang - S - tangent linear sensitivities at the intermediate stages 18376affb6f8SHong Zhang 18386affb6f8SHong Zhang Level: advanced 18396affb6f8SHong Zhang 1840bcf0153eSBarry Smith .seealso: `TS` 18416affb6f8SHong Zhang @*/ 1842d71ae5a4SJacob Faibussowitsch PetscErrorCode TSForwardGetStages(TS ts, PetscInt *ns, Mat **S) 1843d71ae5a4SJacob Faibussowitsch { 18446affb6f8SHong Zhang PetscFunctionBegin; 18456affb6f8SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 18466affb6f8SHong Zhang 18476affb6f8SHong Zhang if (!ts->ops->getstages) *S = NULL; 1848dbbe0bcdSBarry Smith else PetscUseTypeMethod(ts, forwardgetstages, ns, S); 18493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 18506affb6f8SHong Zhang } 1851cd4cee2dSHong Zhang 1852cd4cee2dSHong Zhang /*@ 1853bcf0153eSBarry Smith TSCreateQuadratureTS - Create a sub-`TS` that evaluates integrals over time 1854cd4cee2dSHong Zhang 1855d8d19677SJose E. Roman Input Parameters: 1856bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1857cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1858cd4cee2dSHong Zhang 18592fe279fdSBarry Smith Output Parameter: 1860bcf0153eSBarry Smith . quadts - the child `TS` context 1861cd4cee2dSHong Zhang 1862cd4cee2dSHong Zhang Level: intermediate 1863cd4cee2dSHong Zhang 1864*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSGetQuadratureTS()` 1865cd4cee2dSHong Zhang @*/ 1866d71ae5a4SJacob Faibussowitsch PetscErrorCode TSCreateQuadratureTS(TS ts, PetscBool fwd, TS *quadts) 1867d71ae5a4SJacob Faibussowitsch { 1868cd4cee2dSHong Zhang char prefix[128]; 1869cd4cee2dSHong Zhang 1870cd4cee2dSHong Zhang PetscFunctionBegin; 1871cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1872064a246eSJacob Faibussowitsch PetscValidPointer(quadts, 3); 18739566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts->quadraturets)); 18749566063dSJacob Faibussowitsch PetscCall(TSCreate(PetscObjectComm((PetscObject)ts), &ts->quadraturets)); 18759566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets, (PetscObject)ts, 1)); 18769566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix, sizeof(prefix), "%squad_", ((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "")); 18779566063dSJacob Faibussowitsch PetscCall(TSSetOptionsPrefix(ts->quadraturets, prefix)); 1878cd4cee2dSHong Zhang *quadts = ts->quadraturets; 1879cd4cee2dSHong Zhang 1880cd4cee2dSHong Zhang if (ts->numcost) { 18819566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, ts->numcost, &(*quadts)->vec_sol)); 1882cd4cee2dSHong Zhang } else { 18839566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF, 1, &(*quadts)->vec_sol)); 1884cd4cee2dSHong Zhang } 1885cd4cee2dSHong Zhang ts->costintegralfwd = fwd; 18863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1887cd4cee2dSHong Zhang } 1888cd4cee2dSHong Zhang 1889cd4cee2dSHong Zhang /*@ 1890bcf0153eSBarry Smith TSGetQuadratureTS - Return the sub-`TS` that evaluates integrals over time 1891cd4cee2dSHong Zhang 1892cd4cee2dSHong Zhang Input Parameter: 1893bcf0153eSBarry Smith . ts - the `TS` context obtained from `TSCreate()` 1894cd4cee2dSHong Zhang 1895cd4cee2dSHong Zhang Output Parameters: 1896cd4cee2dSHong Zhang + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1897bcf0153eSBarry Smith - quadts - the child `TS` context 1898cd4cee2dSHong Zhang 1899cd4cee2dSHong Zhang Level: intermediate 1900cd4cee2dSHong Zhang 1901*1cc06b55SBarry Smith .seealso: [](ch_ts), `TSCreateQuadratureTS()` 1902cd4cee2dSHong Zhang @*/ 1903d71ae5a4SJacob Faibussowitsch PetscErrorCode TSGetQuadratureTS(TS ts, PetscBool *fwd, TS *quadts) 1904d71ae5a4SJacob Faibussowitsch { 1905cd4cee2dSHong Zhang PetscFunctionBegin; 1906cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID, 1); 1907cd4cee2dSHong Zhang if (fwd) *fwd = ts->costintegralfwd; 1908cd4cee2dSHong Zhang if (quadts) *quadts = ts->quadraturets; 19093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1910cd4cee2dSHong Zhang } 1911ba0675f6SHong Zhang 1912ba0675f6SHong Zhang /*@ 1913bcf0153eSBarry Smith TSComputeSNESJacobian - Compute the Jacobian needed for the `SNESSolve()` in `TS` 1914bcf0153eSBarry Smith 1915c3339decSBarry Smith Collective 1916ba0675f6SHong Zhang 1917ba0675f6SHong Zhang Input Parameters: 1918bcf0153eSBarry Smith + ts - the `TS` context obtained from `TSCreate()` 1919ba0675f6SHong Zhang - x - state vector 1920ba0675f6SHong Zhang 1921ba0675f6SHong Zhang Output Parameters: 1922ba0675f6SHong Zhang + J - Jacobian matrix 1923ba0675f6SHong Zhang - Jpre - preconditioning matrix for J (may be same as J) 1924ba0675f6SHong Zhang 1925ba0675f6SHong Zhang Level: developer 1926ba0675f6SHong Zhang 1927bcf0153eSBarry Smith Note: 1928bcf0153eSBarry Smith Uses finite differencing when `TS` Jacobian is not available. 1929bcf0153eSBarry Smith 1930bcf0153eSBarry Smith .seealso: `SNES`, `TS`, `SNESSetJacobian()`, TSSetRHSJacobian()`, `TSSetIJacobian()` 1931ba0675f6SHong Zhang @*/ 1932d71ae5a4SJacob Faibussowitsch PetscErrorCode TSComputeSNESJacobian(TS ts, Vec x, Mat J, Mat Jpre) 1933d71ae5a4SJacob Faibussowitsch { 1934ba0675f6SHong Zhang SNES snes = ts->snes; 1935ba0675f6SHong Zhang PetscErrorCode (*jac)(SNES, Vec, Mat, Mat, void *) = NULL; 1936ba0675f6SHong Zhang 1937ba0675f6SHong Zhang PetscFunctionBegin; 1938ba0675f6SHong Zhang /* 1939ba0675f6SHong Zhang Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object 1940ba0675f6SHong Zhang because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for 1941da81f932SPierre Jolivet explicit methods. Instead, we check the Jacobian compute function directly to determine if FD 1942ba0675f6SHong Zhang coloring is used. 1943ba0675f6SHong Zhang */ 19449566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes, NULL, NULL, &jac, NULL)); 1945ba0675f6SHong Zhang if (jac == SNESComputeJacobianDefaultColor) { 1946ba0675f6SHong Zhang Vec f; 19479566063dSJacob Faibussowitsch PetscCall(SNESSetSolution(snes, x)); 19489566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes, &f, NULL, NULL)); 1949ba0675f6SHong Zhang /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */ 19509566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, x, f)); 1951ba0675f6SHong Zhang } 19529566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes, x, J, Jpre)); 19533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1954ba0675f6SHong Zhang } 1955