1814e85d6SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2814e85d6SHong Zhang #include <petscdraw.h> 3814e85d6SHong Zhang 433f72c5dSHong Zhang PetscLogEvent TS_AdjointStep,TS_ForwardStep,TS_JacobianPEval; 5814e85d6SHong Zhang 67f59fb53SHong Zhang /* #define TSADJOINT_STAGE */ 77f59fb53SHong Zhang 8814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/ 9814e85d6SHong Zhang 10814e85d6SHong Zhang /*@C 11b5b4017aSHong Zhang TSSetRHSJacobianP - Sets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix. 12814e85d6SHong Zhang 13814e85d6SHong Zhang Logically Collective on TS 14814e85d6SHong Zhang 15814e85d6SHong Zhang Input Parameters: 16b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 17b5b4017aSHong Zhang . Amat - JacobianP matrix 18b5b4017aSHong Zhang . func - function 19b5b4017aSHong Zhang - ctx - [optional] user-defined function context 20814e85d6SHong Zhang 21814e85d6SHong Zhang Calling sequence of func: 22814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 23814e85d6SHong Zhang + t - current timestep 24c9ad9de0SHong Zhang . U - input vector (current ODE solution) 25814e85d6SHong Zhang . A - output matrix 26814e85d6SHong Zhang - ctx - [optional] user-defined function context 27814e85d6SHong Zhang 28814e85d6SHong Zhang Level: intermediate 29814e85d6SHong Zhang 30814e85d6SHong Zhang Notes: 31814e85d6SHong Zhang Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 32814e85d6SHong Zhang 33db781477SPatrick Sanan .seealso: `TSGetRHSJacobianP()` 34814e85d6SHong Zhang @*/ 35814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 36814e85d6SHong Zhang { 37814e85d6SHong Zhang PetscFunctionBegin; 38814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 39814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 40814e85d6SHong Zhang 41814e85d6SHong Zhang ts->rhsjacobianp = func; 42814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 43814e85d6SHong Zhang if (Amat) { 449566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 459566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacprhs)); 4633f72c5dSHong Zhang ts->Jacprhs = Amat; 47814e85d6SHong Zhang } 48814e85d6SHong Zhang PetscFunctionReturn(0); 49814e85d6SHong Zhang } 50814e85d6SHong Zhang 51814e85d6SHong Zhang /*@C 52cd4cee2dSHong Zhang TSGetRHSJacobianP - Gets the function that computes the Jacobian of G w.r.t. the parameters P where U_t = G(U,P,t), as well as the location to store the matrix. 53cd4cee2dSHong Zhang 54cd4cee2dSHong Zhang Logically Collective on TS 55cd4cee2dSHong Zhang 56f899ff85SJose E. Roman Input Parameter: 57cd4cee2dSHong Zhang . ts - TS context obtained from TSCreate() 58cd4cee2dSHong Zhang 59cd4cee2dSHong Zhang Output Parameters: 60cd4cee2dSHong Zhang + Amat - JacobianP matrix 61cd4cee2dSHong Zhang . func - function 62cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 63cd4cee2dSHong Zhang 64cd4cee2dSHong Zhang Calling sequence of func: 65cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 66cd4cee2dSHong Zhang + t - current timestep 67cd4cee2dSHong Zhang . U - input vector (current ODE solution) 68cd4cee2dSHong Zhang . A - output matrix 69cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 70cd4cee2dSHong Zhang 71cd4cee2dSHong Zhang Level: intermediate 72cd4cee2dSHong Zhang 73cd4cee2dSHong Zhang Notes: 74cd4cee2dSHong Zhang Amat has the same number of rows and the same row parallel layout as u, Amat has the same number of columns and parallel layout as p 75cd4cee2dSHong Zhang 76db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()` 77cd4cee2dSHong Zhang @*/ 78cd4cee2dSHong Zhang PetscErrorCode TSGetRHSJacobianP(TS ts,Mat *Amat,PetscErrorCode (**func)(TS,PetscReal,Vec,Mat,void*),void **ctx) 79cd4cee2dSHong Zhang { 80cd4cee2dSHong Zhang PetscFunctionBegin; 81cd4cee2dSHong Zhang if (func) *func = ts->rhsjacobianp; 82cd4cee2dSHong Zhang if (ctx) *ctx = ts->rhsjacobianpctx; 83cd4cee2dSHong Zhang if (Amat) *Amat = ts->Jacprhs; 84cd4cee2dSHong Zhang PetscFunctionReturn(0); 85cd4cee2dSHong Zhang } 86cd4cee2dSHong Zhang 87cd4cee2dSHong Zhang /*@C 88814e85d6SHong Zhang TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 89814e85d6SHong Zhang 90814e85d6SHong Zhang Collective on TS 91814e85d6SHong Zhang 92814e85d6SHong Zhang Input Parameters: 93814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 94814e85d6SHong Zhang 95814e85d6SHong Zhang Level: developer 96814e85d6SHong Zhang 97db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()` 98814e85d6SHong Zhang @*/ 99c9ad9de0SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat) 100814e85d6SHong Zhang { 101814e85d6SHong Zhang PetscFunctionBegin; 10233f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 103814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 104c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 105814e85d6SHong Zhang 106814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 1079566063dSJacob Faibussowitsch PetscCall((*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx)); 108814e85d6SHong Zhang PetscStackPop; 109814e85d6SHong Zhang PetscFunctionReturn(0); 110814e85d6SHong Zhang } 111814e85d6SHong Zhang 112814e85d6SHong Zhang /*@C 11333f72c5dSHong 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. 11433f72c5dSHong Zhang 11533f72c5dSHong Zhang Logically Collective on TS 11633f72c5dSHong Zhang 11733f72c5dSHong Zhang Input Parameters: 11833f72c5dSHong Zhang + ts - TS context obtained from TSCreate() 11933f72c5dSHong Zhang . Amat - JacobianP matrix 12033f72c5dSHong Zhang . func - function 12133f72c5dSHong Zhang - ctx - [optional] user-defined function context 12233f72c5dSHong Zhang 12333f72c5dSHong Zhang Calling sequence of func: 12433f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 12533f72c5dSHong Zhang + t - current timestep 12633f72c5dSHong Zhang . U - input vector (current ODE solution) 12733f72c5dSHong Zhang . Udot - time derivative of state vector 12833f72c5dSHong Zhang . shift - shift to apply, see note below 12933f72c5dSHong Zhang . A - output matrix 13033f72c5dSHong Zhang - ctx - [optional] user-defined function context 13133f72c5dSHong Zhang 13233f72c5dSHong Zhang Level: intermediate 13333f72c5dSHong Zhang 13433f72c5dSHong Zhang Notes: 13533f72c5dSHong 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 13633f72c5dSHong Zhang 13733f72c5dSHong Zhang .seealso: 13833f72c5dSHong Zhang @*/ 13933f72c5dSHong Zhang PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx) 14033f72c5dSHong Zhang { 14133f72c5dSHong Zhang PetscFunctionBegin; 14233f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 14333f72c5dSHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 14433f72c5dSHong Zhang 14533f72c5dSHong Zhang ts->ijacobianp = func; 14633f72c5dSHong Zhang ts->ijacobianpctx = ctx; 14733f72c5dSHong Zhang if (Amat) { 1489566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 1499566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 15033f72c5dSHong Zhang ts->Jacp = Amat; 15133f72c5dSHong Zhang } 15233f72c5dSHong Zhang PetscFunctionReturn(0); 15333f72c5dSHong Zhang } 15433f72c5dSHong Zhang 15533f72c5dSHong Zhang /*@C 15633f72c5dSHong Zhang TSComputeIJacobianP - Runs the user-defined IJacobianP function. 15733f72c5dSHong Zhang 15833f72c5dSHong Zhang Collective on TS 15933f72c5dSHong Zhang 16033f72c5dSHong Zhang Input Parameters: 16133f72c5dSHong Zhang + ts - the TS context 16233f72c5dSHong Zhang . t - current timestep 16333f72c5dSHong Zhang . U - state vector 16433f72c5dSHong Zhang . Udot - time derivative of state vector 16533f72c5dSHong Zhang . shift - shift to apply, see note below 16633f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 16733f72c5dSHong Zhang 16833f72c5dSHong Zhang Output Parameters: 16933f72c5dSHong Zhang . A - Jacobian matrix 17033f72c5dSHong Zhang 17133f72c5dSHong Zhang Level: developer 17233f72c5dSHong Zhang 173db781477SPatrick Sanan .seealso: `TSSetIJacobianP()` 17433f72c5dSHong Zhang @*/ 17533f72c5dSHong Zhang PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex) 17633f72c5dSHong Zhang { 17733f72c5dSHong Zhang PetscFunctionBegin; 17833f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 17933f72c5dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 18033f72c5dSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 18133f72c5dSHong Zhang PetscValidHeaderSpecific(Udot,VEC_CLASSID,4); 18233f72c5dSHong Zhang 1839566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0)); 18433f72c5dSHong Zhang if (ts->ijacobianp) { 18533f72c5dSHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 1869566063dSJacob Faibussowitsch PetscCall((*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx)); 18733f72c5dSHong Zhang PetscStackPop; 18833f72c5dSHong Zhang } 18933f72c5dSHong Zhang if (imex) { 19033f72c5dSHong Zhang if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */ 19133f72c5dSHong Zhang PetscBool assembled; 1929566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 1939566063dSJacob Faibussowitsch PetscCall(MatAssembled(Amat,&assembled)); 19433f72c5dSHong Zhang if (!assembled) { 1959566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY)); 1969566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY)); 19733f72c5dSHong Zhang } 19833f72c5dSHong Zhang } 19933f72c5dSHong Zhang } else { 20033f72c5dSHong Zhang if (ts->rhsjacobianp) { 2019566063dSJacob Faibussowitsch PetscCall(TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs)); 20233f72c5dSHong Zhang } 20333f72c5dSHong Zhang if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */ 2049566063dSJacob Faibussowitsch PetscCall(MatScale(Amat,-1)); 20533f72c5dSHong Zhang } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */ 20633f72c5dSHong Zhang MatStructure axpy = DIFFERENT_NONZERO_PATTERN; 20733f72c5dSHong Zhang if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */ 2089566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Amat)); 20933f72c5dSHong Zhang } 2109566063dSJacob Faibussowitsch PetscCall(MatAXPY(Amat,-1,ts->Jacprhs,axpy)); 21133f72c5dSHong Zhang } 21233f72c5dSHong Zhang } 2139566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0)); 21433f72c5dSHong Zhang PetscFunctionReturn(0); 21533f72c5dSHong Zhang } 21633f72c5dSHong Zhang 21733f72c5dSHong Zhang /*@C 218814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 219814e85d6SHong Zhang 220814e85d6SHong Zhang Logically Collective on TS 221814e85d6SHong Zhang 222814e85d6SHong Zhang Input Parameters: 223814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 224814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 225814e85d6SHong Zhang . costintegral - vector that stores the integral values 226814e85d6SHong Zhang . rf - routine for evaluating the integrand function 227c9ad9de0SHong Zhang . drduf - function that computes the gradients of the r's with respect to u 228814e85d6SHong Zhang . drdpf - function that computes the gradients of the r's with respect to p, can be NULL if parametric sensitivity is not desired (mu=NULL) 229814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 230814e85d6SHong Zhang - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 231814e85d6SHong Zhang 232814e85d6SHong Zhang Calling sequence of rf: 233c9ad9de0SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx); 234814e85d6SHong Zhang 235c9ad9de0SHong Zhang Calling sequence of drduf: 236c9ad9de0SHong Zhang $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx); 237814e85d6SHong Zhang 238814e85d6SHong Zhang Calling sequence of drdpf: 239c9ad9de0SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx); 240814e85d6SHong Zhang 241cd4cee2dSHong Zhang Level: deprecated 242814e85d6SHong Zhang 243814e85d6SHong Zhang Notes: 244814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 245814e85d6SHong Zhang 246db781477SPatrick Sanan .seealso: `TSSetRHSJacobianP()`, `TSGetCostGradients()`, `TSSetCostGradients()` 247814e85d6SHong Zhang @*/ 248814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*), 249c9ad9de0SHong Zhang PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*), 250814e85d6SHong Zhang PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*), 251814e85d6SHong Zhang PetscBool fwd,void *ctx) 252814e85d6SHong Zhang { 253814e85d6SHong Zhang PetscFunctionBegin; 254814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 255814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3); 2563c633725SBarry 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()"); 257814e85d6SHong Zhang if (!ts->numcost) ts->numcost=numcost; 258814e85d6SHong Zhang 259814e85d6SHong Zhang if (costintegral) { 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)costintegral)); 2619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_costintegral)); 262814e85d6SHong Zhang ts->vec_costintegral = costintegral; 263814e85d6SHong Zhang } else { 264814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 2659566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral)); 266814e85d6SHong Zhang } else { 2679566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegral,0.0)); 268814e85d6SHong Zhang } 269814e85d6SHong Zhang } 270814e85d6SHong Zhang if (!ts->vec_costintegrand) { 2719566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand)); 272814e85d6SHong Zhang } else { 2739566063dSJacob Faibussowitsch PetscCall(VecSet(ts->vec_costintegrand,0.0)); 274814e85d6SHong Zhang } 275814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 276814e85d6SHong Zhang ts->costintegrand = rf; 277814e85d6SHong Zhang ts->costintegrandctx = ctx; 278c9ad9de0SHong Zhang ts->drdufunction = drduf; 279814e85d6SHong Zhang ts->drdpfunction = drdpf; 280814e85d6SHong Zhang PetscFunctionReturn(0); 281814e85d6SHong Zhang } 282814e85d6SHong Zhang 283b5b4017aSHong Zhang /*@C 284814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 285814e85d6SHong Zhang It is valid to call the routine after a backward run. 286814e85d6SHong Zhang 287814e85d6SHong Zhang Not Collective 288814e85d6SHong Zhang 289814e85d6SHong Zhang Input Parameter: 290814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 291814e85d6SHong Zhang 292814e85d6SHong Zhang Output Parameter: 293814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 294814e85d6SHong Zhang 295814e85d6SHong Zhang Level: intermediate 296814e85d6SHong Zhang 297db781477SPatrick Sanan .seealso: `TSSetCostIntegrand()` 298814e85d6SHong Zhang 299814e85d6SHong Zhang @*/ 300814e85d6SHong Zhang PetscErrorCode TSGetCostIntegral(TS ts,Vec *v) 301814e85d6SHong Zhang { 302cd4cee2dSHong Zhang TS quadts; 303cd4cee2dSHong Zhang 304814e85d6SHong Zhang PetscFunctionBegin; 305814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 306814e85d6SHong Zhang PetscValidPointer(v,2); 3079566063dSJacob Faibussowitsch PetscCall(TSGetQuadratureTS(ts,NULL,&quadts)); 308cd4cee2dSHong Zhang *v = quadts->vec_sol; 309814e85d6SHong Zhang PetscFunctionReturn(0); 310814e85d6SHong Zhang } 311814e85d6SHong Zhang 312b5b4017aSHong Zhang /*@C 313814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 314814e85d6SHong Zhang 315814e85d6SHong Zhang Input Parameters: 316814e85d6SHong Zhang + ts - the TS context 317814e85d6SHong Zhang . t - current time 318c9ad9de0SHong Zhang - U - state vector, i.e. current solution 319814e85d6SHong Zhang 320814e85d6SHong Zhang Output Parameter: 321c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs 322814e85d6SHong Zhang 323369cf35fSHong Zhang Notes: 324814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 325814e85d6SHong Zhang is used internally within the sensitivity analysis context. 326814e85d6SHong Zhang 327cd4cee2dSHong Zhang Level: deprecated 328814e85d6SHong Zhang 329db781477SPatrick Sanan .seealso: `TSSetCostIntegrand()` 330814e85d6SHong Zhang @*/ 331c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q) 332814e85d6SHong Zhang { 333814e85d6SHong Zhang PetscFunctionBegin; 334814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 335c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 336c9ad9de0SHong Zhang PetscValidHeaderSpecific(Q,VEC_CLASSID,4); 337814e85d6SHong Zhang 3389566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0)); 339814e85d6SHong Zhang if (ts->costintegrand) { 340814e85d6SHong Zhang PetscStackPush("TS user integrand in the cost function"); 3419566063dSJacob Faibussowitsch PetscCall((*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx)); 342814e85d6SHong Zhang PetscStackPop; 343814e85d6SHong Zhang } else { 3449566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(Q)); 345814e85d6SHong Zhang } 346814e85d6SHong Zhang 3479566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0)); 348814e85d6SHong Zhang PetscFunctionReturn(0); 349814e85d6SHong Zhang } 350814e85d6SHong Zhang 351b5b4017aSHong Zhang /*@C 352d76be551SHong Zhang TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 353814e85d6SHong Zhang 354d76be551SHong Zhang Level: deprecated 355814e85d6SHong Zhang 356814e85d6SHong Zhang @*/ 357c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 358814e85d6SHong Zhang { 359814e85d6SHong Zhang PetscFunctionBegin; 36033f72c5dSHong Zhang if (!DRDU) PetscFunctionReturn(0); 361814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 362c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 363814e85d6SHong Zhang 364c9ad9de0SHong Zhang PetscStackPush("TS user DRDU function for sensitivity analysis"); 3659566063dSJacob Faibussowitsch PetscCall((*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx)); 366814e85d6SHong Zhang PetscStackPop; 367814e85d6SHong Zhang PetscFunctionReturn(0); 368814e85d6SHong Zhang } 369814e85d6SHong Zhang 370b5b4017aSHong Zhang /*@C 371d76be551SHong Zhang TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 372814e85d6SHong Zhang 373d76be551SHong Zhang Level: deprecated 374814e85d6SHong Zhang 375814e85d6SHong Zhang @*/ 376c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 377814e85d6SHong Zhang { 378814e85d6SHong Zhang PetscFunctionBegin; 37933f72c5dSHong Zhang if (!DRDP) PetscFunctionReturn(0); 380814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 381c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 382814e85d6SHong Zhang 383814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 3849566063dSJacob Faibussowitsch PetscCall((*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx)); 385814e85d6SHong Zhang PetscStackPop; 386814e85d6SHong Zhang PetscFunctionReturn(0); 387814e85d6SHong Zhang } 388814e85d6SHong Zhang 389b5b4017aSHong Zhang /*@C 390b5b4017aSHong 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. 391b5b4017aSHong Zhang 392b5b4017aSHong Zhang Logically Collective on TS 393b5b4017aSHong Zhang 394b5b4017aSHong Zhang Input Parameters: 395b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 396b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 397b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 398b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 399b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 400b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 401b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 402b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 403f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for F_PP 404b5b4017aSHong Zhang 405b5b4017aSHong Zhang Calling sequence of ihessianproductfunc: 406369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 407b5b4017aSHong Zhang + t - current timestep 408b5b4017aSHong Zhang . U - input vector (current ODE solution) 409369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 410b5b4017aSHong Zhang . Vr - input vector to be right-multiplied with the Hessian 411369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 412b5b4017aSHong Zhang - ctx - [optional] user-defined function context 413b5b4017aSHong Zhang 414b5b4017aSHong Zhang Level: intermediate 415b5b4017aSHong Zhang 416369cf35fSHong Zhang Notes: 417369cf35fSHong Zhang The first Hessian function and the working array are required. 418369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 419369cf35fSHong Zhang $ Vl_n^T*F_UP*Vr 420369cf35fSHong 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. 421369cf35fSHong Zhang Each entry of F_UP corresponds to the derivative 422369cf35fSHong Zhang $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}. 423369cf35fSHong 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 424369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]} 425369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 426b5b4017aSHong Zhang 427b5b4017aSHong Zhang .seealso: 428b5b4017aSHong Zhang @*/ 429b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 430b5b4017aSHong Zhang Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 431b5b4017aSHong Zhang Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 432b5b4017aSHong Zhang Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 433b5b4017aSHong Zhang void *ctx) 434b5b4017aSHong Zhang { 435b5b4017aSHong Zhang PetscFunctionBegin; 436b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 437b5b4017aSHong Zhang PetscValidPointer(ihp1,2); 438b5b4017aSHong Zhang 439b5b4017aSHong Zhang ts->ihessianproductctx = ctx; 440b5b4017aSHong Zhang if (ihp1) ts->vecs_fuu = ihp1; 441b5b4017aSHong Zhang if (ihp2) ts->vecs_fup = ihp2; 442b5b4017aSHong Zhang if (ihp3) ts->vecs_fpu = ihp3; 443b5b4017aSHong Zhang if (ihp4) ts->vecs_fpp = ihp4; 444b5b4017aSHong Zhang ts->ihessianproduct_fuu = ihessianproductfunc1; 445b5b4017aSHong Zhang ts->ihessianproduct_fup = ihessianproductfunc2; 446b5b4017aSHong Zhang ts->ihessianproduct_fpu = ihessianproductfunc3; 447b5b4017aSHong Zhang ts->ihessianproduct_fpp = ihessianproductfunc4; 448b5b4017aSHong Zhang PetscFunctionReturn(0); 449b5b4017aSHong Zhang } 450b5b4017aSHong Zhang 451b5b4017aSHong Zhang /*@C 452b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu. 453b5b4017aSHong Zhang 454b5b4017aSHong Zhang Collective on TS 455b5b4017aSHong Zhang 456b5b4017aSHong Zhang Input Parameters: 457b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 458b5b4017aSHong Zhang 459b5b4017aSHong Zhang Notes: 460b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation, 461b5b4017aSHong Zhang so most users would not generally call this routine themselves. 462b5b4017aSHong Zhang 463b5b4017aSHong Zhang Level: developer 464b5b4017aSHong Zhang 465db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 466b5b4017aSHong Zhang @*/ 467b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 468b5b4017aSHong Zhang { 469b5b4017aSHong Zhang PetscFunctionBegin; 47033f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 471b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 472b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 473b5b4017aSHong Zhang 47467633408SHong Zhang if (ts->ihessianproduct_fuu) { 475b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis"); 4769566063dSJacob Faibussowitsch PetscCall((*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx)); 477b5b4017aSHong Zhang PetscStackPop; 47867633408SHong Zhang } 47967633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 48067633408SHong Zhang if (ts->rhshessianproduct_guu) { 48167633408SHong Zhang PetscInt nadj; 4829566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUU(ts,t,U,Vl,Vr,VHV)); 48367633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 4849566063dSJacob Faibussowitsch PetscCall(VecScale(VHV[nadj],-1)); 48567633408SHong Zhang } 48667633408SHong Zhang } 487b5b4017aSHong Zhang PetscFunctionReturn(0); 488b5b4017aSHong Zhang } 489b5b4017aSHong Zhang 490b5b4017aSHong Zhang /*@C 491b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup. 492b5b4017aSHong Zhang 493b5b4017aSHong Zhang Collective on TS 494b5b4017aSHong Zhang 495b5b4017aSHong Zhang Input Parameters: 496b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 497b5b4017aSHong Zhang 498b5b4017aSHong Zhang Notes: 499b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation, 500b5b4017aSHong Zhang so most users would not generally call this routine themselves. 501b5b4017aSHong Zhang 502b5b4017aSHong Zhang Level: developer 503b5b4017aSHong Zhang 504db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 505b5b4017aSHong Zhang @*/ 506b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 507b5b4017aSHong Zhang { 508b5b4017aSHong Zhang PetscFunctionBegin; 50933f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 510b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 511b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 512b5b4017aSHong Zhang 51367633408SHong Zhang if (ts->ihessianproduct_fup) { 514b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis"); 5159566063dSJacob Faibussowitsch PetscCall((*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx)); 516b5b4017aSHong Zhang PetscStackPop; 51767633408SHong Zhang } 51867633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 51967633408SHong Zhang if (ts->rhshessianproduct_gup) { 52067633408SHong Zhang PetscInt nadj; 5219566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionUP(ts,t,U,Vl,Vr,VHV)); 52267633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 5239566063dSJacob Faibussowitsch PetscCall(VecScale(VHV[nadj],-1)); 52467633408SHong Zhang } 52567633408SHong Zhang } 526b5b4017aSHong Zhang PetscFunctionReturn(0); 527b5b4017aSHong Zhang } 528b5b4017aSHong Zhang 529b5b4017aSHong Zhang /*@C 530b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu. 531b5b4017aSHong Zhang 532b5b4017aSHong Zhang Collective on TS 533b5b4017aSHong Zhang 534b5b4017aSHong Zhang Input Parameters: 535b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 536b5b4017aSHong Zhang 537b5b4017aSHong Zhang Notes: 538b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation, 539b5b4017aSHong Zhang so most users would not generally call this routine themselves. 540b5b4017aSHong Zhang 541b5b4017aSHong Zhang Level: developer 542b5b4017aSHong Zhang 543db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 544b5b4017aSHong Zhang @*/ 545b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 546b5b4017aSHong Zhang { 547b5b4017aSHong Zhang PetscFunctionBegin; 54833f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 549b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 550b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 551b5b4017aSHong Zhang 55267633408SHong Zhang if (ts->ihessianproduct_fpu) { 553b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 5549566063dSJacob Faibussowitsch PetscCall((*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx)); 555b5b4017aSHong Zhang PetscStackPop; 55667633408SHong Zhang } 55767633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 55867633408SHong Zhang if (ts->rhshessianproduct_gpu) { 55967633408SHong Zhang PetscInt nadj; 5609566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPU(ts,t,U,Vl,Vr,VHV)); 56167633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 5629566063dSJacob Faibussowitsch PetscCall(VecScale(VHV[nadj],-1)); 56367633408SHong Zhang } 56467633408SHong Zhang } 565b5b4017aSHong Zhang PetscFunctionReturn(0); 566b5b4017aSHong Zhang } 567b5b4017aSHong Zhang 568b5b4017aSHong Zhang /*@C 569b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp. 570b5b4017aSHong Zhang 571b5b4017aSHong Zhang Collective on TS 572b5b4017aSHong Zhang 573b5b4017aSHong Zhang Input Parameters: 574b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 575b5b4017aSHong Zhang 576b5b4017aSHong Zhang Notes: 577b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation, 578b5b4017aSHong Zhang so most users would not generally call this routine themselves. 579b5b4017aSHong Zhang 580b5b4017aSHong Zhang Level: developer 581b5b4017aSHong Zhang 582db781477SPatrick Sanan .seealso: `TSSetIHessianProduct()` 583b5b4017aSHong Zhang @*/ 584b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 585b5b4017aSHong Zhang { 586b5b4017aSHong Zhang PetscFunctionBegin; 58733f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 588b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 589b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 590b5b4017aSHong Zhang 59167633408SHong Zhang if (ts->ihessianproduct_fpp) { 592b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 5939566063dSJacob Faibussowitsch PetscCall((*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx)); 594b5b4017aSHong Zhang PetscStackPop; 59567633408SHong Zhang } 59667633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 59767633408SHong Zhang if (ts->rhshessianproduct_gpp) { 59867633408SHong Zhang PetscInt nadj; 5999566063dSJacob Faibussowitsch PetscCall(TSComputeRHSHessianProductFunctionPP(ts,t,U,Vl,Vr,VHV)); 60067633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 6019566063dSJacob Faibussowitsch PetscCall(VecScale(VHV[nadj],-1)); 60267633408SHong Zhang } 60367633408SHong Zhang } 604b5b4017aSHong Zhang PetscFunctionReturn(0); 605b5b4017aSHong Zhang } 606b5b4017aSHong Zhang 60713af1a74SHong Zhang /*@C 6084c500f23SPierre 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. 60913af1a74SHong Zhang 61013af1a74SHong Zhang Logically Collective on TS 61113af1a74SHong Zhang 61213af1a74SHong Zhang Input Parameters: 61313af1a74SHong Zhang + ts - TS context obtained from TSCreate() 61413af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 61513af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 61613af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 61713af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 61813af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 61913af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 62013af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 621f0fc11ceSJed Brown - hessianproductfunc4 - vector-Hessian-vector product function for G_PP 62213af1a74SHong Zhang 62313af1a74SHong Zhang Calling sequence of ihessianproductfunc: 624369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 62513af1a74SHong Zhang + t - current timestep 62613af1a74SHong Zhang . U - input vector (current ODE solution) 627369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 62813af1a74SHong Zhang . Vr - input vector to be right-multiplied with the Hessian 629369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 63013af1a74SHong Zhang - ctx - [optional] user-defined function context 63113af1a74SHong Zhang 63213af1a74SHong Zhang Level: intermediate 63313af1a74SHong Zhang 634369cf35fSHong Zhang Notes: 635369cf35fSHong Zhang The first Hessian function and the working array are required. 636369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 637369cf35fSHong Zhang $ Vl_n^T*G_UP*Vr 638369cf35fSHong 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. 639369cf35fSHong Zhang Each entry of G_UP corresponds to the derivative 640369cf35fSHong Zhang $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}. 641369cf35fSHong 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 642369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]} 643369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 64413af1a74SHong Zhang 64513af1a74SHong Zhang .seealso: 64613af1a74SHong Zhang @*/ 64713af1a74SHong Zhang PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 64813af1a74SHong Zhang Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 64913af1a74SHong Zhang Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 65013af1a74SHong Zhang Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 65113af1a74SHong Zhang void *ctx) 65213af1a74SHong Zhang { 65313af1a74SHong Zhang PetscFunctionBegin; 65413af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 65513af1a74SHong Zhang PetscValidPointer(rhshp1,2); 65613af1a74SHong Zhang 65713af1a74SHong Zhang ts->rhshessianproductctx = ctx; 65813af1a74SHong Zhang if (rhshp1) ts->vecs_guu = rhshp1; 65913af1a74SHong Zhang if (rhshp2) ts->vecs_gup = rhshp2; 66013af1a74SHong Zhang if (rhshp3) ts->vecs_gpu = rhshp3; 66113af1a74SHong Zhang if (rhshp4) ts->vecs_gpp = rhshp4; 66213af1a74SHong Zhang ts->rhshessianproduct_guu = rhshessianproductfunc1; 66313af1a74SHong Zhang ts->rhshessianproduct_gup = rhshessianproductfunc2; 66413af1a74SHong Zhang ts->rhshessianproduct_gpu = rhshessianproductfunc3; 66513af1a74SHong Zhang ts->rhshessianproduct_gpp = rhshessianproductfunc4; 66613af1a74SHong Zhang PetscFunctionReturn(0); 66713af1a74SHong Zhang } 66813af1a74SHong Zhang 66913af1a74SHong Zhang /*@C 670b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu. 67113af1a74SHong Zhang 67213af1a74SHong Zhang Collective on TS 67313af1a74SHong Zhang 67413af1a74SHong Zhang Input Parameters: 67513af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 67613af1a74SHong Zhang 67713af1a74SHong Zhang Notes: 678b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation, 67913af1a74SHong Zhang so most users would not generally call this routine themselves. 68013af1a74SHong Zhang 68113af1a74SHong Zhang Level: developer 68213af1a74SHong Zhang 683db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 68413af1a74SHong Zhang @*/ 685b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 68613af1a74SHong Zhang { 68713af1a74SHong Zhang PetscFunctionBegin; 68813af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 68913af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 69013af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 69113af1a74SHong Zhang 69213af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis"); 6939566063dSJacob Faibussowitsch PetscCall((*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx)); 69413af1a74SHong Zhang PetscStackPop; 69513af1a74SHong Zhang PetscFunctionReturn(0); 69613af1a74SHong Zhang } 69713af1a74SHong Zhang 69813af1a74SHong Zhang /*@C 699b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup. 70013af1a74SHong Zhang 70113af1a74SHong Zhang Collective on TS 70213af1a74SHong Zhang 70313af1a74SHong Zhang Input Parameters: 70413af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 70513af1a74SHong Zhang 70613af1a74SHong Zhang Notes: 707b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation, 70813af1a74SHong Zhang so most users would not generally call this routine themselves. 70913af1a74SHong Zhang 71013af1a74SHong Zhang Level: developer 71113af1a74SHong Zhang 712db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 71313af1a74SHong Zhang @*/ 714b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 71513af1a74SHong Zhang { 71613af1a74SHong Zhang PetscFunctionBegin; 71713af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 71813af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 71913af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 72013af1a74SHong Zhang 72113af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis"); 7229566063dSJacob Faibussowitsch PetscCall((*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx)); 72313af1a74SHong Zhang PetscStackPop; 72413af1a74SHong Zhang PetscFunctionReturn(0); 72513af1a74SHong Zhang } 72613af1a74SHong Zhang 72713af1a74SHong Zhang /*@C 728b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu. 72913af1a74SHong Zhang 73013af1a74SHong Zhang Collective on TS 73113af1a74SHong Zhang 73213af1a74SHong Zhang Input Parameters: 73313af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 73413af1a74SHong Zhang 73513af1a74SHong Zhang Notes: 736b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation, 73713af1a74SHong Zhang so most users would not generally call this routine themselves. 73813af1a74SHong Zhang 73913af1a74SHong Zhang Level: developer 74013af1a74SHong Zhang 741db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 74213af1a74SHong Zhang @*/ 743b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 74413af1a74SHong Zhang { 74513af1a74SHong Zhang PetscFunctionBegin; 74613af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 74713af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 74813af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 74913af1a74SHong Zhang 75013af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis"); 7519566063dSJacob Faibussowitsch PetscCall((*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx)); 75213af1a74SHong Zhang PetscStackPop; 75313af1a74SHong Zhang PetscFunctionReturn(0); 75413af1a74SHong Zhang } 75513af1a74SHong Zhang 75613af1a74SHong Zhang /*@C 757b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp. 75813af1a74SHong Zhang 75913af1a74SHong Zhang Collective on TS 76013af1a74SHong Zhang 76113af1a74SHong Zhang Input Parameters: 76213af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 76313af1a74SHong Zhang 76413af1a74SHong Zhang Notes: 765b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation, 76613af1a74SHong Zhang so most users would not generally call this routine themselves. 76713af1a74SHong Zhang 76813af1a74SHong Zhang Level: developer 76913af1a74SHong Zhang 770db781477SPatrick Sanan .seealso: `TSSetRHSHessianProduct()` 77113af1a74SHong Zhang @*/ 772b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 77313af1a74SHong Zhang { 77413af1a74SHong Zhang PetscFunctionBegin; 77513af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 77613af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 77713af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 77813af1a74SHong Zhang 77913af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis"); 7809566063dSJacob Faibussowitsch PetscCall((*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx)); 78113af1a74SHong Zhang PetscStackPop; 78213af1a74SHong Zhang PetscFunctionReturn(0); 78313af1a74SHong Zhang } 78413af1a74SHong Zhang 785814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 786814e85d6SHong Zhang 787814e85d6SHong Zhang /*@ 788814e85d6SHong 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 789814e85d6SHong Zhang for use by the TSAdjoint routines. 790814e85d6SHong Zhang 791d083f849SBarry Smith Logically Collective on TS 792814e85d6SHong Zhang 793814e85d6SHong Zhang Input Parameters: 794814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 795d2fd7bfcSHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 796814e85d6SHong 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 797814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 798814e85d6SHong Zhang 799814e85d6SHong Zhang Level: beginner 800814e85d6SHong Zhang 801814e85d6SHong Zhang Notes: 802814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 803814e85d6SHong Zhang 804814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 805814e85d6SHong Zhang 806db781477SPatrick Sanan .seealso `TSGetCostGradients()` 807814e85d6SHong Zhang @*/ 808814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu) 809814e85d6SHong Zhang { 810814e85d6SHong Zhang PetscFunctionBegin; 811814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 812064a246eSJacob Faibussowitsch PetscValidPointer(lambda,3); 813814e85d6SHong Zhang ts->vecs_sensi = lambda; 814814e85d6SHong Zhang ts->vecs_sensip = mu; 8153c633725SBarry 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"); 816814e85d6SHong Zhang ts->numcost = numcost; 817814e85d6SHong Zhang PetscFunctionReturn(0); 818814e85d6SHong Zhang } 819814e85d6SHong Zhang 820814e85d6SHong Zhang /*@ 821814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 822814e85d6SHong Zhang 823814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 824814e85d6SHong Zhang 825814e85d6SHong Zhang Input Parameter: 826814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 827814e85d6SHong Zhang 828d8d19677SJose E. Roman Output Parameters: 8296b867d5aSJose E. Roman + numcost - size of returned arrays 8306b867d5aSJose E. Roman . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 831814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 832814e85d6SHong Zhang 833814e85d6SHong Zhang Level: intermediate 834814e85d6SHong Zhang 835db781477SPatrick Sanan .seealso: `TSSetCostGradients()` 836814e85d6SHong Zhang @*/ 837814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu) 838814e85d6SHong Zhang { 839814e85d6SHong Zhang PetscFunctionBegin; 840814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 841814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 842814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 843814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 844814e85d6SHong Zhang PetscFunctionReturn(0); 845814e85d6SHong Zhang } 846814e85d6SHong Zhang 847814e85d6SHong Zhang /*@ 848b5b4017aSHong 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 849b5b4017aSHong Zhang for use by the TSAdjoint routines. 850b5b4017aSHong Zhang 851d083f849SBarry Smith Logically Collective on TS 852b5b4017aSHong Zhang 853b5b4017aSHong Zhang Input Parameters: 854b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 855b5b4017aSHong Zhang . numcost - number of cost functions 856b5b4017aSHong 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 857b5b4017aSHong 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 858b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 859b5b4017aSHong Zhang 860b5b4017aSHong Zhang Level: beginner 861b5b4017aSHong Zhang 862b5b4017aSHong Zhang Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system 863b5b4017aSHong Zhang 864ecf68647SHong Zhang For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve(). 865b5b4017aSHong Zhang 866b5b4017aSHong Zhang After TSAdjointSolve() is called, the lamba2 and the mu2 will contain the computed second-order adjoint sensitivities, and can be used to produce Hessian-vector product (not the full Hessian matrix). Users must provide a direction vector; it is usually generated by an optimization solver. 867b5b4017aSHong Zhang 8683fca9621SHong Zhang Passing NULL for lambda2 disables the second-order calculation. 869db781477SPatrick Sanan .seealso: `TSAdjointSetForward()` 870b5b4017aSHong Zhang @*/ 871b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir) 872b5b4017aSHong Zhang { 873b5b4017aSHong Zhang PetscFunctionBegin; 874b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 8753c633725SBarry 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"); 876b5b4017aSHong Zhang ts->numcost = numcost; 877b5b4017aSHong Zhang ts->vecs_sensi2 = lambda2; 87833f72c5dSHong Zhang ts->vecs_sensi2p = mu2; 879b5b4017aSHong Zhang ts->vec_dir = dir; 880b5b4017aSHong Zhang PetscFunctionReturn(0); 881b5b4017aSHong Zhang } 882b5b4017aSHong Zhang 883b5b4017aSHong Zhang /*@ 884b5b4017aSHong Zhang TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve() 885b5b4017aSHong Zhang 886b5b4017aSHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 887b5b4017aSHong Zhang 888b5b4017aSHong Zhang Input Parameter: 889b5b4017aSHong Zhang . ts - the TS context obtained from TSCreate() 890b5b4017aSHong Zhang 891d8d19677SJose E. Roman Output Parameters: 892b5b4017aSHong Zhang + numcost - number of cost functions 893b5b4017aSHong 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 894b5b4017aSHong 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 895b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 896b5b4017aSHong Zhang 897b5b4017aSHong Zhang Level: intermediate 898b5b4017aSHong Zhang 899db781477SPatrick Sanan .seealso: `TSSetCostHessianProducts()` 900b5b4017aSHong Zhang @*/ 901b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir) 902b5b4017aSHong Zhang { 903b5b4017aSHong Zhang PetscFunctionBegin; 904b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 905b5b4017aSHong Zhang if (numcost) *numcost = ts->numcost; 906b5b4017aSHong Zhang if (lambda2) *lambda2 = ts->vecs_sensi2; 90733f72c5dSHong Zhang if (mu2) *mu2 = ts->vecs_sensi2p; 908b5b4017aSHong Zhang if (dir) *dir = ts->vec_dir; 909b5b4017aSHong Zhang PetscFunctionReturn(0); 910b5b4017aSHong Zhang } 911b5b4017aSHong Zhang 912b5b4017aSHong Zhang /*@ 913ecf68647SHong Zhang TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities 914b5b4017aSHong Zhang 915d083f849SBarry Smith Logically Collective on TS 916b5b4017aSHong Zhang 917b5b4017aSHong Zhang Input Parameters: 918b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 919b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters 920b5b4017aSHong Zhang 921b5b4017aSHong Zhang Level: intermediate 922b5b4017aSHong Zhang 923ecf68647SHong Zhang Notes: When computing sensitivies w.r.t. initial condition, set didp to NULL so that the solver will take it as an identity matrix mathematically. TSAdjoint does not reset the tangent linear solver automatically, TSAdjointResetForward() should be called to reset the tangent linear solver. 924b5b4017aSHong Zhang 925db781477SPatrick Sanan .seealso: `TSSetCostHessianProducts()`, `TSAdjointResetForward()` 926b5b4017aSHong Zhang @*/ 927ecf68647SHong Zhang PetscErrorCode TSAdjointSetForward(TS ts,Mat didp) 928b5b4017aSHong Zhang { 92933f72c5dSHong Zhang Mat A; 93033f72c5dSHong Zhang Vec sp; 93133f72c5dSHong Zhang PetscScalar *xarr; 93233f72c5dSHong Zhang PetscInt lsize; 933b5b4017aSHong Zhang 934b5b4017aSHong Zhang PetscFunctionBegin; 935b5b4017aSHong Zhang ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 9363c633725SBarry Smith PetscCheck(ts->vecs_sensi2,PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first"); 9373c633725SBarry Smith PetscCheck(ts->vec_dir,PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 93833f72c5dSHong Zhang /* create a single-column dense matrix */ 9399566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(ts->vec_sol,&lsize)); 9409566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)ts),lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A)); 94133f72c5dSHong Zhang 9429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&sp)); 9439566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumn(A,0,&xarr)); 9449566063dSJacob Faibussowitsch PetscCall(VecPlaceArray(sp,xarr)); 945ecf68647SHong Zhang if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */ 94633f72c5dSHong Zhang if (didp) { 9479566063dSJacob Faibussowitsch PetscCall(MatMult(didp,ts->vec_dir,sp)); 9489566063dSJacob Faibussowitsch PetscCall(VecScale(sp,2.)); 94933f72c5dSHong Zhang } else { 9509566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(sp)); 95133f72c5dSHong Zhang } 952ecf68647SHong Zhang } else { /* tangent linear variable initialized as dir */ 9539566063dSJacob Faibussowitsch PetscCall(VecCopy(ts->vec_dir,sp)); 95433f72c5dSHong Zhang } 9559566063dSJacob Faibussowitsch PetscCall(VecResetArray(sp)); 9569566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumn(A,&xarr)); 9579566063dSJacob Faibussowitsch PetscCall(VecDestroy(&sp)); 95833f72c5dSHong Zhang 9599566063dSJacob Faibussowitsch PetscCall(TSForwardSetInitialSensitivities(ts,A)); /* if didp is NULL, identity matrix is assumed */ 96033f72c5dSHong Zhang 9619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 962b5b4017aSHong Zhang PetscFunctionReturn(0); 963b5b4017aSHong Zhang } 964b5b4017aSHong Zhang 965b5b4017aSHong Zhang /*@ 966ecf68647SHong Zhang TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context 967ecf68647SHong Zhang 968d083f849SBarry Smith Logically Collective on TS 969ecf68647SHong Zhang 970ecf68647SHong Zhang Input Parameters: 971ecf68647SHong Zhang . ts - the TS context obtained from TSCreate() 972ecf68647SHong Zhang 973ecf68647SHong Zhang Level: intermediate 974ecf68647SHong Zhang 975db781477SPatrick Sanan .seealso: `TSAdjointSetForward()` 976ecf68647SHong Zhang @*/ 977ecf68647SHong Zhang PetscErrorCode TSAdjointResetForward(TS ts) 978ecf68647SHong Zhang { 979ecf68647SHong Zhang PetscFunctionBegin; 980ecf68647SHong Zhang ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */ 9819566063dSJacob Faibussowitsch PetscCall(TSForwardReset(ts)); 982ecf68647SHong Zhang PetscFunctionReturn(0); 983ecf68647SHong Zhang } 984ecf68647SHong Zhang 985ecf68647SHong Zhang /*@ 986814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 987814e85d6SHong Zhang of an adjoint solver 988814e85d6SHong Zhang 989814e85d6SHong Zhang Collective on TS 990814e85d6SHong Zhang 991814e85d6SHong Zhang Input Parameter: 992814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 993814e85d6SHong Zhang 994814e85d6SHong Zhang Level: advanced 995814e85d6SHong Zhang 996db781477SPatrick Sanan .seealso: `TSCreate()`, `TSAdjointStep()`, `TSSetCostGradients()` 997814e85d6SHong Zhang @*/ 998814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts) 999814e85d6SHong Zhang { 1000881c1a9bSHong Zhang TSTrajectory tj; 1001881c1a9bSHong Zhang PetscBool match; 1002814e85d6SHong Zhang 1003814e85d6SHong Zhang PetscFunctionBegin; 1004814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1005814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 10063c633725SBarry Smith PetscCheck(ts->vecs_sensi,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 10073c633725SBarry Smith PetscCheck(!ts->vecs_sensip || ts->Jacp || ts->Jacprhs,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 10089566063dSJacob Faibussowitsch PetscCall(TSGetTrajectory(ts,&tj)); 10099566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match)); 1010881c1a9bSHong Zhang if (match) { 1011881c1a9bSHong Zhang PetscBool solution_only; 10129566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGetSolutionOnly(tj,&solution_only)); 10133c633725SBarry 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"); 1014881c1a9bSHong Zhang } 10159566063dSJacob Faibussowitsch PetscCall(TSTrajectorySetUseHistory(tj,PETSC_FALSE)); /* not use TSHistory */ 101633f72c5dSHong Zhang 1017cd4cee2dSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 10189566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col)); 1019814e85d6SHong Zhang if (ts->vecs_sensip) { 10209566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col)); 1021814e85d6SHong Zhang } 1022814e85d6SHong Zhang } 1023814e85d6SHong Zhang 1024814e85d6SHong Zhang if (ts->ops->adjointsetup) { 10259566063dSJacob Faibussowitsch PetscCall((*ts->ops->adjointsetup)(ts)); 1026814e85d6SHong Zhang } 1027814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 1028814e85d6SHong Zhang PetscFunctionReturn(0); 1029814e85d6SHong Zhang } 1030814e85d6SHong Zhang 1031814e85d6SHong Zhang /*@ 1032ece46509SHong Zhang TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 1033ece46509SHong Zhang 1034ece46509SHong Zhang Collective on TS 1035ece46509SHong Zhang 1036ece46509SHong Zhang Input Parameter: 1037ece46509SHong Zhang . ts - the TS context obtained from TSCreate() 1038ece46509SHong Zhang 1039ece46509SHong Zhang Level: beginner 1040ece46509SHong Zhang 1041db781477SPatrick Sanan .seealso: `TSCreate()`, `TSAdjointSetUp()`, `TSADestroy()` 1042ece46509SHong Zhang @*/ 1043ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts) 1044ece46509SHong Zhang { 1045ece46509SHong Zhang PetscFunctionBegin; 1046ece46509SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1047ece46509SHong Zhang if (ts->ops->adjointreset) { 10489566063dSJacob Faibussowitsch PetscCall((*ts->ops->adjointreset)(ts)); 1049ece46509SHong Zhang } 10507207e0fdSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 10519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_drdu_col)); 10527207e0fdSHong Zhang if (ts->vecs_sensip) { 10539566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_drdp_col)); 10547207e0fdSHong Zhang } 10557207e0fdSHong Zhang } 1056ece46509SHong Zhang ts->vecs_sensi = NULL; 1057ece46509SHong Zhang ts->vecs_sensip = NULL; 1058ece46509SHong Zhang ts->vecs_sensi2 = NULL; 105933f72c5dSHong Zhang ts->vecs_sensi2p = NULL; 1060ece46509SHong Zhang ts->vec_dir = NULL; 1061ece46509SHong Zhang ts->adjointsetupcalled = PETSC_FALSE; 1062ece46509SHong Zhang PetscFunctionReturn(0); 1063ece46509SHong Zhang } 1064ece46509SHong Zhang 1065ece46509SHong Zhang /*@ 1066814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1067814e85d6SHong Zhang 1068814e85d6SHong Zhang Logically Collective on TS 1069814e85d6SHong Zhang 1070814e85d6SHong Zhang Input Parameters: 1071814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1072a2b725a8SWilliam Gropp - steps - number of steps to use 1073814e85d6SHong Zhang 1074814e85d6SHong Zhang Level: intermediate 1075814e85d6SHong Zhang 1076814e85d6SHong Zhang Notes: 1077814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 1078814e85d6SHong Zhang so as to integrate back to less than the original timestep 1079814e85d6SHong Zhang 1080db781477SPatrick Sanan .seealso: `TSSetExactFinalTime()` 1081814e85d6SHong Zhang @*/ 1082814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 1083814e85d6SHong Zhang { 1084814e85d6SHong Zhang PetscFunctionBegin; 1085814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1086814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts,steps,2); 10873c633725SBarry Smith PetscCheck(steps >= 0,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 10883c633725SBarry Smith PetscCheck(steps <= ts->steps,PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 1089814e85d6SHong Zhang ts->adjoint_max_steps = steps; 1090814e85d6SHong Zhang PetscFunctionReturn(0); 1091814e85d6SHong Zhang } 1092814e85d6SHong Zhang 1093814e85d6SHong Zhang /*@C 1094814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1095814e85d6SHong Zhang 1096814e85d6SHong Zhang Level: deprecated 1097814e85d6SHong Zhang 1098814e85d6SHong Zhang @*/ 1099814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 1100814e85d6SHong Zhang { 1101814e85d6SHong Zhang PetscFunctionBegin; 1102814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1103814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 1104814e85d6SHong Zhang 1105814e85d6SHong Zhang ts->rhsjacobianp = func; 1106814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 1107814e85d6SHong Zhang if (Amat) { 11089566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Amat)); 11099566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->Jacp)); 1110814e85d6SHong Zhang ts->Jacp = Amat; 1111814e85d6SHong Zhang } 1112814e85d6SHong Zhang PetscFunctionReturn(0); 1113814e85d6SHong Zhang } 1114814e85d6SHong Zhang 1115814e85d6SHong Zhang /*@C 1116814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1117814e85d6SHong Zhang 1118814e85d6SHong Zhang Level: deprecated 1119814e85d6SHong Zhang 1120814e85d6SHong Zhang @*/ 1121c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat) 1122814e85d6SHong Zhang { 1123814e85d6SHong Zhang PetscFunctionBegin; 1124814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1125c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1126814e85d6SHong Zhang PetscValidPointer(Amat,4); 1127814e85d6SHong Zhang 1128814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 11299566063dSJacob Faibussowitsch PetscCall((*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx)); 1130814e85d6SHong Zhang PetscStackPop; 1131814e85d6SHong Zhang PetscFunctionReturn(0); 1132814e85d6SHong Zhang } 1133814e85d6SHong Zhang 1134814e85d6SHong Zhang /*@ 1135d76be551SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 1136814e85d6SHong Zhang 1137814e85d6SHong Zhang Level: deprecated 1138814e85d6SHong Zhang 1139814e85d6SHong Zhang @*/ 1140c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 1141814e85d6SHong Zhang { 1142814e85d6SHong Zhang PetscFunctionBegin; 1143814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1144c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1145814e85d6SHong Zhang 1146814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 11479566063dSJacob Faibussowitsch PetscCall((*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx)); 1148814e85d6SHong Zhang PetscStackPop; 1149814e85d6SHong Zhang PetscFunctionReturn(0); 1150814e85d6SHong Zhang } 1151814e85d6SHong Zhang 1152814e85d6SHong Zhang /*@ 1153d76be551SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 1154814e85d6SHong Zhang 1155814e85d6SHong Zhang Level: deprecated 1156814e85d6SHong Zhang 1157814e85d6SHong Zhang @*/ 1158c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 1159814e85d6SHong Zhang { 1160814e85d6SHong Zhang PetscFunctionBegin; 1161814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1162c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1163814e85d6SHong Zhang 1164814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 11659566063dSJacob Faibussowitsch PetscCall((*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx)); 1166814e85d6SHong Zhang PetscStackPop; 1167814e85d6SHong Zhang PetscFunctionReturn(0); 1168814e85d6SHong Zhang } 1169814e85d6SHong Zhang 1170814e85d6SHong Zhang /*@C 1171814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 1172814e85d6SHong Zhang 1173814e85d6SHong Zhang Level: intermediate 1174814e85d6SHong Zhang 1175db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1176814e85d6SHong Zhang @*/ 1177814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1178814e85d6SHong Zhang { 1179814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1180814e85d6SHong Zhang 1181814e85d6SHong Zhang PetscFunctionBegin; 1182064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,8); 11839566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,vf->format)); 11849566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0],viewer)); 11859566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1186814e85d6SHong Zhang PetscFunctionReturn(0); 1187814e85d6SHong Zhang } 1188814e85d6SHong Zhang 1189814e85d6SHong Zhang /*@C 1190814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1191814e85d6SHong Zhang 1192814e85d6SHong Zhang Collective on TS 1193814e85d6SHong Zhang 1194814e85d6SHong Zhang Input Parameters: 1195814e85d6SHong Zhang + ts - TS object you wish to monitor 1196814e85d6SHong Zhang . name - the monitor type one is seeking 1197814e85d6SHong Zhang . help - message indicating what monitoring is done 1198814e85d6SHong Zhang . manual - manual page for the monitor 1199814e85d6SHong Zhang . monitor - the monitor function 1200814e85d6SHong Zhang - monitorsetup - a function that is called once ONLY if the user selected this monitor that may set additional features of the TS or PetscViewer objects 1201814e85d6SHong Zhang 1202814e85d6SHong Zhang Level: developer 1203814e85d6SHong Zhang 1204db781477SPatrick Sanan .seealso: `PetscOptionsGetViewer()`, `PetscOptionsGetReal()`, `PetscOptionsHasName()`, `PetscOptionsGetString()`, 1205db781477SPatrick Sanan `PetscOptionsGetIntArray()`, `PetscOptionsGetRealArray()`, `PetscOptionsBool()` 1206db781477SPatrick Sanan `PetscOptionsInt()`, `PetscOptionsString()`, `PetscOptionsReal()`, `PetscOptionsBool()`, 1207db781477SPatrick Sanan `PetscOptionsName()`, `PetscOptionsBegin()`, `PetscOptionsEnd()`, `PetscOptionsHeadBegin()`, 1208*c2e3fba1SPatrick Sanan `PetscOptionsStringArray()`, `PetscOptionsRealArray()`, `PetscOptionsScalar()`, 1209db781477SPatrick Sanan `PetscOptionsBoolGroupBegin()`, `PetscOptionsBoolGroup()`, `PetscOptionsBoolGroupEnd()`, 1210db781477SPatrick Sanan `PetscOptionsFList()`, `PetscOptionsEList()` 1211814e85d6SHong Zhang @*/ 1212814e85d6SHong Zhang 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*)) 1213814e85d6SHong Zhang { 1214814e85d6SHong Zhang PetscViewer viewer; 1215814e85d6SHong Zhang PetscViewerFormat format; 1216814e85d6SHong Zhang PetscBool flg; 1217814e85d6SHong Zhang 1218814e85d6SHong Zhang PetscFunctionBegin; 12199566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg)); 1220814e85d6SHong Zhang if (flg) { 1221814e85d6SHong Zhang PetscViewerAndFormat *vf; 12229566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer,format,&vf)); 12239566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 1224814e85d6SHong Zhang if (monitorsetup) { 12259566063dSJacob Faibussowitsch PetscCall((*monitorsetup)(ts,vf)); 1226814e85d6SHong Zhang } 12279566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy)); 1228814e85d6SHong Zhang } 1229814e85d6SHong Zhang PetscFunctionReturn(0); 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 1236814e85d6SHong Zhang Logically Collective on TS 1237814e85d6SHong Zhang 1238814e85d6SHong Zhang Input Parameters: 1239814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1240814e85d6SHong Zhang . adjointmonitor - monitoring routine 1241814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 1242814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 1243814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 1244814e85d6SHong Zhang (may be NULL) 1245814e85d6SHong Zhang 1246814e85d6SHong Zhang Calling sequence of monitor: 1247814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1248814e85d6SHong Zhang 1249814e85d6SHong Zhang + ts - the TS context 1250814e85d6SHong 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 1251814e85d6SHong Zhang been interpolated to) 1252814e85d6SHong Zhang . time - current time 1253814e85d6SHong Zhang . u - current iterate 1254814e85d6SHong Zhang . numcost - number of cost functionos 1255814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1256814e85d6SHong Zhang . mu - sensitivities to parameters 1257814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 1258814e85d6SHong Zhang 1259814e85d6SHong Zhang Notes: 1260814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 1261814e85d6SHong Zhang already has been loaded. 1262814e85d6SHong Zhang 1263814e85d6SHong Zhang Fortran Notes: 1264814e85d6SHong Zhang Only a single monitor function can be set for each TS object 1265814e85d6SHong Zhang 1266814e85d6SHong Zhang Level: intermediate 1267814e85d6SHong Zhang 1268db781477SPatrick Sanan .seealso: `TSAdjointMonitorCancel()` 1269814e85d6SHong Zhang @*/ 1270814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 1271814e85d6SHong Zhang { 1272814e85d6SHong Zhang PetscInt i; 1273814e85d6SHong Zhang PetscBool identical; 1274814e85d6SHong Zhang 1275814e85d6SHong Zhang PetscFunctionBegin; 1276814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1277814e85d6SHong Zhang for (i=0; i<ts->numbermonitors;i++) { 12789566063dSJacob Faibussowitsch PetscCall(PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical)); 1279814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 1280814e85d6SHong Zhang } 12813c633725SBarry Smith PetscCheck(ts->numberadjointmonitors < MAXTSMONITORS,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 1282814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1283814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1284814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 1285814e85d6SHong Zhang PetscFunctionReturn(0); 1286814e85d6SHong Zhang } 1287814e85d6SHong Zhang 1288814e85d6SHong Zhang /*@C 1289814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1290814e85d6SHong Zhang 1291814e85d6SHong Zhang Logically Collective on TS 1292814e85d6SHong Zhang 1293814e85d6SHong Zhang Input Parameters: 1294814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1295814e85d6SHong Zhang 1296814e85d6SHong Zhang Notes: 1297814e85d6SHong Zhang There is no way to remove a single, specific monitor. 1298814e85d6SHong Zhang 1299814e85d6SHong Zhang Level: intermediate 1300814e85d6SHong Zhang 1301db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1302814e85d6SHong Zhang @*/ 1303814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts) 1304814e85d6SHong Zhang { 1305814e85d6SHong Zhang PetscInt i; 1306814e85d6SHong Zhang 1307814e85d6SHong Zhang PetscFunctionBegin; 1308814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1309814e85d6SHong Zhang for (i=0; i<ts->numberadjointmonitors; i++) { 1310814e85d6SHong Zhang if (ts->adjointmonitordestroy[i]) { 13119566063dSJacob Faibussowitsch PetscCall((*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i])); 1312814e85d6SHong Zhang } 1313814e85d6SHong Zhang } 1314814e85d6SHong Zhang ts->numberadjointmonitors = 0; 1315814e85d6SHong Zhang PetscFunctionReturn(0); 1316814e85d6SHong Zhang } 1317814e85d6SHong Zhang 1318814e85d6SHong Zhang /*@C 1319814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 1320814e85d6SHong Zhang 1321814e85d6SHong Zhang Level: intermediate 1322814e85d6SHong Zhang 1323db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()` 1324814e85d6SHong Zhang @*/ 1325814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1326814e85d6SHong Zhang { 1327814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1328814e85d6SHong Zhang 1329814e85d6SHong Zhang PetscFunctionBegin; 1330064a246eSJacob Faibussowitsch PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,8); 13319566063dSJacob Faibussowitsch PetscCall(PetscViewerPushFormat(viewer,vf->format)); 13329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel)); 133363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"%" PetscInt_FMT " TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n")); 13349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel)); 13359566063dSJacob Faibussowitsch PetscCall(PetscViewerPopFormat(viewer)); 1336814e85d6SHong Zhang PetscFunctionReturn(0); 1337814e85d6SHong Zhang } 1338814e85d6SHong Zhang 1339814e85d6SHong Zhang /*@C 1340814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1341814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 1342814e85d6SHong Zhang 1343814e85d6SHong Zhang Collective on TS 1344814e85d6SHong Zhang 1345814e85d6SHong Zhang Input Parameters: 1346814e85d6SHong Zhang + ts - the TS context 1347814e85d6SHong Zhang . step - current time-step 1348814e85d6SHong Zhang . ptime - current time 1349814e85d6SHong Zhang . u - current state 1350814e85d6SHong Zhang . numcost - number of cost functions 1351814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1352814e85d6SHong Zhang . mu - sensitivities to parameters 1353814e85d6SHong Zhang - dummy - either a viewer or NULL 1354814e85d6SHong Zhang 1355814e85d6SHong Zhang Level: intermediate 1356814e85d6SHong Zhang 1357db781477SPatrick Sanan .seealso: `TSAdjointMonitorSet()`, `TSAdjointMonitorDefault()`, `VecView()` 1358814e85d6SHong Zhang @*/ 1359814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 1360814e85d6SHong Zhang { 1361814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1362814e85d6SHong Zhang PetscDraw draw; 1363814e85d6SHong Zhang PetscReal xl,yl,xr,yr,h; 1364814e85d6SHong Zhang char time[32]; 1365814e85d6SHong Zhang 1366814e85d6SHong Zhang PetscFunctionBegin; 1367814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1368814e85d6SHong Zhang 13699566063dSJacob Faibussowitsch PetscCall(VecView(lambda[0],ictx->viewer)); 13709566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(ictx->viewer,0,&draw)); 13719566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime)); 13729566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr)); 1373814e85d6SHong Zhang h = yl + .95*(yr - yl); 13749566063dSJacob Faibussowitsch PetscCall(PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time)); 13759566063dSJacob Faibussowitsch PetscCall(PetscDrawFlush(draw)); 1376814e85d6SHong Zhang PetscFunctionReturn(0); 1377814e85d6SHong Zhang } 1378814e85d6SHong Zhang 1379814e85d6SHong Zhang /* 1380814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1381814e85d6SHong Zhang 1382814e85d6SHong Zhang Collective on TSAdjoint 1383814e85d6SHong Zhang 1384814e85d6SHong Zhang Input Parameter: 1385814e85d6SHong Zhang . ts - the TS context 1386814e85d6SHong Zhang 1387814e85d6SHong Zhang Options Database Keys: 1388814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1389814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 1390814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1391814e85d6SHong Zhang 1392814e85d6SHong Zhang Level: developer 1393814e85d6SHong Zhang 1394814e85d6SHong Zhang Notes: 1395814e85d6SHong Zhang This is not normally called directly by users 1396814e85d6SHong Zhang 1397db781477SPatrick Sanan .seealso: `TSSetSaveTrajectory()`, `TSTrajectorySetUp()` 1398814e85d6SHong Zhang */ 1399814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 1400814e85d6SHong Zhang { 1401814e85d6SHong Zhang PetscBool tflg,opt; 1402814e85d6SHong Zhang 1403814e85d6SHong Zhang PetscFunctionBegin; 1404814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,2); 1405d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"TS Adjoint options"); 1406814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 14079566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt)); 1408814e85d6SHong Zhang if (opt) { 14099566063dSJacob Faibussowitsch PetscCall(TSSetSaveTrajectory(ts)); 1410814e85d6SHong Zhang ts->adjoint_solve = tflg; 1411814e85d6SHong Zhang } 14129566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL)); 14139566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL)); 1414814e85d6SHong Zhang opt = PETSC_FALSE; 14159566063dSJacob Faibussowitsch PetscCall(PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt)); 1416814e85d6SHong Zhang if (opt) { 1417814e85d6SHong Zhang TSMonitorDrawCtx ctx; 1418814e85d6SHong Zhang PetscInt howoften = 1; 1419814e85d6SHong Zhang 14209566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL)); 14219566063dSJacob Faibussowitsch PetscCall(TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),NULL,NULL,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx)); 14229566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy)); 1423814e85d6SHong Zhang } 1424814e85d6SHong Zhang PetscFunctionReturn(0); 1425814e85d6SHong Zhang } 1426814e85d6SHong Zhang 1427814e85d6SHong Zhang /*@ 1428814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 1429814e85d6SHong Zhang 1430814e85d6SHong Zhang Collective on TS 1431814e85d6SHong Zhang 1432814e85d6SHong Zhang Input Parameter: 1433814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1434814e85d6SHong Zhang 1435814e85d6SHong Zhang Level: intermediate 1436814e85d6SHong Zhang 1437db781477SPatrick Sanan .seealso: `TSAdjointSetUp()`, `TSAdjointSolve()` 1438814e85d6SHong Zhang @*/ 1439814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts) 1440814e85d6SHong Zhang { 1441814e85d6SHong Zhang DM dm; 1442814e85d6SHong Zhang 1443814e85d6SHong Zhang PetscFunctionBegin; 1444814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 14459566063dSJacob Faibussowitsch PetscCall(TSGetDM(ts,&dm)); 14469566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 14477b0e2f17SHong Zhang ts->steps--; /* must decrease the step index before the adjoint step is taken. */ 1448814e85d6SHong Zhang 1449814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1450814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 14513c633725SBarry Smith PetscCheck(ts->ops->adjointstep,PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed because the adjoint of %s has not been implemented, try other time stepping methods for adjoint sensitivity analysis",((PetscObject)ts)->type_name); 14529566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_AdjointStep,ts,0,0,0)); 14539566063dSJacob Faibussowitsch PetscCall((*ts->ops->adjointstep)(ts)); 14549566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_AdjointStep,ts,0,0,0)); 14557b0e2f17SHong Zhang ts->adjoint_steps++; 1456814e85d6SHong Zhang 1457814e85d6SHong Zhang if (ts->reason < 0) { 14583c633725SBarry Smith PetscCheck(!ts->errorifstepfailed,PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]); 1459814e85d6SHong Zhang } else if (!ts->reason) { 1460814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1461814e85d6SHong Zhang } 1462814e85d6SHong Zhang PetscFunctionReturn(0); 1463814e85d6SHong Zhang } 1464814e85d6SHong Zhang 1465814e85d6SHong Zhang /*@ 1466814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1467814e85d6SHong Zhang 1468814e85d6SHong Zhang Collective on TS 1469814e85d6SHong Zhang 1470814e85d6SHong Zhang Input Parameter: 1471814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1472814e85d6SHong Zhang 1473814e85d6SHong Zhang Options Database: 1474814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1475814e85d6SHong Zhang 1476814e85d6SHong Zhang Level: intermediate 1477814e85d6SHong Zhang 1478814e85d6SHong Zhang Notes: 1479814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 1480814e85d6SHong Zhang 1481814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1482814e85d6SHong Zhang 1483db781477SPatrick Sanan .seealso: `TSCreate()`, `TSSetCostGradients()`, `TSSetSolution()`, `TSAdjointStep()` 1484814e85d6SHong Zhang @*/ 1485814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts) 1486814e85d6SHong Zhang { 1487f1d62c27SHong Zhang static PetscBool cite = PETSC_FALSE; 14887f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 14897f59fb53SHong Zhang PetscLogStage adjoint_stage; 14907f59fb53SHong Zhang #endif 1491814e85d6SHong Zhang 1492814e85d6SHong Zhang PetscFunctionBegin; 1493814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 14949566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister("@article{tsadjointpaper,\n" 1495f1d62c27SHong Zhang " title = {{PETSc TSAdjoint: a discrete adjoint ODE solver for first-order and second-order sensitivity analysis}},\n" 1496f1d62c27SHong Zhang " author = {Zhang, Hong and Constantinescu, Emil M. and Smith, Barry F.},\n" 1497f1d62c27SHong Zhang " journal = {arXiv e-preprints},\n" 1498f1d62c27SHong Zhang " eprint = {1912.07696},\n" 1499f1d62c27SHong Zhang " archivePrefix = {arXiv},\n" 1500b122ec5aSJacob Faibussowitsch " year = {2019}\n}\n",&cite)); 15017f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 15029566063dSJacob Faibussowitsch PetscCall(PetscLogStageRegister("TSAdjoint",&adjoint_stage)); 15039566063dSJacob Faibussowitsch PetscCall(PetscLogStagePush(adjoint_stage)); 15047f59fb53SHong Zhang #endif 15059566063dSJacob Faibussowitsch PetscCall(TSAdjointSetUp(ts)); 1506814e85d6SHong Zhang 1507814e85d6SHong Zhang /* reset time step and iteration counters */ 1508814e85d6SHong Zhang ts->adjoint_steps = 0; 1509814e85d6SHong Zhang ts->ksp_its = 0; 1510814e85d6SHong Zhang ts->snes_its = 0; 1511814e85d6SHong Zhang ts->num_snes_failures = 0; 1512814e85d6SHong Zhang ts->reject = 0; 1513814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1514814e85d6SHong Zhang 1515814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1516814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1517814e85d6SHong Zhang 1518814e85d6SHong Zhang while (!ts->reason) { 15199566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime)); 15209566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip)); 15219566063dSJacob Faibussowitsch PetscCall(TSAdjointEventHandler(ts)); 15229566063dSJacob Faibussowitsch PetscCall(TSAdjointStep(ts)); 1523cd4cee2dSHong Zhang if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) { 15249566063dSJacob Faibussowitsch PetscCall(TSAdjointCostIntegral(ts)); 1525814e85d6SHong Zhang } 1526814e85d6SHong Zhang } 1527bdbff044SHong Zhang if (!ts->steps) { 15289566063dSJacob Faibussowitsch PetscCall(TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime)); 15299566063dSJacob Faibussowitsch PetscCall(TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip)); 1530bdbff044SHong Zhang } 1531814e85d6SHong Zhang ts->solvetime = ts->ptime; 15329566063dSJacob Faibussowitsch PetscCall(TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view")); 15339566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution")); 1534814e85d6SHong Zhang ts->adjoint_max_steps = 0; 15357f59fb53SHong Zhang #if defined(TSADJOINT_STAGE) 15369566063dSJacob Faibussowitsch PetscCall(PetscLogStagePop()); 15377f59fb53SHong Zhang #endif 1538814e85d6SHong Zhang PetscFunctionReturn(0); 1539814e85d6SHong Zhang } 1540814e85d6SHong Zhang 1541814e85d6SHong Zhang /*@C 1542814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1543814e85d6SHong Zhang 1544814e85d6SHong Zhang Collective on TS 1545814e85d6SHong Zhang 1546814e85d6SHong Zhang Input Parameters: 1547814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 1548814e85d6SHong Zhang . step - step number that has just completed 1549814e85d6SHong Zhang . ptime - model time of the state 1550814e85d6SHong Zhang . u - state at the current model time 1551814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 1552814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1553814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1554814e85d6SHong Zhang 1555814e85d6SHong Zhang Notes: 1556814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1557814e85d6SHong Zhang Users would almost never call this routine directly. 1558814e85d6SHong Zhang 1559814e85d6SHong Zhang Level: developer 1560814e85d6SHong Zhang 1561814e85d6SHong Zhang @*/ 1562814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 1563814e85d6SHong Zhang { 1564814e85d6SHong Zhang PetscInt i,n = ts->numberadjointmonitors; 1565814e85d6SHong Zhang 1566814e85d6SHong Zhang PetscFunctionBegin; 1567814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1568814e85d6SHong Zhang PetscValidHeaderSpecific(u,VEC_CLASSID,4); 15699566063dSJacob Faibussowitsch PetscCall(VecLockReadPush(u)); 1570814e85d6SHong Zhang for (i=0; i<n; i++) { 15719566063dSJacob Faibussowitsch PetscCall((*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i])); 1572814e85d6SHong Zhang } 15739566063dSJacob Faibussowitsch PetscCall(VecLockReadPop(u)); 1574814e85d6SHong Zhang PetscFunctionReturn(0); 1575814e85d6SHong Zhang } 1576814e85d6SHong Zhang 1577814e85d6SHong Zhang /*@ 1578814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1579814e85d6SHong Zhang 1580814e85d6SHong Zhang Collective on TS 1581814e85d6SHong Zhang 15824165533cSJose E. Roman Input Parameter: 1583814e85d6SHong Zhang . ts - time stepping context 1584814e85d6SHong Zhang 1585814e85d6SHong Zhang Level: advanced 1586814e85d6SHong Zhang 1587814e85d6SHong Zhang Notes: 1588814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 1589814e85d6SHong Zhang 1590db781477SPatrick Sanan .seealso: `TSAdjointSolve()`, `TSAdjointStep` 1591814e85d6SHong Zhang @*/ 1592814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts) 1593814e85d6SHong Zhang { 1594362febeeSStefano Zampini PetscFunctionBegin; 1595814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 15963c633725SBarry Smith PetscCheck(ts->ops->adjointintegral,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the adjoint run",((PetscObject)ts)->type_name); 15979566063dSJacob Faibussowitsch PetscCall((*ts->ops->adjointintegral)(ts)); 1598814e85d6SHong Zhang PetscFunctionReturn(0); 1599814e85d6SHong Zhang } 1600814e85d6SHong Zhang 1601814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1602814e85d6SHong Zhang 1603814e85d6SHong Zhang /*@ 1604814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 1605814e85d6SHong Zhang of forward sensitivity analysis 1606814e85d6SHong Zhang 1607814e85d6SHong Zhang Collective on TS 1608814e85d6SHong Zhang 1609814e85d6SHong Zhang Input Parameter: 1610814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1611814e85d6SHong Zhang 1612814e85d6SHong Zhang Level: advanced 1613814e85d6SHong Zhang 1614db781477SPatrick Sanan .seealso: `TSCreate()`, `TSDestroy()`, `TSSetUp()` 1615814e85d6SHong Zhang @*/ 1616814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts) 1617814e85d6SHong Zhang { 1618814e85d6SHong Zhang PetscFunctionBegin; 1619814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1620814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1621814e85d6SHong Zhang if (ts->ops->forwardsetup) { 16229566063dSJacob Faibussowitsch PetscCall((*ts->ops->forwardsetup)(ts)); 1623814e85d6SHong Zhang } 16249566063dSJacob Faibussowitsch PetscCall(VecDuplicate(ts->vec_sol,&ts->vec_sensip_col)); 1625814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 1626814e85d6SHong Zhang PetscFunctionReturn(0); 1627814e85d6SHong Zhang } 1628814e85d6SHong Zhang 1629814e85d6SHong Zhang /*@ 16307adebddeSHong Zhang TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 16317adebddeSHong Zhang 16327adebddeSHong Zhang Collective on TS 16337adebddeSHong Zhang 16347adebddeSHong Zhang Input Parameter: 16357adebddeSHong Zhang . ts - the TS context obtained from TSCreate() 16367adebddeSHong Zhang 16377adebddeSHong Zhang Level: advanced 16387adebddeSHong Zhang 1639db781477SPatrick Sanan .seealso: `TSCreate()`, `TSDestroy()`, `TSForwardSetUp()` 16407adebddeSHong Zhang @*/ 16417adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts) 16427adebddeSHong Zhang { 16437207e0fdSHong Zhang TS quadts = ts->quadraturets; 16447adebddeSHong Zhang 16457adebddeSHong Zhang PetscFunctionBegin; 16467adebddeSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 16477adebddeSHong Zhang if (ts->ops->forwardreset) { 16489566063dSJacob Faibussowitsch PetscCall((*ts->ops->forwardreset)(ts)); 16497adebddeSHong Zhang } 16509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 16517207e0fdSHong Zhang if (quadts) { 16529566063dSJacob Faibussowitsch PetscCall(MatDestroy(&quadts->mat_sensip)); 16537207e0fdSHong Zhang } 16549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&ts->vec_sensip_col)); 16557adebddeSHong Zhang ts->forward_solve = PETSC_FALSE; 16567adebddeSHong Zhang ts->forwardsetupcalled = PETSC_FALSE; 16577adebddeSHong Zhang PetscFunctionReturn(0); 16587adebddeSHong Zhang } 16597adebddeSHong Zhang 16607adebddeSHong Zhang /*@ 1661814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1662814e85d6SHong Zhang 1663d8d19677SJose E. Roman Input Parameters: 1664a2b725a8SWilliam Gropp + ts - the TS context obtained from TSCreate() 1665814e85d6SHong Zhang . numfwdint - number of integrals 166667b8a455SSatish Balay - vp - the vectors containing the gradients for each integral w.r.t. parameters 1667814e85d6SHong Zhang 16687207e0fdSHong Zhang Level: deprecated 1669814e85d6SHong Zhang 1670db781477SPatrick Sanan .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1671814e85d6SHong Zhang @*/ 1672814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 1673814e85d6SHong Zhang { 1674814e85d6SHong Zhang PetscFunctionBegin; 1675814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 16763c633725SBarry 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()"); 1677814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1678814e85d6SHong Zhang 1679814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 1680814e85d6SHong Zhang PetscFunctionReturn(0); 1681814e85d6SHong Zhang } 1682814e85d6SHong Zhang 1683814e85d6SHong Zhang /*@ 1684814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1685814e85d6SHong Zhang 1686814e85d6SHong Zhang Input Parameter: 1687814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1688814e85d6SHong Zhang 1689814e85d6SHong Zhang Output Parameter: 169067b8a455SSatish Balay . vp - the vectors containing the gradients for each integral w.r.t. parameters 1691814e85d6SHong Zhang 16927207e0fdSHong Zhang Level: deprecated 1693814e85d6SHong Zhang 1694db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1695814e85d6SHong Zhang @*/ 1696814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1697814e85d6SHong Zhang { 1698814e85d6SHong Zhang PetscFunctionBegin; 1699814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1700814e85d6SHong Zhang PetscValidPointer(vp,3); 1701814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1702814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1703814e85d6SHong Zhang PetscFunctionReturn(0); 1704814e85d6SHong Zhang } 1705814e85d6SHong Zhang 1706814e85d6SHong Zhang /*@ 1707814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1708814e85d6SHong Zhang 1709814e85d6SHong Zhang Collective on TS 1710814e85d6SHong Zhang 17114165533cSJose E. Roman Input Parameter: 1712814e85d6SHong Zhang . ts - time stepping context 1713814e85d6SHong Zhang 1714814e85d6SHong Zhang Level: advanced 1715814e85d6SHong Zhang 1716814e85d6SHong Zhang Notes: 1717814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1718814e85d6SHong Zhang 1719db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardSetUp()` 1720814e85d6SHong Zhang @*/ 1721814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts) 1722814e85d6SHong Zhang { 1723362febeeSStefano Zampini PetscFunctionBegin; 1724814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 17253c633725SBarry Smith PetscCheck(ts->ops->forwardstep,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 17269566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(TS_ForwardStep,ts,0,0,0)); 17279566063dSJacob Faibussowitsch PetscCall((*ts->ops->forwardstep)(ts)); 17289566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(TS_ForwardStep,ts,0,0,0)); 17293c633725SBarry Smith PetscCheck(ts->reason >= 0 || !ts->errorifstepfailed,PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSFowardStep has failed due to %s",TSConvergedReasons[ts->reason]); 1730814e85d6SHong Zhang PetscFunctionReturn(0); 1731814e85d6SHong Zhang } 1732814e85d6SHong Zhang 1733814e85d6SHong Zhang /*@ 1734814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1735814e85d6SHong Zhang 1736d083f849SBarry Smith Logically Collective on TS 1737814e85d6SHong Zhang 1738814e85d6SHong Zhang Input Parameters: 1739814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1740814e85d6SHong Zhang . nump - number of parameters 1741814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1742814e85d6SHong Zhang 1743814e85d6SHong Zhang Level: beginner 1744814e85d6SHong Zhang 1745814e85d6SHong Zhang Notes: 1746814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1747814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1748814e85d6SHong Zhang You must call this function before TSSolve(). 1749814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1750814e85d6SHong Zhang 1751db781477SPatrick Sanan .seealso: `TSForwardGetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1752814e85d6SHong Zhang @*/ 1753814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1754814e85d6SHong Zhang { 1755814e85d6SHong Zhang PetscFunctionBegin; 1756814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1757814e85d6SHong Zhang PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1758814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1759b5b4017aSHong Zhang if (nump == PETSC_DEFAULT) { 17609566063dSJacob Faibussowitsch PetscCall(MatGetSize(Smat,NULL,&ts->num_parameters)); 1761b5b4017aSHong Zhang } else ts->num_parameters = nump; 17629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)Smat)); 17639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&ts->mat_sensip)); 1764814e85d6SHong Zhang ts->mat_sensip = Smat; 1765814e85d6SHong Zhang PetscFunctionReturn(0); 1766814e85d6SHong Zhang } 1767814e85d6SHong Zhang 1768814e85d6SHong Zhang /*@ 1769814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1770814e85d6SHong Zhang 1771814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1772814e85d6SHong Zhang 1773d8d19677SJose E. Roman Output Parameters: 1774814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1775814e85d6SHong Zhang . nump - number of parameters 1776814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1777814e85d6SHong Zhang 1778814e85d6SHong Zhang Level: intermediate 1779814e85d6SHong Zhang 1780db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()`, `TSForwardSetIntegralGradients()`, `TSForwardGetIntegralGradients()`, `TSForwardStep()` 1781814e85d6SHong Zhang @*/ 1782814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1783814e85d6SHong Zhang { 1784814e85d6SHong Zhang PetscFunctionBegin; 1785814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1786814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1787814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1788814e85d6SHong Zhang PetscFunctionReturn(0); 1789814e85d6SHong Zhang } 1790814e85d6SHong Zhang 1791814e85d6SHong Zhang /*@ 1792814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1793814e85d6SHong Zhang 1794814e85d6SHong Zhang Collective on TS 1795814e85d6SHong Zhang 17964165533cSJose E. Roman Input Parameter: 1797814e85d6SHong Zhang . ts - time stepping context 1798814e85d6SHong Zhang 1799814e85d6SHong Zhang Level: advanced 1800814e85d6SHong Zhang 1801814e85d6SHong Zhang Notes: 1802814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1803814e85d6SHong Zhang 1804db781477SPatrick Sanan .seealso: `TSSolve()`, `TSAdjointCostIntegral()` 1805814e85d6SHong Zhang @*/ 1806814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts) 1807814e85d6SHong Zhang { 1808362febeeSStefano Zampini PetscFunctionBegin; 1809814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 18103c633725SBarry Smith PetscCheck(ts->ops->forwardintegral,PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the forward run",((PetscObject)ts)->type_name); 18119566063dSJacob Faibussowitsch PetscCall((*ts->ops->forwardintegral)(ts)); 1812814e85d6SHong Zhang PetscFunctionReturn(0); 1813814e85d6SHong Zhang } 1814b5b4017aSHong Zhang 1815b5b4017aSHong Zhang /*@ 1816b5b4017aSHong Zhang TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1817b5b4017aSHong Zhang 1818d083f849SBarry Smith Collective on TS 1819b5b4017aSHong Zhang 1820d8d19677SJose E. Roman Input Parameters: 1821b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 1822b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition 1823b5b4017aSHong Zhang 1824b5b4017aSHong Zhang Level: intermediate 1825b5b4017aSHong Zhang 1826b5b4017aSHong Zhang Notes: TSSolve() allows users to pass the initial solution directly to TS. But the tangent linear variables cannot be initialized in this way. This function is used to set initial values for tangent linear variables. 1827b5b4017aSHong Zhang 1828db781477SPatrick Sanan .seealso: `TSForwardSetSensitivities()` 1829b5b4017aSHong Zhang @*/ 1830b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp) 1831b5b4017aSHong Zhang { 1832362febeeSStefano Zampini PetscFunctionBegin; 1833b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1834b5b4017aSHong Zhang PetscValidHeaderSpecific(didp,MAT_CLASSID,2); 1835b5b4017aSHong Zhang if (!ts->mat_sensip) { 18369566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp)); 1837b5b4017aSHong Zhang } 1838b5b4017aSHong Zhang PetscFunctionReturn(0); 1839b5b4017aSHong Zhang } 18406affb6f8SHong Zhang 18416affb6f8SHong Zhang /*@ 18426affb6f8SHong Zhang TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 18436affb6f8SHong Zhang 18446affb6f8SHong Zhang Input Parameter: 18456affb6f8SHong Zhang . ts - the TS context obtained from TSCreate() 18466affb6f8SHong Zhang 18476affb6f8SHong Zhang Output Parameters: 1848cd4cee2dSHong Zhang + ns - number of stages 18496affb6f8SHong Zhang - S - tangent linear sensitivities at the intermediate stages 18506affb6f8SHong Zhang 18516affb6f8SHong Zhang Level: advanced 18526affb6f8SHong Zhang 18536affb6f8SHong Zhang @*/ 18546affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S) 18556affb6f8SHong Zhang { 18566affb6f8SHong Zhang PetscFunctionBegin; 18576affb6f8SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 18586affb6f8SHong Zhang 18596affb6f8SHong Zhang if (!ts->ops->getstages) *S=NULL; 18606affb6f8SHong Zhang else { 18619566063dSJacob Faibussowitsch PetscCall((*ts->ops->forwardgetstages)(ts,ns,S)); 18626affb6f8SHong Zhang } 18636affb6f8SHong Zhang PetscFunctionReturn(0); 18646affb6f8SHong Zhang } 1865cd4cee2dSHong Zhang 1866cd4cee2dSHong Zhang /*@ 1867cd4cee2dSHong Zhang TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 1868cd4cee2dSHong Zhang 1869d8d19677SJose E. Roman Input Parameters: 1870cd4cee2dSHong Zhang + ts - the TS context obtained from TSCreate() 1871cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1872cd4cee2dSHong Zhang 1873cd4cee2dSHong Zhang Output Parameters: 1874cd4cee2dSHong Zhang . quadts - the child TS context 1875cd4cee2dSHong Zhang 1876cd4cee2dSHong Zhang Level: intermediate 1877cd4cee2dSHong Zhang 1878db781477SPatrick Sanan .seealso: `TSGetQuadratureTS()` 1879cd4cee2dSHong Zhang @*/ 1880cd4cee2dSHong Zhang PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts) 1881cd4cee2dSHong Zhang { 1882cd4cee2dSHong Zhang char prefix[128]; 1883cd4cee2dSHong Zhang 1884cd4cee2dSHong Zhang PetscFunctionBegin; 1885cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1886064a246eSJacob Faibussowitsch PetscValidPointer(quadts,3); 18879566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts->quadraturets)); 18889566063dSJacob Faibussowitsch PetscCall(TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets)); 18899566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1)); 18909566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets)); 18919566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "")); 18929566063dSJacob Faibussowitsch PetscCall(TSSetOptionsPrefix(ts->quadraturets,prefix)); 1893cd4cee2dSHong Zhang *quadts = ts->quadraturets; 1894cd4cee2dSHong Zhang 1895cd4cee2dSHong Zhang if (ts->numcost) { 18969566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol)); 1897cd4cee2dSHong Zhang } else { 18989566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol)); 1899cd4cee2dSHong Zhang } 1900cd4cee2dSHong Zhang ts->costintegralfwd = fwd; 1901cd4cee2dSHong Zhang PetscFunctionReturn(0); 1902cd4cee2dSHong Zhang } 1903cd4cee2dSHong Zhang 1904cd4cee2dSHong Zhang /*@ 1905cd4cee2dSHong Zhang TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 1906cd4cee2dSHong Zhang 1907cd4cee2dSHong Zhang Input Parameter: 1908cd4cee2dSHong Zhang . ts - the TS context obtained from TSCreate() 1909cd4cee2dSHong Zhang 1910cd4cee2dSHong Zhang Output Parameters: 1911cd4cee2dSHong Zhang + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1912cd4cee2dSHong Zhang - quadts - the child TS context 1913cd4cee2dSHong Zhang 1914cd4cee2dSHong Zhang Level: intermediate 1915cd4cee2dSHong Zhang 1916db781477SPatrick Sanan .seealso: `TSCreateQuadratureTS()` 1917cd4cee2dSHong Zhang @*/ 1918cd4cee2dSHong Zhang PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts) 1919cd4cee2dSHong Zhang { 1920cd4cee2dSHong Zhang PetscFunctionBegin; 1921cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1922cd4cee2dSHong Zhang if (fwd) *fwd = ts->costintegralfwd; 1923cd4cee2dSHong Zhang if (quadts) *quadts = ts->quadraturets; 1924cd4cee2dSHong Zhang PetscFunctionReturn(0); 1925cd4cee2dSHong Zhang } 1926ba0675f6SHong Zhang 1927ba0675f6SHong Zhang /*@ 1928ba0675f6SHong Zhang TSComputeSNESJacobian - Compute the SNESJacobian 1929ba0675f6SHong Zhang 1930ba0675f6SHong Zhang Input Parameters: 1931ba0675f6SHong Zhang + ts - the TS context obtained from TSCreate() 1932ba0675f6SHong Zhang - x - state vector 1933ba0675f6SHong Zhang 1934ba0675f6SHong Zhang Output Parameters: 1935ba0675f6SHong Zhang + J - Jacobian matrix 1936ba0675f6SHong Zhang - Jpre - preconditioning matrix for J (may be same as J) 1937ba0675f6SHong Zhang 1938ba0675f6SHong Zhang Level: developer 1939ba0675f6SHong Zhang 1940ba0675f6SHong Zhang Notes: 1941ba0675f6SHong Zhang Using SNES to compute the Jacobian enables finite differencing when TS Jacobian is not available. 1942ba0675f6SHong Zhang @*/ 1943ba0675f6SHong Zhang PetscErrorCode TSComputeSNESJacobian(TS ts,Vec x,Mat J,Mat Jpre) 1944ba0675f6SHong Zhang { 1945ba0675f6SHong Zhang SNES snes = ts->snes; 1946ba0675f6SHong Zhang PetscErrorCode (*jac)(SNES,Vec,Mat,Mat,void*) = NULL; 1947ba0675f6SHong Zhang 1948ba0675f6SHong Zhang PetscFunctionBegin; 1949ba0675f6SHong Zhang /* 1950ba0675f6SHong Zhang Unlike implicit methods, explicit methods do not have SNESMatFDColoring in the snes object 1951ba0675f6SHong Zhang because SNESSolve() has not been called yet; so querying SNESMatFDColoring does not work for 1952ba0675f6SHong Zhang explicit methods. Instead, we check the Jacobian compute function directly to determin if FD 1953ba0675f6SHong Zhang coloring is used. 1954ba0675f6SHong Zhang */ 19559566063dSJacob Faibussowitsch PetscCall(SNESGetJacobian(snes,NULL,NULL,&jac,NULL)); 1956ba0675f6SHong Zhang if (jac == SNESComputeJacobianDefaultColor) { 1957ba0675f6SHong Zhang Vec f; 19589566063dSJacob Faibussowitsch PetscCall(SNESSetSolution(snes,x)); 19599566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(snes,&f,NULL,NULL)); 1960ba0675f6SHong Zhang /* Force MatFDColoringApply to evaluate the SNES residual function for the base vector */ 19619566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,x,f)); 1962ba0675f6SHong Zhang } 19639566063dSJacob Faibussowitsch PetscCall(SNESComputeJacobian(snes,x,J,Jpre)); 1964ba0675f6SHong Zhang PetscFunctionReturn(0); 1965ba0675f6SHong Zhang } 1966