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 31814e85d6SHong Zhang .keywords: TS, sensitivity 32cd4cee2dSHong Zhang .seealso: TSGetRHSJacobianP() 33814e85d6SHong Zhang @*/ 34814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 35814e85d6SHong Zhang { 36814e85d6SHong Zhang PetscErrorCode ierr; 37814e85d6SHong Zhang 38814e85d6SHong Zhang PetscFunctionBegin; 39814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 40814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 41814e85d6SHong Zhang 42814e85d6SHong Zhang ts->rhsjacobianp = func; 43814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 44814e85d6SHong Zhang if(Amat) { 45814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 4633f72c5dSHong Zhang ierr = MatDestroy(&ts->Jacprhs);CHKERRQ(ierr); 4733f72c5dSHong Zhang ts->Jacprhs = Amat; 48814e85d6SHong Zhang } 49814e85d6SHong Zhang PetscFunctionReturn(0); 50814e85d6SHong Zhang } 51814e85d6SHong Zhang 52814e85d6SHong Zhang /*@C 53cd4cee2dSHong 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. 54cd4cee2dSHong Zhang 55cd4cee2dSHong Zhang Logically Collective on TS 56cd4cee2dSHong Zhang 57cd4cee2dSHong Zhang Input Parameters: 58cd4cee2dSHong Zhang . ts - TS context obtained from TSCreate() 59cd4cee2dSHong Zhang 60cd4cee2dSHong Zhang Output Parameters: 61cd4cee2dSHong Zhang + Amat - JacobianP matrix 62cd4cee2dSHong Zhang . func - function 63cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 64cd4cee2dSHong Zhang 65cd4cee2dSHong Zhang Calling sequence of func: 66cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 67cd4cee2dSHong Zhang + t - current timestep 68cd4cee2dSHong Zhang . U - input vector (current ODE solution) 69cd4cee2dSHong Zhang . A - output matrix 70cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 71cd4cee2dSHong Zhang 72cd4cee2dSHong Zhang Level: intermediate 73cd4cee2dSHong Zhang 74cd4cee2dSHong Zhang Notes: 75cd4cee2dSHong 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 76cd4cee2dSHong Zhang 77cd4cee2dSHong Zhang .keywords: TS, sensitivity 78cd4cee2dSHong Zhang .seealso: TSSetRHSJacobianP() 79cd4cee2dSHong Zhang @*/ 80cd4cee2dSHong Zhang PetscErrorCode TSGetRHSJacobianP(TS ts,Mat *Amat,PetscErrorCode (**func)(TS,PetscReal,Vec,Mat,void*),void **ctx) 81cd4cee2dSHong Zhang { 82cd4cee2dSHong Zhang PetscFunctionBegin; 83cd4cee2dSHong Zhang if (func) *func = ts->rhsjacobianp; 84cd4cee2dSHong Zhang if (ctx) *ctx = ts->rhsjacobianpctx; 85cd4cee2dSHong Zhang if (Amat) *Amat = ts->Jacprhs; 86cd4cee2dSHong Zhang PetscFunctionReturn(0); 87cd4cee2dSHong Zhang } 88cd4cee2dSHong Zhang 89cd4cee2dSHong Zhang /*@C 90814e85d6SHong Zhang TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 91814e85d6SHong Zhang 92814e85d6SHong Zhang Collective on TS 93814e85d6SHong Zhang 94814e85d6SHong Zhang Input Parameters: 95814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 96814e85d6SHong Zhang 97814e85d6SHong Zhang Level: developer 98814e85d6SHong Zhang 99814e85d6SHong Zhang .keywords: TS, sensitivity 100814e85d6SHong Zhang .seealso: TSSetRHSJacobianP() 101814e85d6SHong Zhang @*/ 102c9ad9de0SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat) 103814e85d6SHong Zhang { 104814e85d6SHong Zhang PetscErrorCode ierr; 105814e85d6SHong Zhang 106814e85d6SHong Zhang PetscFunctionBegin; 10733f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 108814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 109c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 110814e85d6SHong Zhang 111814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 112c9ad9de0SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr); 113814e85d6SHong Zhang PetscStackPop; 114814e85d6SHong Zhang PetscFunctionReturn(0); 115814e85d6SHong Zhang } 116814e85d6SHong Zhang 117814e85d6SHong Zhang /*@C 11833f72c5dSHong 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. 11933f72c5dSHong Zhang 12033f72c5dSHong Zhang Logically Collective on TS 12133f72c5dSHong Zhang 12233f72c5dSHong Zhang Input Parameters: 12333f72c5dSHong Zhang + ts - TS context obtained from TSCreate() 12433f72c5dSHong Zhang . Amat - JacobianP matrix 12533f72c5dSHong Zhang . func - function 12633f72c5dSHong Zhang - ctx - [optional] user-defined function context 12733f72c5dSHong Zhang 12833f72c5dSHong Zhang Calling sequence of func: 12933f72c5dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 13033f72c5dSHong Zhang + t - current timestep 13133f72c5dSHong Zhang . U - input vector (current ODE solution) 13233f72c5dSHong Zhang . Udot - time derivative of state vector 13333f72c5dSHong Zhang . shift - shift to apply, see note below 13433f72c5dSHong Zhang . A - output matrix 13533f72c5dSHong Zhang - ctx - [optional] user-defined function context 13633f72c5dSHong Zhang 13733f72c5dSHong Zhang Level: intermediate 13833f72c5dSHong Zhang 13933f72c5dSHong Zhang Notes: 14033f72c5dSHong 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 14133f72c5dSHong Zhang 14233f72c5dSHong Zhang .keywords: TS, sensitivity 14333f72c5dSHong Zhang .seealso: 14433f72c5dSHong Zhang @*/ 14533f72c5dSHong Zhang PetscErrorCode TSSetIJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Vec,PetscReal,Mat,void*),void *ctx) 14633f72c5dSHong Zhang { 14733f72c5dSHong Zhang PetscErrorCode ierr; 14833f72c5dSHong Zhang 14933f72c5dSHong Zhang PetscFunctionBegin; 15033f72c5dSHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 15133f72c5dSHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 15233f72c5dSHong Zhang 15333f72c5dSHong Zhang ts->ijacobianp = func; 15433f72c5dSHong Zhang ts->ijacobianpctx = ctx; 15533f72c5dSHong Zhang if(Amat) { 15633f72c5dSHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 15733f72c5dSHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 15833f72c5dSHong Zhang ts->Jacp = Amat; 15933f72c5dSHong Zhang } 16033f72c5dSHong Zhang PetscFunctionReturn(0); 16133f72c5dSHong Zhang } 16233f72c5dSHong Zhang 16333f72c5dSHong Zhang /*@C 16433f72c5dSHong Zhang TSComputeIJacobianP - Runs the user-defined IJacobianP function. 16533f72c5dSHong Zhang 16633f72c5dSHong Zhang Collective on TS 16733f72c5dSHong Zhang 16833f72c5dSHong Zhang Input Parameters: 16933f72c5dSHong Zhang + ts - the TS context 17033f72c5dSHong Zhang . t - current timestep 17133f72c5dSHong Zhang . U - state vector 17233f72c5dSHong Zhang . Udot - time derivative of state vector 17333f72c5dSHong Zhang . shift - shift to apply, see note below 17433f72c5dSHong Zhang - imex - flag indicates if the method is IMEX so that the RHSJacobian should be kept separate 17533f72c5dSHong Zhang 17633f72c5dSHong Zhang Output Parameters: 17733f72c5dSHong Zhang . A - Jacobian matrix 17833f72c5dSHong Zhang 17933f72c5dSHong Zhang Level: developer 18033f72c5dSHong Zhang 18133f72c5dSHong Zhang .keywords: TS, sensitivity 18233f72c5dSHong Zhang .seealso: TSSetIJacobianP() 18333f72c5dSHong Zhang @*/ 18433f72c5dSHong Zhang PetscErrorCode TSComputeIJacobianP(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal shift,Mat Amat,PetscBool imex) 18533f72c5dSHong Zhang { 18633f72c5dSHong Zhang PetscErrorCode ierr; 18733f72c5dSHong Zhang 18833f72c5dSHong Zhang PetscFunctionBegin; 18933f72c5dSHong Zhang if (!Amat) PetscFunctionReturn(0); 19033f72c5dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 19133f72c5dSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 19233f72c5dSHong Zhang PetscValidHeaderSpecific(Udot,VEC_CLASSID,4); 19333f72c5dSHong Zhang 19433f72c5dSHong Zhang ierr = PetscLogEventBegin(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr); 19533f72c5dSHong Zhang if (ts->ijacobianp) { 19633f72c5dSHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 19733f72c5dSHong Zhang ierr = (*ts->ijacobianp)(ts,t,U,Udot,shift,Amat,ts->ijacobianpctx);CHKERRQ(ierr); 19833f72c5dSHong Zhang PetscStackPop; 19933f72c5dSHong Zhang } 20033f72c5dSHong Zhang if (imex) { 20133f72c5dSHong Zhang if (!ts->ijacobianp) { /* system was written as Udot = G(t,U) */ 20233f72c5dSHong Zhang PetscBool assembled; 20333f72c5dSHong Zhang ierr = MatZeroEntries(Amat);CHKERRQ(ierr); 20433f72c5dSHong Zhang ierr = MatAssembled(Amat,&assembled);CHKERRQ(ierr); 20533f72c5dSHong Zhang if (!assembled) { 20633f72c5dSHong Zhang ierr = MatAssemblyBegin(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20733f72c5dSHong Zhang ierr = MatAssemblyEnd(Amat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 20833f72c5dSHong Zhang } 20933f72c5dSHong Zhang } 21033f72c5dSHong Zhang } else { 21133f72c5dSHong Zhang if (ts->rhsjacobianp) { 21233f72c5dSHong Zhang ierr = TSComputeRHSJacobianP(ts,t,U,ts->Jacprhs);CHKERRQ(ierr); 21333f72c5dSHong Zhang } 21433f72c5dSHong Zhang if (ts->Jacprhs == Amat) { /* No IJacobian, so we only have the RHS matrix */ 21533f72c5dSHong Zhang ierr = MatScale(Amat,-1);CHKERRQ(ierr); 21633f72c5dSHong Zhang } else if (ts->Jacprhs) { /* Both IJacobian and RHSJacobian */ 21733f72c5dSHong Zhang MatStructure axpy = DIFFERENT_NONZERO_PATTERN; 21833f72c5dSHong Zhang if (!ts->ijacobianp) { /* No IJacobianp provided, but we have a separate RHS matrix */ 21933f72c5dSHong Zhang ierr = MatZeroEntries(Amat);CHKERRQ(ierr); 22033f72c5dSHong Zhang } 22133f72c5dSHong Zhang ierr = MatAXPY(Amat,-1,ts->Jacprhs,axpy);CHKERRQ(ierr); 22233f72c5dSHong Zhang } 22333f72c5dSHong Zhang } 22433f72c5dSHong Zhang ierr = PetscLogEventEnd(TS_JacobianPEval,ts,U,Amat,0);CHKERRQ(ierr); 22533f72c5dSHong Zhang PetscFunctionReturn(0); 22633f72c5dSHong Zhang } 22733f72c5dSHong Zhang 22833f72c5dSHong Zhang /*@C 229814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 230814e85d6SHong Zhang 231814e85d6SHong Zhang Logically Collective on TS 232814e85d6SHong Zhang 233814e85d6SHong Zhang Input Parameters: 234814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 235814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 236814e85d6SHong Zhang . costintegral - vector that stores the integral values 237814e85d6SHong Zhang . rf - routine for evaluating the integrand function 238c9ad9de0SHong Zhang . drduf - function that computes the gradients of the r's with respect to u 239814e85d6SHong 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) 240814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 241814e85d6SHong Zhang - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 242814e85d6SHong Zhang 243814e85d6SHong Zhang Calling sequence of rf: 244c9ad9de0SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx); 245814e85d6SHong Zhang 246c9ad9de0SHong Zhang Calling sequence of drduf: 247c9ad9de0SHong Zhang $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx); 248814e85d6SHong Zhang 249814e85d6SHong Zhang Calling sequence of drdpf: 250c9ad9de0SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx); 251814e85d6SHong Zhang 252cd4cee2dSHong Zhang Level: deprecated 253814e85d6SHong Zhang 254814e85d6SHong Zhang Notes: 255814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 256814e85d6SHong Zhang 257814e85d6SHong Zhang .keywords: TS, sensitivity analysis, timestep, set, quadrature, function 258814e85d6SHong Zhang 259814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients() 260814e85d6SHong Zhang @*/ 261814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*), 262c9ad9de0SHong Zhang PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*), 263814e85d6SHong Zhang PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*), 264814e85d6SHong Zhang PetscBool fwd,void *ctx) 265814e85d6SHong Zhang { 266814e85d6SHong Zhang PetscErrorCode ierr; 267814e85d6SHong Zhang 268814e85d6SHong Zhang PetscFunctionBegin; 269814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 270814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3); 271814e85d6SHong 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()"); 272814e85d6SHong Zhang if (!ts->numcost) ts->numcost=numcost; 273814e85d6SHong Zhang 274814e85d6SHong Zhang if (costintegral) { 275814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr); 276814e85d6SHong Zhang ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr); 277814e85d6SHong Zhang ts->vec_costintegral = costintegral; 278814e85d6SHong Zhang } else { 279814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 280814e85d6SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr); 281814e85d6SHong Zhang } else { 282814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr); 283814e85d6SHong Zhang } 284814e85d6SHong Zhang } 285814e85d6SHong Zhang if (!ts->vec_costintegrand) { 286814e85d6SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr); 287814e85d6SHong Zhang } else { 288814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr); 289814e85d6SHong Zhang } 290814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 291814e85d6SHong Zhang ts->costintegrand = rf; 292814e85d6SHong Zhang ts->costintegrandctx = ctx; 293c9ad9de0SHong Zhang ts->drdufunction = drduf; 294814e85d6SHong Zhang ts->drdpfunction = drdpf; 295814e85d6SHong Zhang PetscFunctionReturn(0); 296814e85d6SHong Zhang } 297814e85d6SHong Zhang 298b5b4017aSHong Zhang /*@C 299814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 300814e85d6SHong Zhang It is valid to call the routine after a backward run. 301814e85d6SHong Zhang 302814e85d6SHong Zhang Not Collective 303814e85d6SHong Zhang 304814e85d6SHong Zhang Input Parameter: 305814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 306814e85d6SHong Zhang 307814e85d6SHong Zhang Output Parameter: 308814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 309814e85d6SHong Zhang 310814e85d6SHong Zhang Level: intermediate 311814e85d6SHong Zhang 312814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 313814e85d6SHong Zhang 314814e85d6SHong Zhang .keywords: TS, sensitivity analysis 315814e85d6SHong Zhang @*/ 316814e85d6SHong Zhang PetscErrorCode TSGetCostIntegral(TS ts,Vec *v) 317814e85d6SHong Zhang { 318cd4cee2dSHong Zhang TS quadts; 319cd4cee2dSHong Zhang PetscErrorCode ierr; 320cd4cee2dSHong Zhang 321814e85d6SHong Zhang PetscFunctionBegin; 322814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 323814e85d6SHong Zhang PetscValidPointer(v,2); 324cd4cee2dSHong Zhang ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr); 325cd4cee2dSHong Zhang *v = quadts->vec_sol; 326814e85d6SHong Zhang PetscFunctionReturn(0); 327814e85d6SHong Zhang } 328814e85d6SHong Zhang 329b5b4017aSHong Zhang /*@C 330814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 331814e85d6SHong Zhang 332814e85d6SHong Zhang Input Parameters: 333814e85d6SHong Zhang + ts - the TS context 334814e85d6SHong Zhang . t - current time 335c9ad9de0SHong Zhang - U - state vector, i.e. current solution 336814e85d6SHong Zhang 337814e85d6SHong Zhang Output Parameter: 338c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs 339814e85d6SHong Zhang 340369cf35fSHong Zhang Notes: 341814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 342814e85d6SHong Zhang is used internally within the sensitivity analysis context. 343814e85d6SHong Zhang 344cd4cee2dSHong Zhang Level: deprecated 345814e85d6SHong Zhang 346814e85d6SHong Zhang .keywords: TS, compute 347814e85d6SHong Zhang 348814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 349814e85d6SHong Zhang @*/ 350c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q) 351814e85d6SHong Zhang { 352814e85d6SHong Zhang PetscErrorCode ierr; 353814e85d6SHong Zhang 354814e85d6SHong Zhang PetscFunctionBegin; 355814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 356c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 357c9ad9de0SHong Zhang PetscValidHeaderSpecific(Q,VEC_CLASSID,4); 358814e85d6SHong Zhang 359c9ad9de0SHong Zhang ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr); 360814e85d6SHong Zhang if (ts->costintegrand) { 361814e85d6SHong Zhang PetscStackPush("TS user integrand in the cost function"); 362c9ad9de0SHong Zhang ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr); 363814e85d6SHong Zhang PetscStackPop; 364814e85d6SHong Zhang } else { 365c9ad9de0SHong Zhang ierr = VecZeroEntries(Q);CHKERRQ(ierr); 366814e85d6SHong Zhang } 367814e85d6SHong Zhang 368c9ad9de0SHong Zhang ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr); 369814e85d6SHong Zhang PetscFunctionReturn(0); 370814e85d6SHong Zhang } 371814e85d6SHong Zhang 372b5b4017aSHong Zhang /*@C 373*d76be551SHong Zhang TSComputeDRDUFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 374814e85d6SHong Zhang 375*d76be551SHong Zhang Level: deprecated 376814e85d6SHong Zhang 377814e85d6SHong Zhang @*/ 378c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 379814e85d6SHong Zhang { 380814e85d6SHong Zhang PetscErrorCode ierr; 381814e85d6SHong Zhang 382814e85d6SHong Zhang PetscFunctionBegin; 38333f72c5dSHong Zhang if (!DRDU) PetscFunctionReturn(0); 384814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 385c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 386814e85d6SHong Zhang 387c9ad9de0SHong Zhang PetscStackPush("TS user DRDU function for sensitivity analysis"); 388c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 389814e85d6SHong Zhang PetscStackPop; 390814e85d6SHong Zhang PetscFunctionReturn(0); 391814e85d6SHong Zhang } 392814e85d6SHong Zhang 393b5b4017aSHong Zhang /*@C 394*d76be551SHong Zhang TSComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 395814e85d6SHong Zhang 396*d76be551SHong Zhang Level: deprecated 397814e85d6SHong Zhang 398814e85d6SHong Zhang @*/ 399c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 400814e85d6SHong Zhang { 401814e85d6SHong Zhang PetscErrorCode ierr; 402814e85d6SHong Zhang 403814e85d6SHong Zhang PetscFunctionBegin; 40433f72c5dSHong Zhang if (!DRDP) PetscFunctionReturn(0); 405814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 406c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 407814e85d6SHong Zhang 408814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 409c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 410814e85d6SHong Zhang PetscStackPop; 411814e85d6SHong Zhang PetscFunctionReturn(0); 412814e85d6SHong Zhang } 413814e85d6SHong Zhang 414b5b4017aSHong Zhang /*@C 415b5b4017aSHong 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. 416b5b4017aSHong Zhang 417b5b4017aSHong Zhang Logically Collective on TS 418b5b4017aSHong Zhang 419b5b4017aSHong Zhang Input Parameters: 420b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 421b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 422b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 423b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 424b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 425b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 426b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 427b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 428b5b4017aSHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for F_PP 429b5b4017aSHong Zhang 430b5b4017aSHong Zhang Calling sequence of ihessianproductfunc: 431369cf35fSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 432b5b4017aSHong Zhang + t - current timestep 433b5b4017aSHong Zhang . U - input vector (current ODE solution) 434369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 435b5b4017aSHong Zhang . Vr - input vector to be right-multiplied with the Hessian 436369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 437b5b4017aSHong Zhang - ctx - [optional] user-defined function context 438b5b4017aSHong Zhang 439b5b4017aSHong Zhang Level: intermediate 440b5b4017aSHong Zhang 441369cf35fSHong Zhang Notes: 442369cf35fSHong Zhang The first Hessian function and the working array are required. 443369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 444369cf35fSHong Zhang $ Vl_n^T*F_UP*Vr 445369cf35fSHong 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. 446369cf35fSHong Zhang Each entry of F_UP corresponds to the derivative 447369cf35fSHong Zhang $ F_UP[i][j][k] = \frac{\partial^2 F[i]}{\partial U[j] \partial P[k]}. 448369cf35fSHong 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 449369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * F_UP[i][j][k] * Vr[k]} 450369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 451b5b4017aSHong Zhang 452b5b4017aSHong Zhang .keywords: TS, sensitivity 453b5b4017aSHong Zhang 454b5b4017aSHong Zhang .seealso: 455b5b4017aSHong Zhang @*/ 456b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 457b5b4017aSHong Zhang Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 458b5b4017aSHong Zhang Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 459b5b4017aSHong Zhang Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 460b5b4017aSHong Zhang void *ctx) 461b5b4017aSHong Zhang { 462b5b4017aSHong Zhang PetscFunctionBegin; 463b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 464b5b4017aSHong Zhang PetscValidPointer(ihp1,2); 465b5b4017aSHong Zhang 466b5b4017aSHong Zhang ts->ihessianproductctx = ctx; 467b5b4017aSHong Zhang if (ihp1) ts->vecs_fuu = ihp1; 468b5b4017aSHong Zhang if (ihp2) ts->vecs_fup = ihp2; 469b5b4017aSHong Zhang if (ihp3) ts->vecs_fpu = ihp3; 470b5b4017aSHong Zhang if (ihp4) ts->vecs_fpp = ihp4; 471b5b4017aSHong Zhang ts->ihessianproduct_fuu = ihessianproductfunc1; 472b5b4017aSHong Zhang ts->ihessianproduct_fup = ihessianproductfunc2; 473b5b4017aSHong Zhang ts->ihessianproduct_fpu = ihessianproductfunc3; 474b5b4017aSHong Zhang ts->ihessianproduct_fpp = ihessianproductfunc4; 475b5b4017aSHong Zhang PetscFunctionReturn(0); 476b5b4017aSHong Zhang } 477b5b4017aSHong Zhang 478b5b4017aSHong Zhang /*@C 479b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Fuu. 480b5b4017aSHong Zhang 481b5b4017aSHong Zhang Collective on TS 482b5b4017aSHong Zhang 483b5b4017aSHong Zhang Input Parameters: 484b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 485b5b4017aSHong Zhang 486b5b4017aSHong Zhang Notes: 487b1e111ebSHong Zhang TSComputeIHessianProductFunctionUU() is typically used for sensitivity implementation, 488b5b4017aSHong Zhang so most users would not generally call this routine themselves. 489b5b4017aSHong Zhang 490b5b4017aSHong Zhang Level: developer 491b5b4017aSHong Zhang 492b5b4017aSHong Zhang .keywords: TS, sensitivity 493b5b4017aSHong Zhang 494b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 495b5b4017aSHong Zhang @*/ 496b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 497b5b4017aSHong Zhang { 498b5b4017aSHong Zhang PetscErrorCode ierr; 499b5b4017aSHong Zhang 500b5b4017aSHong Zhang PetscFunctionBegin; 50133f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 502b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 503b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 504b5b4017aSHong Zhang 50567633408SHong Zhang if (ts->ihessianproduct_fuu) { 506b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis"); 507b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 508b5b4017aSHong Zhang PetscStackPop; 50967633408SHong Zhang } 51067633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 51167633408SHong Zhang if (ts->rhshessianproduct_guu) { 51267633408SHong Zhang PetscInt nadj; 513b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 51467633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 51567633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 51667633408SHong Zhang } 51767633408SHong Zhang } 518b5b4017aSHong Zhang PetscFunctionReturn(0); 519b5b4017aSHong Zhang } 520b5b4017aSHong Zhang 521b5b4017aSHong Zhang /*@C 522b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Fup. 523b5b4017aSHong Zhang 524b5b4017aSHong Zhang Collective on TS 525b5b4017aSHong Zhang 526b5b4017aSHong Zhang Input Parameters: 527b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 528b5b4017aSHong Zhang 529b5b4017aSHong Zhang Notes: 530b1e111ebSHong Zhang TSComputeIHessianProductFunctionUP() is typically used for sensitivity implementation, 531b5b4017aSHong Zhang so most users would not generally call this routine themselves. 532b5b4017aSHong Zhang 533b5b4017aSHong Zhang Level: developer 534b5b4017aSHong Zhang 535b5b4017aSHong Zhang .keywords: TS, sensitivity 536b5b4017aSHong Zhang 537b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 538b5b4017aSHong Zhang @*/ 539b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 540b5b4017aSHong Zhang { 541b5b4017aSHong Zhang PetscErrorCode ierr; 542b5b4017aSHong Zhang 543b5b4017aSHong Zhang PetscFunctionBegin; 54433f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 545b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 546b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 547b5b4017aSHong Zhang 54867633408SHong Zhang if (ts->ihessianproduct_fup) { 549b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis"); 550b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 551b5b4017aSHong Zhang PetscStackPop; 55267633408SHong Zhang } 55367633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 55467633408SHong Zhang if (ts->rhshessianproduct_gup) { 55567633408SHong Zhang PetscInt nadj; 556b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionUP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 55767633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 55867633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 55967633408SHong Zhang } 56067633408SHong Zhang } 561b5b4017aSHong Zhang PetscFunctionReturn(0); 562b5b4017aSHong Zhang } 563b5b4017aSHong Zhang 564b5b4017aSHong Zhang /*@C 565b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Fpu. 566b5b4017aSHong Zhang 567b5b4017aSHong Zhang Collective on TS 568b5b4017aSHong Zhang 569b5b4017aSHong Zhang Input Parameters: 570b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 571b5b4017aSHong Zhang 572b5b4017aSHong Zhang Notes: 573b1e111ebSHong Zhang TSComputeIHessianProductFunctionPU() is typically used for sensitivity implementation, 574b5b4017aSHong Zhang so most users would not generally call this routine themselves. 575b5b4017aSHong Zhang 576b5b4017aSHong Zhang Level: developer 577b5b4017aSHong Zhang 578b5b4017aSHong Zhang .keywords: TS, sensitivity 579b5b4017aSHong Zhang 580b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 581b5b4017aSHong Zhang @*/ 582b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 583b5b4017aSHong Zhang { 584b5b4017aSHong Zhang PetscErrorCode ierr; 585b5b4017aSHong Zhang 586b5b4017aSHong Zhang PetscFunctionBegin; 58733f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 588b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 589b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 590b5b4017aSHong Zhang 59167633408SHong Zhang if (ts->ihessianproduct_fpu) { 592b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 593b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 594b5b4017aSHong Zhang PetscStackPop; 59567633408SHong Zhang } 59667633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 59767633408SHong Zhang if (ts->rhshessianproduct_gpu) { 59867633408SHong Zhang PetscInt nadj; 599b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPU(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 60067633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 60167633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 60267633408SHong Zhang } 60367633408SHong Zhang } 604b5b4017aSHong Zhang PetscFunctionReturn(0); 605b5b4017aSHong Zhang } 606b5b4017aSHong Zhang 607b5b4017aSHong Zhang /*@C 608b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Fpp. 609b5b4017aSHong Zhang 610b5b4017aSHong Zhang Collective on TS 611b5b4017aSHong Zhang 612b5b4017aSHong Zhang Input Parameters: 613b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 614b5b4017aSHong Zhang 615b5b4017aSHong Zhang Notes: 616b1e111ebSHong Zhang TSComputeIHessianProductFunctionPP() is typically used for sensitivity implementation, 617b5b4017aSHong Zhang so most users would not generally call this routine themselves. 618b5b4017aSHong Zhang 619b5b4017aSHong Zhang Level: developer 620b5b4017aSHong Zhang 621b5b4017aSHong Zhang .keywords: TS, sensitivity 622b5b4017aSHong Zhang 623b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 624b5b4017aSHong Zhang @*/ 625b1e111ebSHong Zhang PetscErrorCode TSComputeIHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 626b5b4017aSHong Zhang { 627b5b4017aSHong Zhang PetscErrorCode ierr; 628b5b4017aSHong Zhang 629b5b4017aSHong Zhang PetscFunctionBegin; 63033f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 631b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 632b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 633b5b4017aSHong Zhang 63467633408SHong Zhang if (ts->ihessianproduct_fpp) { 635b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 636b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 637b5b4017aSHong Zhang PetscStackPop; 63867633408SHong Zhang } 63967633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 64067633408SHong Zhang if (ts->rhshessianproduct_gpp) { 64167633408SHong Zhang PetscInt nadj; 642b1e111ebSHong Zhang ierr = TSComputeRHSHessianProductFunctionPP(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 64367633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 64467633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 64567633408SHong Zhang } 64667633408SHong Zhang } 647b5b4017aSHong Zhang PetscFunctionReturn(0); 648b5b4017aSHong Zhang } 649b5b4017aSHong Zhang 65013af1a74SHong Zhang /*@C 65113af1a74SHong Zhang TSSetRHSHessianProduct - Sets the function that computes the vecotr-Hessian-vector product. The Hessian is the second-order derivative of G (RHSFunction) w.r.t. the state variable. 65213af1a74SHong Zhang 65313af1a74SHong Zhang Logically Collective on TS 65413af1a74SHong Zhang 65513af1a74SHong Zhang Input Parameters: 65613af1a74SHong Zhang + ts - TS context obtained from TSCreate() 65713af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 65813af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 65913af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 66013af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 66113af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 66213af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 66313af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 66413af1a74SHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for G_PP 66513af1a74SHong Zhang 66613af1a74SHong Zhang Calling sequence of ihessianproductfunc: 667369cf35fSHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV,void *ctx); 66813af1a74SHong Zhang + t - current timestep 66913af1a74SHong Zhang . U - input vector (current ODE solution) 670369cf35fSHong Zhang . Vl - an array of input vectors to be left-multiplied with the Hessian 67113af1a74SHong Zhang . Vr - input vector to be right-multiplied with the Hessian 672369cf35fSHong Zhang . VHV - an array of output vectors for vector-Hessian-vector product 67313af1a74SHong Zhang - ctx - [optional] user-defined function context 67413af1a74SHong Zhang 67513af1a74SHong Zhang Level: intermediate 67613af1a74SHong Zhang 677369cf35fSHong Zhang Notes: 678369cf35fSHong Zhang The first Hessian function and the working array are required. 679369cf35fSHong Zhang As an example to implement the callback functions, the second callback function calculates the vector-Hessian-vector product 680369cf35fSHong Zhang $ Vl_n^T*G_UP*Vr 681369cf35fSHong 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. 682369cf35fSHong Zhang Each entry of G_UP corresponds to the derivative 683369cf35fSHong Zhang $ G_UP[i][j][k] = \frac{\partial^2 G[i]}{\partial U[j] \partial P[k]}. 684369cf35fSHong 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 685369cf35fSHong Zhang $ VHV_n[j] = \sum_i \sum_k {Vl_n[i] * G_UP[i][j][k] * Vr[k]} 686369cf35fSHong Zhang If the cost function is a scalar, there will be only one vector in Vl and VHV. 68713af1a74SHong Zhang 68813af1a74SHong Zhang .keywords: TS, sensitivity 68913af1a74SHong Zhang 69013af1a74SHong Zhang .seealso: 69113af1a74SHong Zhang @*/ 69213af1a74SHong Zhang PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 69313af1a74SHong Zhang Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 69413af1a74SHong Zhang Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 69513af1a74SHong Zhang Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 69613af1a74SHong Zhang void *ctx) 69713af1a74SHong Zhang { 69813af1a74SHong Zhang PetscFunctionBegin; 69913af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 70013af1a74SHong Zhang PetscValidPointer(rhshp1,2); 70113af1a74SHong Zhang 70213af1a74SHong Zhang ts->rhshessianproductctx = ctx; 70313af1a74SHong Zhang if (rhshp1) ts->vecs_guu = rhshp1; 70413af1a74SHong Zhang if (rhshp2) ts->vecs_gup = rhshp2; 70513af1a74SHong Zhang if (rhshp3) ts->vecs_gpu = rhshp3; 70613af1a74SHong Zhang if (rhshp4) ts->vecs_gpp = rhshp4; 70713af1a74SHong Zhang ts->rhshessianproduct_guu = rhshessianproductfunc1; 70813af1a74SHong Zhang ts->rhshessianproduct_gup = rhshessianproductfunc2; 70913af1a74SHong Zhang ts->rhshessianproduct_gpu = rhshessianproductfunc3; 71013af1a74SHong Zhang ts->rhshessianproduct_gpp = rhshessianproductfunc4; 71113af1a74SHong Zhang PetscFunctionReturn(0); 71213af1a74SHong Zhang } 71313af1a74SHong Zhang 71413af1a74SHong Zhang /*@C 715b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU - Runs the user-defined vector-Hessian-vector product function for Guu. 71613af1a74SHong Zhang 71713af1a74SHong Zhang Collective on TS 71813af1a74SHong Zhang 71913af1a74SHong Zhang Input Parameters: 72013af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 72113af1a74SHong Zhang 72213af1a74SHong Zhang Notes: 723b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUU() is typically used for sensitivity implementation, 72413af1a74SHong Zhang so most users would not generally call this routine themselves. 72513af1a74SHong Zhang 72613af1a74SHong Zhang Level: developer 72713af1a74SHong Zhang 72813af1a74SHong Zhang .keywords: TS, sensitivity 72913af1a74SHong Zhang 73013af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 73113af1a74SHong Zhang @*/ 732b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 73313af1a74SHong Zhang { 73413af1a74SHong Zhang PetscErrorCode ierr; 73513af1a74SHong Zhang 73613af1a74SHong Zhang PetscFunctionBegin; 73713af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 73813af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 73913af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 74013af1a74SHong Zhang 74113af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis"); 74213af1a74SHong Zhang ierr = (*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 74313af1a74SHong Zhang PetscStackPop; 74413af1a74SHong Zhang PetscFunctionReturn(0); 74513af1a74SHong Zhang } 74613af1a74SHong Zhang 74713af1a74SHong Zhang /*@C 748b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP - Runs the user-defined vector-Hessian-vector product function for Gup. 74913af1a74SHong Zhang 75013af1a74SHong Zhang Collective on TS 75113af1a74SHong Zhang 75213af1a74SHong Zhang Input Parameters: 75313af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 75413af1a74SHong Zhang 75513af1a74SHong Zhang Notes: 756b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionUP() is typically used for sensitivity implementation, 75713af1a74SHong Zhang so most users would not generally call this routine themselves. 75813af1a74SHong Zhang 75913af1a74SHong Zhang Level: developer 76013af1a74SHong Zhang 76113af1a74SHong Zhang .keywords: TS, sensitivity 76213af1a74SHong Zhang 76313af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 76413af1a74SHong Zhang @*/ 765b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionUP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 76613af1a74SHong Zhang { 76713af1a74SHong Zhang PetscErrorCode ierr; 76813af1a74SHong Zhang 76913af1a74SHong Zhang PetscFunctionBegin; 77013af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 77113af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 77213af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 77313af1a74SHong Zhang 77413af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis"); 77513af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 77613af1a74SHong Zhang PetscStackPop; 77713af1a74SHong Zhang PetscFunctionReturn(0); 77813af1a74SHong Zhang } 77913af1a74SHong Zhang 78013af1a74SHong Zhang /*@C 781b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU - Runs the user-defined vector-Hessian-vector product function for Gpu. 78213af1a74SHong Zhang 78313af1a74SHong Zhang Collective on TS 78413af1a74SHong Zhang 78513af1a74SHong Zhang Input Parameters: 78613af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 78713af1a74SHong Zhang 78813af1a74SHong Zhang Notes: 789b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPU() is typically used for sensitivity implementation, 79013af1a74SHong Zhang so most users would not generally call this routine themselves. 79113af1a74SHong Zhang 79213af1a74SHong Zhang Level: developer 79313af1a74SHong Zhang 79413af1a74SHong Zhang .keywords: TS, sensitivity 79513af1a74SHong Zhang 79613af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 79713af1a74SHong Zhang @*/ 798b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPU(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 79913af1a74SHong Zhang { 80013af1a74SHong Zhang PetscErrorCode ierr; 80113af1a74SHong Zhang 80213af1a74SHong Zhang PetscFunctionBegin; 80313af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 80413af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 80513af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 80613af1a74SHong Zhang 80713af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis"); 80813af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 80913af1a74SHong Zhang PetscStackPop; 81013af1a74SHong Zhang PetscFunctionReturn(0); 81113af1a74SHong Zhang } 81213af1a74SHong Zhang 81313af1a74SHong Zhang /*@C 814b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP - Runs the user-defined vector-Hessian-vector product function for Gpp. 81513af1a74SHong Zhang 81613af1a74SHong Zhang Collective on TS 81713af1a74SHong Zhang 81813af1a74SHong Zhang Input Parameters: 81913af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 82013af1a74SHong Zhang 82113af1a74SHong Zhang Notes: 822b1e111ebSHong Zhang TSComputeRHSHessianProductFunctionPP() is typically used for sensitivity implementation, 82313af1a74SHong Zhang so most users would not generally call this routine themselves. 82413af1a74SHong Zhang 82513af1a74SHong Zhang Level: developer 82613af1a74SHong Zhang 82713af1a74SHong Zhang .keywords: TS, sensitivity 82813af1a74SHong Zhang 82913af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 83013af1a74SHong Zhang @*/ 831b1e111ebSHong Zhang PetscErrorCode TSComputeRHSHessianProductFunctionPP(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 83213af1a74SHong Zhang { 83313af1a74SHong Zhang PetscErrorCode ierr; 83413af1a74SHong Zhang 83513af1a74SHong Zhang PetscFunctionBegin; 83613af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 83713af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 83813af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 83913af1a74SHong Zhang 84013af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis"); 84113af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 84213af1a74SHong Zhang PetscStackPop; 84313af1a74SHong Zhang PetscFunctionReturn(0); 84413af1a74SHong Zhang } 84513af1a74SHong Zhang 846814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 847814e85d6SHong Zhang 848814e85d6SHong Zhang /*@ 849814e85d6SHong 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 850814e85d6SHong Zhang for use by the TSAdjoint routines. 851814e85d6SHong Zhang 852814e85d6SHong Zhang Logically Collective on TS and Vec 853814e85d6SHong Zhang 854814e85d6SHong Zhang Input Parameters: 855814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 856814e85d6SHong 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 857814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 858814e85d6SHong Zhang 859814e85d6SHong Zhang Level: beginner 860814e85d6SHong Zhang 861814e85d6SHong Zhang Notes: 862814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 863814e85d6SHong Zhang 864814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 865814e85d6SHong Zhang 866814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values 867b5b4017aSHong Zhang 868b5b4017aSHong Zhang .seealso TSGetCostGradients() 869814e85d6SHong Zhang @*/ 870814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu) 871814e85d6SHong Zhang { 872814e85d6SHong Zhang PetscFunctionBegin; 873814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 874814e85d6SHong Zhang PetscValidPointer(lambda,2); 875814e85d6SHong Zhang ts->vecs_sensi = lambda; 876814e85d6SHong Zhang ts->vecs_sensip = mu; 877814e85d6SHong 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"); 878814e85d6SHong Zhang ts->numcost = numcost; 879814e85d6SHong Zhang PetscFunctionReturn(0); 880814e85d6SHong Zhang } 881814e85d6SHong Zhang 882814e85d6SHong Zhang /*@ 883814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 884814e85d6SHong Zhang 885814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 886814e85d6SHong Zhang 887814e85d6SHong Zhang Input Parameter: 888814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 889814e85d6SHong Zhang 890814e85d6SHong Zhang Output Parameter: 891814e85d6SHong Zhang + lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 892814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 893814e85d6SHong Zhang 894814e85d6SHong Zhang Level: intermediate 895814e85d6SHong Zhang 896814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity 897b5b4017aSHong Zhang 898b5b4017aSHong Zhang .seealso: TSSetCostGradients() 899814e85d6SHong Zhang @*/ 900814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu) 901814e85d6SHong Zhang { 902814e85d6SHong Zhang PetscFunctionBegin; 903814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 904814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 905814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 906814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 907814e85d6SHong Zhang PetscFunctionReturn(0); 908814e85d6SHong Zhang } 909814e85d6SHong Zhang 910814e85d6SHong Zhang /*@ 911b5b4017aSHong 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 912b5b4017aSHong Zhang for use by the TSAdjoint routines. 913b5b4017aSHong Zhang 914b5b4017aSHong Zhang Logically Collective on TS and Vec 915b5b4017aSHong Zhang 916b5b4017aSHong Zhang Input Parameters: 917b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 918b5b4017aSHong Zhang . numcost - number of cost functions 919b5b4017aSHong 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 920b5b4017aSHong 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 921b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 922b5b4017aSHong Zhang 923b5b4017aSHong Zhang Level: beginner 924b5b4017aSHong Zhang 925b5b4017aSHong Zhang Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system 926b5b4017aSHong Zhang 927ecf68647SHong Zhang For second-order adjoint, one needs to call this function and then TSAdjointSetForward() before TSSolve(). 928b5b4017aSHong Zhang 929b5b4017aSHong 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. 930b5b4017aSHong Zhang 9313fca9621SHong Zhang Passing NULL for lambda2 disables the second-order calculation. 932b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint 933b5b4017aSHong Zhang 934ecf68647SHong Zhang .seealso: TSAdjointSetForward() 935b5b4017aSHong Zhang @*/ 936b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir) 937b5b4017aSHong Zhang { 938b5b4017aSHong Zhang PetscFunctionBegin; 939b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 940b5b4017aSHong 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"); 941b5b4017aSHong Zhang ts->numcost = numcost; 942b5b4017aSHong Zhang ts->vecs_sensi2 = lambda2; 94333f72c5dSHong Zhang ts->vecs_sensi2p = mu2; 944b5b4017aSHong Zhang ts->vec_dir = dir; 945b5b4017aSHong Zhang PetscFunctionReturn(0); 946b5b4017aSHong Zhang } 947b5b4017aSHong Zhang 948b5b4017aSHong Zhang /*@ 949b5b4017aSHong Zhang TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve() 950b5b4017aSHong Zhang 951b5b4017aSHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 952b5b4017aSHong Zhang 953b5b4017aSHong Zhang Input Parameter: 954b5b4017aSHong Zhang . ts - the TS context obtained from TSCreate() 955b5b4017aSHong Zhang 956b5b4017aSHong Zhang Output Parameter: 957b5b4017aSHong Zhang + numcost - number of cost functions 958b5b4017aSHong 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 959b5b4017aSHong 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 960b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 961b5b4017aSHong Zhang 962b5b4017aSHong Zhang Level: intermediate 963b5b4017aSHong Zhang 964b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint 965b5b4017aSHong Zhang 966b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts() 967b5b4017aSHong Zhang @*/ 968b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir) 969b5b4017aSHong Zhang { 970b5b4017aSHong Zhang PetscFunctionBegin; 971b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 972b5b4017aSHong Zhang if (numcost) *numcost = ts->numcost; 973b5b4017aSHong Zhang if (lambda2) *lambda2 = ts->vecs_sensi2; 97433f72c5dSHong Zhang if (mu2) *mu2 = ts->vecs_sensi2p; 975b5b4017aSHong Zhang if (dir) *dir = ts->vec_dir; 976b5b4017aSHong Zhang PetscFunctionReturn(0); 977b5b4017aSHong Zhang } 978b5b4017aSHong Zhang 979b5b4017aSHong Zhang /*@ 980ecf68647SHong Zhang TSAdjointSetForward - Trigger the tangent linear solver and initialize the forward sensitivities 981b5b4017aSHong Zhang 982b5b4017aSHong Zhang Logically Collective on TS and Mat 983b5b4017aSHong Zhang 984b5b4017aSHong Zhang Input Parameters: 985b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 986b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters 987b5b4017aSHong Zhang 988b5b4017aSHong Zhang Level: intermediate 989b5b4017aSHong Zhang 990ecf68647SHong 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. 991b5b4017aSHong Zhang 992b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint 993b5b4017aSHong Zhang 994ecf68647SHong Zhang .seealso: TSSetCostHessianProducts(), TSAdjointResetForward() 995b5b4017aSHong Zhang @*/ 996ecf68647SHong Zhang PetscErrorCode TSAdjointSetForward(TS ts,Mat didp) 997b5b4017aSHong Zhang { 99833f72c5dSHong Zhang Mat A; 99933f72c5dSHong Zhang Vec sp; 100033f72c5dSHong Zhang PetscScalar *xarr; 100133f72c5dSHong Zhang PetscInt lsize; 1002b5b4017aSHong Zhang PetscErrorCode ierr; 1003b5b4017aSHong Zhang 1004b5b4017aSHong Zhang PetscFunctionBegin; 1005b5b4017aSHong Zhang ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 1006b5b4017aSHong Zhang if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first"); 100733f72c5dSHong Zhang if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 100833f72c5dSHong Zhang /* create a single-column dense matrix */ 100933f72c5dSHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&lsize);CHKERRQ(ierr); 1010f63bf25fSHong Zhang ierr = MatCreateDense(PetscObjectComm((PetscObject)ts),lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr); 101133f72c5dSHong Zhang 101233f72c5dSHong Zhang ierr = VecDuplicate(ts->vec_sol,&sp);CHKERRQ(ierr); 101333f72c5dSHong Zhang ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr); 101433f72c5dSHong Zhang ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr); 1015ecf68647SHong Zhang if (ts->vecs_sensi2p) { /* tangent linear variable initialized as 2*dIdP*dir */ 101633f72c5dSHong Zhang if (didp) { 101733f72c5dSHong Zhang ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr); 101833f72c5dSHong Zhang ierr = VecScale(sp,2.);CHKERRQ(ierr); 101933f72c5dSHong Zhang } else { 102033f72c5dSHong Zhang ierr = VecZeroEntries(sp);CHKERRQ(ierr); 102133f72c5dSHong Zhang } 1022ecf68647SHong Zhang } else { /* tangent linear variable initialized as dir */ 102333f72c5dSHong Zhang ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr); 102433f72c5dSHong Zhang } 102533f72c5dSHong Zhang ierr = VecResetArray(sp);CHKERRQ(ierr); 102633f72c5dSHong Zhang ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr); 102733f72c5dSHong Zhang ierr = VecDestroy(&sp);CHKERRQ(ierr); 102833f72c5dSHong Zhang 102933f72c5dSHong Zhang ierr = TSForwardSetInitialSensitivities(ts,A);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */ 103033f72c5dSHong Zhang 103133f72c5dSHong Zhang ierr = MatDestroy(&A);CHKERRQ(ierr); 1032b5b4017aSHong Zhang PetscFunctionReturn(0); 1033b5b4017aSHong Zhang } 1034b5b4017aSHong Zhang 1035b5b4017aSHong Zhang /*@ 1036ecf68647SHong Zhang TSAdjointResetForward - Reset the tangent linear solver and destroy the tangent linear context 1037ecf68647SHong Zhang 1038ecf68647SHong Zhang Logically Collective on TS and Mat 1039ecf68647SHong Zhang 1040ecf68647SHong Zhang Input Parameters: 1041ecf68647SHong Zhang . ts - the TS context obtained from TSCreate() 1042ecf68647SHong Zhang 1043ecf68647SHong Zhang Level: intermediate 1044ecf68647SHong Zhang 1045ecf68647SHong Zhang .keywords: TS, sensitivity, second-order adjoint 1046ecf68647SHong Zhang 1047ecf68647SHong Zhang .seealso: TSAdjointSetForward() 1048ecf68647SHong Zhang @*/ 1049ecf68647SHong Zhang PetscErrorCode TSAdjointResetForward(TS ts) 1050ecf68647SHong Zhang { 1051ecf68647SHong Zhang PetscErrorCode ierr; 1052ecf68647SHong Zhang 1053ecf68647SHong Zhang PetscFunctionBegin; 1054ecf68647SHong Zhang ts->forward_solve = PETSC_FALSE; /* turn off tangent linear mode */ 1055ecf68647SHong Zhang ierr = TSForwardReset(ts);CHKERRQ(ierr); 1056ecf68647SHong Zhang PetscFunctionReturn(0); 1057ecf68647SHong Zhang } 1058ecf68647SHong Zhang 1059ecf68647SHong Zhang /*@ 1060814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 1061814e85d6SHong Zhang of an adjoint solver 1062814e85d6SHong Zhang 1063814e85d6SHong Zhang Collective on TS 1064814e85d6SHong Zhang 1065814e85d6SHong Zhang Input Parameter: 1066814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1067814e85d6SHong Zhang 1068814e85d6SHong Zhang Level: advanced 1069814e85d6SHong Zhang 1070814e85d6SHong Zhang .keywords: TS, timestep, setup 1071814e85d6SHong Zhang 1072814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients() 1073814e85d6SHong Zhang @*/ 1074814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts) 1075814e85d6SHong Zhang { 1076881c1a9bSHong Zhang TSTrajectory tj; 1077881c1a9bSHong Zhang PetscBool match; 1078814e85d6SHong Zhang PetscErrorCode ierr; 1079814e85d6SHong Zhang 1080814e85d6SHong Zhang PetscFunctionBegin; 1081814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1082814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 1083814e85d6SHong Zhang if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 108433f72c5dSHong Zhang if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 1085881c1a9bSHong Zhang ierr = TSGetTrajectory(ts,&tj);CHKERRQ(ierr); 10868e224c67SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)tj,TSTRAJECTORYBASIC,&match);CHKERRQ(ierr); 1087881c1a9bSHong Zhang if (match) { 1088881c1a9bSHong Zhang PetscBool solution_only; 1089881c1a9bSHong Zhang ierr = TSTrajectoryGetSolutionOnly(tj,&solution_only);CHKERRQ(ierr); 1090881c1a9bSHong 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"); 1091881c1a9bSHong Zhang } 10929ffb3502SHong Zhang ierr = TSTrajectorySetUseHistory(tj,PETSC_FALSE);CHKERRQ(ierr); /* not use TSHistory */ 109333f72c5dSHong Zhang 1094cd4cee2dSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 1095cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);CHKERRQ(ierr); 1096814e85d6SHong Zhang if (ts->vecs_sensip){ 1097cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr); 1098814e85d6SHong Zhang } 1099814e85d6SHong Zhang } 1100814e85d6SHong Zhang 1101814e85d6SHong Zhang if (ts->ops->adjointsetup) { 1102814e85d6SHong Zhang ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr); 1103814e85d6SHong Zhang } 1104814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 1105814e85d6SHong Zhang PetscFunctionReturn(0); 1106814e85d6SHong Zhang } 1107814e85d6SHong Zhang 1108814e85d6SHong Zhang /*@ 1109ece46509SHong Zhang TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 1110ece46509SHong Zhang 1111ece46509SHong Zhang Collective on TS 1112ece46509SHong Zhang 1113ece46509SHong Zhang Input Parameter: 1114ece46509SHong Zhang . ts - the TS context obtained from TSCreate() 1115ece46509SHong Zhang 1116ece46509SHong Zhang Level: beginner 1117ece46509SHong Zhang 1118ece46509SHong Zhang .keywords: TS, timestep, reset 1119ece46509SHong Zhang 11209ffb3502SHong Zhang .seealso: TSCreate(), TSAdjointSetUp(), TSADestroy() 1121ece46509SHong Zhang @*/ 1122ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts) 1123ece46509SHong Zhang { 1124ece46509SHong Zhang PetscErrorCode ierr; 1125ece46509SHong Zhang 1126ece46509SHong Zhang PetscFunctionBegin; 1127ece46509SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1128ece46509SHong Zhang if (ts->ops->adjointreset) { 1129ece46509SHong Zhang ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr); 1130ece46509SHong Zhang } 11317207e0fdSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 11327207e0fdSHong Zhang ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr); 11337207e0fdSHong Zhang if (ts->vecs_sensip){ 11347207e0fdSHong Zhang ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr); 11357207e0fdSHong Zhang } 11367207e0fdSHong Zhang } 1137ece46509SHong Zhang ts->vecs_sensi = NULL; 1138ece46509SHong Zhang ts->vecs_sensip = NULL; 1139ece46509SHong Zhang ts->vecs_sensi2 = NULL; 114033f72c5dSHong Zhang ts->vecs_sensi2p = NULL; 1141ece46509SHong Zhang ts->vec_dir = NULL; 1142ece46509SHong Zhang ts->adjointsetupcalled = PETSC_FALSE; 1143ece46509SHong Zhang PetscFunctionReturn(0); 1144ece46509SHong Zhang } 1145ece46509SHong Zhang 1146ece46509SHong Zhang /*@ 1147814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1148814e85d6SHong Zhang 1149814e85d6SHong Zhang Logically Collective on TS 1150814e85d6SHong Zhang 1151814e85d6SHong Zhang Input Parameters: 1152814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1153814e85d6SHong Zhang . steps - number of steps to use 1154814e85d6SHong Zhang 1155814e85d6SHong Zhang Level: intermediate 1156814e85d6SHong Zhang 1157814e85d6SHong Zhang Notes: 1158814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 1159814e85d6SHong Zhang so as to integrate back to less than the original timestep 1160814e85d6SHong Zhang 1161814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations 1162814e85d6SHong Zhang 1163814e85d6SHong Zhang .seealso: TSSetExactFinalTime() 1164814e85d6SHong Zhang @*/ 1165814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 1166814e85d6SHong Zhang { 1167814e85d6SHong Zhang PetscFunctionBegin; 1168814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1169814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts,steps,2); 1170814e85d6SHong Zhang if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 1171814e85d6SHong Zhang if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 1172814e85d6SHong Zhang ts->adjoint_max_steps = steps; 1173814e85d6SHong Zhang PetscFunctionReturn(0); 1174814e85d6SHong Zhang } 1175814e85d6SHong Zhang 1176814e85d6SHong Zhang /*@C 1177814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1178814e85d6SHong Zhang 1179814e85d6SHong Zhang Level: deprecated 1180814e85d6SHong Zhang 1181814e85d6SHong Zhang @*/ 1182814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 1183814e85d6SHong Zhang { 1184814e85d6SHong Zhang PetscErrorCode ierr; 1185814e85d6SHong Zhang 1186814e85d6SHong Zhang PetscFunctionBegin; 1187814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1188814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 1189814e85d6SHong Zhang 1190814e85d6SHong Zhang ts->rhsjacobianp = func; 1191814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 1192814e85d6SHong Zhang if(Amat) { 1193814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 1194814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 1195814e85d6SHong Zhang ts->Jacp = Amat; 1196814e85d6SHong Zhang } 1197814e85d6SHong Zhang PetscFunctionReturn(0); 1198814e85d6SHong Zhang } 1199814e85d6SHong Zhang 1200814e85d6SHong Zhang /*@C 1201814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1202814e85d6SHong Zhang 1203814e85d6SHong Zhang Level: deprecated 1204814e85d6SHong Zhang 1205814e85d6SHong Zhang @*/ 1206c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat) 1207814e85d6SHong Zhang { 1208814e85d6SHong Zhang PetscErrorCode ierr; 1209814e85d6SHong Zhang 1210814e85d6SHong Zhang PetscFunctionBegin; 1211814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1212c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1213814e85d6SHong Zhang PetscValidPointer(Amat,4); 1214814e85d6SHong Zhang 1215814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 1216c9ad9de0SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 1217814e85d6SHong Zhang PetscStackPop; 1218814e85d6SHong Zhang PetscFunctionReturn(0); 1219814e85d6SHong Zhang } 1220814e85d6SHong Zhang 1221814e85d6SHong Zhang /*@ 1222*d76be551SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobian() 1223814e85d6SHong Zhang 1224814e85d6SHong Zhang Level: deprecated 1225814e85d6SHong Zhang 1226814e85d6SHong Zhang @*/ 1227c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 1228814e85d6SHong Zhang { 1229814e85d6SHong Zhang PetscErrorCode ierr; 1230814e85d6SHong Zhang 1231814e85d6SHong Zhang PetscFunctionBegin; 1232814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1233c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1234814e85d6SHong Zhang 1235814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 1236c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 1237814e85d6SHong Zhang PetscStackPop; 1238814e85d6SHong Zhang PetscFunctionReturn(0); 1239814e85d6SHong Zhang } 1240814e85d6SHong Zhang 1241814e85d6SHong Zhang /*@ 1242*d76be551SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSGetQuadratureTS() then TSComputeRHSJacobianP() 1243814e85d6SHong Zhang 1244814e85d6SHong Zhang Level: deprecated 1245814e85d6SHong Zhang 1246814e85d6SHong Zhang @*/ 1247c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 1248814e85d6SHong Zhang { 1249814e85d6SHong Zhang PetscErrorCode ierr; 1250814e85d6SHong Zhang 1251814e85d6SHong Zhang PetscFunctionBegin; 1252814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1253c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1254814e85d6SHong Zhang 1255814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 1256c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 1257814e85d6SHong Zhang PetscStackPop; 1258814e85d6SHong Zhang PetscFunctionReturn(0); 1259814e85d6SHong Zhang } 1260814e85d6SHong Zhang 1261814e85d6SHong Zhang /*@C 1262814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 1263814e85d6SHong Zhang 1264814e85d6SHong Zhang Level: intermediate 1265814e85d6SHong Zhang 1266814e85d6SHong Zhang .keywords: TS, set, monitor 1267814e85d6SHong Zhang 1268814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1269814e85d6SHong Zhang @*/ 1270814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1271814e85d6SHong Zhang { 1272814e85d6SHong Zhang PetscErrorCode ierr; 1273814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1274814e85d6SHong Zhang 1275814e85d6SHong Zhang PetscFunctionBegin; 1276814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1277814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1278814e85d6SHong Zhang ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 1279814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1280814e85d6SHong Zhang PetscFunctionReturn(0); 1281814e85d6SHong Zhang } 1282814e85d6SHong Zhang 1283814e85d6SHong Zhang /*@C 1284814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1285814e85d6SHong Zhang 1286814e85d6SHong Zhang Collective on TS 1287814e85d6SHong Zhang 1288814e85d6SHong Zhang Input Parameters: 1289814e85d6SHong Zhang + ts - TS object you wish to monitor 1290814e85d6SHong Zhang . name - the monitor type one is seeking 1291814e85d6SHong Zhang . help - message indicating what monitoring is done 1292814e85d6SHong Zhang . manual - manual page for the monitor 1293814e85d6SHong Zhang . monitor - the monitor function 1294814e85d6SHong 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 1295814e85d6SHong Zhang 1296814e85d6SHong Zhang Level: developer 1297814e85d6SHong Zhang 1298814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 1299814e85d6SHong Zhang PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1300814e85d6SHong Zhang PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 1301814e85d6SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1302814e85d6SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1303814e85d6SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1304814e85d6SHong Zhang PetscOptionsFList(), PetscOptionsEList() 1305814e85d6SHong Zhang @*/ 1306814e85d6SHong 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*)) 1307814e85d6SHong Zhang { 1308814e85d6SHong Zhang PetscErrorCode ierr; 1309814e85d6SHong Zhang PetscViewer viewer; 1310814e85d6SHong Zhang PetscViewerFormat format; 1311814e85d6SHong Zhang PetscBool flg; 1312814e85d6SHong Zhang 1313814e85d6SHong Zhang PetscFunctionBegin; 131416413a6aSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 1315814e85d6SHong Zhang if (flg) { 1316814e85d6SHong Zhang PetscViewerAndFormat *vf; 1317814e85d6SHong Zhang ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 1318814e85d6SHong Zhang ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 1319814e85d6SHong Zhang if (monitorsetup) { 1320814e85d6SHong Zhang ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 1321814e85d6SHong Zhang } 1322814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 1323814e85d6SHong Zhang } 1324814e85d6SHong Zhang PetscFunctionReturn(0); 1325814e85d6SHong Zhang } 1326814e85d6SHong Zhang 1327814e85d6SHong Zhang /*@C 1328814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1329814e85d6SHong Zhang timestep to display the iteration's progress. 1330814e85d6SHong Zhang 1331814e85d6SHong Zhang Logically Collective on TS 1332814e85d6SHong Zhang 1333814e85d6SHong Zhang Input Parameters: 1334814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1335814e85d6SHong Zhang . adjointmonitor - monitoring routine 1336814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 1337814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 1338814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 1339814e85d6SHong Zhang (may be NULL) 1340814e85d6SHong Zhang 1341814e85d6SHong Zhang Calling sequence of monitor: 1342814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1343814e85d6SHong Zhang 1344814e85d6SHong Zhang + ts - the TS context 1345814e85d6SHong 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 1346814e85d6SHong Zhang been interpolated to) 1347814e85d6SHong Zhang . time - current time 1348814e85d6SHong Zhang . u - current iterate 1349814e85d6SHong Zhang . numcost - number of cost functionos 1350814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1351814e85d6SHong Zhang . mu - sensitivities to parameters 1352814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 1353814e85d6SHong Zhang 1354814e85d6SHong Zhang Notes: 1355814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 1356814e85d6SHong Zhang already has been loaded. 1357814e85d6SHong Zhang 1358814e85d6SHong Zhang Fortran Notes: 1359814e85d6SHong Zhang Only a single monitor function can be set for each TS object 1360814e85d6SHong Zhang 1361814e85d6SHong Zhang Level: intermediate 1362814e85d6SHong Zhang 1363814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 1364814e85d6SHong Zhang 1365814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel() 1366814e85d6SHong Zhang @*/ 1367814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 1368814e85d6SHong Zhang { 1369814e85d6SHong Zhang PetscErrorCode ierr; 1370814e85d6SHong Zhang PetscInt i; 1371814e85d6SHong Zhang PetscBool identical; 1372814e85d6SHong Zhang 1373814e85d6SHong Zhang PetscFunctionBegin; 1374814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1375814e85d6SHong Zhang for (i=0; i<ts->numbermonitors;i++) { 1376814e85d6SHong Zhang ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 1377814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 1378814e85d6SHong Zhang } 1379814e85d6SHong Zhang if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 1380814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1381814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1382814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 1383814e85d6SHong Zhang PetscFunctionReturn(0); 1384814e85d6SHong Zhang } 1385814e85d6SHong Zhang 1386814e85d6SHong Zhang /*@C 1387814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1388814e85d6SHong Zhang 1389814e85d6SHong Zhang Logically Collective on TS 1390814e85d6SHong Zhang 1391814e85d6SHong Zhang Input Parameters: 1392814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1393814e85d6SHong Zhang 1394814e85d6SHong Zhang Notes: 1395814e85d6SHong Zhang There is no way to remove a single, specific monitor. 1396814e85d6SHong Zhang 1397814e85d6SHong Zhang Level: intermediate 1398814e85d6SHong Zhang 1399814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 1400814e85d6SHong Zhang 1401814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1402814e85d6SHong Zhang @*/ 1403814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts) 1404814e85d6SHong Zhang { 1405814e85d6SHong Zhang PetscErrorCode ierr; 1406814e85d6SHong Zhang PetscInt i; 1407814e85d6SHong Zhang 1408814e85d6SHong Zhang PetscFunctionBegin; 1409814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1410814e85d6SHong Zhang for (i=0; i<ts->numberadjointmonitors; i++) { 1411814e85d6SHong Zhang if (ts->adjointmonitordestroy[i]) { 1412814e85d6SHong Zhang ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1413814e85d6SHong Zhang } 1414814e85d6SHong Zhang } 1415814e85d6SHong Zhang ts->numberadjointmonitors = 0; 1416814e85d6SHong Zhang PetscFunctionReturn(0); 1417814e85d6SHong Zhang } 1418814e85d6SHong Zhang 1419814e85d6SHong Zhang /*@C 1420814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 1421814e85d6SHong Zhang 1422814e85d6SHong Zhang Level: intermediate 1423814e85d6SHong Zhang 1424814e85d6SHong Zhang .keywords: TS, set, monitor 1425814e85d6SHong Zhang 1426814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1427814e85d6SHong Zhang @*/ 1428814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1429814e85d6SHong Zhang { 1430814e85d6SHong Zhang PetscErrorCode ierr; 1431814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1432814e85d6SHong Zhang 1433814e85d6SHong Zhang PetscFunctionBegin; 1434814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1435814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1436814e85d6SHong Zhang ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1437814e85d6SHong 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); 1438814e85d6SHong Zhang ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1439814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1440814e85d6SHong Zhang PetscFunctionReturn(0); 1441814e85d6SHong Zhang } 1442814e85d6SHong Zhang 1443814e85d6SHong Zhang /*@C 1444814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1445814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 1446814e85d6SHong Zhang 1447814e85d6SHong Zhang Collective on TS 1448814e85d6SHong Zhang 1449814e85d6SHong Zhang Input Parameters: 1450814e85d6SHong Zhang + ts - the TS context 1451814e85d6SHong Zhang . step - current time-step 1452814e85d6SHong Zhang . ptime - current time 1453814e85d6SHong Zhang . u - current state 1454814e85d6SHong Zhang . numcost - number of cost functions 1455814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1456814e85d6SHong Zhang . mu - sensitivities to parameters 1457814e85d6SHong Zhang - dummy - either a viewer or NULL 1458814e85d6SHong Zhang 1459814e85d6SHong Zhang Level: intermediate 1460814e85d6SHong Zhang 1461814e85d6SHong Zhang .keywords: TS, vector, adjoint, monitor, view 1462814e85d6SHong Zhang 1463814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 1464814e85d6SHong Zhang @*/ 1465814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 1466814e85d6SHong Zhang { 1467814e85d6SHong Zhang PetscErrorCode ierr; 1468814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1469814e85d6SHong Zhang PetscDraw draw; 1470814e85d6SHong Zhang PetscReal xl,yl,xr,yr,h; 1471814e85d6SHong Zhang char time[32]; 1472814e85d6SHong Zhang 1473814e85d6SHong Zhang PetscFunctionBegin; 1474814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1475814e85d6SHong Zhang 1476814e85d6SHong Zhang ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 1477814e85d6SHong Zhang ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 1478814e85d6SHong Zhang ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 1479814e85d6SHong Zhang ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1480814e85d6SHong Zhang h = yl + .95*(yr - yl); 1481814e85d6SHong Zhang ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 1482814e85d6SHong Zhang ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 1483814e85d6SHong Zhang PetscFunctionReturn(0); 1484814e85d6SHong Zhang } 1485814e85d6SHong Zhang 1486814e85d6SHong Zhang /* 1487814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1488814e85d6SHong Zhang 1489814e85d6SHong Zhang Collective on TSAdjoint 1490814e85d6SHong Zhang 1491814e85d6SHong Zhang Input Parameter: 1492814e85d6SHong Zhang . ts - the TS context 1493814e85d6SHong Zhang 1494814e85d6SHong Zhang Options Database Keys: 1495814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1496814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 1497814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1498814e85d6SHong Zhang 1499814e85d6SHong Zhang Level: developer 1500814e85d6SHong Zhang 1501814e85d6SHong Zhang Notes: 1502814e85d6SHong Zhang This is not normally called directly by users 1503814e85d6SHong Zhang 1504814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database 1505814e85d6SHong Zhang 1506814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 1507814e85d6SHong Zhang */ 1508814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 1509814e85d6SHong Zhang { 1510814e85d6SHong Zhang PetscBool tflg,opt; 1511814e85d6SHong Zhang PetscErrorCode ierr; 1512814e85d6SHong Zhang 1513814e85d6SHong Zhang PetscFunctionBegin; 1514814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,2); 1515814e85d6SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 1516814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 1517814e85d6SHong Zhang ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 1518814e85d6SHong Zhang if (opt) { 1519814e85d6SHong Zhang ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 1520814e85d6SHong Zhang ts->adjoint_solve = tflg; 1521814e85d6SHong Zhang } 1522814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 1523814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 1524814e85d6SHong Zhang opt = PETSC_FALSE; 1525814e85d6SHong Zhang ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 1526814e85d6SHong Zhang if (opt) { 1527814e85d6SHong Zhang TSMonitorDrawCtx ctx; 1528814e85d6SHong Zhang PetscInt howoften = 1; 1529814e85d6SHong Zhang 1530814e85d6SHong Zhang ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 1531814e85d6SHong Zhang ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 1532814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 1533814e85d6SHong Zhang } 1534814e85d6SHong Zhang PetscFunctionReturn(0); 1535814e85d6SHong Zhang } 1536814e85d6SHong Zhang 1537814e85d6SHong Zhang /*@ 1538814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 1539814e85d6SHong Zhang 1540814e85d6SHong Zhang Collective on TS 1541814e85d6SHong Zhang 1542814e85d6SHong Zhang Input Parameter: 1543814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1544814e85d6SHong Zhang 1545814e85d6SHong Zhang Level: intermediate 1546814e85d6SHong Zhang 1547814e85d6SHong Zhang .keywords: TS, adjoint, step 1548814e85d6SHong Zhang 1549814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve() 1550814e85d6SHong Zhang @*/ 1551814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts) 1552814e85d6SHong Zhang { 1553814e85d6SHong Zhang DM dm; 1554814e85d6SHong Zhang PetscErrorCode ierr; 1555814e85d6SHong Zhang 1556814e85d6SHong Zhang PetscFunctionBegin; 1557814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1558814e85d6SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1559814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1560814e85d6SHong Zhang 1561814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1562814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 1563814e85d6SHong 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); 1564814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1565814e85d6SHong Zhang ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 1566814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1567814e85d6SHong Zhang ts->adjoint_steps++; ts->steps--; 1568814e85d6SHong Zhang 1569814e85d6SHong Zhang if (ts->reason < 0) { 1570814e85d6SHong Zhang if (ts->errorifstepfailed) { 157105c9950eSHong 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]); 157205c9950eSHong Zhang else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]); 1573814e85d6SHong Zhang } 1574814e85d6SHong Zhang } else if (!ts->reason) { 1575814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1576814e85d6SHong Zhang } 1577814e85d6SHong Zhang PetscFunctionReturn(0); 1578814e85d6SHong Zhang } 1579814e85d6SHong Zhang 1580814e85d6SHong Zhang /*@ 1581814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1582814e85d6SHong Zhang 1583814e85d6SHong Zhang Collective on TS 1584814e85d6SHong Zhang 1585814e85d6SHong Zhang Input Parameter: 1586814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1587814e85d6SHong Zhang 1588814e85d6SHong Zhang Options Database: 1589814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1590814e85d6SHong Zhang 1591814e85d6SHong Zhang Level: intermediate 1592814e85d6SHong Zhang 1593814e85d6SHong Zhang Notes: 1594814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 1595814e85d6SHong Zhang 1596814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1597814e85d6SHong Zhang 1598814e85d6SHong Zhang .keywords: TS, timestep, solve 1599814e85d6SHong Zhang 1600814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 1601814e85d6SHong Zhang @*/ 1602814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts) 1603814e85d6SHong Zhang { 1604814e85d6SHong Zhang PetscErrorCode ierr; 1605814e85d6SHong Zhang 1606814e85d6SHong Zhang PetscFunctionBegin; 1607814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1608814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1609814e85d6SHong Zhang 1610814e85d6SHong Zhang /* reset time step and iteration counters */ 1611814e85d6SHong Zhang ts->adjoint_steps = 0; 1612814e85d6SHong Zhang ts->ksp_its = 0; 1613814e85d6SHong Zhang ts->snes_its = 0; 1614814e85d6SHong Zhang ts->num_snes_failures = 0; 1615814e85d6SHong Zhang ts->reject = 0; 1616814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1617814e85d6SHong Zhang 1618814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1619814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1620814e85d6SHong Zhang 1621814e85d6SHong Zhang while (!ts->reason) { 1622814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1623814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1624814e85d6SHong Zhang ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 1625814e85d6SHong Zhang ierr = TSAdjointStep(ts);CHKERRQ(ierr); 1626cd4cee2dSHong Zhang if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) { 1627814e85d6SHong Zhang ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 1628814e85d6SHong Zhang } 1629814e85d6SHong Zhang } 1630814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1631814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1632814e85d6SHong Zhang ts->solvetime = ts->ptime; 1633814e85d6SHong Zhang ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 1634814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 1635814e85d6SHong Zhang ts->adjoint_max_steps = 0; 1636814e85d6SHong Zhang PetscFunctionReturn(0); 1637814e85d6SHong Zhang } 1638814e85d6SHong Zhang 1639814e85d6SHong Zhang /*@C 1640814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1641814e85d6SHong Zhang 1642814e85d6SHong Zhang Collective on TS 1643814e85d6SHong Zhang 1644814e85d6SHong Zhang Input Parameters: 1645814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 1646814e85d6SHong Zhang . step - step number that has just completed 1647814e85d6SHong Zhang . ptime - model time of the state 1648814e85d6SHong Zhang . u - state at the current model time 1649814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 1650814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1651814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1652814e85d6SHong Zhang 1653814e85d6SHong Zhang Notes: 1654814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1655814e85d6SHong Zhang Users would almost never call this routine directly. 1656814e85d6SHong Zhang 1657814e85d6SHong Zhang Level: developer 1658814e85d6SHong Zhang 1659814e85d6SHong Zhang .keywords: TS, timestep 1660814e85d6SHong Zhang @*/ 1661814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 1662814e85d6SHong Zhang { 1663814e85d6SHong Zhang PetscErrorCode ierr; 1664814e85d6SHong Zhang PetscInt i,n = ts->numberadjointmonitors; 1665814e85d6SHong Zhang 1666814e85d6SHong Zhang PetscFunctionBegin; 1667814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1668814e85d6SHong Zhang PetscValidHeaderSpecific(u,VEC_CLASSID,4); 16698860a134SJunchao Zhang ierr = VecLockReadPush(u);CHKERRQ(ierr); 1670814e85d6SHong Zhang for (i=0; i<n; i++) { 1671814e85d6SHong Zhang ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1672814e85d6SHong Zhang } 16738860a134SJunchao Zhang ierr = VecLockReadPop(u);CHKERRQ(ierr); 1674814e85d6SHong Zhang PetscFunctionReturn(0); 1675814e85d6SHong Zhang } 1676814e85d6SHong Zhang 1677814e85d6SHong Zhang /*@ 1678814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1679814e85d6SHong Zhang 1680814e85d6SHong Zhang Collective on TS 1681814e85d6SHong Zhang 1682814e85d6SHong Zhang Input Arguments: 1683814e85d6SHong Zhang . ts - time stepping context 1684814e85d6SHong Zhang 1685814e85d6SHong Zhang Level: advanced 1686814e85d6SHong Zhang 1687814e85d6SHong Zhang Notes: 1688814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 1689814e85d6SHong Zhang 1690814e85d6SHong Zhang .seealso: TSAdjointSolve(), TSAdjointStep 1691814e85d6SHong Zhang @*/ 1692814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts) 1693814e85d6SHong Zhang { 1694814e85d6SHong Zhang PetscErrorCode ierr; 1695814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1696814e85d6SHong 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); 1697814e85d6SHong Zhang ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 1698814e85d6SHong Zhang PetscFunctionReturn(0); 1699814e85d6SHong Zhang } 1700814e85d6SHong Zhang 1701814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1702814e85d6SHong Zhang 1703814e85d6SHong Zhang /*@ 1704814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 1705814e85d6SHong Zhang of forward sensitivity analysis 1706814e85d6SHong Zhang 1707814e85d6SHong Zhang Collective on TS 1708814e85d6SHong Zhang 1709814e85d6SHong Zhang Input Parameter: 1710814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1711814e85d6SHong Zhang 1712814e85d6SHong Zhang Level: advanced 1713814e85d6SHong Zhang 1714814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup 1715814e85d6SHong Zhang 1716814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp() 1717814e85d6SHong Zhang @*/ 1718814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts) 1719814e85d6SHong Zhang { 1720814e85d6SHong Zhang PetscErrorCode ierr; 1721814e85d6SHong Zhang 1722814e85d6SHong Zhang PetscFunctionBegin; 1723814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1724814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1725814e85d6SHong Zhang if (ts->ops->forwardsetup) { 1726814e85d6SHong Zhang ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 1727814e85d6SHong Zhang } 17287207e0fdSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ts->vec_sensip_col);CHKERRQ(ierr); 1729814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 1730814e85d6SHong Zhang PetscFunctionReturn(0); 1731814e85d6SHong Zhang } 1732814e85d6SHong Zhang 1733814e85d6SHong Zhang /*@ 17347adebddeSHong Zhang TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 17357adebddeSHong Zhang 17367adebddeSHong Zhang Collective on TS 17377adebddeSHong Zhang 17387adebddeSHong Zhang Input Parameter: 17397adebddeSHong Zhang . ts - the TS context obtained from TSCreate() 17407adebddeSHong Zhang 17417adebddeSHong Zhang Level: advanced 17427adebddeSHong Zhang 17437adebddeSHong Zhang .keywords: TS, forward sensitivity, reset 17447adebddeSHong Zhang 17457adebddeSHong Zhang .seealso: TSCreate(), TSDestroy(), TSForwardSetUp() 17467adebddeSHong Zhang @*/ 17477adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts) 17487adebddeSHong Zhang { 17497207e0fdSHong Zhang TS quadts = ts->quadraturets; 17507adebddeSHong Zhang PetscErrorCode ierr; 17517adebddeSHong Zhang 17527adebddeSHong Zhang PetscFunctionBegin; 17537adebddeSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 17547adebddeSHong Zhang if (ts->ops->forwardreset) { 17557adebddeSHong Zhang ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr); 17567adebddeSHong Zhang } 17577adebddeSHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 17587207e0fdSHong Zhang if (quadts) { 17597207e0fdSHong Zhang ierr = MatDestroy(&quadts->mat_sensip);CHKERRQ(ierr); 17607207e0fdSHong Zhang } 17617207e0fdSHong Zhang ierr = VecDestroy(&ts->vec_sensip_col);CHKERRQ(ierr); 17627adebddeSHong Zhang ts->forward_solve = PETSC_FALSE; 17637adebddeSHong Zhang ts->forwardsetupcalled = PETSC_FALSE; 17647adebddeSHong Zhang PetscFunctionReturn(0); 17657adebddeSHong Zhang } 17667adebddeSHong Zhang 17677adebddeSHong Zhang /*@ 1768814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1769814e85d6SHong Zhang 1770814e85d6SHong Zhang Input Parameter: 1771814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1772814e85d6SHong Zhang . numfwdint- number of integrals 1773814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1774814e85d6SHong Zhang 17757207e0fdSHong Zhang Level: deprecated 1776814e85d6SHong Zhang 1777814e85d6SHong Zhang .keywords: TS, forward sensitivity 1778814e85d6SHong Zhang 1779814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1780814e85d6SHong Zhang @*/ 1781814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 1782814e85d6SHong Zhang { 1783814e85d6SHong Zhang PetscFunctionBegin; 1784814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1785814e85d6SHong 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()"); 1786814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1787814e85d6SHong Zhang 1788814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 1789814e85d6SHong Zhang PetscFunctionReturn(0); 1790814e85d6SHong Zhang } 1791814e85d6SHong Zhang 1792814e85d6SHong Zhang /*@ 1793814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1794814e85d6SHong Zhang 1795814e85d6SHong Zhang Input Parameter: 1796814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1797814e85d6SHong Zhang 1798814e85d6SHong Zhang Output Parameter: 1799814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1800814e85d6SHong Zhang 18017207e0fdSHong Zhang Level: deprecated 1802814e85d6SHong Zhang 1803814e85d6SHong Zhang .keywords: TS, forward sensitivity 1804814e85d6SHong Zhang 1805814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1806814e85d6SHong Zhang @*/ 1807814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1808814e85d6SHong Zhang { 1809814e85d6SHong Zhang PetscFunctionBegin; 1810814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1811814e85d6SHong Zhang PetscValidPointer(vp,3); 1812814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1813814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1814814e85d6SHong Zhang PetscFunctionReturn(0); 1815814e85d6SHong Zhang } 1816814e85d6SHong Zhang 1817814e85d6SHong Zhang /*@ 1818814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1819814e85d6SHong Zhang 1820814e85d6SHong Zhang Collective on TS 1821814e85d6SHong Zhang 1822814e85d6SHong Zhang Input Arguments: 1823814e85d6SHong Zhang . ts - time stepping context 1824814e85d6SHong Zhang 1825814e85d6SHong Zhang Level: advanced 1826814e85d6SHong Zhang 1827814e85d6SHong Zhang Notes: 1828814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1829814e85d6SHong Zhang 1830814e85d6SHong Zhang .keywords: TS, forward sensitivity 1831814e85d6SHong Zhang 1832814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1833814e85d6SHong Zhang @*/ 1834814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts) 1835814e85d6SHong Zhang { 1836814e85d6SHong Zhang PetscErrorCode ierr; 1837814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1838814e85d6SHong Zhang if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1839814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1840814e85d6SHong Zhang ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1841814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1842814e85d6SHong Zhang PetscFunctionReturn(0); 1843814e85d6SHong Zhang } 1844814e85d6SHong Zhang 1845814e85d6SHong Zhang /*@ 1846814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1847814e85d6SHong Zhang 1848814e85d6SHong Zhang Logically Collective on TS and Vec 1849814e85d6SHong Zhang 1850814e85d6SHong Zhang Input Parameters: 1851814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1852814e85d6SHong Zhang . nump - number of parameters 1853814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1854814e85d6SHong Zhang 1855814e85d6SHong Zhang Level: beginner 1856814e85d6SHong Zhang 1857814e85d6SHong Zhang Notes: 1858814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1859814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1860814e85d6SHong Zhang You must call this function before TSSolve(). 1861814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1862814e85d6SHong Zhang 1863814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values 1864814e85d6SHong Zhang 1865814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1866814e85d6SHong Zhang @*/ 1867814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1868814e85d6SHong Zhang { 1869814e85d6SHong Zhang PetscErrorCode ierr; 1870814e85d6SHong Zhang 1871814e85d6SHong Zhang PetscFunctionBegin; 1872814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1873814e85d6SHong Zhang PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1874814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1875b5b4017aSHong Zhang if (nump == PETSC_DEFAULT) { 1876b5b4017aSHong Zhang ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr); 1877b5b4017aSHong Zhang } else ts->num_parameters = nump; 1878814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1879814e85d6SHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1880814e85d6SHong Zhang ts->mat_sensip = Smat; 1881814e85d6SHong Zhang PetscFunctionReturn(0); 1882814e85d6SHong Zhang } 1883814e85d6SHong Zhang 1884814e85d6SHong Zhang /*@ 1885814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1886814e85d6SHong Zhang 1887814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1888814e85d6SHong Zhang 1889814e85d6SHong Zhang Output Parameter: 1890814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1891814e85d6SHong Zhang . nump - number of parameters 1892814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1893814e85d6SHong Zhang 1894814e85d6SHong Zhang Level: intermediate 1895814e85d6SHong Zhang 1896814e85d6SHong Zhang .keywords: TS, forward sensitivity 1897814e85d6SHong Zhang 1898814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1899814e85d6SHong Zhang @*/ 1900814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1901814e85d6SHong Zhang { 1902814e85d6SHong Zhang PetscFunctionBegin; 1903814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1904814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1905814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1906814e85d6SHong Zhang PetscFunctionReturn(0); 1907814e85d6SHong Zhang } 1908814e85d6SHong Zhang 1909814e85d6SHong Zhang /*@ 1910814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1911814e85d6SHong Zhang 1912814e85d6SHong Zhang Collective on TS 1913814e85d6SHong Zhang 1914814e85d6SHong Zhang Input Arguments: 1915814e85d6SHong Zhang . ts - time stepping context 1916814e85d6SHong Zhang 1917814e85d6SHong Zhang Level: advanced 1918814e85d6SHong Zhang 1919814e85d6SHong Zhang Notes: 1920814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1921814e85d6SHong Zhang 1922814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral() 1923814e85d6SHong Zhang @*/ 1924814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts) 1925814e85d6SHong Zhang { 1926814e85d6SHong Zhang PetscErrorCode ierr; 1927814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1928814e85d6SHong 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); 1929814e85d6SHong Zhang ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1930814e85d6SHong Zhang PetscFunctionReturn(0); 1931814e85d6SHong Zhang } 1932b5b4017aSHong Zhang 1933b5b4017aSHong Zhang /*@ 1934b5b4017aSHong Zhang TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1935b5b4017aSHong Zhang 1936b5b4017aSHong Zhang Collective on TS and Mat 1937b5b4017aSHong Zhang 1938b5b4017aSHong Zhang Input Parameter 1939b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 1940b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition 1941b5b4017aSHong Zhang 1942b5b4017aSHong Zhang Level: intermediate 1943b5b4017aSHong Zhang 1944b5b4017aSHong 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. 1945b5b4017aSHong Zhang 1946b5b4017aSHong Zhang .seealso: TSForwardSetSensitivities() 1947b5b4017aSHong Zhang @*/ 1948b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp) 1949b5b4017aSHong Zhang { 1950b5b4017aSHong Zhang PetscErrorCode ierr; 1951b5b4017aSHong Zhang 1952b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1953b5b4017aSHong Zhang PetscValidHeaderSpecific(didp,MAT_CLASSID,2); 1954b5b4017aSHong Zhang if (!ts->mat_sensip) { 1955b5b4017aSHong Zhang ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr); 1956b5b4017aSHong Zhang } 1957b5b4017aSHong Zhang PetscFunctionReturn(0); 1958b5b4017aSHong Zhang } 19596affb6f8SHong Zhang 19606affb6f8SHong Zhang /*@ 19616affb6f8SHong Zhang TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 19626affb6f8SHong Zhang 19636affb6f8SHong Zhang Input Parameter: 19646affb6f8SHong Zhang . ts - the TS context obtained from TSCreate() 19656affb6f8SHong Zhang 19666affb6f8SHong Zhang Output Parameters: 1967cd4cee2dSHong Zhang + ns - number of stages 19686affb6f8SHong Zhang - S - tangent linear sensitivities at the intermediate stages 19696affb6f8SHong Zhang 19706affb6f8SHong Zhang Level: advanced 19716affb6f8SHong Zhang 19726affb6f8SHong Zhang .keywords: TS, second-order adjoint, forward sensitivity 19736affb6f8SHong Zhang @*/ 19746affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S) 19756affb6f8SHong Zhang { 19766affb6f8SHong Zhang PetscErrorCode ierr; 19776affb6f8SHong Zhang 19786affb6f8SHong Zhang PetscFunctionBegin; 19796affb6f8SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 19806affb6f8SHong Zhang 19816affb6f8SHong Zhang if (!ts->ops->getstages) *S=NULL; 19826affb6f8SHong Zhang else { 19836affb6f8SHong Zhang ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr); 19846affb6f8SHong Zhang } 19856affb6f8SHong Zhang PetscFunctionReturn(0); 19866affb6f8SHong Zhang } 1987cd4cee2dSHong Zhang 1988cd4cee2dSHong Zhang /*@ 1989cd4cee2dSHong Zhang TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 1990cd4cee2dSHong Zhang 1991cd4cee2dSHong Zhang Input Parameter: 1992cd4cee2dSHong Zhang + ts - the TS context obtained from TSCreate() 1993cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1994cd4cee2dSHong Zhang 1995cd4cee2dSHong Zhang Output Parameters: 1996cd4cee2dSHong Zhang . quadts - the child TS context 1997cd4cee2dSHong Zhang 1998cd4cee2dSHong Zhang Level: intermediate 1999cd4cee2dSHong Zhang 2000cd4cee2dSHong Zhang .keywords: TS, quadrature evaluation 2001cd4cee2dSHong Zhang 2002cd4cee2dSHong Zhang .seealso: TSGetQuadratureTS() 2003cd4cee2dSHong Zhang @*/ 2004cd4cee2dSHong Zhang PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts) 2005cd4cee2dSHong Zhang { 2006cd4cee2dSHong Zhang char prefix[128]; 2007cd4cee2dSHong Zhang PetscErrorCode ierr; 2008cd4cee2dSHong Zhang 2009cd4cee2dSHong Zhang PetscFunctionBegin; 2010cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2011cd4cee2dSHong Zhang PetscValidPointer(quadts,2); 2012cd4cee2dSHong Zhang ierr = TSDestroy(&ts->quadraturets);CHKERRQ(ierr); 2013cd4cee2dSHong Zhang ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);CHKERRQ(ierr); 2014cd4cee2dSHong Zhang ierr = PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);CHKERRQ(ierr); 2015cd4cee2dSHong Zhang ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);CHKERRQ(ierr); 2016cd4cee2dSHong Zhang ierr = PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");CHKERRQ(ierr); 2017cd4cee2dSHong Zhang ierr = TSSetOptionsPrefix(ts->quadraturets,prefix);CHKERRQ(ierr); 2018cd4cee2dSHong Zhang *quadts = ts->quadraturets; 2019cd4cee2dSHong Zhang 2020cd4cee2dSHong Zhang if (ts->numcost) { 2021cd4cee2dSHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);CHKERRQ(ierr); 2022cd4cee2dSHong Zhang } else { 2023cd4cee2dSHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);CHKERRQ(ierr); 2024cd4cee2dSHong Zhang } 2025cd4cee2dSHong Zhang ts->costintegralfwd = fwd; 2026cd4cee2dSHong Zhang PetscFunctionReturn(0); 2027cd4cee2dSHong Zhang } 2028cd4cee2dSHong Zhang 2029cd4cee2dSHong Zhang /*@ 2030cd4cee2dSHong Zhang TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 2031cd4cee2dSHong Zhang 2032cd4cee2dSHong Zhang Input Parameter: 2033cd4cee2dSHong Zhang . ts - the TS context obtained from TSCreate() 2034cd4cee2dSHong Zhang 2035cd4cee2dSHong Zhang Output Parameters: 2036cd4cee2dSHong Zhang + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 2037cd4cee2dSHong Zhang - quadts - the child TS context 2038cd4cee2dSHong Zhang 2039cd4cee2dSHong Zhang Level: intermediate 2040cd4cee2dSHong Zhang 2041cd4cee2dSHong Zhang .keywords: TS, quadrature evaluation 2042cd4cee2dSHong Zhang 2043cd4cee2dSHong Zhang .seealso: TSCreateQuadratureTS() 2044cd4cee2dSHong Zhang @*/ 2045cd4cee2dSHong Zhang PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts) 2046cd4cee2dSHong Zhang { 2047cd4cee2dSHong Zhang PetscFunctionBegin; 2048cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2049cd4cee2dSHong Zhang if (fwd) *fwd = ts->costintegralfwd; 2050cd4cee2dSHong Zhang if (quadts) *quadts = ts->quadraturets; 2051cd4cee2dSHong Zhang PetscFunctionReturn(0); 2052cd4cee2dSHong Zhang } 2053