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 32*cd4cee2dSHong 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 53*cd4cee2dSHong 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. 54*cd4cee2dSHong Zhang 55*cd4cee2dSHong Zhang Logically Collective on TS 56*cd4cee2dSHong Zhang 57*cd4cee2dSHong Zhang Input Parameters: 58*cd4cee2dSHong Zhang . ts - TS context obtained from TSCreate() 59*cd4cee2dSHong Zhang 60*cd4cee2dSHong Zhang Output Parameters: 61*cd4cee2dSHong Zhang + Amat - JacobianP matrix 62*cd4cee2dSHong Zhang . func - function 63*cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 64*cd4cee2dSHong Zhang 65*cd4cee2dSHong Zhang Calling sequence of func: 66*cd4cee2dSHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 67*cd4cee2dSHong Zhang + t - current timestep 68*cd4cee2dSHong Zhang . U - input vector (current ODE solution) 69*cd4cee2dSHong Zhang . A - output matrix 70*cd4cee2dSHong Zhang - ctx - [optional] user-defined function context 71*cd4cee2dSHong Zhang 72*cd4cee2dSHong Zhang Level: intermediate 73*cd4cee2dSHong Zhang 74*cd4cee2dSHong Zhang Notes: 75*cd4cee2dSHong 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 76*cd4cee2dSHong Zhang 77*cd4cee2dSHong Zhang .keywords: TS, sensitivity 78*cd4cee2dSHong Zhang .seealso: TSSetRHSJacobianP() 79*cd4cee2dSHong Zhang @*/ 80*cd4cee2dSHong Zhang PetscErrorCode TSGetRHSJacobianP(TS ts,Mat *Amat,PetscErrorCode (**func)(TS,PetscReal,Vec,Mat,void*),void **ctx) 81*cd4cee2dSHong Zhang { 82*cd4cee2dSHong Zhang PetscFunctionBegin; 83*cd4cee2dSHong Zhang if (func) *func = ts->rhsjacobianp; 84*cd4cee2dSHong Zhang if (ctx) *ctx = ts->rhsjacobianpctx; 85*cd4cee2dSHong Zhang if (Amat) *Amat = ts->Jacprhs; 86*cd4cee2dSHong Zhang PetscFunctionReturn(0); 87*cd4cee2dSHong Zhang } 88*cd4cee2dSHong Zhang 89*cd4cee2dSHong 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 252*cd4cee2dSHong 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 { 318*cd4cee2dSHong Zhang TS quadts; 319*cd4cee2dSHong Zhang PetscErrorCode ierr; 320*cd4cee2dSHong Zhang 321814e85d6SHong Zhang PetscFunctionBegin; 322814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 323814e85d6SHong Zhang PetscValidPointer(v,2); 324*cd4cee2dSHong Zhang ierr = TSGetQuadratureTS(ts,NULL,&quadts);CHKERRQ(ierr); 325*cd4cee2dSHong 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 340814e85d6SHong Zhang Note: 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 344*cd4cee2dSHong 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 373c9ad9de0SHong Zhang TSComputeDRDUFunction - Runs the user-defined DRDU function. 374814e85d6SHong Zhang 375814e85d6SHong Zhang Collective on TS 376814e85d6SHong Zhang 377814e85d6SHong Zhang Input Parameters: 378c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate() 379c9ad9de0SHong Zhang . t - current time 380c9ad9de0SHong Zhang - U - stata vector 381c9ad9de0SHong Zhang 382c9ad9de0SHong Zhang Output Parameters: 383b5b4017aSHong Zhang . DRDU - vector array to hold the outputs 384814e85d6SHong Zhang 385814e85d6SHong Zhang Notes: 386c9ad9de0SHong Zhang TSComputeDRDUFunction() is typically used for sensitivity implementation, 387814e85d6SHong Zhang so most users would not generally call this routine themselves. 388814e85d6SHong Zhang 389814e85d6SHong Zhang Level: developer 390814e85d6SHong Zhang 391814e85d6SHong Zhang .keywords: TS, sensitivity 392814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 393814e85d6SHong Zhang @*/ 394c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 395814e85d6SHong Zhang { 396814e85d6SHong Zhang PetscErrorCode ierr; 397814e85d6SHong Zhang 398814e85d6SHong Zhang PetscFunctionBegin; 39933f72c5dSHong Zhang if (!DRDU) PetscFunctionReturn(0); 400814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 401c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 402814e85d6SHong Zhang 403c9ad9de0SHong Zhang PetscStackPush("TS user DRDU function for sensitivity analysis"); 404c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 405814e85d6SHong Zhang PetscStackPop; 406814e85d6SHong Zhang PetscFunctionReturn(0); 407814e85d6SHong Zhang } 408814e85d6SHong Zhang 409b5b4017aSHong Zhang /*@C 410814e85d6SHong Zhang TSComputeDRDPFunction - Runs the user-defined DRDP function. 411814e85d6SHong Zhang 412814e85d6SHong Zhang Collective on TS 413814e85d6SHong Zhang 414814e85d6SHong Zhang Input Parameters: 415c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate() 416c9ad9de0SHong Zhang . t - current time 417c9ad9de0SHong Zhang - U - stata vector 418c9ad9de0SHong Zhang 419c9ad9de0SHong Zhang Output Parameters: 420b5b4017aSHong Zhang . DRDP - vector array to hold the outputs 421814e85d6SHong Zhang 422814e85d6SHong Zhang Notes: 423c9ad9de0SHong Zhang TSComputeDRDPFunction() is typically used for sensitivity implementation, 424814e85d6SHong Zhang so most users would not generally call this routine themselves. 425814e85d6SHong Zhang 426814e85d6SHong Zhang Level: developer 427814e85d6SHong Zhang 428814e85d6SHong Zhang .keywords: TS, sensitivity 429814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 430814e85d6SHong Zhang @*/ 431c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 432814e85d6SHong Zhang { 433814e85d6SHong Zhang PetscErrorCode ierr; 434814e85d6SHong Zhang 435814e85d6SHong Zhang PetscFunctionBegin; 43633f72c5dSHong Zhang if (!DRDP) PetscFunctionReturn(0); 437814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 438c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 439814e85d6SHong Zhang 440814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 441c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 442814e85d6SHong Zhang PetscStackPop; 443814e85d6SHong Zhang PetscFunctionReturn(0); 444814e85d6SHong Zhang } 445814e85d6SHong Zhang 446b5b4017aSHong Zhang /*@C 447b5b4017aSHong 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. 448b5b4017aSHong Zhang 449b5b4017aSHong Zhang Logically Collective on TS 450b5b4017aSHong Zhang 451b5b4017aSHong Zhang Input Parameters: 452b5b4017aSHong Zhang + ts - TS context obtained from TSCreate() 453b5b4017aSHong Zhang . ihp1 - an array of vectors storing the result of vector-Hessian-vector product for F_UU 454b5b4017aSHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for F_UU 455b5b4017aSHong Zhang . ihp2 - an array of vectors storing the result of vector-Hessian-vector product for F_UP 456b5b4017aSHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for F_UP 457b5b4017aSHong Zhang . ihp3 - an array of vectors storing the result of vector-Hessian-vector product for F_PU 458b5b4017aSHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for F_PU 459b5b4017aSHong Zhang . ihp4 - an array of vectors storing the result of vector-Hessian-vector product for F_PP 460b5b4017aSHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for F_PP 461b5b4017aSHong Zhang 462b5b4017aSHong Zhang Calling sequence of ihessianproductfunc: 463b5b4017aSHong Zhang $ ihessianproductfunc (TS ts,PetscReal t,Vec U,Vec Vl,Vec Vr,Vec VHV,void *ctx); 464b5b4017aSHong Zhang + t - current timestep 465b5b4017aSHong Zhang . U - input vector (current ODE solution) 466b5b4017aSHong Zhang . Vl - input vector to be left-multiplied with the Hessian 467b5b4017aSHong Zhang . Vr - input vector to be right-multiplied with the Hessian 468b5b4017aSHong Zhang . VHV - output vector for vector-Hessian-vector product 469b5b4017aSHong Zhang - ctx - [optional] user-defined function context 470b5b4017aSHong Zhang 471b5b4017aSHong Zhang Level: intermediate 472b5b4017aSHong Zhang 473b5b4017aSHong Zhang Note: The first Hessian function and the working array are required. 474b5b4017aSHong Zhang 475b5b4017aSHong Zhang .keywords: TS, sensitivity 476b5b4017aSHong Zhang 477b5b4017aSHong Zhang .seealso: 478b5b4017aSHong Zhang @*/ 479b5b4017aSHong Zhang PetscErrorCode TSSetIHessianProduct(TS ts,Vec *ihp1,PetscErrorCode (*ihessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 480b5b4017aSHong Zhang Vec *ihp2,PetscErrorCode (*ihessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 481b5b4017aSHong Zhang Vec *ihp3,PetscErrorCode (*ihessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 482b5b4017aSHong Zhang Vec *ihp4,PetscErrorCode (*ihessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 483b5b4017aSHong Zhang void *ctx) 484b5b4017aSHong Zhang { 485b5b4017aSHong Zhang PetscFunctionBegin; 486b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 487b5b4017aSHong Zhang PetscValidPointer(ihp1,2); 488b5b4017aSHong Zhang 489b5b4017aSHong Zhang ts->ihessianproductctx = ctx; 490b5b4017aSHong Zhang if (ihp1) ts->vecs_fuu = ihp1; 491b5b4017aSHong Zhang if (ihp2) ts->vecs_fup = ihp2; 492b5b4017aSHong Zhang if (ihp3) ts->vecs_fpu = ihp3; 493b5b4017aSHong Zhang if (ihp4) ts->vecs_fpp = ihp4; 494b5b4017aSHong Zhang ts->ihessianproduct_fuu = ihessianproductfunc1; 495b5b4017aSHong Zhang ts->ihessianproduct_fup = ihessianproductfunc2; 496b5b4017aSHong Zhang ts->ihessianproduct_fpu = ihessianproductfunc3; 497b5b4017aSHong Zhang ts->ihessianproduct_fpp = ihessianproductfunc4; 498b5b4017aSHong Zhang PetscFunctionReturn(0); 499b5b4017aSHong Zhang } 500b5b4017aSHong Zhang 501b5b4017aSHong Zhang /*@C 502b5b4017aSHong Zhang TSComputeIHessianProductFunction1 - Runs the user-defined vector-Hessian-vector product function for Fuu. 503b5b4017aSHong Zhang 504b5b4017aSHong Zhang Collective on TS 505b5b4017aSHong Zhang 506b5b4017aSHong Zhang Input Parameters: 507b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 508b5b4017aSHong Zhang 509b5b4017aSHong Zhang Notes: 510b5b4017aSHong Zhang TSComputeIHessianProductFunction1() is typically used for sensitivity implementation, 511b5b4017aSHong Zhang so most users would not generally call this routine themselves. 512b5b4017aSHong Zhang 513b5b4017aSHong Zhang Level: developer 514b5b4017aSHong Zhang 515b5b4017aSHong Zhang .keywords: TS, sensitivity 516b5b4017aSHong Zhang 517b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 518b5b4017aSHong Zhang @*/ 519b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction1(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 520b5b4017aSHong Zhang { 521b5b4017aSHong Zhang PetscErrorCode ierr; 522b5b4017aSHong Zhang 523b5b4017aSHong Zhang PetscFunctionBegin; 52433f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 525b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 526b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 527b5b4017aSHong Zhang 52867633408SHong Zhang if (ts->ihessianproduct_fuu) { 529b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 1 for sensitivity analysis"); 530b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fuu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 531b5b4017aSHong Zhang PetscStackPop; 53267633408SHong Zhang } 53367633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 53467633408SHong Zhang if (ts->rhshessianproduct_guu) { 53567633408SHong Zhang PetscInt nadj; 53667633408SHong Zhang ierr = TSComputeRHSHessianProductFunction1(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 53767633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 53867633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 53967633408SHong Zhang } 54067633408SHong Zhang } 541b5b4017aSHong Zhang PetscFunctionReturn(0); 542b5b4017aSHong Zhang } 543b5b4017aSHong Zhang 544b5b4017aSHong Zhang /*@C 545b5b4017aSHong Zhang TSComputeIHessianProductFunction2 - Runs the user-defined vector-Hessian-vector product function for Fup. 546b5b4017aSHong Zhang 547b5b4017aSHong Zhang Collective on TS 548b5b4017aSHong Zhang 549b5b4017aSHong Zhang Input Parameters: 550b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 551b5b4017aSHong Zhang 552b5b4017aSHong Zhang Notes: 553b5b4017aSHong Zhang TSComputeIHessianProductFunction2() is typically used for sensitivity implementation, 554b5b4017aSHong Zhang so most users would not generally call this routine themselves. 555b5b4017aSHong Zhang 556b5b4017aSHong Zhang Level: developer 557b5b4017aSHong Zhang 558b5b4017aSHong Zhang .keywords: TS, sensitivity 559b5b4017aSHong Zhang 560b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 561b5b4017aSHong Zhang @*/ 562b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction2(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 563b5b4017aSHong Zhang { 564b5b4017aSHong Zhang PetscErrorCode ierr; 565b5b4017aSHong Zhang 566b5b4017aSHong Zhang PetscFunctionBegin; 56733f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 568b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 569b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 570b5b4017aSHong Zhang 57167633408SHong Zhang if (ts->ihessianproduct_fup) { 572b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 2 for sensitivity analysis"); 573b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fup)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 574b5b4017aSHong Zhang PetscStackPop; 57567633408SHong Zhang } 57667633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 57767633408SHong Zhang if (ts->rhshessianproduct_gup) { 57867633408SHong Zhang PetscInt nadj; 57967633408SHong Zhang ierr = TSComputeRHSHessianProductFunction2(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 58067633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 58167633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 58267633408SHong Zhang } 58367633408SHong Zhang } 584b5b4017aSHong Zhang PetscFunctionReturn(0); 585b5b4017aSHong Zhang } 586b5b4017aSHong Zhang 587b5b4017aSHong Zhang /*@C 588b5b4017aSHong Zhang TSComputeIHessianProductFunction3 - Runs the user-defined vector-Hessian-vector product function for Fpu. 589b5b4017aSHong Zhang 590b5b4017aSHong Zhang Collective on TS 591b5b4017aSHong Zhang 592b5b4017aSHong Zhang Input Parameters: 593b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 594b5b4017aSHong Zhang 595b5b4017aSHong Zhang Notes: 596b5b4017aSHong Zhang TSComputeIHessianProductFunction3() is typically used for sensitivity implementation, 597b5b4017aSHong Zhang so most users would not generally call this routine themselves. 598b5b4017aSHong Zhang 599b5b4017aSHong Zhang Level: developer 600b5b4017aSHong Zhang 601b5b4017aSHong Zhang .keywords: TS, sensitivity 602b5b4017aSHong Zhang 603b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 604b5b4017aSHong Zhang @*/ 605b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction3(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 606b5b4017aSHong Zhang { 607b5b4017aSHong Zhang PetscErrorCode ierr; 608b5b4017aSHong Zhang 609b5b4017aSHong Zhang PetscFunctionBegin; 61033f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 611b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 612b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 613b5b4017aSHong Zhang 61467633408SHong Zhang if (ts->ihessianproduct_fpu) { 615b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 616b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fpu)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 617b5b4017aSHong Zhang PetscStackPop; 61867633408SHong Zhang } 61967633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 62067633408SHong Zhang if (ts->rhshessianproduct_gpu) { 62167633408SHong Zhang PetscInt nadj; 62267633408SHong Zhang ierr = TSComputeRHSHessianProductFunction3(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 62367633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 62467633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 62567633408SHong Zhang } 62667633408SHong Zhang } 627b5b4017aSHong Zhang PetscFunctionReturn(0); 628b5b4017aSHong Zhang } 629b5b4017aSHong Zhang 630b5b4017aSHong Zhang /*@C 631b5b4017aSHong Zhang TSComputeIHessianProductFunction4 - Runs the user-defined vector-Hessian-vector product function for Fpp. 632b5b4017aSHong Zhang 633b5b4017aSHong Zhang Collective on TS 634b5b4017aSHong Zhang 635b5b4017aSHong Zhang Input Parameters: 636b5b4017aSHong Zhang . ts - The TS context obtained from TSCreate() 637b5b4017aSHong Zhang 638b5b4017aSHong Zhang Notes: 639b5b4017aSHong Zhang TSComputeIHessianProductFunction4() is typically used for sensitivity implementation, 640b5b4017aSHong Zhang so most users would not generally call this routine themselves. 641b5b4017aSHong Zhang 642b5b4017aSHong Zhang Level: developer 643b5b4017aSHong Zhang 644b5b4017aSHong Zhang .keywords: TS, sensitivity 645b5b4017aSHong Zhang 646b5b4017aSHong Zhang .seealso: TSSetIHessianProduct() 647b5b4017aSHong Zhang @*/ 648b5b4017aSHong Zhang PetscErrorCode TSComputeIHessianProductFunction4(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 649b5b4017aSHong Zhang { 650b5b4017aSHong Zhang PetscErrorCode ierr; 651b5b4017aSHong Zhang 652b5b4017aSHong Zhang PetscFunctionBegin; 65333f72c5dSHong Zhang if (!VHV) PetscFunctionReturn(0); 654b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 655b5b4017aSHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 656b5b4017aSHong Zhang 65767633408SHong Zhang if (ts->ihessianproduct_fpp) { 658b5b4017aSHong Zhang PetscStackPush("TS user IHessianProduct function 3 for sensitivity analysis"); 659b5b4017aSHong Zhang ierr = (*ts->ihessianproduct_fpp)(ts,t,U,Vl,Vr,VHV,ts->ihessianproductctx);CHKERRQ(ierr); 660b5b4017aSHong Zhang PetscStackPop; 66167633408SHong Zhang } 66267633408SHong Zhang /* does not consider IMEX for now, so either IHessian or RHSHessian will be calculated, using the same output VHV */ 66367633408SHong Zhang if (ts->rhshessianproduct_gpp) { 66467633408SHong Zhang PetscInt nadj; 66567633408SHong Zhang ierr = TSComputeRHSHessianProductFunction4(ts,t,U,Vl,Vr,VHV);CHKERRQ(ierr); 66667633408SHong Zhang for (nadj=0; nadj<ts->numcost; nadj++) { 66767633408SHong Zhang ierr = VecScale(VHV[nadj],-1);CHKERRQ(ierr); 66867633408SHong Zhang } 66967633408SHong Zhang } 670b5b4017aSHong Zhang PetscFunctionReturn(0); 671b5b4017aSHong Zhang } 672b5b4017aSHong Zhang 67313af1a74SHong Zhang /*@C 67413af1a74SHong 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. 67513af1a74SHong Zhang 67613af1a74SHong Zhang Logically Collective on TS 67713af1a74SHong Zhang 67813af1a74SHong Zhang Input Parameters: 67913af1a74SHong Zhang + ts - TS context obtained from TSCreate() 68013af1a74SHong Zhang . rhshp1 - an array of vectors storing the result of vector-Hessian-vector product for G_UU 68113af1a74SHong Zhang . hessianproductfunc1 - vector-Hessian-vector product function for G_UU 68213af1a74SHong Zhang . rhshp2 - an array of vectors storing the result of vector-Hessian-vector product for G_UP 68313af1a74SHong Zhang . hessianproductfunc2 - vector-Hessian-vector product function for G_UP 68413af1a74SHong Zhang . rhshp3 - an array of vectors storing the result of vector-Hessian-vector product for G_PU 68513af1a74SHong Zhang . hessianproductfunc3 - vector-Hessian-vector product function for G_PU 68613af1a74SHong Zhang . rhshp4 - an array of vectors storing the result of vector-Hessian-vector product for G_PP 68713af1a74SHong Zhang . hessianproductfunc4 - vector-Hessian-vector product function for G_PP 68813af1a74SHong Zhang 68913af1a74SHong Zhang Calling sequence of ihessianproductfunc: 69013af1a74SHong Zhang $ rhshessianproductfunc (TS ts,PetscReal t,Vec U,Vec Vl,Vec Vr,Vec VHV,void *ctx); 69113af1a74SHong Zhang + t - current timestep 69213af1a74SHong Zhang . U - input vector (current ODE solution) 69313af1a74SHong Zhang . Vl - input vector to be left-multiplied with the Hessian 69413af1a74SHong Zhang . Vr - input vector to be right-multiplied with the Hessian 69513af1a74SHong Zhang . VHV - output vector for vector-Hessian-vector product 69613af1a74SHong Zhang - ctx - [optional] user-defined function context 69713af1a74SHong Zhang 69813af1a74SHong Zhang Level: intermediate 69913af1a74SHong Zhang 70013af1a74SHong Zhang Note: The first Hessian function and the working array are required. 70113af1a74SHong Zhang 70213af1a74SHong Zhang .keywords: TS, sensitivity 70313af1a74SHong Zhang 70413af1a74SHong Zhang .seealso: 70513af1a74SHong Zhang @*/ 70613af1a74SHong Zhang PetscErrorCode TSSetRHSHessianProduct(TS ts,Vec *rhshp1,PetscErrorCode (*rhshessianproductfunc1)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 70713af1a74SHong Zhang Vec *rhshp2,PetscErrorCode (*rhshessianproductfunc2)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 70813af1a74SHong Zhang Vec *rhshp3,PetscErrorCode (*rhshessianproductfunc3)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 70913af1a74SHong Zhang Vec *rhshp4,PetscErrorCode (*rhshessianproductfunc4)(TS,PetscReal,Vec,Vec*,Vec,Vec*,void*), 71013af1a74SHong Zhang void *ctx) 71113af1a74SHong Zhang { 71213af1a74SHong Zhang PetscFunctionBegin; 71313af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 71413af1a74SHong Zhang PetscValidPointer(rhshp1,2); 71513af1a74SHong Zhang 71613af1a74SHong Zhang ts->rhshessianproductctx = ctx; 71713af1a74SHong Zhang if (rhshp1) ts->vecs_guu = rhshp1; 71813af1a74SHong Zhang if (rhshp2) ts->vecs_gup = rhshp2; 71913af1a74SHong Zhang if (rhshp3) ts->vecs_gpu = rhshp3; 72013af1a74SHong Zhang if (rhshp4) ts->vecs_gpp = rhshp4; 72113af1a74SHong Zhang ts->rhshessianproduct_guu = rhshessianproductfunc1; 72213af1a74SHong Zhang ts->rhshessianproduct_gup = rhshessianproductfunc2; 72313af1a74SHong Zhang ts->rhshessianproduct_gpu = rhshessianproductfunc3; 72413af1a74SHong Zhang ts->rhshessianproduct_gpp = rhshessianproductfunc4; 72513af1a74SHong Zhang PetscFunctionReturn(0); 72613af1a74SHong Zhang } 72713af1a74SHong Zhang 72813af1a74SHong Zhang /*@C 72913af1a74SHong Zhang TSComputeRHSHessianProductFunction1 - Runs the user-defined vector-Hessian-vector product function for Guu. 73013af1a74SHong Zhang 73113af1a74SHong Zhang Collective on TS 73213af1a74SHong Zhang 73313af1a74SHong Zhang Input Parameters: 73413af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 73513af1a74SHong Zhang 73613af1a74SHong Zhang Notes: 73713af1a74SHong Zhang TSComputeRHSHessianProductFunction1() is typically used for sensitivity implementation, 73813af1a74SHong Zhang so most users would not generally call this routine themselves. 73913af1a74SHong Zhang 74013af1a74SHong Zhang Level: developer 74113af1a74SHong Zhang 74213af1a74SHong Zhang .keywords: TS, sensitivity 74313af1a74SHong Zhang 74413af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 74513af1a74SHong Zhang @*/ 74613af1a74SHong Zhang PetscErrorCode TSComputeRHSHessianProductFunction1(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 74713af1a74SHong Zhang { 74813af1a74SHong Zhang PetscErrorCode ierr; 74913af1a74SHong Zhang 75013af1a74SHong Zhang PetscFunctionBegin; 75113af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 75213af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 75313af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 75413af1a74SHong Zhang 75513af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 1 for sensitivity analysis"); 75613af1a74SHong Zhang ierr = (*ts->rhshessianproduct_guu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 75713af1a74SHong Zhang PetscStackPop; 75813af1a74SHong Zhang PetscFunctionReturn(0); 75913af1a74SHong Zhang } 76013af1a74SHong Zhang 76113af1a74SHong Zhang /*@C 76213af1a74SHong Zhang TSComputeRHSHessianProductFunction2 - Runs the user-defined vector-Hessian-vector product function for Gup. 76313af1a74SHong Zhang 76413af1a74SHong Zhang Collective on TS 76513af1a74SHong Zhang 76613af1a74SHong Zhang Input Parameters: 76713af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 76813af1a74SHong Zhang 76913af1a74SHong Zhang Notes: 77013af1a74SHong Zhang TSComputeRHSHessianProductFunction2() is typically used for sensitivity implementation, 77113af1a74SHong Zhang so most users would not generally call this routine themselves. 77213af1a74SHong Zhang 77313af1a74SHong Zhang Level: developer 77413af1a74SHong Zhang 77513af1a74SHong Zhang .keywords: TS, sensitivity 77613af1a74SHong Zhang 77713af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 77813af1a74SHong Zhang @*/ 77913af1a74SHong Zhang PetscErrorCode TSComputeRHSHessianProductFunction2(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 78013af1a74SHong Zhang { 78113af1a74SHong Zhang PetscErrorCode ierr; 78213af1a74SHong Zhang 78313af1a74SHong Zhang PetscFunctionBegin; 78413af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 78513af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 78613af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 78713af1a74SHong Zhang 78813af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 2 for sensitivity analysis"); 78913af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gup)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 79013af1a74SHong Zhang PetscStackPop; 79113af1a74SHong Zhang PetscFunctionReturn(0); 79213af1a74SHong Zhang } 79313af1a74SHong Zhang 79413af1a74SHong Zhang /*@C 79513af1a74SHong Zhang TSComputeRHSHessianProductFunction3 - Runs the user-defined vector-Hessian-vector product function for Gpu. 79613af1a74SHong Zhang 79713af1a74SHong Zhang Collective on TS 79813af1a74SHong Zhang 79913af1a74SHong Zhang Input Parameters: 80013af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 80113af1a74SHong Zhang 80213af1a74SHong Zhang Notes: 80313af1a74SHong Zhang TSComputeRHSHessianProductFunction3() is typically used for sensitivity implementation, 80413af1a74SHong Zhang so most users would not generally call this routine themselves. 80513af1a74SHong Zhang 80613af1a74SHong Zhang Level: developer 80713af1a74SHong Zhang 80813af1a74SHong Zhang .keywords: TS, sensitivity 80913af1a74SHong Zhang 81013af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 81113af1a74SHong Zhang @*/ 81213af1a74SHong Zhang PetscErrorCode TSComputeRHSHessianProductFunction3(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 81313af1a74SHong Zhang { 81413af1a74SHong Zhang PetscErrorCode ierr; 81513af1a74SHong Zhang 81613af1a74SHong Zhang PetscFunctionBegin; 81713af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 81813af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 81913af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 82013af1a74SHong Zhang 82113af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis"); 82213af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gpu)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 82313af1a74SHong Zhang PetscStackPop; 82413af1a74SHong Zhang PetscFunctionReturn(0); 82513af1a74SHong Zhang } 82613af1a74SHong Zhang 82713af1a74SHong Zhang /*@C 82813af1a74SHong Zhang TSComputeRHSHessianProductFunction4 - Runs the user-defined vector-Hessian-vector product function for Gpp. 82913af1a74SHong Zhang 83013af1a74SHong Zhang Collective on TS 83113af1a74SHong Zhang 83213af1a74SHong Zhang Input Parameters: 83313af1a74SHong Zhang . ts - The TS context obtained from TSCreate() 83413af1a74SHong Zhang 83513af1a74SHong Zhang Notes: 83613af1a74SHong Zhang TSComputeRHSHessianProductFunction4() is typically used for sensitivity implementation, 83713af1a74SHong Zhang so most users would not generally call this routine themselves. 83813af1a74SHong Zhang 83913af1a74SHong Zhang Level: developer 84013af1a74SHong Zhang 84113af1a74SHong Zhang .keywords: TS, sensitivity 84213af1a74SHong Zhang 84313af1a74SHong Zhang .seealso: TSSetRHSHessianProduct() 84413af1a74SHong Zhang @*/ 84513af1a74SHong Zhang PetscErrorCode TSComputeRHSHessianProductFunction4(TS ts,PetscReal t,Vec U,Vec *Vl,Vec Vr,Vec *VHV) 84613af1a74SHong Zhang { 84713af1a74SHong Zhang PetscErrorCode ierr; 84813af1a74SHong Zhang 84913af1a74SHong Zhang PetscFunctionBegin; 85013af1a74SHong Zhang if (!VHV) PetscFunctionReturn(0); 85113af1a74SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 85213af1a74SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 85313af1a74SHong Zhang 85413af1a74SHong Zhang PetscStackPush("TS user RHSHessianProduct function 3 for sensitivity analysis"); 85513af1a74SHong Zhang ierr = (*ts->rhshessianproduct_gpp)(ts,t,U,Vl,Vr,VHV,ts->rhshessianproductctx);CHKERRQ(ierr); 85613af1a74SHong Zhang PetscStackPop; 85713af1a74SHong Zhang PetscFunctionReturn(0); 85813af1a74SHong Zhang } 85913af1a74SHong Zhang 860814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 861814e85d6SHong Zhang 862814e85d6SHong Zhang /*@ 863814e85d6SHong 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 864814e85d6SHong Zhang for use by the TSAdjoint routines. 865814e85d6SHong Zhang 866814e85d6SHong Zhang Logically Collective on TS and Vec 867814e85d6SHong Zhang 868814e85d6SHong Zhang Input Parameters: 869814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 870814e85d6SHong 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 871814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 872814e85d6SHong Zhang 873814e85d6SHong Zhang Level: beginner 874814e85d6SHong Zhang 875814e85d6SHong Zhang Notes: 876814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 877814e85d6SHong Zhang 878814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 879814e85d6SHong Zhang 880814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values 881b5b4017aSHong Zhang 882b5b4017aSHong Zhang .seealso TSGetCostGradients() 883814e85d6SHong Zhang @*/ 884814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu) 885814e85d6SHong Zhang { 886814e85d6SHong Zhang PetscFunctionBegin; 887814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 888814e85d6SHong Zhang PetscValidPointer(lambda,2); 889814e85d6SHong Zhang ts->vecs_sensi = lambda; 890814e85d6SHong Zhang ts->vecs_sensip = mu; 891814e85d6SHong 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"); 892814e85d6SHong Zhang ts->numcost = numcost; 893814e85d6SHong Zhang PetscFunctionReturn(0); 894814e85d6SHong Zhang } 895814e85d6SHong Zhang 896814e85d6SHong Zhang /*@ 897814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 898814e85d6SHong Zhang 899814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 900814e85d6SHong Zhang 901814e85d6SHong Zhang Input Parameter: 902814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 903814e85d6SHong Zhang 904814e85d6SHong Zhang Output Parameter: 905814e85d6SHong Zhang + lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 906814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 907814e85d6SHong Zhang 908814e85d6SHong Zhang Level: intermediate 909814e85d6SHong Zhang 910814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity 911b5b4017aSHong Zhang 912b5b4017aSHong Zhang .seealso: TSSetCostGradients() 913814e85d6SHong Zhang @*/ 914814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu) 915814e85d6SHong Zhang { 916814e85d6SHong Zhang PetscFunctionBegin; 917814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 918814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 919814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 920814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 921814e85d6SHong Zhang PetscFunctionReturn(0); 922814e85d6SHong Zhang } 923814e85d6SHong Zhang 924814e85d6SHong Zhang /*@ 925b5b4017aSHong 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 926b5b4017aSHong Zhang for use by the TSAdjoint routines. 927b5b4017aSHong Zhang 928b5b4017aSHong Zhang Logically Collective on TS and Vec 929b5b4017aSHong Zhang 930b5b4017aSHong Zhang Input Parameters: 931b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 932b5b4017aSHong Zhang . numcost - number of cost functions 933b5b4017aSHong 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 934b5b4017aSHong 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 935b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 936b5b4017aSHong Zhang 937b5b4017aSHong Zhang Level: beginner 938b5b4017aSHong Zhang 939b5b4017aSHong Zhang Notes: Hessian of the cost function is completely different from Hessian of the ODE/DAE system 940b5b4017aSHong Zhang 941b5b4017aSHong Zhang For second-order adjoint, one needs to call this function and then TSAdjointInitializeForward() before TSSolve(). 942b5b4017aSHong Zhang 943b5b4017aSHong 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. 944b5b4017aSHong Zhang 9453fca9621SHong Zhang Passing NULL for lambda2 disables the second-order calculation. 946b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint 947b5b4017aSHong Zhang 948b5b4017aSHong Zhang .seealso: TSAdjointInitializeForward() 949b5b4017aSHong Zhang @*/ 950b5b4017aSHong Zhang PetscErrorCode TSSetCostHessianProducts(TS ts,PetscInt numcost,Vec *lambda2,Vec *mu2,Vec dir) 951b5b4017aSHong Zhang { 952b5b4017aSHong Zhang PetscFunctionBegin; 953b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 954b5b4017aSHong 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"); 955b5b4017aSHong Zhang ts->numcost = numcost; 956b5b4017aSHong Zhang ts->vecs_sensi2 = lambda2; 95733f72c5dSHong Zhang ts->vecs_sensi2p = mu2; 958b5b4017aSHong Zhang ts->vec_dir = dir; 959b5b4017aSHong Zhang PetscFunctionReturn(0); 960b5b4017aSHong Zhang } 961b5b4017aSHong Zhang 962b5b4017aSHong Zhang /*@ 963b5b4017aSHong Zhang TSGetCostHessianProducts - Returns the gradients from the TSAdjointSolve() 964b5b4017aSHong Zhang 965b5b4017aSHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 966b5b4017aSHong Zhang 967b5b4017aSHong Zhang Input Parameter: 968b5b4017aSHong Zhang . ts - the TS context obtained from TSCreate() 969b5b4017aSHong Zhang 970b5b4017aSHong Zhang Output Parameter: 971b5b4017aSHong Zhang + numcost - number of cost functions 972b5b4017aSHong 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 973b5b4017aSHong 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 974b5b4017aSHong Zhang - dir - the direction vector that are multiplied with the Hessian of the cost functions 975b5b4017aSHong Zhang 976b5b4017aSHong Zhang Level: intermediate 977b5b4017aSHong Zhang 978b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint 979b5b4017aSHong Zhang 980b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts() 981b5b4017aSHong Zhang @*/ 982b5b4017aSHong Zhang PetscErrorCode TSGetCostHessianProducts(TS ts,PetscInt *numcost,Vec **lambda2,Vec **mu2, Vec *dir) 983b5b4017aSHong Zhang { 984b5b4017aSHong Zhang PetscFunctionBegin; 985b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 986b5b4017aSHong Zhang if (numcost) *numcost = ts->numcost; 987b5b4017aSHong Zhang if (lambda2) *lambda2 = ts->vecs_sensi2; 98833f72c5dSHong Zhang if (mu2) *mu2 = ts->vecs_sensi2p; 989b5b4017aSHong Zhang if (dir) *dir = ts->vec_dir; 990b5b4017aSHong Zhang PetscFunctionReturn(0); 991b5b4017aSHong Zhang } 992b5b4017aSHong Zhang 993b5b4017aSHong Zhang /*@ 994b5b4017aSHong Zhang TSAdjointInitializeForward - Trigger the tangent linear solver and initialize the forward sensitivities 995b5b4017aSHong Zhang 996b5b4017aSHong Zhang Logically Collective on TS and Mat 997b5b4017aSHong Zhang 998b5b4017aSHong Zhang Input Parameters: 999b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 1000b5b4017aSHong Zhang - didp - the derivative of initial values w.r.t. parameters 1001b5b4017aSHong Zhang 1002b5b4017aSHong Zhang Level: intermediate 1003b5b4017aSHong Zhang 1004b5b4017aSHong 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. 1005b5b4017aSHong Zhang 1006b5b4017aSHong Zhang .keywords: TS, sensitivity, second-order adjoint 1007b5b4017aSHong Zhang 1008b5b4017aSHong Zhang .seealso: TSSetCostHessianProducts() 1009b5b4017aSHong Zhang @*/ 1010b5b4017aSHong Zhang PetscErrorCode TSAdjointInitializeForward(TS ts,Mat didp) 1011b5b4017aSHong Zhang { 101233f72c5dSHong Zhang Mat A; 101333f72c5dSHong Zhang Vec sp; 101433f72c5dSHong Zhang PetscScalar *xarr; 101533f72c5dSHong Zhang PetscInt lsize; 1016b5b4017aSHong Zhang PetscErrorCode ierr; 1017b5b4017aSHong Zhang 1018b5b4017aSHong Zhang PetscFunctionBegin; 1019b5b4017aSHong Zhang ts->forward_solve = PETSC_TRUE; /* turn on tangent linear mode */ 1020b5b4017aSHong Zhang if (!ts->vecs_sensi2) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Must call TSSetCostHessianProducts() first"); 102133f72c5dSHong Zhang if (!ts->vec_dir) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_USER,"Directional vector is missing. Call TSSetCostHessianProducts() to set it."); 102233f72c5dSHong Zhang /* create a single-column dense matrix */ 102333f72c5dSHong Zhang ierr = VecGetLocalSize(ts->vec_sol,&lsize);CHKERRQ(ierr); 102433f72c5dSHong Zhang ierr = MatCreateDense(PETSC_COMM_WORLD,lsize,PETSC_DECIDE,PETSC_DECIDE,1,NULL,&A);CHKERRQ(ierr); 102533f72c5dSHong Zhang 102633f72c5dSHong Zhang ierr = VecDuplicate(ts->vec_sol,&sp);CHKERRQ(ierr); 102733f72c5dSHong Zhang ierr = MatDenseGetColumn(A,0,&xarr);CHKERRQ(ierr); 102833f72c5dSHong Zhang ierr = VecPlaceArray(sp,xarr);CHKERRQ(ierr); 102933f72c5dSHong Zhang if (ts->vecs_sensi2p) { /* TLM variable initialized as 2*dIdP*dir */ 103033f72c5dSHong Zhang if (didp) { 103133f72c5dSHong Zhang ierr = MatMult(didp,ts->vec_dir,sp);CHKERRQ(ierr); 103233f72c5dSHong Zhang ierr = VecScale(sp,2.);CHKERRQ(ierr); 103333f72c5dSHong Zhang } else { 103433f72c5dSHong Zhang ierr = VecZeroEntries(sp);CHKERRQ(ierr); 103533f72c5dSHong Zhang } 103633f72c5dSHong Zhang } else { /* TLM variable initialized as dir */ 103733f72c5dSHong Zhang ierr = VecCopy(ts->vec_dir,sp);CHKERRQ(ierr); 103833f72c5dSHong Zhang } 103933f72c5dSHong Zhang ierr = VecResetArray(sp);CHKERRQ(ierr); 104033f72c5dSHong Zhang ierr = MatDenseRestoreColumn(A,&xarr);CHKERRQ(ierr); 104133f72c5dSHong Zhang ierr = VecDestroy(&sp);CHKERRQ(ierr); 104233f72c5dSHong Zhang 104333f72c5dSHong Zhang ierr = TSForwardSetInitialSensitivities(ts,A);CHKERRQ(ierr); /* if didp is NULL, identity matrix is assumed */ 104433f72c5dSHong Zhang 104533f72c5dSHong Zhang ierr = MatDestroy(&A);CHKERRQ(ierr); 1046b5b4017aSHong Zhang PetscFunctionReturn(0); 1047b5b4017aSHong Zhang } 1048b5b4017aSHong Zhang 1049b5b4017aSHong Zhang /*@ 1050814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 1051814e85d6SHong Zhang of an adjoint solver 1052814e85d6SHong Zhang 1053814e85d6SHong Zhang Collective on TS 1054814e85d6SHong Zhang 1055814e85d6SHong Zhang Input Parameter: 1056814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1057814e85d6SHong Zhang 1058814e85d6SHong Zhang Level: advanced 1059814e85d6SHong Zhang 1060814e85d6SHong Zhang .keywords: TS, timestep, setup 1061814e85d6SHong Zhang 1062814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients() 1063814e85d6SHong Zhang @*/ 1064814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts) 1065814e85d6SHong Zhang { 1066814e85d6SHong Zhang PetscErrorCode ierr; 1067814e85d6SHong Zhang 1068814e85d6SHong Zhang PetscFunctionBegin; 1069814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1070814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 1071814e85d6SHong Zhang if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 107233f72c5dSHong Zhang if (ts->vecs_sensip && !ts->Jacp && !ts->Jacprhs) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetRHSJacobianP() or TSSetIJacobianP() first"); 107333f72c5dSHong Zhang 107433f72c5dSHong Zhang if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs; 1075814e85d6SHong Zhang 1076*cd4cee2dSHong Zhang if (ts->quadraturets) { /* if there is integral in the cost function */ 1077*cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vecs_sensi[0],&ts->vec_drdu_col);CHKERRQ(ierr); 1078814e85d6SHong Zhang if (ts->vecs_sensip){ 1079*cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vecs_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr); 1080814e85d6SHong Zhang } 1081814e85d6SHong Zhang } 1082814e85d6SHong Zhang 1083814e85d6SHong Zhang if (ts->ops->adjointsetup) { 1084814e85d6SHong Zhang ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr); 1085814e85d6SHong Zhang } 1086814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 1087814e85d6SHong Zhang PetscFunctionReturn(0); 1088814e85d6SHong Zhang } 1089814e85d6SHong Zhang 1090814e85d6SHong Zhang /*@ 1091ece46509SHong Zhang TSAdjointReset - Resets a TSAdjoint context and removes any allocated Vecs and Mats. 1092ece46509SHong Zhang 1093ece46509SHong Zhang Collective on TS 1094ece46509SHong Zhang 1095ece46509SHong Zhang Input Parameter: 1096ece46509SHong Zhang . ts - the TS context obtained from TSCreate() 1097ece46509SHong Zhang 1098ece46509SHong Zhang Level: beginner 1099ece46509SHong Zhang 1100ece46509SHong Zhang .keywords: TS, timestep, reset 1101ece46509SHong Zhang 1102ece46509SHong Zhang .seealso: TSCreate(), TSAdjointSetup(), TSADestroy() 1103ece46509SHong Zhang @*/ 1104ece46509SHong Zhang PetscErrorCode TSAdjointReset(TS ts) 1105ece46509SHong Zhang { 1106ece46509SHong Zhang PetscErrorCode ierr; 1107ece46509SHong Zhang 1108ece46509SHong Zhang PetscFunctionBegin; 1109ece46509SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1110*cd4cee2dSHong Zhang if (ts->vecs_integral_sensip) { 1111*cd4cee2dSHong Zhang ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr); 1112*cd4cee2dSHong Zhang ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr); 1113*cd4cee2dSHong Zhang } 1114ece46509SHong Zhang if (ts->ops->adjointreset) { 1115ece46509SHong Zhang ierr = (*ts->ops->adjointreset)(ts);CHKERRQ(ierr); 1116ece46509SHong Zhang } 1117cda2db4bSHong Zhang if (ts->vec_dir) { /* second-order adjoint */ 1118cda2db4bSHong Zhang ierr = TSForwardReset(ts);CHKERRQ(ierr); 1119cda2db4bSHong Zhang } 1120ece46509SHong Zhang ts->vecs_sensi = NULL; 1121ece46509SHong Zhang ts->vecs_sensip = NULL; 1122ece46509SHong Zhang ts->vecs_sensi2 = NULL; 112333f72c5dSHong Zhang ts->vecs_sensi2p = NULL; 1124ece46509SHong Zhang ts->vec_dir = NULL; 1125ece46509SHong Zhang ts->adjointsetupcalled = PETSC_FALSE; 1126ece46509SHong Zhang PetscFunctionReturn(0); 1127ece46509SHong Zhang } 1128ece46509SHong Zhang 1129ece46509SHong Zhang /*@ 1130814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 1131814e85d6SHong Zhang 1132814e85d6SHong Zhang Logically Collective on TS 1133814e85d6SHong Zhang 1134814e85d6SHong Zhang Input Parameters: 1135814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1136814e85d6SHong Zhang . steps - number of steps to use 1137814e85d6SHong Zhang 1138814e85d6SHong Zhang Level: intermediate 1139814e85d6SHong Zhang 1140814e85d6SHong Zhang Notes: 1141814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 1142814e85d6SHong Zhang so as to integrate back to less than the original timestep 1143814e85d6SHong Zhang 1144814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations 1145814e85d6SHong Zhang 1146814e85d6SHong Zhang .seealso: TSSetExactFinalTime() 1147814e85d6SHong Zhang @*/ 1148814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 1149814e85d6SHong Zhang { 1150814e85d6SHong Zhang PetscFunctionBegin; 1151814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1152814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts,steps,2); 1153814e85d6SHong Zhang if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 1154814e85d6SHong Zhang if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 1155814e85d6SHong Zhang ts->adjoint_max_steps = steps; 1156814e85d6SHong Zhang PetscFunctionReturn(0); 1157814e85d6SHong Zhang } 1158814e85d6SHong Zhang 1159814e85d6SHong Zhang /*@C 1160814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 1161814e85d6SHong Zhang 1162814e85d6SHong Zhang Level: deprecated 1163814e85d6SHong Zhang 1164814e85d6SHong Zhang @*/ 1165814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 1166814e85d6SHong Zhang { 1167814e85d6SHong Zhang PetscErrorCode ierr; 1168814e85d6SHong Zhang 1169814e85d6SHong Zhang PetscFunctionBegin; 1170814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 1171814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 1172814e85d6SHong Zhang 1173814e85d6SHong Zhang ts->rhsjacobianp = func; 1174814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 1175814e85d6SHong Zhang if(Amat) { 1176814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 1177814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 1178814e85d6SHong Zhang ts->Jacp = Amat; 1179814e85d6SHong Zhang } 1180814e85d6SHong Zhang PetscFunctionReturn(0); 1181814e85d6SHong Zhang } 1182814e85d6SHong Zhang 1183814e85d6SHong Zhang /*@C 1184814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 1185814e85d6SHong Zhang 1186814e85d6SHong Zhang Level: deprecated 1187814e85d6SHong Zhang 1188814e85d6SHong Zhang @*/ 1189c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat) 1190814e85d6SHong Zhang { 1191814e85d6SHong Zhang PetscErrorCode ierr; 1192814e85d6SHong Zhang 1193814e85d6SHong Zhang PetscFunctionBegin; 1194814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1195c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1196814e85d6SHong Zhang PetscValidPointer(Amat,4); 1197814e85d6SHong Zhang 1198814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 1199c9ad9de0SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 1200814e85d6SHong Zhang PetscStackPop; 1201814e85d6SHong Zhang PetscFunctionReturn(0); 1202814e85d6SHong Zhang } 1203814e85d6SHong Zhang 1204814e85d6SHong Zhang /*@ 1205c9ad9de0SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction() 1206814e85d6SHong Zhang 1207814e85d6SHong Zhang Level: deprecated 1208814e85d6SHong Zhang 1209814e85d6SHong Zhang @*/ 1210c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 1211814e85d6SHong Zhang { 1212814e85d6SHong Zhang PetscErrorCode ierr; 1213814e85d6SHong Zhang 1214814e85d6SHong Zhang PetscFunctionBegin; 1215814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1216c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1217814e85d6SHong Zhang 1218814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 1219c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 1220814e85d6SHong Zhang PetscStackPop; 1221814e85d6SHong Zhang PetscFunctionReturn(0); 1222814e85d6SHong Zhang } 1223814e85d6SHong Zhang 1224814e85d6SHong Zhang /*@ 1225814e85d6SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction() 1226814e85d6SHong Zhang 1227814e85d6SHong Zhang Level: deprecated 1228814e85d6SHong Zhang 1229814e85d6SHong Zhang @*/ 1230c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 1231814e85d6SHong Zhang { 1232814e85d6SHong Zhang PetscErrorCode ierr; 1233814e85d6SHong Zhang 1234814e85d6SHong Zhang PetscFunctionBegin; 1235814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1236c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 1237814e85d6SHong Zhang 1238814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 1239c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 1240814e85d6SHong Zhang PetscStackPop; 1241814e85d6SHong Zhang PetscFunctionReturn(0); 1242814e85d6SHong Zhang } 1243814e85d6SHong Zhang 1244814e85d6SHong Zhang /*@C 1245814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 1246814e85d6SHong Zhang 1247814e85d6SHong Zhang Level: intermediate 1248814e85d6SHong Zhang 1249814e85d6SHong Zhang .keywords: TS, set, monitor 1250814e85d6SHong Zhang 1251814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1252814e85d6SHong Zhang @*/ 1253814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1254814e85d6SHong Zhang { 1255814e85d6SHong Zhang PetscErrorCode ierr; 1256814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1257814e85d6SHong Zhang 1258814e85d6SHong Zhang PetscFunctionBegin; 1259814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1260814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1261814e85d6SHong Zhang ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 1262814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1263814e85d6SHong Zhang PetscFunctionReturn(0); 1264814e85d6SHong Zhang } 1265814e85d6SHong Zhang 1266814e85d6SHong Zhang /*@C 1267814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 1268814e85d6SHong Zhang 1269814e85d6SHong Zhang Collective on TS 1270814e85d6SHong Zhang 1271814e85d6SHong Zhang Input Parameters: 1272814e85d6SHong Zhang + ts - TS object you wish to monitor 1273814e85d6SHong Zhang . name - the monitor type one is seeking 1274814e85d6SHong Zhang . help - message indicating what monitoring is done 1275814e85d6SHong Zhang . manual - manual page for the monitor 1276814e85d6SHong Zhang . monitor - the monitor function 1277814e85d6SHong 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 1278814e85d6SHong Zhang 1279814e85d6SHong Zhang Level: developer 1280814e85d6SHong Zhang 1281814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 1282814e85d6SHong Zhang PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 1283814e85d6SHong Zhang PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 1284814e85d6SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 1285814e85d6SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 1286814e85d6SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 1287814e85d6SHong Zhang PetscOptionsFList(), PetscOptionsEList() 1288814e85d6SHong Zhang @*/ 1289814e85d6SHong 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*)) 1290814e85d6SHong Zhang { 1291814e85d6SHong Zhang PetscErrorCode ierr; 1292814e85d6SHong Zhang PetscViewer viewer; 1293814e85d6SHong Zhang PetscViewerFormat format; 1294814e85d6SHong Zhang PetscBool flg; 1295814e85d6SHong Zhang 1296814e85d6SHong Zhang PetscFunctionBegin; 129716413a6aSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 1298814e85d6SHong Zhang if (flg) { 1299814e85d6SHong Zhang PetscViewerAndFormat *vf; 1300814e85d6SHong Zhang ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 1301814e85d6SHong Zhang ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 1302814e85d6SHong Zhang if (monitorsetup) { 1303814e85d6SHong Zhang ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 1304814e85d6SHong Zhang } 1305814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 1306814e85d6SHong Zhang } 1307814e85d6SHong Zhang PetscFunctionReturn(0); 1308814e85d6SHong Zhang } 1309814e85d6SHong Zhang 1310814e85d6SHong Zhang /*@C 1311814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 1312814e85d6SHong Zhang timestep to display the iteration's progress. 1313814e85d6SHong Zhang 1314814e85d6SHong Zhang Logically Collective on TS 1315814e85d6SHong Zhang 1316814e85d6SHong Zhang Input Parameters: 1317814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1318814e85d6SHong Zhang . adjointmonitor - monitoring routine 1319814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 1320814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 1321814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 1322814e85d6SHong Zhang (may be NULL) 1323814e85d6SHong Zhang 1324814e85d6SHong Zhang Calling sequence of monitor: 1325814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 1326814e85d6SHong Zhang 1327814e85d6SHong Zhang + ts - the TS context 1328814e85d6SHong 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 1329814e85d6SHong Zhang been interpolated to) 1330814e85d6SHong Zhang . time - current time 1331814e85d6SHong Zhang . u - current iterate 1332814e85d6SHong Zhang . numcost - number of cost functionos 1333814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1334814e85d6SHong Zhang . mu - sensitivities to parameters 1335814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 1336814e85d6SHong Zhang 1337814e85d6SHong Zhang Notes: 1338814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 1339814e85d6SHong Zhang already has been loaded. 1340814e85d6SHong Zhang 1341814e85d6SHong Zhang Fortran Notes: 1342814e85d6SHong Zhang Only a single monitor function can be set for each TS object 1343814e85d6SHong Zhang 1344814e85d6SHong Zhang Level: intermediate 1345814e85d6SHong Zhang 1346814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 1347814e85d6SHong Zhang 1348814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel() 1349814e85d6SHong Zhang @*/ 1350814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 1351814e85d6SHong Zhang { 1352814e85d6SHong Zhang PetscErrorCode ierr; 1353814e85d6SHong Zhang PetscInt i; 1354814e85d6SHong Zhang PetscBool identical; 1355814e85d6SHong Zhang 1356814e85d6SHong Zhang PetscFunctionBegin; 1357814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1358814e85d6SHong Zhang for (i=0; i<ts->numbermonitors;i++) { 1359814e85d6SHong Zhang ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 1360814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 1361814e85d6SHong Zhang } 1362814e85d6SHong Zhang if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 1363814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 1364814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 1365814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 1366814e85d6SHong Zhang PetscFunctionReturn(0); 1367814e85d6SHong Zhang } 1368814e85d6SHong Zhang 1369814e85d6SHong Zhang /*@C 1370814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 1371814e85d6SHong Zhang 1372814e85d6SHong Zhang Logically Collective on TS 1373814e85d6SHong Zhang 1374814e85d6SHong Zhang Input Parameters: 1375814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1376814e85d6SHong Zhang 1377814e85d6SHong Zhang Notes: 1378814e85d6SHong Zhang There is no way to remove a single, specific monitor. 1379814e85d6SHong Zhang 1380814e85d6SHong Zhang Level: intermediate 1381814e85d6SHong Zhang 1382814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 1383814e85d6SHong Zhang 1384814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1385814e85d6SHong Zhang @*/ 1386814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts) 1387814e85d6SHong Zhang { 1388814e85d6SHong Zhang PetscErrorCode ierr; 1389814e85d6SHong Zhang PetscInt i; 1390814e85d6SHong Zhang 1391814e85d6SHong Zhang PetscFunctionBegin; 1392814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1393814e85d6SHong Zhang for (i=0; i<ts->numberadjointmonitors; i++) { 1394814e85d6SHong Zhang if (ts->adjointmonitordestroy[i]) { 1395814e85d6SHong Zhang ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1396814e85d6SHong Zhang } 1397814e85d6SHong Zhang } 1398814e85d6SHong Zhang ts->numberadjointmonitors = 0; 1399814e85d6SHong Zhang PetscFunctionReturn(0); 1400814e85d6SHong Zhang } 1401814e85d6SHong Zhang 1402814e85d6SHong Zhang /*@C 1403814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 1404814e85d6SHong Zhang 1405814e85d6SHong Zhang Level: intermediate 1406814e85d6SHong Zhang 1407814e85d6SHong Zhang .keywords: TS, set, monitor 1408814e85d6SHong Zhang 1409814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 1410814e85d6SHong Zhang @*/ 1411814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 1412814e85d6SHong Zhang { 1413814e85d6SHong Zhang PetscErrorCode ierr; 1414814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 1415814e85d6SHong Zhang 1416814e85d6SHong Zhang PetscFunctionBegin; 1417814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 1418814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 1419814e85d6SHong Zhang ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1420814e85d6SHong 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); 1421814e85d6SHong Zhang ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 1422814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 1423814e85d6SHong Zhang PetscFunctionReturn(0); 1424814e85d6SHong Zhang } 1425814e85d6SHong Zhang 1426814e85d6SHong Zhang /*@C 1427814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 1428814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 1429814e85d6SHong Zhang 1430814e85d6SHong Zhang Collective on TS 1431814e85d6SHong Zhang 1432814e85d6SHong Zhang Input Parameters: 1433814e85d6SHong Zhang + ts - the TS context 1434814e85d6SHong Zhang . step - current time-step 1435814e85d6SHong Zhang . ptime - current time 1436814e85d6SHong Zhang . u - current state 1437814e85d6SHong Zhang . numcost - number of cost functions 1438814e85d6SHong Zhang . lambda - sensitivities to initial conditions 1439814e85d6SHong Zhang . mu - sensitivities to parameters 1440814e85d6SHong Zhang - dummy - either a viewer or NULL 1441814e85d6SHong Zhang 1442814e85d6SHong Zhang Level: intermediate 1443814e85d6SHong Zhang 1444814e85d6SHong Zhang .keywords: TS, vector, adjoint, monitor, view 1445814e85d6SHong Zhang 1446814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 1447814e85d6SHong Zhang @*/ 1448814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 1449814e85d6SHong Zhang { 1450814e85d6SHong Zhang PetscErrorCode ierr; 1451814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 1452814e85d6SHong Zhang PetscDraw draw; 1453814e85d6SHong Zhang PetscReal xl,yl,xr,yr,h; 1454814e85d6SHong Zhang char time[32]; 1455814e85d6SHong Zhang 1456814e85d6SHong Zhang PetscFunctionBegin; 1457814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 1458814e85d6SHong Zhang 1459814e85d6SHong Zhang ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 1460814e85d6SHong Zhang ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 1461814e85d6SHong Zhang ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 1462814e85d6SHong Zhang ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 1463814e85d6SHong Zhang h = yl + .95*(yr - yl); 1464814e85d6SHong Zhang ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 1465814e85d6SHong Zhang ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 1466814e85d6SHong Zhang PetscFunctionReturn(0); 1467814e85d6SHong Zhang } 1468814e85d6SHong Zhang 1469814e85d6SHong Zhang /* 1470814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 1471814e85d6SHong Zhang 1472814e85d6SHong Zhang Collective on TSAdjoint 1473814e85d6SHong Zhang 1474814e85d6SHong Zhang Input Parameter: 1475814e85d6SHong Zhang . ts - the TS context 1476814e85d6SHong Zhang 1477814e85d6SHong Zhang Options Database Keys: 1478814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 1479814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 1480814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 1481814e85d6SHong Zhang 1482814e85d6SHong Zhang Level: developer 1483814e85d6SHong Zhang 1484814e85d6SHong Zhang Notes: 1485814e85d6SHong Zhang This is not normally called directly by users 1486814e85d6SHong Zhang 1487814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database 1488814e85d6SHong Zhang 1489814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 1490814e85d6SHong Zhang */ 1491814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 1492814e85d6SHong Zhang { 1493814e85d6SHong Zhang PetscBool tflg,opt; 1494814e85d6SHong Zhang PetscErrorCode ierr; 1495814e85d6SHong Zhang 1496814e85d6SHong Zhang PetscFunctionBegin; 1497814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,2); 1498814e85d6SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 1499814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 1500814e85d6SHong Zhang ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 1501814e85d6SHong Zhang if (opt) { 1502814e85d6SHong Zhang ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 1503814e85d6SHong Zhang ts->adjoint_solve = tflg; 1504814e85d6SHong Zhang } 1505814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 1506814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 1507814e85d6SHong Zhang opt = PETSC_FALSE; 1508814e85d6SHong Zhang ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 1509814e85d6SHong Zhang if (opt) { 1510814e85d6SHong Zhang TSMonitorDrawCtx ctx; 1511814e85d6SHong Zhang PetscInt howoften = 1; 1512814e85d6SHong Zhang 1513814e85d6SHong Zhang ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 1514814e85d6SHong Zhang ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 1515814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 1516814e85d6SHong Zhang } 1517814e85d6SHong Zhang PetscFunctionReturn(0); 1518814e85d6SHong Zhang } 1519814e85d6SHong Zhang 1520814e85d6SHong Zhang /*@ 1521814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 1522814e85d6SHong Zhang 1523814e85d6SHong Zhang Collective on TS 1524814e85d6SHong Zhang 1525814e85d6SHong Zhang Input Parameter: 1526814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1527814e85d6SHong Zhang 1528814e85d6SHong Zhang Level: intermediate 1529814e85d6SHong Zhang 1530814e85d6SHong Zhang .keywords: TS, adjoint, step 1531814e85d6SHong Zhang 1532814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve() 1533814e85d6SHong Zhang @*/ 1534814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts) 1535814e85d6SHong Zhang { 1536814e85d6SHong Zhang DM dm; 1537814e85d6SHong Zhang PetscErrorCode ierr; 1538814e85d6SHong Zhang 1539814e85d6SHong Zhang PetscFunctionBegin; 1540814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1541814e85d6SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 1542814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1543814e85d6SHong Zhang 1544814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1545814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 1546814e85d6SHong 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); 1547814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1548814e85d6SHong Zhang ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 1549814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 1550814e85d6SHong Zhang ts->adjoint_steps++; ts->steps--; 1551814e85d6SHong Zhang 1552814e85d6SHong Zhang if (ts->reason < 0) { 1553814e85d6SHong Zhang if (ts->errorifstepfailed) { 155405c9950eSHong 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]); 155505c9950eSHong Zhang else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSAdjointStep has failed due to %s",TSConvergedReasons[ts->reason]); 1556814e85d6SHong Zhang } 1557814e85d6SHong Zhang } else if (!ts->reason) { 1558814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1559814e85d6SHong Zhang } 1560814e85d6SHong Zhang PetscFunctionReturn(0); 1561814e85d6SHong Zhang } 1562814e85d6SHong Zhang 1563814e85d6SHong Zhang /*@ 1564814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 1565814e85d6SHong Zhang 1566814e85d6SHong Zhang Collective on TS 1567814e85d6SHong Zhang 1568814e85d6SHong Zhang Input Parameter: 1569814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1570814e85d6SHong Zhang 1571814e85d6SHong Zhang Options Database: 1572814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 1573814e85d6SHong Zhang 1574814e85d6SHong Zhang Level: intermediate 1575814e85d6SHong Zhang 1576814e85d6SHong Zhang Notes: 1577814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 1578814e85d6SHong Zhang 1579814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 1580814e85d6SHong Zhang 1581814e85d6SHong Zhang .keywords: TS, timestep, solve 1582814e85d6SHong Zhang 1583814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 1584814e85d6SHong Zhang @*/ 1585814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts) 1586814e85d6SHong Zhang { 1587814e85d6SHong Zhang PetscErrorCode ierr; 1588814e85d6SHong Zhang 1589814e85d6SHong Zhang PetscFunctionBegin; 1590814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1591814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 1592814e85d6SHong Zhang 1593814e85d6SHong Zhang /* reset time step and iteration counters */ 1594814e85d6SHong Zhang ts->adjoint_steps = 0; 1595814e85d6SHong Zhang ts->ksp_its = 0; 1596814e85d6SHong Zhang ts->snes_its = 0; 1597814e85d6SHong Zhang ts->num_snes_failures = 0; 1598814e85d6SHong Zhang ts->reject = 0; 1599814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 1600814e85d6SHong Zhang 1601814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 1602814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 1603814e85d6SHong Zhang 1604814e85d6SHong Zhang while (!ts->reason) { 1605814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1606814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1607814e85d6SHong Zhang ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 1608814e85d6SHong Zhang ierr = TSAdjointStep(ts);CHKERRQ(ierr); 1609*cd4cee2dSHong Zhang if ((ts->vec_costintegral || ts->quadraturets) && !ts->costintegralfwd) { 1610814e85d6SHong Zhang ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 1611814e85d6SHong Zhang } 1612814e85d6SHong Zhang } 1613814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 1614814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 1615814e85d6SHong Zhang ts->solvetime = ts->ptime; 1616814e85d6SHong Zhang ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 1617814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 1618814e85d6SHong Zhang ts->adjoint_max_steps = 0; 1619814e85d6SHong Zhang PetscFunctionReturn(0); 1620814e85d6SHong Zhang } 1621814e85d6SHong Zhang 1622814e85d6SHong Zhang /*@C 1623814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 1624814e85d6SHong Zhang 1625814e85d6SHong Zhang Collective on TS 1626814e85d6SHong Zhang 1627814e85d6SHong Zhang Input Parameters: 1628814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 1629814e85d6SHong Zhang . step - step number that has just completed 1630814e85d6SHong Zhang . ptime - model time of the state 1631814e85d6SHong Zhang . u - state at the current model time 1632814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 1633814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 1634814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 1635814e85d6SHong Zhang 1636814e85d6SHong Zhang Notes: 1637814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 1638814e85d6SHong Zhang Users would almost never call this routine directly. 1639814e85d6SHong Zhang 1640814e85d6SHong Zhang Level: developer 1641814e85d6SHong Zhang 1642814e85d6SHong Zhang .keywords: TS, timestep 1643814e85d6SHong Zhang @*/ 1644814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 1645814e85d6SHong Zhang { 1646814e85d6SHong Zhang PetscErrorCode ierr; 1647814e85d6SHong Zhang PetscInt i,n = ts->numberadjointmonitors; 1648814e85d6SHong Zhang 1649814e85d6SHong Zhang PetscFunctionBegin; 1650814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1651814e85d6SHong Zhang PetscValidHeaderSpecific(u,VEC_CLASSID,4); 16528860a134SJunchao Zhang ierr = VecLockReadPush(u);CHKERRQ(ierr); 1653814e85d6SHong Zhang for (i=0; i<n; i++) { 1654814e85d6SHong Zhang ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 1655814e85d6SHong Zhang } 16568860a134SJunchao Zhang ierr = VecLockReadPop(u);CHKERRQ(ierr); 1657814e85d6SHong Zhang PetscFunctionReturn(0); 1658814e85d6SHong Zhang } 1659814e85d6SHong Zhang 1660814e85d6SHong Zhang /*@ 1661814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 1662814e85d6SHong Zhang 1663814e85d6SHong Zhang Collective on TS 1664814e85d6SHong Zhang 1665814e85d6SHong Zhang Input Arguments: 1666814e85d6SHong Zhang . ts - time stepping context 1667814e85d6SHong Zhang 1668814e85d6SHong Zhang Level: advanced 1669814e85d6SHong Zhang 1670814e85d6SHong Zhang Notes: 1671814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 1672814e85d6SHong Zhang 1673814e85d6SHong Zhang .seealso: TSAdjointSolve(), TSAdjointStep 1674814e85d6SHong Zhang @*/ 1675814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts) 1676814e85d6SHong Zhang { 1677814e85d6SHong Zhang PetscErrorCode ierr; 1678814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1679814e85d6SHong 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); 1680814e85d6SHong Zhang ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 1681814e85d6SHong Zhang PetscFunctionReturn(0); 1682814e85d6SHong Zhang } 1683814e85d6SHong Zhang 1684814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 1685814e85d6SHong Zhang 1686814e85d6SHong Zhang /*@ 1687814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 1688814e85d6SHong Zhang of forward sensitivity analysis 1689814e85d6SHong Zhang 1690814e85d6SHong Zhang Collective on TS 1691814e85d6SHong Zhang 1692814e85d6SHong Zhang Input Parameter: 1693814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 1694814e85d6SHong Zhang 1695814e85d6SHong Zhang Level: advanced 1696814e85d6SHong Zhang 1697814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup 1698814e85d6SHong Zhang 1699814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp() 1700814e85d6SHong Zhang @*/ 1701814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts) 1702814e85d6SHong Zhang { 1703814e85d6SHong Zhang PetscErrorCode ierr; 1704814e85d6SHong Zhang 1705814e85d6SHong Zhang PetscFunctionBegin; 1706814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1707814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 1708814e85d6SHong Zhang if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()"); 1709814e85d6SHong Zhang if (ts->vecs_integral_sensip) { 1710*cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vec_sol,&ts->vec_drdu_col);CHKERRQ(ierr); 1711*cd4cee2dSHong Zhang ierr = VecDuplicate(ts->vecs_integral_sensip[0],&ts->vec_drdp_col);CHKERRQ(ierr); 1712814e85d6SHong Zhang } 171333f72c5dSHong Zhang if (!ts->Jacp && ts->Jacprhs) ts->Jacp = ts->Jacprhs; 1714814e85d6SHong Zhang 1715814e85d6SHong Zhang if (ts->ops->forwardsetup) { 1716814e85d6SHong Zhang ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 1717814e85d6SHong Zhang } 1718814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 1719814e85d6SHong Zhang PetscFunctionReturn(0); 1720814e85d6SHong Zhang } 1721814e85d6SHong Zhang 1722814e85d6SHong Zhang /*@ 17237adebddeSHong Zhang TSForwardReset - Reset the internal data structures used by forward sensitivity analysis 17247adebddeSHong Zhang 17257adebddeSHong Zhang Collective on TS 17267adebddeSHong Zhang 17277adebddeSHong Zhang Input Parameter: 17287adebddeSHong Zhang . ts - the TS context obtained from TSCreate() 17297adebddeSHong Zhang 17307adebddeSHong Zhang Level: advanced 17317adebddeSHong Zhang 17327adebddeSHong Zhang .keywords: TS, forward sensitivity, reset 17337adebddeSHong Zhang 17347adebddeSHong Zhang .seealso: TSCreate(), TSDestroy(), TSForwardSetUp() 17357adebddeSHong Zhang @*/ 17367adebddeSHong Zhang PetscErrorCode TSForwardReset(TS ts) 17377adebddeSHong Zhang { 17387adebddeSHong Zhang PetscErrorCode ierr; 17397adebddeSHong Zhang 17407adebddeSHong Zhang PetscFunctionBegin; 17417adebddeSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1742*cd4cee2dSHong Zhang if (ts->vecs_integral_sensip) { 1743*cd4cee2dSHong Zhang ierr = VecDestroy(&ts->vec_drdu_col);CHKERRQ(ierr); 1744*cd4cee2dSHong Zhang ierr = VecDestroy(&ts->vec_drdp_col);CHKERRQ(ierr); 1745*cd4cee2dSHong Zhang } 17467adebddeSHong Zhang if (ts->ops->forwardreset) { 17477adebddeSHong Zhang ierr = (*ts->ops->forwardreset)(ts);CHKERRQ(ierr); 17487adebddeSHong Zhang } 17497adebddeSHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 17507adebddeSHong Zhang ts->vecs_integral_sensip = NULL; 17517adebddeSHong Zhang ts->forward_solve = PETSC_FALSE; 17527adebddeSHong Zhang ts->forwardsetupcalled = PETSC_FALSE; 17537adebddeSHong Zhang PetscFunctionReturn(0); 17547adebddeSHong Zhang } 17557adebddeSHong Zhang 17567adebddeSHong Zhang /*@ 1757814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 1758814e85d6SHong Zhang 1759814e85d6SHong Zhang Input Parameter: 1760814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1761814e85d6SHong Zhang . numfwdint- number of integrals 1762814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1763814e85d6SHong Zhang 1764814e85d6SHong Zhang Level: intermediate 1765814e85d6SHong Zhang 1766814e85d6SHong Zhang .keywords: TS, forward sensitivity 1767814e85d6SHong Zhang 1768814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1769814e85d6SHong Zhang @*/ 1770814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 1771814e85d6SHong Zhang { 1772814e85d6SHong Zhang PetscFunctionBegin; 1773814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1774814e85d6SHong 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()"); 1775814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1776814e85d6SHong Zhang 1777814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 1778814e85d6SHong Zhang PetscFunctionReturn(0); 1779814e85d6SHong Zhang } 1780814e85d6SHong Zhang 1781814e85d6SHong Zhang /*@ 1782814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1783814e85d6SHong Zhang 1784814e85d6SHong Zhang Input Parameter: 1785814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1786814e85d6SHong Zhang 1787814e85d6SHong Zhang Output Parameter: 1788814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1789814e85d6SHong Zhang 1790814e85d6SHong Zhang Level: intermediate 1791814e85d6SHong Zhang 1792814e85d6SHong Zhang .keywords: TS, forward sensitivity 1793814e85d6SHong Zhang 1794814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1795814e85d6SHong Zhang @*/ 1796814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1797814e85d6SHong Zhang { 1798814e85d6SHong Zhang PetscFunctionBegin; 1799814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1800814e85d6SHong Zhang PetscValidPointer(vp,3); 1801814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1802814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1803814e85d6SHong Zhang PetscFunctionReturn(0); 1804814e85d6SHong Zhang } 1805814e85d6SHong Zhang 1806814e85d6SHong Zhang /*@ 1807814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1808814e85d6SHong Zhang 1809814e85d6SHong Zhang Collective on TS 1810814e85d6SHong Zhang 1811814e85d6SHong Zhang Input Arguments: 1812814e85d6SHong Zhang . ts - time stepping context 1813814e85d6SHong Zhang 1814814e85d6SHong Zhang Level: advanced 1815814e85d6SHong Zhang 1816814e85d6SHong Zhang Notes: 1817814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1818814e85d6SHong Zhang 1819814e85d6SHong Zhang .keywords: TS, forward sensitivity 1820814e85d6SHong Zhang 1821814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1822814e85d6SHong Zhang @*/ 1823814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts) 1824814e85d6SHong Zhang { 1825814e85d6SHong Zhang PetscErrorCode ierr; 1826814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1827814e85d6SHong Zhang if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1828814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1829814e85d6SHong Zhang ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1830814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1831814e85d6SHong Zhang PetscFunctionReturn(0); 1832814e85d6SHong Zhang } 1833814e85d6SHong Zhang 1834814e85d6SHong Zhang /*@ 1835814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1836814e85d6SHong Zhang 1837814e85d6SHong Zhang Logically Collective on TS and Vec 1838814e85d6SHong Zhang 1839814e85d6SHong Zhang Input Parameters: 1840814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1841814e85d6SHong Zhang . nump - number of parameters 1842814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1843814e85d6SHong Zhang 1844814e85d6SHong Zhang Level: beginner 1845814e85d6SHong Zhang 1846814e85d6SHong Zhang Notes: 1847814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1848814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1849814e85d6SHong Zhang You must call this function before TSSolve(). 1850814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1851814e85d6SHong Zhang 1852814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values 1853814e85d6SHong Zhang 1854814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1855814e85d6SHong Zhang @*/ 1856814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1857814e85d6SHong Zhang { 1858814e85d6SHong Zhang PetscErrorCode ierr; 1859814e85d6SHong Zhang 1860814e85d6SHong Zhang PetscFunctionBegin; 1861814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1862814e85d6SHong Zhang PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1863814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1864b5b4017aSHong Zhang if (nump == PETSC_DEFAULT) { 1865b5b4017aSHong Zhang ierr = MatGetSize(Smat,NULL,&ts->num_parameters);CHKERRQ(ierr); 1866b5b4017aSHong Zhang } else ts->num_parameters = nump; 1867814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1868814e85d6SHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1869814e85d6SHong Zhang ts->mat_sensip = Smat; 1870814e85d6SHong Zhang PetscFunctionReturn(0); 1871814e85d6SHong Zhang } 1872814e85d6SHong Zhang 1873814e85d6SHong Zhang /*@ 1874814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1875814e85d6SHong Zhang 1876814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1877814e85d6SHong Zhang 1878814e85d6SHong Zhang Output Parameter: 1879814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1880814e85d6SHong Zhang . nump - number of parameters 1881814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1882814e85d6SHong Zhang 1883814e85d6SHong Zhang Level: intermediate 1884814e85d6SHong Zhang 1885814e85d6SHong Zhang .keywords: TS, forward sensitivity 1886814e85d6SHong Zhang 1887814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1888814e85d6SHong Zhang @*/ 1889814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1890814e85d6SHong Zhang { 1891814e85d6SHong Zhang PetscFunctionBegin; 1892814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1893814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1894814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1895814e85d6SHong Zhang PetscFunctionReturn(0); 1896814e85d6SHong Zhang } 1897814e85d6SHong Zhang 1898814e85d6SHong Zhang /*@ 1899814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1900814e85d6SHong Zhang 1901814e85d6SHong Zhang Collective on TS 1902814e85d6SHong Zhang 1903814e85d6SHong Zhang Input Arguments: 1904814e85d6SHong Zhang . ts - time stepping context 1905814e85d6SHong Zhang 1906814e85d6SHong Zhang Level: advanced 1907814e85d6SHong Zhang 1908814e85d6SHong Zhang Notes: 1909814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1910814e85d6SHong Zhang 1911814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral() 1912814e85d6SHong Zhang @*/ 1913814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts) 1914814e85d6SHong Zhang { 1915814e85d6SHong Zhang PetscErrorCode ierr; 1916814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1917814e85d6SHong 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); 1918814e85d6SHong Zhang ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1919814e85d6SHong Zhang PetscFunctionReturn(0); 1920814e85d6SHong Zhang } 1921b5b4017aSHong Zhang 1922b5b4017aSHong Zhang /*@ 1923b5b4017aSHong Zhang TSForwardSetInitialSensitivities - Set initial values for tangent linear sensitivities 1924b5b4017aSHong Zhang 1925b5b4017aSHong Zhang Collective on TS and Mat 1926b5b4017aSHong Zhang 1927b5b4017aSHong Zhang Input Parameter 1928b5b4017aSHong Zhang + ts - the TS context obtained from TSCreate() 1929b5b4017aSHong Zhang - didp - parametric sensitivities of the initial condition 1930b5b4017aSHong Zhang 1931b5b4017aSHong Zhang Level: intermediate 1932b5b4017aSHong Zhang 1933b5b4017aSHong 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. 1934b5b4017aSHong Zhang 1935b5b4017aSHong Zhang .seealso: TSForwardSetSensitivities() 1936b5b4017aSHong Zhang @*/ 1937b5b4017aSHong Zhang PetscErrorCode TSForwardSetInitialSensitivities(TS ts,Mat didp) 1938b5b4017aSHong Zhang { 1939b5b4017aSHong Zhang PetscErrorCode ierr; 1940b5b4017aSHong Zhang 1941b5b4017aSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1942b5b4017aSHong Zhang PetscValidHeaderSpecific(didp,MAT_CLASSID,2); 1943b5b4017aSHong Zhang if (!ts->mat_sensip) { 1944b5b4017aSHong Zhang ierr = TSForwardSetSensitivities(ts,PETSC_DEFAULT,didp);CHKERRQ(ierr); 1945b5b4017aSHong Zhang } 1946b5b4017aSHong Zhang PetscFunctionReturn(0); 1947b5b4017aSHong Zhang } 19486affb6f8SHong Zhang 19496affb6f8SHong Zhang /*@ 19506affb6f8SHong Zhang TSForwardGetStages - Get the number of stages and the tangent linear sensitivities at the intermediate stages 19516affb6f8SHong Zhang 19526affb6f8SHong Zhang Input Parameter: 19536affb6f8SHong Zhang . ts - the TS context obtained from TSCreate() 19546affb6f8SHong Zhang 19556affb6f8SHong Zhang Output Parameters: 1956*cd4cee2dSHong Zhang + ns - number of stages 19576affb6f8SHong Zhang - S - tangent linear sensitivities at the intermediate stages 19586affb6f8SHong Zhang 19596affb6f8SHong Zhang Level: advanced 19606affb6f8SHong Zhang 19616affb6f8SHong Zhang .keywords: TS, second-order adjoint, forward sensitivity 19626affb6f8SHong Zhang @*/ 19636affb6f8SHong Zhang PetscErrorCode TSForwardGetStages(TS ts,PetscInt *ns,Mat **S) 19646affb6f8SHong Zhang { 19656affb6f8SHong Zhang PetscErrorCode ierr; 19666affb6f8SHong Zhang 19676affb6f8SHong Zhang PetscFunctionBegin; 19686affb6f8SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 19696affb6f8SHong Zhang 19706affb6f8SHong Zhang if (!ts->ops->getstages) *S=NULL; 19716affb6f8SHong Zhang else { 19726affb6f8SHong Zhang ierr = (*ts->ops->forwardgetstages)(ts,ns,S);CHKERRQ(ierr); 19736affb6f8SHong Zhang } 19746affb6f8SHong Zhang PetscFunctionReturn(0); 19756affb6f8SHong Zhang } 1976*cd4cee2dSHong Zhang 1977*cd4cee2dSHong Zhang /*@ 1978*cd4cee2dSHong Zhang TSCreateQuadratureTS - Create a sub-TS that evaluates integrals over time 1979*cd4cee2dSHong Zhang 1980*cd4cee2dSHong Zhang Input Parameter: 1981*cd4cee2dSHong Zhang + ts - the TS context obtained from TSCreate() 1982*cd4cee2dSHong Zhang - fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 1983*cd4cee2dSHong Zhang 1984*cd4cee2dSHong Zhang Output Parameters: 1985*cd4cee2dSHong Zhang . quadts - the child TS context 1986*cd4cee2dSHong Zhang 1987*cd4cee2dSHong Zhang Level: intermediate 1988*cd4cee2dSHong Zhang 1989*cd4cee2dSHong Zhang .keywords: TS, quadrature evaluation 1990*cd4cee2dSHong Zhang 1991*cd4cee2dSHong Zhang .seealso: TSGetQuadratureTS() 1992*cd4cee2dSHong Zhang @*/ 1993*cd4cee2dSHong Zhang PetscErrorCode TSCreateQuadratureTS(TS ts,PetscBool fwd,TS *quadts) 1994*cd4cee2dSHong Zhang { 1995*cd4cee2dSHong Zhang char prefix[128]; 1996*cd4cee2dSHong Zhang PetscErrorCode ierr; 1997*cd4cee2dSHong Zhang 1998*cd4cee2dSHong Zhang PetscFunctionBegin; 1999*cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2000*cd4cee2dSHong Zhang PetscValidPointer(quadts,2); 2001*cd4cee2dSHong Zhang ierr = TSDestroy(&ts->quadraturets);CHKERRQ(ierr); 2002*cd4cee2dSHong Zhang ierr = TSCreate(PetscObjectComm((PetscObject)ts),&ts->quadraturets);CHKERRQ(ierr); 2003*cd4cee2dSHong Zhang ierr = PetscObjectIncrementTabLevel((PetscObject)ts->quadraturets,(PetscObject)ts,1);CHKERRQ(ierr); 2004*cd4cee2dSHong Zhang ierr = PetscLogObjectParent((PetscObject)ts,(PetscObject)ts->quadraturets);CHKERRQ(ierr); 2005*cd4cee2dSHong Zhang ierr = PetscSNPrintf(prefix,sizeof(prefix),"%squad_",((PetscObject)ts)->prefix ? ((PetscObject)ts)->prefix : "");CHKERRQ(ierr); 2006*cd4cee2dSHong Zhang ierr = TSSetOptionsPrefix(ts->quadraturets,prefix);CHKERRQ(ierr); 2007*cd4cee2dSHong Zhang *quadts = ts->quadraturets; 2008*cd4cee2dSHong Zhang 2009*cd4cee2dSHong Zhang if (ts->numcost) { 2010*cd4cee2dSHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,ts->numcost,&(*quadts)->vec_sol);CHKERRQ(ierr); 2011*cd4cee2dSHong Zhang } else { 2012*cd4cee2dSHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,1,&(*quadts)->vec_sol);CHKERRQ(ierr); 2013*cd4cee2dSHong Zhang } 2014*cd4cee2dSHong Zhang ierr = VecDuplicate((*quadts)->vec_sol,&ts->vec_costintegrand);CHKERRQ(ierr); 2015*cd4cee2dSHong Zhang ts->costintegralfwd = fwd; 2016*cd4cee2dSHong Zhang PetscFunctionReturn(0); 2017*cd4cee2dSHong Zhang } 2018*cd4cee2dSHong Zhang 2019*cd4cee2dSHong Zhang /*@ 2020*cd4cee2dSHong Zhang TSGetQuadratureTS - Return the sub-TS that evaluates integrals over time 2021*cd4cee2dSHong Zhang 2022*cd4cee2dSHong Zhang Input Parameter: 2023*cd4cee2dSHong Zhang . ts - the TS context obtained from TSCreate() 2024*cd4cee2dSHong Zhang 2025*cd4cee2dSHong Zhang Output Parameters: 2026*cd4cee2dSHong Zhang + fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 2027*cd4cee2dSHong Zhang - quadts - the child TS context 2028*cd4cee2dSHong Zhang 2029*cd4cee2dSHong Zhang Level: intermediate 2030*cd4cee2dSHong Zhang 2031*cd4cee2dSHong Zhang .keywords: TS, quadrature evaluation 2032*cd4cee2dSHong Zhang 2033*cd4cee2dSHong Zhang .seealso: TSCreateQuadratureTS() 2034*cd4cee2dSHong Zhang @*/ 2035*cd4cee2dSHong Zhang PetscErrorCode TSGetQuadratureTS(TS ts,PetscBool *fwd,TS *quadts) 2036*cd4cee2dSHong Zhang { 2037*cd4cee2dSHong Zhang PetscFunctionBegin; 2038*cd4cee2dSHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 2039*cd4cee2dSHong Zhang if (fwd) *fwd = ts->costintegralfwd; 2040*cd4cee2dSHong Zhang if (quadts) *quadts = ts->quadraturets; 2041*cd4cee2dSHong Zhang PetscFunctionReturn(0); 2042*cd4cee2dSHong Zhang } 2043