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 6814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/ 7814e85d6SHong Zhang 8814e85d6SHong Zhang /*@C 9b5b4017aSHong 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. 10814e85d6SHong Zhang 11814e85d6SHong Zhang Logically Collective on TS 12814e85d6SHong Zhang 13814e85d6SHong Zhang Input Parameters: 14b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 15b5b4017aSHong Zhang . Amat - JacobianP matrix 16b5b4017aSHong Zhang . func - function 17b5b4017aSHong Zhang - ctx - [optional] user-defined function context 18814e85d6SHong Zhang 19814e85d6SHong Zhang Calling sequence of func: 20814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 21814e85d6SHong Zhang + t - current timestep 22c9ad9de0SHong Zhang . U - input vector (current ODE solution) 23814e85d6SHong Zhang . A - output matrix 24814e85d6SHong Zhang - ctx - [optional] user-defined function context 25814e85d6SHong Zhang 26814e85d6SHong Zhang Level: intermediate 27814e85d6SHong Zhang 28814e85d6SHong Zhang Notes: 29814e85d6SHong 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 30814e85d6SHong Zhang 31cd4cee2dSHong Zhang .seealso: TSGetRHSJacobianP() 32814e85d6SHong Zhang @*/ 33814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 34814e85d6SHong Zhang { 35814e85d6SHong Zhang PetscErrorCode ierr; 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) { 44814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 4533f72c5dSHong Zhang ierr = MatDestroy(&ts->Jacprhs);CHKERRQ(ierr); 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 56cd4cee2dSHong Zhang Input Parameters: 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 76cd4cee2dSHong Zhang .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 97814e85d6SHong Zhang .seealso: TSSetRHSJacobianP() 98814e85d6SHong Zhang @*/ 99c9ad9de0SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat) 100814e85d6SHong Zhang { 101814e85d6SHong Zhang PetscErrorCode ierr; 102814e85d6SHong Zhang 103814e85d6SHong Zhang PetscFunctionBegin; 10433f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 105814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 106c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 107814e85d6SHong Zhang 108814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 109c9ad9de0SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr); 110814e85d6SHong Zhang PetscStackPop; 111814e85d6SHong Zhang PetscFunctionReturn(0); 112814e85d6SHong Zhang } 113814e85d6SHong Zhang 114814e85d6SHong Zhang /*@C 11533f72c5dSHong 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. 11633f72c5dSHong Zhang 11733f72c5dSHong Zhang Logically Collective on TS 11833f72c5dSHong Zhang 11933f72c5dSHong Zhang Input Parameters: 12033f72c5dSHong Zhang + ts - TS context obtained from TSCreate() 12133f72c5dSHong Zhang . Amat - JacobianP matrix 12233f72c5dSHong Zhang . func - function 12333f72c5dSHong Zhang - ctx - [optional] user-defined function context 12433f72c5dSHong Zhang 12533f72c5dSHong Zhang Calling sequence of func: 12633f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 12733f72c5dSHong Zhang + t - current timestep 12833f72c5dSHong Zhang . U - input vector (current ODE solution) 12933f72c5dSHong Zhang . Udot - time derivative of state vector 13033f72c5dSHong Zhang . shift - shift to apply, see note below 13133f72c5dSHong Zhang . A - output matrix 13233f72c5dSHong Zhang - ctx - [optional] user-defined function context 13333f72c5dSHong Zhang 13433f72c5dSHong Zhang Level: intermediate 13533f72c5dSHong Zhang 13633f72c5dSHong Zhang Notes: 13733f72c5dSHong 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 13833f72c5dSHong Zhang 13933f72c5dSHong Zhang .seealso: 14033f72c5dSHong Zhang @*/ 14133f72c5dSHong Zhang PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx) 14233f72c5dSHong Zhang { 14333f72c5dSHong Zhang PetscErrorCode ierr; 14433f72c5dSHong Zhang 14533f72c5dSHong Zhang PetscFunctionBegin; 14633f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 14733f72c5dSHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 14833f72c5dSHong Zhang 14933f72c5dSHong Zhang ts->ijacobianp = func; 15033f72c5dSHong Zhang ts->ijacobianpctx = ctx; 15133f72c5dSHong Zhang if(Amat) { 15233f72c5dSHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 15333f72c5dSHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 15433f72c5dSHong Zhang ts->Jacp = Amat; 15533f72c5dSHong Zhang } 15633f72c5dSHong Zhang PetscFunctionReturn(0); 15733f72c5dSHong Zhang } 15833f72c5dSHong Zhang 15933f72c5dSHong Zhang /*@C 16033f72c5dSHong Zhang TSComputeIJacobianP - Runs the user-defined IJacobianP function. 16133f72c5dSHong Zhang 16233f72c5dSHong Zhang Collective on TS 16333f72c5dSHong Zhang 16433f72c5dSHong Zhang Input Parameters: 16533f72c5dSHong Zhang + ts - the TS context 16633f72c5dSHong Zhang . t - current timestep 16733f72c5dSHong Zhang . U - state vector 16833f72c5dSHong Zhang . Udot - time derivative of state vector 16933f72c5dSHong Zhang . shift - shift to apply, see note below 17033f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 17133f72c5dSHong Zhang 17233f72c5dSHong Zhang Output Parameters: 17333f72c5dSHong Zhang . A - Jacobian matrix 17433f72c5dSHong Zhang 17533f72c5dSHong Zhang Level: developer 17633f72c5dSHong Zhang 17733f72c5dSHong Zhang .seealso: TSSetIJacobianP() 17833f72c5dSHong Zhang @*/ 17933f72c5dSHong Zhang PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex) 18033f72c5dSHong Zhang { 18133f72c5dSHong Zhang PetscErrorCode ierr; 18233f72c5dSHong Zhang 18333f72c5dSHong Zhang PetscFunctionBegin; 18433f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 18533f72c5dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 18633f72c5dSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 18733f72c5dSHong Zhang PetscValidHeaderSpecific(Udot,VEC_CLASSID,4); 18833f72c5dSHong Zhang 18933f72c5dSHong Zhang ierr = PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr); 19033f72c5dSHong Zhang if (ts->ijacobianp) { 19133f72c5dSHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 19233f72c5dSHong Zhang ierr = (*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx);CHKERRQ(ierr); 19333f72c5dSHong Zhang PetscStackPop; 19433f72c5dSHong Zhang } 19533f72c5dSHong Zhang if (imex) { 19633f72c5dSHong Zhang if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */ 19733f72c5dSHong Zhang PetscBool assembled; 19833f72c5dSHong Zhang ierr = MatZeroEntries(Amat);CHKERRQ(ierr); 19933f72c5dSHong Zhang ierr = MatAssembled(Amat,&assembled);CHKERRQ(ierr); 20033f72c5dSHong Zhang if (!assembled) { 20133f72c5dSHong Zhang ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20233f72c5dSHong Zhang ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20333f72c5dSHong Zhang } 20433f72c5dSHong Zhang } 20533f72c5dSHong Zhang } else { 20633f72c5dSHong Zhang if (ts->rhsjacobianp) { 20733f72c5dSHong Zhang ierr = TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs);CHKERRQ(ierr); 20833f72c5dSHong Zhang } 20933f72c5dSHong Zhang if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */ 21033f72c5dSHong Zhang ierr = MatScale(Amat,-1);CHKERRQ(ierr); 21133f72c5dSHong Zhang } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */ 21233f72c5dSHong Zhang MatStructure axpy = DIFFERENT_NONZERO_PATTERN; 21333f72c5dSHong Zhang if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */ 21433f72c5dSHong Zhang ierr = MatZeroEntries(Amat);CHKERRQ(ierr); 21533f72c5dSHong Zhang } 21633f72c5dSHong Zhang ierr = MatAXPY(Amat,-1,ts->Jacprhs,axpy);CHKERRQ(ierr); 21733f72c5dSHong Zhang } 21833f72c5dSHong Zhang } 21933f72c5dSHong Zhang ierr = PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr); 22033f72c5dSHong Zhang PetscFunctionReturn(0); 22133f72c5dSHong Zhang } 22233f72c5dSHong Zhang 22333f72c5dSHong Zhang /*@C 224814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 225814e85d6SHong Zhang 226814e85d6SHong Zhang Logically Collective on TS 227814e85d6SHong Zhang 228814e85d6SHong Zhang Input Parameters: 229814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 230814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 231814e85d6SHong Zhang . costintegral - vector that stores the integral values 232814e85d6SHong Zhang . rf - routine for evaluating the integrand function 233c9ad9de0SHong Zhang . drduf - function that computes the gradients of the r's with respect to u 234814e85d6SHong 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) 235814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 236814e85d6SHong Zhang - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 237814e85d6SHong Zhang 238814e85d6SHong Zhang Calling sequence of rf: 239c9ad9de0SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx); 240814e85d6SHong Zhang 241c9ad9de0SHong Zhang Calling sequence of drduf: 242c9ad9de0SHong Zhang $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx); 243814e85d6SHong Zhang 244814e85d6SHong Zhang Calling sequence of drdpf: 245c9ad9de0SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx); 246814e85d6SHong Zhang 247cd4cee2dSHong Zhang Level: deprecated 248814e85d6SHong Zhang 249814e85d6SHong Zhang Notes: 250814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 251814e85d6SHong Zhang 252814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients() 253814e85d6SHong Zhang @*/ 254814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*), 255c9ad9de0SHong Zhang PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*), 256814e85d6SHong Zhang PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*), 257814e85d6SHong Zhang PetscBool fwd,void *ctx) 258814e85d6SHong Zhang { 259814e85d6SHong Zhang PetscErrorCode ierr; 260814e85d6SHong Zhang 261814e85d6SHong Zhang PetscFunctionBegin; 262814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 263814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3); 264814e85d6SHong Zhang if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostGradients() or TSForwardSetIntegralGradients()"); 265814e85d6SHong Zhang if (!ts->numcost) ts->numcost=numcost; 266814e85d6SHong Zhang 267814e85d6SHong Zhang if (costintegral) { 268814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr); 269814e85d6SHong Zhang ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr); 270814e85d6SHong Zhang ts->vec_costintegral = costintegral; 271814e85d6SHong Zhang } else { 272814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 273814e85d6SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr); 274814e85d6SHong Zhang } else { 275814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr); 276814e85d6SHong Zhang } 277814e85d6SHong Zhang } 278814e85d6SHong Zhang if (!ts->vec_costintegrand) { 279814e85d6SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr); 280814e85d6SHong Zhang } else { 281814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr); 282814e85d6SHong Zhang } 283814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 284814e85d6SHong Zhang ts->costintegrand = rf; 285814e85d6SHong Zhang ts->costintegrandctx = ctx; 286c9ad9de0SHong Zhang ts->drdufunction = drduf; 287814e85d6SHong Zhang ts->drdpfunction = drdpf; 288814e85d6SHong Zhang PetscFunctionReturn(0); 289814e85d6SHong Zhang } 290814e85d6SHong Zhang 291b5b4017aSHong Zhang /*@C 292814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 293814e85d6SHong Zhang It is valid to call the routine after a backward run. 294814e85d6SHong Zhang 295814e85d6SHong Zhang Not Collective 296814e85d6SHong Zhang 297814e85d6SHong Zhang Input Parameter: 298814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 299814e85d6SHong Zhang 300814e85d6SHong Zhang Output Parameter: 301814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 302814e85d6SHong Zhang 303814e85d6SHong Zhang Level: intermediate 304814e85d6SHong Zhang 305814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 306814e85d6SHong Zhang 307814e85d6SHong Zhang @*/ 308814e85d6SHong Zhang PetscErrorCode TSGetCostIntegral(TS ts,Vec *v) 309814e85d6SHong Zhang { 310cd4cee2dSHong Zhang TS quadts; 311cd4cee2dSHong Zhang PetscErrorCode ierr; 312cd4cee2dSHong Zhang 313814e85d6SHong Zhang PetscFunctionBegin; 314814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 315814e85d6SHong Zhang PetscValidPointer(v,2); 316cd4cee2dSHong Zhang ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr); 317cd4cee2dSHong Zhang *v = quadts->vec_sol; 318814e85d6SHong Zhang PetscFunctionReturn(0); 319814e85d6SHong Zhang } 320814e85d6SHong Zhang 321b5b4017aSHong Zhang /*@C 322814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 323814e85d6SHong Zhang 324814e85d6SHong Zhang Input Parameters: 325814e85d6SHong Zhang + ts - the TS context 326814e85d6SHong Zhang . t - current time 327c9ad9de0SHong Zhang - U - state vector, i.e. current solution 328814e85d6SHong Zhang 329814e85d6SHong Zhang Output Parameter: 330c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs 331814e85d6SHong Zhang 332369cf35fSHong Zhang Notes: 333814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 334814e85d6SHong Zhang is used internally within the sensitivity analysis context. 335814e85d6SHong Zhang 336cd4cee2dSHong Zhang Level: deprecated 337814e85d6SHong Zhang 338814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 339814e85d6SHong Zhang @*/ 340c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q) 341814e85d6SHong Zhang { 342814e85d6SHong Zhang PetscErrorCode ierr; 343814e85d6SHong Zhang 344814e85d6SHong Zhang PetscFunctionBegin; 345814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 346c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 347c9ad9de0SHong Zhang PetscValidHeaderSpecific(Q,VEC_CLASSID,4); 348814e85d6SHong Zhang 349c9ad9de0SHong Zhang ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr); 350814e85d6SHong Zhang if (ts->costintegrand) { 351814e85d6SHong Zhang PetscStackPush("TS user integrand in the cost function"); 352c9ad9de0SHong Zhang ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr); 353814e85d6SHong Zhang PetscStackPop; 354814e85d6SHong Zhang } else { 355c9ad9de0SHong Zhang ierr = VecZeroEntries(Q);CHKERRQ(ierr); 356814e85d6SHong Zhang } 357814e85d6SHong Zhang 358c9ad9de0SHong Zhang ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr); 359814e85d6SHong Zhang PetscFunctionReturn(0); 360814e85d6SHong Zhang } 361814e85d6SHong Zhang 362b5b4017aSHong Zhang /*@C 363d76be551SHong Zhang TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 364814e85d6SHong Zhang 365d76be551SHong Zhang Level: deprecated 366814e85d6SHong Zhang 367814e85d6SHong Zhang @*/ 368c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 369814e85d6SHong Zhang { 370814e85d6SHong Zhang PetscErrorCode ierr; 371814e85d6SHong Zhang 372814e85d6SHong Zhang PetscFunctionBegin; 37333f72c5dSHong Zhang if (!DRDU) PetscFunctionReturn(0); 374814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 375c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 376814e85d6SHong Zhang 377c9ad9de0SHong Zhang PetscStackPush("TS user DRDU function for sensitivity analysis"); 378c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 379814e85d6SHong Zhang PetscStackPop; 380814e85d6SHong Zhang PetscFunctionReturn(0); 381814e85d6SHong Zhang } 382814e85d6SHong Zhang 383b5b4017aSHong Zhang /*@C 384d76be551SHong Zhang TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 385814e85d6SHong Zhang 386d76be551SHong Zhang Level: deprecated 387814e85d6SHong Zhang 388814e85d6SHong Zhang @*/ 389c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 390814e85d6SHong Zhang { 391814e85d6SHong Zhang PetscErrorCode ierr; 392814e85d6SHong Zhang 393814e85d6SHong Zhang PetscFunctionBegin; 39433f72c5dSHong Zhang if (!DRDP) PetscFunctionReturn(0); 395814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 396c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 397814e85d6SHong Zhang 398814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 399c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 400814e85d6SHong Zhang PetscStackPop; 401814e85d6SHong Zhang PetscFunctionReturn(0); 402814e85d6SHong Zhang } 403814e85d6SHong Zhang 404b5b4017aSHong Zhang /*@C 405b5b4017aSHong 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. 406b5b4017aSHong Zhang 407b5b4017aSHong Zhang Logically Collective on TS 408b5b4017aSHong Zhang 409b5b4017aSHong Zhang Input Parameters: 410b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 411b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 412b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 413b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 414b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 415b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 416b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 417b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 418b5b4017aSHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for F_PP 419b5b4017aSHong Zhang 420b5b4017aSHong Zhang Calling sequence of ihessianproductfunc: 421369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 422b5b4017aSHong Zhang + t - current timestep 423b5b4017aSHong Zhang . U - input vector (current ODE solution) 424369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 425b5b4017aSHong Zhang . Vr - input vector to be right-multiplied with the Hessian 426369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 427b5b4017aSHong Zhang - ctx - [optional] user-defined function context 428b5b4017aSHong Zhang 429b5b4017aSHong Zhang Level: intermediate 430b5b4017aSHong Zhang 431369cf35fSHong Zhang Notes: 432369cf35fSHong Zhang The first Hessian function and the working array are required. 433369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 434369cf35fSHong Zhang $ Vl_n^T*F_UP*Vr 435369cf35fSHong 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. 436369cf35fSHong Zhang Each entry of F_UP corresponds to the derivative 437369cf35fSHong Zhang $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}. 438369cf35fSHong 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 439369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]} 440369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 441b5b4017aSHong Zhang 442b5b4017aSHong Zhang .seealso: 443b5b4017aSHong Zhang @*/ 444b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 445b5b4017aSHong Zhang Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 446b5b4017aSHong Zhang Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 447b5b4017aSHong Zhang Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 448b5b4017aSHong Zhang void *ctx) 449b5b4017aSHong Zhang { 450b5b4017aSHong Zhang PetscFunctionBegin; 451b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 452b5b4017aSHong Zhang PetscValidPointer(ihp1,2); 453b5b4017aSHong Zhang 454b5b4017aSHong Zhang ts->ihessianproductctx = ctx; 455b5b4017aSHong Zhang if (ihp1) ts->vecs_fuu = ihp1; 456b5b4017aSHong Zhang if (ihp2) ts->vecs_fup = ihp2; 457b5b4017aSHong Zhang if (ihp3) ts->vecs_fpu = ihp3; 458b5b4017aSHong Zhang if (ihp4) ts->vecs_fpp = ihp4; 459b5b4017aSHong Zhang ts->ihessianproduct_fuu = ihessianproductfunc1; 460b5b4017aSHong Zhang ts->ihessianproduct_fup = ihessianproductfunc2; 461b5b4017aSHong Zhang ts->ihessianproduct_fpu = ihessianproductfunc3; 462b5b4017aSHong Zhang ts->ihessianproduct_fpp = ihessianproductfunc4; 463b5b4017aSHong Zhang PetscFunctionReturn(0); 464b5b4017aSHong Zhang } 465b5b4017aSHong Zhang 466b5b4017aSHong Zhang /*@C 467b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu. 468b5b4017aSHong Zhang 469b5b4017aSHong Zhang Collective on TS 470b5b4017aSHong Zhang 471b5b4017aSHong Zhang Input Parameters: 472b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 473b5b4017aSHong Zhang 474b5b4017aSHong Zhang Notes: 475b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation, 476b5b4017aSHong Zhang so most users would not generally call this routine themselves. 477b5b4017aSHong Zhang 478b5b4017aSHong Zhang Level: developer 479b5b4017aSHong Zhang 480b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 481b5b4017aSHong Zhang @*/ 482b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 483b5b4017aSHong Zhang { 484b5b4017aSHong Zhang PetscErrorCode ierr; 485b5b4017aSHong Zhang 486b5b4017aSHong Zhang PetscFunctionBegin; 48733f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 488b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 489b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 490b5b4017aSHong Zhang 49167633408SHong Zhang if (ts->ihessianproduct_fuu) { 492b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis"); 493b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 494b5b4017aSHong Zhang PetscStackPop; 49567633408SHong Zhang } 49667633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 49767633408SHong Zhang if (ts->rhshessianproduct_guu) { 49867633408SHong Zhang PetscInt nadj; 499b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 50067633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 50167633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 50267633408SHong Zhang } 50367633408SHong Zhang } 504b5b4017aSHong Zhang PetscFunctionReturn(0); 505b5b4017aSHong Zhang } 506b5b4017aSHong Zhang 507b5b4017aSHong Zhang /*@C 508b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup. 509b5b4017aSHong Zhang 510b5b4017aSHong Zhang Collective on TS 511b5b4017aSHong Zhang 512b5b4017aSHong Zhang Input Parameters: 513b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 514b5b4017aSHong Zhang 515b5b4017aSHong Zhang Notes: 516b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation, 517b5b4017aSHong Zhang so most users would not generally call this routine themselves. 518b5b4017aSHong Zhang 519b5b4017aSHong Zhang Level: developer 520b5b4017aSHong Zhang 521b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 522b5b4017aSHong Zhang @*/ 523b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 524b5b4017aSHong Zhang { 525b5b4017aSHong Zhang PetscErrorCode ierr; 526b5b4017aSHong Zhang 527b5b4017aSHong Zhang PetscFunctionBegin; 52833f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 529b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 530b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 531b5b4017aSHong Zhang 53267633408SHong Zhang if (ts->ihessianproduct_fup) { 533b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis"); 534b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 535b5b4017aSHong Zhang PetscStackPop; 53667633408SHong Zhang } 53767633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 53867633408SHong Zhang if (ts->rhshessianproduct_gup) { 53967633408SHong Zhang PetscInt nadj; 540b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 54167633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 54267633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 54367633408SHong Zhang } 54467633408SHong Zhang } 545b5b4017aSHong Zhang PetscFunctionReturn(0); 546b5b4017aSHong Zhang } 547b5b4017aSHong Zhang 548b5b4017aSHong Zhang /*@C 549b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu. 550b5b4017aSHong Zhang 551b5b4017aSHong Zhang Collective on TS 552b5b4017aSHong Zhang 553b5b4017aSHong Zhang Input Parameters: 554b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 555b5b4017aSHong Zhang 556b5b4017aSHong Zhang Notes: 557b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation, 558b5b4017aSHong Zhang so most users would not generally call this routine themselves. 559b5b4017aSHong Zhang 560b5b4017aSHong Zhang Level: developer 561b5b4017aSHong Zhang 562b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 563b5b4017aSHong Zhang @*/ 564b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 565b5b4017aSHong Zhang { 566b5b4017aSHong Zhang PetscErrorCode ierr; 567b5b4017aSHong Zhang 568b5b4017aSHong Zhang PetscFunctionBegin; 56933f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 570b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 571b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 572b5b4017aSHong Zhang 57367633408SHong Zhang if (ts->ihessianproduct_fpu) { 574b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 575b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 576b5b4017aSHong Zhang PetscStackPop; 57767633408SHong Zhang } 57867633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 57967633408SHong Zhang if (ts->rhshessianproduct_gpu) { 58067633408SHong Zhang PetscInt nadj; 581b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 58267633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 58367633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 58467633408SHong Zhang } 58567633408SHong Zhang } 586b5b4017aSHong Zhang PetscFunctionReturn(0); 587b5b4017aSHong Zhang } 588b5b4017aSHong Zhang 589b5b4017aSHong Zhang /*@C 590b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp. 591b5b4017aSHong Zhang 592b5b4017aSHong Zhang Collective on TS 593b5b4017aSHong Zhang 594b5b4017aSHong Zhang Input Parameters: 595b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 596b5b4017aSHong Zhang 597b5b4017aSHong Zhang Notes: 598b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation, 599b5b4017aSHong Zhang so most users would not generally call this routine themselves. 600b5b4017aSHong Zhang 601b5b4017aSHong Zhang Level: developer 602b5b4017aSHong Zhang 603b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 604b5b4017aSHong Zhang @*/ 605b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 606b5b4017aSHong Zhang { 607b5b4017aSHong Zhang PetscErrorCode ierr; 608b5b4017aSHong Zhang 609b5b4017aSHong Zhang PetscFunctionBegin; 61033f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 611b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 612b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 613b5b4017aSHong Zhang 61467633408SHong Zhang if (ts->ihessianproduct_fpp) { 615b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 616b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 617b5b4017aSHong Zhang PetscStackPop; 61867633408SHong Zhang } 61967633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 62067633408SHong Zhang if (ts->rhshessianproduct_gpp) { 62167633408SHong Zhang PetscInt nadj; 622b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 62367633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 62467633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 62567633408SHong Zhang } 62667633408SHong Zhang } 627b5b4017aSHong Zhang PetscFunctionReturn(0); 628b5b4017aSHong Zhang } 629b5b4017aSHong Zhang 63013af1a74SHong Zhang /*@C 631*4c500f23SPierre 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. 63213af1a74SHong Zhang 63313af1a74SHong Zhang Logically Collective on TS 63413af1a74SHong Zhang 63513af1a74SHong Zhang Input Parameters: 63613af1a74SHong Zhang + ts - TS context obtained from TSCreate() 63713af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 63813af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 63913af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 64013af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 64113af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 64213af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 64313af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 64413af1a74SHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for G_PP 64513af1a74SHong Zhang 64613af1a74SHong Zhang Calling sequence of ihessianproductfunc: 647369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 64813af1a74SHong Zhang + t - current timestep 64913af1a74SHong Zhang . U - input vector (current ODE solution) 650369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 65113af1a74SHong Zhang . Vr - input vector to be right-multiplied with the Hessian 652369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 65313af1a74SHong Zhang - ctx - [optional] user-defined function context 65413af1a74SHong Zhang 65513af1a74SHong Zhang Level: intermediate 65613af1a74SHong Zhang 657369cf35fSHong Zhang Notes: 658369cf35fSHong Zhang The first Hessian function and the working array are required. 659369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 660369cf35fSHong Zhang $ Vl_n^T*G_UP*Vr 661369cf35fSHong 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. 662369cf35fSHong Zhang Each entry of G_UP corresponds to the derivative 663369cf35fSHong Zhang $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}. 664369cf35fSHong 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 665369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]} 666369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 66713af1a74SHong Zhang 66813af1a74SHong Zhang .seealso: 66913af1a74SHong Zhang @*/ 67013af1a74SHong Zhang PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 67113af1a74SHong Zhang Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 67213af1a74SHong Zhang Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 67313af1a74SHong Zhang Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 67413af1a74SHong Zhang void *ctx) 67513af1a74SHong Zhang { 67613af1a74SHong Zhang PetscFunctionBegin; 67713af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 67813af1a74SHong Zhang PetscValidPointer(rhshp1,2); 67913af1a74SHong Zhang 68013af1a74SHong Zhang ts->rhshessianproductctx = ctx; 68113af1a74SHong Zhang if (rhshp1) ts->vecs_guu = rhshp1; 68213af1a74SHong Zhang if (rhshp2) ts->vecs_gup = rhshp2; 68313af1a74SHong Zhang if (rhshp3) ts->vecs_gpu = rhshp3; 68413af1a74SHong Zhang if (rhshp4) ts->vecs_gpp = rhshp4; 68513af1a74SHong Zhang ts->rhshessianproduct_guu = rhshessianproductfunc1; 68613af1a74SHong Zhang ts->rhshessianproduct_gup = rhshessianproductfunc2; 68713af1a74SHong Zhang ts->rhshessianproduct_gpu = rhshessianproductfunc3; 68813af1a74SHong Zhang ts->rhshessianproduct_gpp = rhshessianproductfunc4; 68913af1a74SHong Zhang PetscFunctionReturn(0); 69013af1a74SHong Zhang } 69113af1a74SHong Zhang 69213af1a74SHong Zhang /*@C 693b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu. 69413af1a74SHong Zhang 69513af1a74SHong Zhang Collective on TS 69613af1a74SHong Zhang 69713af1a74SHong Zhang Input Parameters: 69813af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 69913af1a74SHong Zhang 70013af1a74SHong Zhang Notes: 701b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation, 70213af1a74SHong Zhang so most users would not generally call this routine themselves. 70313af1a74SHong Zhang 70413af1a74SHong Zhang Level: developer 70513af1a74SHong Zhang 70613af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 70713af1a74SHong Zhang @*/ 708b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 70913af1a74SHong Zhang { 71013af1a74SHong Zhang PetscErrorCode ierr; 71113af1a74SHong Zhang 71213af1a74SHong Zhang PetscFunctionBegin; 71313af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 71413af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 71513af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 71613af1a74SHong Zhang 71713af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis"); 71813af1a74SHong Zhang ierr = (*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 71913af1a74SHong Zhang PetscStackPop; 72013af1a74SHong Zhang PetscFunctionReturn(0); 72113af1a74SHong Zhang } 72213af1a74SHong Zhang 72313af1a74SHong Zhang /*@C 724b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup. 72513af1a74SHong Zhang 72613af1a74SHong Zhang Collective on TS 72713af1a74SHong Zhang 72813af1a74SHong Zhang Input Parameters: 72913af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 73013af1a74SHong Zhang 73113af1a74SHong Zhang Notes: 732b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation, 73313af1a74SHong Zhang so most users would not generally call this routine themselves. 73413af1a74SHong Zhang 73513af1a74SHong Zhang Level: developer 73613af1a74SHong Zhang 73713af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 73813af1a74SHong Zhang @*/ 739b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 74013af1a74SHong Zhang { 74113af1a74SHong Zhang PetscErrorCode ierr; 74213af1a74SHong Zhang 74313af1a74SHong Zhang PetscFunctionBegin; 74413af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 74513af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 74613af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 74713af1a74SHong Zhang 74813af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis"); 74913af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 75013af1a74SHong Zhang PetscStackPop; 75113af1a74SHong Zhang PetscFunctionReturn(0); 75213af1a74SHong Zhang } 75313af1a74SHong Zhang 75413af1a74SHong Zhang /*@C 755b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu. 75613af1a74SHong Zhang 75713af1a74SHong Zhang Collective on TS 75813af1a74SHong Zhang 75913af1a74SHong Zhang Input Parameters: 76013af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 76113af1a74SHong Zhang 76213af1a74SHong Zhang Notes: 763b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation, 76413af1a74SHong Zhang so most users would not generally call this routine themselves. 76513af1a74SHong Zhang 76613af1a74SHong Zhang Level: developer 76713af1a74SHong Zhang 76813af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 76913af1a74SHong Zhang @*/ 770b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 77113af1a74SHong Zhang { 77213af1a74SHong Zhang PetscErrorCode ierr; 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"); 78013af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 78113af1a74SHong Zhang PetscStackPop; 78213af1a74SHong Zhang PetscFunctionReturn(0); 78313af1a74SHong Zhang } 78413af1a74SHong Zhang 78513af1a74SHong Zhang /*@C 786b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp. 78713af1a74SHong Zhang 78813af1a74SHong Zhang Collective on TS 78913af1a74SHong Zhang 79013af1a74SHong Zhang Input Parameters: 79113af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 79213af1a74SHong Zhang 79313af1a74SHong Zhang Notes: 794b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation, 79513af1a74SHong Zhang so most users would not generally call this routine themselves. 79613af1a74SHong Zhang 79713af1a74SHong Zhang Level: developer 79813af1a74SHong Zhang 79913af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 80013af1a74SHong Zhang @*/ 801b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 80213af1a74SHong Zhang { 80313af1a74SHong Zhang PetscErrorCode ierr; 80413af1a74SHong Zhang 80513af1a74SHong Zhang PetscFunctionBegin; 80613af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 80713af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 80813af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 80913af1a74SHong Zhang 81013af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis"); 81113af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 81213af1a74SHong Zhang PetscStackPop; 81313af1a74SHong Zhang PetscFunctionReturn(0); 81413af1a74SHong Zhang } 81513af1a74SHong Zhang 816814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 817814e85d6SHong Zhang 818814e85d6SHong Zhang /*@ 819814e85d6SHong 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 820814e85d6SHong Zhang for use by the TSAdjoint routines. 821814e85d6SHong Zhang 822d083f849SBarry Smith Logically Collective on TS 823814e85d6SHong Zhang 824814e85d6SHong Zhang Input Parameters: 825814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 826814e85d6SHong 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 827814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 828814e85d6SHong Zhang 829814e85d6SHong Zhang Level: beginner 830814e85d6SHong Zhang 831814e85d6SHong Zhang Notes: 832814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 833814e85d6SHong Zhang 834814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 835814e85d6SHong Zhang 836b5b4017aSHong Zhang .seealso TSGetCostGradients() 837814e85d6SHong Zhang @*/ 838814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu) 839814e85d6SHong Zhang { 840814e85d6SHong Zhang PetscFunctionBegin; 841814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 842814e85d6SHong Zhang PetscValidPointer(lambda,2); 843814e85d6SHong Zhang ts->vecs_sensi = lambda; 844814e85d6SHong Zhang ts->vecs_sensip = mu; 845814e85d6SHong Zhang if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand"); 846814e85d6SHong Zhang ts->numcost = numcost; 847814e85d6SHong Zhang PetscFunctionReturn(0); 848814e85d6SHong Zhang } 849814e85d6SHong Zhang 850814e85d6SHong Zhang /*@ 851814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 852814e85d6SHong Zhang 853814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 854814e85d6SHong Zhang 855814e85d6SHong Zhang Input Parameter: 856814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 857814e85d6SHong Zhang 858814e85d6SHong Zhang Output Parameter: 859814e85d6SHong Zhang + lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 860814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 861814e85d6SHong Zhang 862814e85d6SHong Zhang Level: intermediate 863814e85d6SHong Zhang 864b5b4017aSHong Zhang .seealso: TSSetCostGradients() 865814e85d6SHong Zhang @*/ 866814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu) 867814e85d6SHong Zhang { 868814e85d6SHong Zhang PetscFunctionBegin; 869814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 870814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 871814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 872814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 873814e85d6SHong Zhang PetscFunctionReturn(0); 874814e85d6SHong Zhang } 875814e85d6SHong Zhang 876814e85d6SHong Zhang /*@ 877b5b4017aSHong 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 878b5b4017aSHong Zhang for use by the TSAdjoint routines. 879b5b4017aSHong Zhang 880d083f849SBarry Smith Logically Collective on TS 881b5b4017aSHong Zhang 882b5b4017aSHong Zhang Input Parameters: 883b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 884b5b4017aSHong Zhang . numcost - number of cost functions 885b5b4017aSHong 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 886b5b4017aSHong 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 887b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 888b5b4017aSHong Zhang 889b5b4017aSHong Zhang Level: beginner 890b5b4017aSHong Zhang 891b5b4017aSHong Zhang Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system 892b5b4017aSHong Zhang 893ecf68647SHong Zhang For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve(). 894b5b4017aSHong Zhang 895b5b4017aSHong 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. 896b5b4017aSHong Zhang 8973fca9621SHong Zhang Passing NULL for lambda2 disables the second-order calculation. 898ecf68647SHong Zhang .seealso: TSAdjointSetForward() 899b5b4017aSHong Zhang @*/ 900b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir) 901b5b4017aSHong Zhang { 902b5b4017aSHong Zhang PetscFunctionBegin; 903b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 904b5b4017aSHong Zhang if (ts->numcost && ts->numcost!=numcost) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand"); 905b5b4017aSHong Zhang ts->numcost = numcost; 906b5b4017aSHong Zhang ts->vecs_sensi2 = lambda2; 90733f72c5dSHong Zhang ts->vecs_sensi2p = mu2; 908b5b4017aSHong Zhang ts->vec_dir = dir; 909b5b4017aSHong Zhang PetscFunctionReturn(0); 910b5b4017aSHong Zhang } 911b5b4017aSHong Zhang 912b5b4017aSHong Zhang /*@ 913b5b4017aSHong Zhang TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve() 914b5b4017aSHong Zhang 915b5b4017aSHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 916b5b4017aSHong Zhang 917b5b4017aSHong Zhang Input Parameter: 918b5b4017aSHong Zhang . ts - the TS context obtained from TSCreate() 919b5b4017aSHong Zhang 920b5b4017aSHong Zhang Output Parameter: 921b5b4017aSHong Zhang + numcost - number of cost functions 922b5b4017aSHong 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 923b5b4017aSHong 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 924b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 925b5b4017aSHong Zhang 926b5b4017aSHong Zhang Level: intermediate 927b5b4017aSHong Zhang 928b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts() 929b5b4017aSHong Zhang @*/ 930b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir) 931b5b4017aSHong Zhang { 932b5b4017aSHong Zhang PetscFunctionBegin; 933b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 934b5b4017aSHong Zhang if (numcost) *numcost = ts->numcost; 935b5b4017aSHong Zhang if (lambda2) *lambda2 = ts->vecs_sensi2; 93633f72c5dSHong Zhang if (mu2) *mu2 = ts->vecs_sensi2p; 937b5b4017aSHong Zhang if (dir) *dir = ts->vec_dir; 938b5b4017aSHong Zhang PetscFunctionReturn(0); 939b5b4017aSHong Zhang } 940b5b4017aSHong Zhang 941b5b4017aSHong Zhang /*@ 942ecf68647SHong Zhang TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities 943b5b4017aSHong Zhang 944d083f849SBarry Smith Logically Collective on TS 945b5b4017aSHong Zhang 946b5b4017aSHong Zhang Input Parameters: 947b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 948b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters 949b5b4017aSHong Zhang 950b5b4017aSHong Zhang Level: intermediate 951b5b4017aSHong Zhang 952ecf68647SHong 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. 953b5b4017aSHong Zhang 954ecf68647SHong Zhang .seealso: TSSetCostHessianProducts(), TSAdjointResetForward() 955b5b4017aSHong Zhang @*/ 956ecf68647SHong Zhang PetscErrorCode TSAdjointSetForward(TS ts,Mat didp) 957b5b4017aSHong Zhang { 95833f72c5dSHong Zhang Mat A; 95933f72c5dSHong Zhang Vec sp; 96033f72c5dSHong Zhang PetscScalar *xarr; 96133f72c5dSHong Zhang PetscInt lsize; 962b5b4017aSHong Zhang PetscErrorCode ierr; 963b5b4017aSHong Zhang 964b5b4017aSHong Zhang PetscFunctionBegin; 965b5b4017aSHong Zhang ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 966b5b4017aSHong Zhang if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first"); 96733f72c5dSHong Zhang if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 96833f72c5dSHong Zhang /* create a single-column dense matrix */ 96933f72c5dSHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&lsize);CHKERRQ(ierr); 970f63bf25fSHong Zhang ierr = MatCreateDense(PetscObjectComm((PetscObject)ts),lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr); 97133f72c5dSHong Zhang 97233f72c5dSHong Zhang ierr = VecDuplicate(ts->vec_sol,&sp);CHKERRQ(ierr); 97333f72c5dSHong Zhang ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr); 97433f72c5dSHong Zhang ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr); 975ecf68647SHong Zhang if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */ 97633f72c5dSHong Zhang if (didp) { 97733f72c5dSHong Zhang ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr); 97833f72c5dSHong Zhang ierr = VecScale(sp,2.);CHKERRQ(ierr); 97933f72c5dSHong Zhang } else { 98033f72c5dSHong Zhang ierr = VecZeroEntries(sp);CHKERRQ(ierr); 98133f72c5dSHong Zhang } 982ecf68647SHong Zhang } else { /* tangent linear variable initialized as dir */ 98333f72c5dSHong Zhang ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr); 98433f72c5dSHong Zhang } 98533f72c5dSHong Zhang ierr = VecResetArray(sp);CHKERRQ(ierr); 98633f72c5dSHong Zhang ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr); 98733f72c5dSHong Zhang ierr = VecDestroy(&sp);CHKERRQ(ierr); 98833f72c5dSHong Zhang 98933f72c5dSHong Zhang ierr = TSForwardSetInitialSensitivities(ts,A);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */ 99033f72c5dSHong Zhang 99133f72c5dSHong Zhang ierr = MatDestroy(&A);CHKERRQ(ierr); 992b5b4017aSHong Zhang PetscFunctionReturn(0); 993b5b4017aSHong Zhang } 994b5b4017aSHong Zhang 995b5b4017aSHong Zhang /*@ 996ecf68647SHong Zhang TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context 997ecf68647SHong Zhang 998d083f849SBarry Smith Logically Collective on TS 999ecf68647SHong Zhang 1000ecf68647SHong Zhang Input Parameters: 1001ecf68647SHong Zhang . ts - the TS context obtained from TSCreate() 1002ecf68647SHong Zhang 1003ecf68647SHong Zhang Level: intermediate 1004ecf68647SHong Zhang 1005ecf68647SHong Zhang .seealso: TSAdjointSetForward() 1006ecf68647SHong Zhang @*/ 1007ecf68647SHong Zhang PetscErrorCode TSAdjointResetForward(TS ts) 1008ecf68647SHong Zhang { 1009ecf68647SHong Zhang PetscErrorCode ierr; 1010ecf68647SHong Zhang 1011ecf68647SHong Zhang PetscFunctionBegin; 1012ecf68647SHong Zhang ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */ 1013ecf68647SHong Zhang ierr = TSForwardReset(ts);CHKERRQ(ierr); 1014ecf68647SHong Zhang PetscFunctionReturn(0); 1015ecf68647SHong Zhang } 1016ecf68647SHong Zhang 1017ecf68647SHong Zhang /*@ 1018814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 1019814e85d6SHong Zhang of an adjoint solver 1020814e85d6SHong Zhang 1021814e85d6SHong Zhang Collective on TS 1022814e85d6SHong Zhang 1023814e85d6SHong Zhang Input Parameter: 1024814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1025814e85d6SHong Zhang 1026814e85d6SHong Zhang Level: advanced 1027814e85d6SHong Zhang 1028814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients() 1029814e85d6SHong Zhang @*/ 1030814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts) 1031814e85d6SHong Zhang { 1032881c1a9bSHong Zhang TSTrajectory tj; 1033881c1a9bSHong Zhang PetscBool match; 1034814e85d6SHong Zhang PetscErrorCode ierr; 1035814e85d6SHong Zhang 1036814e85d6SHong Zhang PetscFunctionBegin; 1037814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1038814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 1039814e85d6SHong Zhang if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 104033f72c5dSHong Zhang if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 1041881c1a9bSHong Zhang ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr); 10428e224c67SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match);CHKERRQ(ierr); 1043881c1a9bSHong Zhang if (match) { 1044881c1a9bSHong Zhang PetscBool solution_only; 1045881c1a9bSHong Zhang ierr = TSTrajectoryGetSolutionOnly(tj,&solution_only);CHKERRQ(ierr); 1046881c1a9bSHong Zhang if (solution_only) SETERRQ(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"); 1047881c1a9bSHong Zhang } 10489ffb3502SHong Zhang ierr = TSTrajectorySetUseHistory(tj,PETSC_FALSE);CHKERRQ(ierr); /* not use TSHistory */ 104933f72c5dSHong Zhang 1050cd4cee2dSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 1051cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);CHKERRQ(ierr); 1052814e85d6SHong Zhang if (ts->vecs_sensip){ 1053cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr); 1054814e85d6SHong Zhang } 1055814e85d6SHong Zhang } 1056814e85d6SHong Zhang 1057814e85d6SHong Zhang if (ts->ops->adjointsetup) { 1058814e85d6SHong Zhang ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr); 1059814e85d6SHong Zhang } 1060814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 1061814e85d6SHong Zhang PetscFunctionReturn(0); 1062814e85d6SHong Zhang } 1063814e85d6SHong Zhang 1064814e85d6SHong Zhang /*@ 1065ece46509SHong Zhang TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 1066ece46509SHong Zhang 1067ece46509SHong Zhang Collective on TS 1068ece46509SHong Zhang 1069ece46509SHong Zhang Input Parameter: 1070ece46509SHong Zhang . ts - the TS context obtained from TSCreate() 1071ece46509SHong Zhang 1072ece46509SHong Zhang Level: beginner 1073ece46509SHong Zhang 10749ffb3502SHong Zhang .seealso: TSCreate(), TSAdjointSetUp(), TSADestroy() 1075ece46509SHong Zhang @*/ 1076ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts) 1077ece46509SHong Zhang { 1078ece46509SHong Zhang PetscErrorCode ierr; 1079ece46509SHong Zhang 1080ece46509SHong Zhang PetscFunctionBegin; 1081ece46509SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1082ece46509SHong Zhang if (ts->ops->adjointreset) { 1083ece46509SHong Zhang ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr); 1084ece46509SHong Zhang } 10857207e0fdSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 10867207e0fdSHong Zhang ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr); 10877207e0fdSHong Zhang if (ts->vecs_sensip){ 10887207e0fdSHong Zhang ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr); 10897207e0fdSHong Zhang } 10907207e0fdSHong Zhang } 1091ece46509SHong Zhang ts->vecs_sensi = NULL; 1092ece46509SHong Zhang ts->vecs_sensip = NULL; 1093ece46509SHong Zhang ts->vecs_sensi2 = NULL; 109433f72c5dSHong Zhang ts->vecs_sensi2p = NULL; 1095ece46509SHong Zhang ts->vec_dir = NULL; 1096ece46509SHong Zhang ts->adjointsetupcalled = PETSC_FALSE; 1097ece46509SHong Zhang PetscFunctionReturn(0); 1098ece46509SHong Zhang } 1099ece46509SHong Zhang 1100ece46509SHong Zhang /*@ 1101814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1102814e85d6SHong Zhang 1103814e85d6SHong Zhang Logically Collective on TS 1104814e85d6SHong Zhang 1105814e85d6SHong Zhang Input Parameters: 1106814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1107a2b725a8SWilliam Gropp - steps - number of steps to use 1108814e85d6SHong Zhang 1109814e85d6SHong Zhang Level: intermediate 1110814e85d6SHong Zhang 1111814e85d6SHong Zhang Notes: 1112814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 1113814e85d6SHong Zhang so as to integrate back to less than the original timestep 1114814e85d6SHong Zhang 1115814e85d6SHong Zhang .seealso: TSSetExactFinalTime() 1116814e85d6SHong Zhang @*/ 1117814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 1118814e85d6SHong Zhang { 1119814e85d6SHong Zhang PetscFunctionBegin; 1120814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1121814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts,steps,2); 1122814e85d6SHong Zhang if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 1123814e85d6SHong Zhang if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 1124814e85d6SHong Zhang ts->adjoint_max_steps = steps; 1125814e85d6SHong Zhang PetscFunctionReturn(0); 1126814e85d6SHong Zhang } 1127814e85d6SHong Zhang 1128814e85d6SHong Zhang /*@C 1129814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1130814e85d6SHong Zhang 1131814e85d6SHong Zhang Level: deprecated 1132814e85d6SHong Zhang 1133814e85d6SHong Zhang @*/ 1134814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 1135814e85d6SHong Zhang { 1136814e85d6SHong Zhang PetscErrorCode ierr; 1137814e85d6SHong Zhang 1138814e85d6SHong Zhang PetscFunctionBegin; 1139814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1140814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 1141814e85d6SHong Zhang 1142814e85d6SHong Zhang ts->rhsjacobianp = func; 1143814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 1144814e85d6SHong Zhang if(Amat) { 1145814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 1146814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 1147814e85d6SHong Zhang ts->Jacp = Amat; 1148814e85d6SHong Zhang } 1149814e85d6SHong Zhang PetscFunctionReturn(0); 1150814e85d6SHong Zhang } 1151814e85d6SHong Zhang 1152814e85d6SHong Zhang /*@C 1153814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1154814e85d6SHong Zhang 1155814e85d6SHong Zhang Level: deprecated 1156814e85d6SHong Zhang 1157814e85d6SHong Zhang @*/ 1158c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat) 1159814e85d6SHong Zhang { 1160814e85d6SHong Zhang PetscErrorCode ierr; 1161814e85d6SHong Zhang 1162814e85d6SHong Zhang PetscFunctionBegin; 1163814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1164c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1165814e85d6SHong Zhang PetscValidPointer(Amat,4); 1166814e85d6SHong Zhang 1167814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 1168c9ad9de0SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 1169814e85d6SHong Zhang PetscStackPop; 1170814e85d6SHong Zhang PetscFunctionReturn(0); 1171814e85d6SHong Zhang } 1172814e85d6SHong Zhang 1173814e85d6SHong Zhang /*@ 1174d76be551SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 1175814e85d6SHong Zhang 1176814e85d6SHong Zhang Level: deprecated 1177814e85d6SHong Zhang 1178814e85d6SHong Zhang @*/ 1179c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 1180814e85d6SHong Zhang { 1181814e85d6SHong Zhang PetscErrorCode ierr; 1182814e85d6SHong Zhang 1183814e85d6SHong Zhang PetscFunctionBegin; 1184814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1185c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1186814e85d6SHong Zhang 1187814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 1188c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 1189814e85d6SHong Zhang PetscStackPop; 1190814e85d6SHong Zhang PetscFunctionReturn(0); 1191814e85d6SHong Zhang } 1192814e85d6SHong Zhang 1193814e85d6SHong Zhang /*@ 1194d76be551SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 1195814e85d6SHong Zhang 1196814e85d6SHong Zhang Level: deprecated 1197814e85d6SHong Zhang 1198814e85d6SHong Zhang @*/ 1199c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 1200814e85d6SHong Zhang { 1201814e85d6SHong Zhang PetscErrorCode ierr; 1202814e85d6SHong Zhang 1203814e85d6SHong Zhang PetscFunctionBegin; 1204814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1205c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1206814e85d6SHong Zhang 1207814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 1208c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 1209814e85d6SHong Zhang PetscStackPop; 1210814e85d6SHong Zhang PetscFunctionReturn(0); 1211814e85d6SHong Zhang } 1212814e85d6SHong Zhang 1213814e85d6SHong Zhang /*@C 1214814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 1215814e85d6SHong Zhang 1216814e85d6SHong Zhang Level: intermediate 1217814e85d6SHong Zhang 1218814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1219814e85d6SHong Zhang @*/ 1220814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1221814e85d6SHong Zhang { 1222814e85d6SHong Zhang PetscErrorCode ierr; 1223814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1224814e85d6SHong Zhang 1225814e85d6SHong Zhang PetscFunctionBegin; 1226814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1227814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1228814e85d6SHong Zhang ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 1229814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1230814e85d6SHong Zhang PetscFunctionReturn(0); 1231814e85d6SHong Zhang } 1232814e85d6SHong Zhang 1233814e85d6SHong Zhang /*@C 1234814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1235814e85d6SHong Zhang 1236814e85d6SHong Zhang Collective on TS 1237814e85d6SHong Zhang 1238814e85d6SHong Zhang Input Parameters: 1239814e85d6SHong Zhang + ts - TS object you wish to monitor 1240814e85d6SHong Zhang . name - the monitor type one is seeking 1241814e85d6SHong Zhang . help - message indicating what monitoring is done 1242814e85d6SHong Zhang . manual - manual page for the monitor 1243814e85d6SHong Zhang . monitor - the monitor function 1244814e85d6SHong 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 1245814e85d6SHong Zhang 1246814e85d6SHong Zhang Level: developer 1247814e85d6SHong Zhang 1248814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 1249814e85d6SHong Zhang PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1250814e85d6SHong Zhang PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 1251814e85d6SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1252814e85d6SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1253814e85d6SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1254814e85d6SHong Zhang PetscOptionsFList(), PetscOptionsEList() 1255814e85d6SHong Zhang @*/ 1256814e85d6SHong 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*)) 1257814e85d6SHong Zhang { 1258814e85d6SHong Zhang PetscErrorCode ierr; 1259814e85d6SHong Zhang PetscViewer viewer; 1260814e85d6SHong Zhang PetscViewerFormat format; 1261814e85d6SHong Zhang PetscBool flg; 1262814e85d6SHong Zhang 1263814e85d6SHong Zhang PetscFunctionBegin; 126416413a6aSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 1265814e85d6SHong Zhang if (flg) { 1266814e85d6SHong Zhang PetscViewerAndFormat *vf; 1267814e85d6SHong Zhang ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 1268814e85d6SHong Zhang ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 1269814e85d6SHong Zhang if (monitorsetup) { 1270814e85d6SHong Zhang ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 1271814e85d6SHong Zhang } 1272814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 1273814e85d6SHong Zhang } 1274814e85d6SHong Zhang PetscFunctionReturn(0); 1275814e85d6SHong Zhang } 1276814e85d6SHong Zhang 1277814e85d6SHong Zhang /*@C 1278814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1279814e85d6SHong Zhang timestep to display the iteration's progress. 1280814e85d6SHong Zhang 1281814e85d6SHong Zhang Logically Collective on TS 1282814e85d6SHong Zhang 1283814e85d6SHong Zhang Input Parameters: 1284814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1285814e85d6SHong Zhang . adjointmonitor - monitoring routine 1286814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 1287814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 1288814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 1289814e85d6SHong Zhang (may be NULL) 1290814e85d6SHong Zhang 1291814e85d6SHong Zhang Calling sequence of monitor: 1292814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1293814e85d6SHong Zhang 1294814e85d6SHong Zhang + ts - the TS context 1295814e85d6SHong 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 1296814e85d6SHong Zhang been interpolated to) 1297814e85d6SHong Zhang . time - current time 1298814e85d6SHong Zhang . u - current iterate 1299814e85d6SHong Zhang . numcost - number of cost functionos 1300814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1301814e85d6SHong Zhang . mu - sensitivities to parameters 1302814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 1303814e85d6SHong Zhang 1304814e85d6SHong Zhang Notes: 1305814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 1306814e85d6SHong Zhang already has been loaded. 1307814e85d6SHong Zhang 1308814e85d6SHong Zhang Fortran Notes: 1309814e85d6SHong Zhang Only a single monitor function can be set for each TS object 1310814e85d6SHong Zhang 1311814e85d6SHong Zhang Level: intermediate 1312814e85d6SHong Zhang 1313814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel() 1314814e85d6SHong Zhang @*/ 1315814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 1316814e85d6SHong Zhang { 1317814e85d6SHong Zhang PetscErrorCode ierr; 1318814e85d6SHong Zhang PetscInt i; 1319814e85d6SHong Zhang PetscBool identical; 1320814e85d6SHong Zhang 1321814e85d6SHong Zhang PetscFunctionBegin; 1322814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1323814e85d6SHong Zhang for (i=0; i<ts->numbermonitors;i++) { 1324814e85d6SHong Zhang ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 1325814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 1326814e85d6SHong Zhang } 1327814e85d6SHong Zhang if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 1328814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1329814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1330814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 1331814e85d6SHong Zhang PetscFunctionReturn(0); 1332814e85d6SHong Zhang } 1333814e85d6SHong Zhang 1334814e85d6SHong Zhang /*@C 1335814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1336814e85d6SHong Zhang 1337814e85d6SHong Zhang Logically Collective on TS 1338814e85d6SHong Zhang 1339814e85d6SHong Zhang Input Parameters: 1340814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1341814e85d6SHong Zhang 1342814e85d6SHong Zhang Notes: 1343814e85d6SHong Zhang There is no way to remove a single, specific monitor. 1344814e85d6SHong Zhang 1345814e85d6SHong Zhang Level: intermediate 1346814e85d6SHong Zhang 1347814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1348814e85d6SHong Zhang @*/ 1349814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts) 1350814e85d6SHong Zhang { 1351814e85d6SHong Zhang PetscErrorCode ierr; 1352814e85d6SHong Zhang PetscInt i; 1353814e85d6SHong Zhang 1354814e85d6SHong Zhang PetscFunctionBegin; 1355814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1356814e85d6SHong Zhang for (i=0; i<ts->numberadjointmonitors; i++) { 1357814e85d6SHong Zhang if (ts->adjointmonitordestroy[i]) { 1358814e85d6SHong Zhang ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1359814e85d6SHong Zhang } 1360814e85d6SHong Zhang } 1361814e85d6SHong Zhang ts->numberadjointmonitors = 0; 1362814e85d6SHong Zhang PetscFunctionReturn(0); 1363814e85d6SHong Zhang } 1364814e85d6SHong Zhang 1365814e85d6SHong Zhang /*@C 1366814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 1367814e85d6SHong Zhang 1368814e85d6SHong Zhang Level: intermediate 1369814e85d6SHong Zhang 1370814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1371814e85d6SHong Zhang @*/ 1372814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1373814e85d6SHong Zhang { 1374814e85d6SHong Zhang PetscErrorCode ierr; 1375814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1376814e85d6SHong Zhang 1377814e85d6SHong Zhang PetscFunctionBegin; 1378814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1379814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1380814e85d6SHong Zhang ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1381814e85d6SHong Zhang ierr = PetscViewerASCIIPrintf(viewer,"%D TS dt %g time %g%s",step,(double)ts->time_step,(double)ptime,ts->steprollback ? " (r)\n" : "\n");CHKERRQ(ierr); 1382814e85d6SHong Zhang ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1383814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1384814e85d6SHong Zhang PetscFunctionReturn(0); 1385814e85d6SHong Zhang } 1386814e85d6SHong Zhang 1387814e85d6SHong Zhang /*@C 1388814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1389814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 1390814e85d6SHong Zhang 1391814e85d6SHong Zhang Collective on TS 1392814e85d6SHong Zhang 1393814e85d6SHong Zhang Input Parameters: 1394814e85d6SHong Zhang + ts - the TS context 1395814e85d6SHong Zhang . step - current time-step 1396814e85d6SHong Zhang . ptime - current time 1397814e85d6SHong Zhang . u - current state 1398814e85d6SHong Zhang . numcost - number of cost functions 1399814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1400814e85d6SHong Zhang . mu - sensitivities to parameters 1401814e85d6SHong Zhang - dummy - either a viewer or NULL 1402814e85d6SHong Zhang 1403814e85d6SHong Zhang Level: intermediate 1404814e85d6SHong Zhang 1405814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 1406814e85d6SHong Zhang @*/ 1407814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 1408814e85d6SHong Zhang { 1409814e85d6SHong Zhang PetscErrorCode ierr; 1410814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1411814e85d6SHong Zhang PetscDraw draw; 1412814e85d6SHong Zhang PetscReal xl,yl,xr,yr,h; 1413814e85d6SHong Zhang char time[32]; 1414814e85d6SHong Zhang 1415814e85d6SHong Zhang PetscFunctionBegin; 1416814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1417814e85d6SHong Zhang 1418814e85d6SHong Zhang ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 1419814e85d6SHong Zhang ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 1420814e85d6SHong Zhang ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 1421814e85d6SHong Zhang ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1422814e85d6SHong Zhang h = yl + .95*(yr - yl); 1423814e85d6SHong Zhang ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 1424814e85d6SHong Zhang ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 1425814e85d6SHong Zhang PetscFunctionReturn(0); 1426814e85d6SHong Zhang } 1427814e85d6SHong Zhang 1428814e85d6SHong Zhang /* 1429814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1430814e85d6SHong Zhang 1431814e85d6SHong Zhang Collective on TSAdjoint 1432814e85d6SHong Zhang 1433814e85d6SHong Zhang Input Parameter: 1434814e85d6SHong Zhang . ts - the TS context 1435814e85d6SHong Zhang 1436814e85d6SHong Zhang Options Database Keys: 1437814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1438814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 1439814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1440814e85d6SHong Zhang 1441814e85d6SHong Zhang Level: developer 1442814e85d6SHong Zhang 1443814e85d6SHong Zhang Notes: 1444814e85d6SHong Zhang This is not normally called directly by users 1445814e85d6SHong Zhang 1446814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 1447814e85d6SHong Zhang */ 1448814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 1449814e85d6SHong Zhang { 1450814e85d6SHong Zhang PetscBool tflg,opt; 1451814e85d6SHong Zhang PetscErrorCode ierr; 1452814e85d6SHong Zhang 1453814e85d6SHong Zhang PetscFunctionBegin; 1454814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,2); 1455814e85d6SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 1456814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 1457814e85d6SHong Zhang ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 1458814e85d6SHong Zhang if (opt) { 1459814e85d6SHong Zhang ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 1460814e85d6SHong Zhang ts->adjoint_solve = tflg; 1461814e85d6SHong Zhang } 1462814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 1463814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 1464814e85d6SHong Zhang opt = PETSC_FALSE; 1465814e85d6SHong Zhang ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 1466814e85d6SHong Zhang if (opt) { 1467814e85d6SHong Zhang TSMonitorDrawCtx ctx; 1468814e85d6SHong Zhang PetscInt howoften = 1; 1469814e85d6SHong Zhang 1470814e85d6SHong Zhang ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 1471814e85d6SHong Zhang ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 1472814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 1473814e85d6SHong Zhang } 1474814e85d6SHong Zhang PetscFunctionReturn(0); 1475814e85d6SHong Zhang } 1476814e85d6SHong Zhang 1477814e85d6SHong Zhang /*@ 1478814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 1479814e85d6SHong Zhang 1480814e85d6SHong Zhang Collective on TS 1481814e85d6SHong Zhang 1482814e85d6SHong Zhang Input Parameter: 1483814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1484814e85d6SHong Zhang 1485814e85d6SHong Zhang Level: intermediate 1486814e85d6SHong Zhang 1487814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve() 1488814e85d6SHong Zhang @*/ 1489814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts) 1490814e85d6SHong Zhang { 1491814e85d6SHong Zhang DM dm; 1492814e85d6SHong Zhang PetscErrorCode ierr; 1493814e85d6SHong Zhang 1494814e85d6SHong Zhang PetscFunctionBegin; 1495814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1496814e85d6SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1497814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1498814e85d6SHong Zhang 1499814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1500814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 1501814e85d6SHong Zhang if (!ts->ops->adjointstep) SETERRQ1(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); 1502814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1503814e85d6SHong Zhang ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 1504814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1505814e85d6SHong Zhang ts->adjoint_steps++; ts->steps--; 1506814e85d6SHong Zhang 1507814e85d6SHong Zhang if (ts->reason < 0) { 1508814e85d6SHong Zhang if (ts->errorifstepfailed) { 150905c9950eSHong Zhang if (ts->reason == TSADJOINT_DIVERGED_LINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]); 151005c9950eSHong Zhang else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]); 1511814e85d6SHong Zhang } 1512814e85d6SHong Zhang } else if (!ts->reason) { 1513814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1514814e85d6SHong Zhang } 1515814e85d6SHong Zhang PetscFunctionReturn(0); 1516814e85d6SHong Zhang } 1517814e85d6SHong Zhang 1518814e85d6SHong Zhang /*@ 1519814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1520814e85d6SHong Zhang 1521814e85d6SHong Zhang Collective on TS 1522814e85d6SHong Zhang 1523814e85d6SHong Zhang Input Parameter: 1524814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1525814e85d6SHong Zhang 1526814e85d6SHong Zhang Options Database: 1527814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1528814e85d6SHong Zhang 1529814e85d6SHong Zhang Level: intermediate 1530814e85d6SHong Zhang 1531814e85d6SHong Zhang Notes: 1532814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 1533814e85d6SHong Zhang 1534814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1535814e85d6SHong Zhang 1536814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 1537814e85d6SHong Zhang @*/ 1538814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts) 1539814e85d6SHong Zhang { 1540814e85d6SHong Zhang PetscErrorCode ierr; 1541814e85d6SHong Zhang 1542814e85d6SHong Zhang PetscFunctionBegin; 1543814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1544814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1545814e85d6SHong Zhang 1546814e85d6SHong Zhang /* reset time step and iteration counters */ 1547814e85d6SHong Zhang ts->adjoint_steps = 0; 1548814e85d6SHong Zhang ts->ksp_its = 0; 1549814e85d6SHong Zhang ts->snes_its = 0; 1550814e85d6SHong Zhang ts->num_snes_failures = 0; 1551814e85d6SHong Zhang ts->reject = 0; 1552814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1553814e85d6SHong Zhang 1554814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1555814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1556814e85d6SHong Zhang 1557814e85d6SHong Zhang while (!ts->reason) { 1558814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1559814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1560814e85d6SHong Zhang ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 1561814e85d6SHong Zhang ierr = TSAdjointStep(ts);CHKERRQ(ierr); 1562cd4cee2dSHong Zhang if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) { 1563814e85d6SHong Zhang ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 1564814e85d6SHong Zhang } 1565814e85d6SHong Zhang } 1566814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1567814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1568814e85d6SHong Zhang ts->solvetime = ts->ptime; 1569814e85d6SHong Zhang ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 1570814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 1571814e85d6SHong Zhang ts->adjoint_max_steps = 0; 1572814e85d6SHong Zhang PetscFunctionReturn(0); 1573814e85d6SHong Zhang } 1574814e85d6SHong Zhang 1575814e85d6SHong Zhang /*@C 1576814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1577814e85d6SHong Zhang 1578814e85d6SHong Zhang Collective on TS 1579814e85d6SHong Zhang 1580814e85d6SHong Zhang Input Parameters: 1581814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 1582814e85d6SHong Zhang . step - step number that has just completed 1583814e85d6SHong Zhang . ptime - model time of the state 1584814e85d6SHong Zhang . u - state at the current model time 1585814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 1586814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1587814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1588814e85d6SHong Zhang 1589814e85d6SHong Zhang Notes: 1590814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1591814e85d6SHong Zhang Users would almost never call this routine directly. 1592814e85d6SHong Zhang 1593814e85d6SHong Zhang Level: developer 1594814e85d6SHong Zhang 1595814e85d6SHong Zhang @*/ 1596814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 1597814e85d6SHong Zhang { 1598814e85d6SHong Zhang PetscErrorCode ierr; 1599814e85d6SHong Zhang PetscInt i,n = ts->numberadjointmonitors; 1600814e85d6SHong Zhang 1601814e85d6SHong Zhang PetscFunctionBegin; 1602814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1603814e85d6SHong Zhang PetscValidHeaderSpecific(u,VEC_CLASSID,4); 16048860a134SJunchao Zhang ierr = VecLockReadPush(u);CHKERRQ(ierr); 1605814e85d6SHong Zhang for (i=0; i<n; i++) { 1606814e85d6SHong Zhang ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1607814e85d6SHong Zhang } 16088860a134SJunchao Zhang ierr = VecLockReadPop(u);CHKERRQ(ierr); 1609814e85d6SHong Zhang PetscFunctionReturn(0); 1610814e85d6SHong Zhang } 1611814e85d6SHong Zhang 1612814e85d6SHong Zhang /*@ 1613814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1614814e85d6SHong Zhang 1615814e85d6SHong Zhang Collective on TS 1616814e85d6SHong Zhang 1617814e85d6SHong Zhang Input Arguments: 1618814e85d6SHong Zhang . ts - time stepping context 1619814e85d6SHong Zhang 1620814e85d6SHong Zhang Level: advanced 1621814e85d6SHong Zhang 1622814e85d6SHong Zhang Notes: 1623814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 1624814e85d6SHong Zhang 1625814e85d6SHong Zhang .seealso: TSAdjointSolve(), TSAdjointStep 1626814e85d6SHong Zhang @*/ 1627814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts) 1628814e85d6SHong Zhang { 1629814e85d6SHong Zhang PetscErrorCode ierr; 1630814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1631814e85d6SHong Zhang if (!ts->ops->adjointintegral) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the adjoint run",((PetscObject)ts)->type_name); 1632814e85d6SHong Zhang ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 1633814e85d6SHong Zhang PetscFunctionReturn(0); 1634814e85d6SHong Zhang } 1635814e85d6SHong Zhang 1636814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1637814e85d6SHong Zhang 1638814e85d6SHong Zhang /*@ 1639814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 1640814e85d6SHong Zhang of forward sensitivity analysis 1641814e85d6SHong Zhang 1642814e85d6SHong Zhang Collective on TS 1643814e85d6SHong Zhang 1644814e85d6SHong Zhang Input Parameter: 1645814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1646814e85d6SHong Zhang 1647814e85d6SHong Zhang Level: advanced 1648814e85d6SHong Zhang 1649814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp() 1650814e85d6SHong Zhang @*/ 1651814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts) 1652814e85d6SHong Zhang { 1653814e85d6SHong Zhang PetscErrorCode ierr; 1654814e85d6SHong Zhang 1655814e85d6SHong Zhang PetscFunctionBegin; 1656814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1657814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1658814e85d6SHong Zhang if (ts->ops->forwardsetup) { 1659814e85d6SHong Zhang ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 1660814e85d6SHong Zhang } 16617207e0fdSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 1662814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 1663814e85d6SHong Zhang PetscFunctionReturn(0); 1664814e85d6SHong Zhang } 1665814e85d6SHong Zhang 1666814e85d6SHong Zhang /*@ 16677adebddeSHong Zhang TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 16687adebddeSHong Zhang 16697adebddeSHong Zhang Collective on TS 16707adebddeSHong Zhang 16717adebddeSHong Zhang Input Parameter: 16727adebddeSHong Zhang . ts - the TS context obtained from TSCreate() 16737adebddeSHong Zhang 16747adebddeSHong Zhang Level: advanced 16757adebddeSHong Zhang 16767adebddeSHong Zhang .seealso: TSCreate(), TSDestroy(), TSForwardSetUp() 16777adebddeSHong Zhang @*/ 16787adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts) 16797adebddeSHong Zhang { 16807207e0fdSHong Zhang TS quadts = ts->quadraturets; 16817adebddeSHong Zhang PetscErrorCode ierr; 16827adebddeSHong Zhang 16837adebddeSHong Zhang PetscFunctionBegin; 16847adebddeSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 16857adebddeSHong Zhang if (ts->ops->forwardreset) { 16867adebddeSHong Zhang ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr); 16877adebddeSHong Zhang } 16887adebddeSHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 16897207e0fdSHong Zhang if (quadts) { 16907207e0fdSHong Zhang ierr = MatDestroy(&quadts->mat_sensip);CHKERRQ(ierr); 16917207e0fdSHong Zhang } 16927207e0fdSHong Zhang ierr = VecDestroy(&ts->vec_sensip_col);CHKERRQ(ierr); 16937adebddeSHong Zhang ts->forward_solve = PETSC_FALSE; 16947adebddeSHong Zhang ts->forwardsetupcalled = PETSC_FALSE; 16957adebddeSHong Zhang PetscFunctionReturn(0); 16967adebddeSHong Zhang } 16977adebddeSHong Zhang 16987adebddeSHong Zhang /*@ 1699814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1700814e85d6SHong Zhang 1701814e85d6SHong Zhang Input Parameter: 1702a2b725a8SWilliam Gropp + ts- the TS context obtained from TSCreate() 1703814e85d6SHong Zhang . numfwdint- number of integrals 1704a2b725a8SWilliam Gropp - vp = the vectors containing the gradients for each integral w.r.t. parameters 1705814e85d6SHong Zhang 17067207e0fdSHong Zhang Level: deprecated 1707814e85d6SHong Zhang 1708814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1709814e85d6SHong Zhang @*/ 1710814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 1711814e85d6SHong Zhang { 1712814e85d6SHong Zhang PetscFunctionBegin; 1713814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1714814e85d6SHong Zhang if (ts->numcost && ts->numcost!=numfwdint) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"The number of cost functions (2rd parameter of TSSetCostIntegrand()) is inconsistent with the one set by TSSetCostIntegrand()"); 1715814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1716814e85d6SHong Zhang 1717814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 1718814e85d6SHong Zhang PetscFunctionReturn(0); 1719814e85d6SHong Zhang } 1720814e85d6SHong Zhang 1721814e85d6SHong Zhang /*@ 1722814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1723814e85d6SHong Zhang 1724814e85d6SHong Zhang Input Parameter: 1725814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1726814e85d6SHong Zhang 1727814e85d6SHong Zhang Output Parameter: 1728814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1729814e85d6SHong Zhang 17307207e0fdSHong Zhang Level: deprecated 1731814e85d6SHong Zhang 1732814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1733814e85d6SHong Zhang @*/ 1734814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1735814e85d6SHong Zhang { 1736814e85d6SHong Zhang PetscFunctionBegin; 1737814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1738814e85d6SHong Zhang PetscValidPointer(vp,3); 1739814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1740814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1741814e85d6SHong Zhang PetscFunctionReturn(0); 1742814e85d6SHong Zhang } 1743814e85d6SHong Zhang 1744814e85d6SHong Zhang /*@ 1745814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1746814e85d6SHong Zhang 1747814e85d6SHong Zhang Collective on TS 1748814e85d6SHong Zhang 1749814e85d6SHong Zhang Input Arguments: 1750814e85d6SHong Zhang . ts - time stepping context 1751814e85d6SHong Zhang 1752814e85d6SHong Zhang Level: advanced 1753814e85d6SHong Zhang 1754814e85d6SHong Zhang Notes: 1755814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1756814e85d6SHong Zhang 1757814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1758814e85d6SHong Zhang @*/ 1759814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts) 1760814e85d6SHong Zhang { 1761814e85d6SHong Zhang PetscErrorCode ierr; 1762814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1763814e85d6SHong Zhang if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1764814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1765814e85d6SHong Zhang ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1766814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1767814e85d6SHong Zhang PetscFunctionReturn(0); 1768814e85d6SHong Zhang } 1769814e85d6SHong Zhang 1770814e85d6SHong Zhang /*@ 1771814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1772814e85d6SHong Zhang 1773d083f849SBarry Smith Logically Collective on TS 1774814e85d6SHong Zhang 1775814e85d6SHong Zhang Input Parameters: 1776814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1777814e85d6SHong Zhang . nump - number of parameters 1778814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1779814e85d6SHong Zhang 1780814e85d6SHong Zhang Level: beginner 1781814e85d6SHong Zhang 1782814e85d6SHong Zhang Notes: 1783814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1784814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1785814e85d6SHong Zhang You must call this function before TSSolve(). 1786814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1787814e85d6SHong Zhang 1788814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1789814e85d6SHong Zhang @*/ 1790814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1791814e85d6SHong Zhang { 1792814e85d6SHong Zhang PetscErrorCode ierr; 1793814e85d6SHong Zhang 1794814e85d6SHong Zhang PetscFunctionBegin; 1795814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1796814e85d6SHong Zhang PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1797814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1798b5b4017aSHong Zhang if (nump == PETSC_DEFAULT) { 1799b5b4017aSHong Zhang ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr); 1800b5b4017aSHong Zhang } else ts->num_parameters = nump; 1801814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1802814e85d6SHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1803814e85d6SHong Zhang ts->mat_sensip = Smat; 1804814e85d6SHong Zhang PetscFunctionReturn(0); 1805814e85d6SHong Zhang } 1806814e85d6SHong Zhang 1807814e85d6SHong Zhang /*@ 1808814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1809814e85d6SHong Zhang 1810814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1811814e85d6SHong Zhang 1812814e85d6SHong Zhang Output Parameter: 1813814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1814814e85d6SHong Zhang . nump - number of parameters 1815814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1816814e85d6SHong Zhang 1817814e85d6SHong Zhang Level: intermediate 1818814e85d6SHong Zhang 1819814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1820814e85d6SHong Zhang @*/ 1821814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1822814e85d6SHong Zhang { 1823814e85d6SHong Zhang PetscFunctionBegin; 1824814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1825814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1826814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1827814e85d6SHong Zhang PetscFunctionReturn(0); 1828814e85d6SHong Zhang } 1829814e85d6SHong Zhang 1830814e85d6SHong Zhang /*@ 1831814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1832814e85d6SHong Zhang 1833814e85d6SHong Zhang Collective on TS 1834814e85d6SHong Zhang 1835814e85d6SHong Zhang Input Arguments: 1836814e85d6SHong Zhang . ts - time stepping context 1837814e85d6SHong Zhang 1838814e85d6SHong Zhang Level: advanced 1839814e85d6SHong Zhang 1840814e85d6SHong Zhang Notes: 1841814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1842814e85d6SHong Zhang 1843814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral() 1844814e85d6SHong Zhang @*/ 1845814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts) 1846814e85d6SHong Zhang { 1847814e85d6SHong Zhang PetscErrorCode ierr; 1848814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1849814e85d6SHong Zhang if (!ts->ops->forwardintegral) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide integral evaluation in the forward run",((PetscObject)ts)->type_name); 1850814e85d6SHong Zhang ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1851814e85d6SHong Zhang PetscFunctionReturn(0); 1852814e85d6SHong Zhang } 1853b5b4017aSHong Zhang 1854b5b4017aSHong Zhang /*@ 1855b5b4017aSHong Zhang TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1856b5b4017aSHong Zhang 1857d083f849SBarry Smith Collective on TS 1858b5b4017aSHong Zhang 1859b5b4017aSHong Zhang Input Parameter 1860b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 1861b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition 1862b5b4017aSHong Zhang 1863b5b4017aSHong Zhang Level: intermediate 1864b5b4017aSHong Zhang 1865b5b4017aSHong 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. 1866b5b4017aSHong Zhang 1867b5b4017aSHong Zhang .seealso: TSForwardSetSensitivities() 1868b5b4017aSHong Zhang @*/ 1869b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp) 1870b5b4017aSHong Zhang { 1871b5b4017aSHong Zhang PetscErrorCode ierr; 1872b5b4017aSHong Zhang 1873b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1874b5b4017aSHong Zhang PetscValidHeaderSpecific(didp,MAT_CLASSID,2); 1875b5b4017aSHong Zhang if (!ts->mat_sensip) { 1876b5b4017aSHong Zhang ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr); 1877b5b4017aSHong Zhang } 1878b5b4017aSHong Zhang PetscFunctionReturn(0); 1879b5b4017aSHong Zhang } 18806affb6f8SHong Zhang 18816affb6f8SHong Zhang /*@ 18826affb6f8SHong Zhang TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 18836affb6f8SHong Zhang 18846affb6f8SHong Zhang Input Parameter: 18856affb6f8SHong Zhang . ts - the TS context obtained from TSCreate() 18866affb6f8SHong Zhang 18876affb6f8SHong Zhang Output Parameters: 1888cd4cee2dSHong Zhang + ns - number of stages 18896affb6f8SHong Zhang - S - tangent linear sensitivities at the intermediate stages 18906affb6f8SHong Zhang 18916affb6f8SHong Zhang Level: advanced 18926affb6f8SHong Zhang 18936affb6f8SHong Zhang @*/ 18946affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S) 18956affb6f8SHong Zhang { 18966affb6f8SHong Zhang PetscErrorCode ierr; 18976affb6f8SHong Zhang 18986affb6f8SHong Zhang PetscFunctionBegin; 18996affb6f8SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 19006affb6f8SHong Zhang 19016affb6f8SHong Zhang if (!ts->ops->getstages) *S=NULL; 19026affb6f8SHong Zhang else { 19036affb6f8SHong Zhang ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr); 19046affb6f8SHong Zhang } 19056affb6f8SHong Zhang PetscFunctionReturn(0); 19066affb6f8SHong Zhang } 1907cd4cee2dSHong Zhang 1908cd4cee2dSHong Zhang /*@ 1909cd4cee2dSHong Zhang TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 1910cd4cee2dSHong Zhang 1911cd4cee2dSHong Zhang Input Parameter: 1912cd4cee2dSHong Zhang + ts - the TS context obtained from TSCreate() 1913cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1914cd4cee2dSHong Zhang 1915cd4cee2dSHong Zhang Output Parameters: 1916cd4cee2dSHong Zhang . quadts - the child TS context 1917cd4cee2dSHong Zhang 1918cd4cee2dSHong Zhang Level: intermediate 1919cd4cee2dSHong Zhang 1920cd4cee2dSHong Zhang .seealso: TSGetQuadratureTS() 1921cd4cee2dSHong Zhang @*/ 1922cd4cee2dSHong Zhang PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts) 1923cd4cee2dSHong Zhang { 1924cd4cee2dSHong Zhang char prefix[128]; 1925cd4cee2dSHong Zhang PetscErrorCode ierr; 1926cd4cee2dSHong Zhang 1927cd4cee2dSHong Zhang PetscFunctionBegin; 1928cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1929cd4cee2dSHong Zhang PetscValidPointer(quadts,2); 1930cd4cee2dSHong Zhang ierr = TSDestroy(&ts->quadraturets);CHKERRQ(ierr); 1931cd4cee2dSHong Zhang ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);CHKERRQ(ierr); 1932cd4cee2dSHong Zhang ierr = PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);CHKERRQ(ierr); 1933cd4cee2dSHong Zhang ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);CHKERRQ(ierr); 1934cd4cee2dSHong Zhang ierr = PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");CHKERRQ(ierr); 1935cd4cee2dSHong Zhang ierr = TSSetOptionsPrefix(ts->quadraturets,prefix);CHKERRQ(ierr); 1936cd4cee2dSHong Zhang *quadts = ts->quadraturets; 1937cd4cee2dSHong Zhang 1938cd4cee2dSHong Zhang if (ts->numcost) { 1939cd4cee2dSHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);CHKERRQ(ierr); 1940cd4cee2dSHong Zhang } else { 1941cd4cee2dSHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);CHKERRQ(ierr); 1942cd4cee2dSHong Zhang } 1943cd4cee2dSHong Zhang ts->costintegralfwd = fwd; 1944cd4cee2dSHong Zhang PetscFunctionReturn(0); 1945cd4cee2dSHong Zhang } 1946cd4cee2dSHong Zhang 1947cd4cee2dSHong Zhang /*@ 1948cd4cee2dSHong Zhang TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 1949cd4cee2dSHong Zhang 1950cd4cee2dSHong Zhang Input Parameter: 1951cd4cee2dSHong Zhang . ts - the TS context obtained from TSCreate() 1952cd4cee2dSHong Zhang 1953cd4cee2dSHong Zhang Output Parameters: 1954cd4cee2dSHong Zhang + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1955cd4cee2dSHong Zhang - quadts - the child TS context 1956cd4cee2dSHong Zhang 1957cd4cee2dSHong Zhang Level: intermediate 1958cd4cee2dSHong Zhang 1959cd4cee2dSHong Zhang .seealso: TSCreateQuadratureTS() 1960cd4cee2dSHong Zhang @*/ 1961cd4cee2dSHong Zhang PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts) 1962cd4cee2dSHong Zhang { 1963cd4cee2dSHong Zhang PetscFunctionBegin; 1964cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1965cd4cee2dSHong Zhang if (fwd) *fwd = ts->costintegralfwd; 1966cd4cee2dSHong Zhang if (quadts) *quadts = ts->quadraturets; 1967cd4cee2dSHong Zhang PetscFunctionReturn(0); 1968cd4cee2dSHong Zhang } 1969