1814e85d6SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2814e85d6SHong Zhang #include <petscdraw.h> 3814e85d6SHong Zhang 4814e85d6SHong Zhang PetscLogEvent TS_AdjointStep, TS_ForwardStep; 5814e85d6SHong Zhang 6814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/ 7814e85d6SHong Zhang 8814e85d6SHong Zhang /*@C 9*c9ad9de0SHong 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: 14814e85d6SHong Zhang + ts - The TS context obtained from TSCreate() 15814e85d6SHong Zhang - func - The function 16814e85d6SHong Zhang 17814e85d6SHong Zhang Calling sequence of func: 18814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 19814e85d6SHong Zhang + t - current timestep 20*c9ad9de0SHong Zhang . U - input vector (current ODE solution) 21814e85d6SHong Zhang . A - output matrix 22814e85d6SHong Zhang - ctx - [optional] user-defined function context 23814e85d6SHong Zhang 24814e85d6SHong Zhang Level: intermediate 25814e85d6SHong Zhang 26814e85d6SHong Zhang Notes: 27814e85d6SHong 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 28814e85d6SHong Zhang 29814e85d6SHong Zhang .keywords: TS, sensitivity 30814e85d6SHong Zhang .seealso: 31814e85d6SHong Zhang @*/ 32814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 33814e85d6SHong Zhang { 34814e85d6SHong Zhang PetscErrorCode ierr; 35814e85d6SHong Zhang 36814e85d6SHong Zhang PetscFunctionBegin; 37814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 38814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 39814e85d6SHong Zhang 40814e85d6SHong Zhang ts->rhsjacobianp = func; 41814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 42814e85d6SHong Zhang if(Amat) { 43814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 44814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 45814e85d6SHong Zhang ts->Jacp = Amat; 46814e85d6SHong Zhang } 47814e85d6SHong Zhang PetscFunctionReturn(0); 48814e85d6SHong Zhang } 49814e85d6SHong Zhang 50814e85d6SHong Zhang /*@C 51814e85d6SHong Zhang TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 52814e85d6SHong Zhang 53814e85d6SHong Zhang Collective on TS 54814e85d6SHong Zhang 55814e85d6SHong Zhang Input Parameters: 56814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 57814e85d6SHong Zhang 58814e85d6SHong Zhang Level: developer 59814e85d6SHong Zhang 60814e85d6SHong Zhang .keywords: TS, sensitivity 61814e85d6SHong Zhang .seealso: TSSetRHSJacobianP() 62814e85d6SHong Zhang @*/ 63*c9ad9de0SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec U,Mat Amat) 64814e85d6SHong Zhang { 65814e85d6SHong Zhang PetscErrorCode ierr; 66814e85d6SHong Zhang 67814e85d6SHong Zhang PetscFunctionBegin; 68814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 69*c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 70814e85d6SHong Zhang PetscValidPointer(Amat,4); 71814e85d6SHong Zhang 72814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 73*c9ad9de0SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx);CHKERRQ(ierr); 74814e85d6SHong Zhang PetscStackPop; 75814e85d6SHong Zhang PetscFunctionReturn(0); 76814e85d6SHong Zhang } 77814e85d6SHong Zhang 78814e85d6SHong Zhang /*@C 79814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 80814e85d6SHong Zhang 81814e85d6SHong Zhang Logically Collective on TS 82814e85d6SHong Zhang 83814e85d6SHong Zhang Input Parameters: 84814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 85814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 86814e85d6SHong Zhang . costintegral - vector that stores the integral values 87814e85d6SHong Zhang . rf - routine for evaluating the integrand function 88*c9ad9de0SHong Zhang . drduf - function that computes the gradients of the r's with respect to u 89814e85d6SHong 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) 90814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 91814e85d6SHong Zhang - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 92814e85d6SHong Zhang 93814e85d6SHong Zhang Calling sequence of rf: 94*c9ad9de0SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec U,Vec F,void *ctx); 95814e85d6SHong Zhang 96*c9ad9de0SHong Zhang Calling sequence of drduf: 97*c9ad9de0SHong Zhang $ PetscErroCode drduf(TS ts,PetscReal t,Vec U,Vec *dRdU,void *ctx); 98814e85d6SHong Zhang 99814e85d6SHong Zhang Calling sequence of drdpf: 100*c9ad9de0SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec U,Vec *dRdP,void *ctx); 101814e85d6SHong Zhang 102814e85d6SHong Zhang Level: intermediate 103814e85d6SHong Zhang 104814e85d6SHong Zhang Notes: 105814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 106814e85d6SHong Zhang 107814e85d6SHong Zhang .keywords: TS, sensitivity analysis, timestep, set, quadrature, function 108814e85d6SHong Zhang 109814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients() 110814e85d6SHong Zhang @*/ 111814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*), 112*c9ad9de0SHong Zhang PetscErrorCode (*drduf)(TS,PetscReal,Vec,Vec*,void*), 113814e85d6SHong Zhang PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*), 114814e85d6SHong Zhang PetscBool fwd,void *ctx) 115814e85d6SHong Zhang { 116814e85d6SHong Zhang PetscErrorCode ierr; 117814e85d6SHong Zhang 118814e85d6SHong Zhang PetscFunctionBegin; 119814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 120814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3); 121814e85d6SHong 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()"); 122814e85d6SHong Zhang if (!ts->numcost) ts->numcost=numcost; 123814e85d6SHong Zhang 124814e85d6SHong Zhang if (costintegral) { 125814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr); 126814e85d6SHong Zhang ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr); 127814e85d6SHong Zhang ts->vec_costintegral = costintegral; 128814e85d6SHong Zhang } else { 129814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 130814e85d6SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr); 131814e85d6SHong Zhang } else { 132814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr); 133814e85d6SHong Zhang } 134814e85d6SHong Zhang } 135814e85d6SHong Zhang if (!ts->vec_costintegrand) { 136814e85d6SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr); 137814e85d6SHong Zhang } else { 138814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr); 139814e85d6SHong Zhang } 140814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 141814e85d6SHong Zhang ts->costintegrand = rf; 142814e85d6SHong Zhang ts->costintegrandctx = ctx; 143*c9ad9de0SHong Zhang ts->drdufunction = drduf; 144814e85d6SHong Zhang ts->drdpfunction = drdpf; 145814e85d6SHong Zhang PetscFunctionReturn(0); 146814e85d6SHong Zhang } 147814e85d6SHong Zhang 148814e85d6SHong Zhang /*@ 149814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 150814e85d6SHong Zhang It is valid to call the routine after a backward run. 151814e85d6SHong Zhang 152814e85d6SHong Zhang Not Collective 153814e85d6SHong Zhang 154814e85d6SHong Zhang Input Parameter: 155814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 156814e85d6SHong Zhang 157814e85d6SHong Zhang Output Parameter: 158814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 159814e85d6SHong Zhang 160814e85d6SHong Zhang Level: intermediate 161814e85d6SHong Zhang 162814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 163814e85d6SHong Zhang 164814e85d6SHong Zhang .keywords: TS, sensitivity analysis 165814e85d6SHong Zhang @*/ 166814e85d6SHong Zhang PetscErrorCode TSGetCostIntegral(TS ts,Vec *v) 167814e85d6SHong Zhang { 168814e85d6SHong Zhang PetscFunctionBegin; 169814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 170814e85d6SHong Zhang PetscValidPointer(v,2); 171814e85d6SHong Zhang *v = ts->vec_costintegral; 172814e85d6SHong Zhang PetscFunctionReturn(0); 173814e85d6SHong Zhang } 174814e85d6SHong Zhang 175814e85d6SHong Zhang /*@ 176814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 177814e85d6SHong Zhang 178814e85d6SHong Zhang Input Parameters: 179814e85d6SHong Zhang + ts - the TS context 180814e85d6SHong Zhang . t - current time 181*c9ad9de0SHong Zhang - U - state vector, i.e. current solution 182814e85d6SHong Zhang 183814e85d6SHong Zhang Output Parameter: 184*c9ad9de0SHong Zhang . Q - vector of size numcost to hold the outputs 185814e85d6SHong Zhang 186814e85d6SHong Zhang Note: 187814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 188814e85d6SHong Zhang is used internally within the sensitivity analysis context. 189814e85d6SHong Zhang 190814e85d6SHong Zhang Level: developer 191814e85d6SHong Zhang 192814e85d6SHong Zhang .keywords: TS, compute 193814e85d6SHong Zhang 194814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 195814e85d6SHong Zhang @*/ 196*c9ad9de0SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec U,Vec Q) 197814e85d6SHong Zhang { 198814e85d6SHong Zhang PetscErrorCode ierr; 199814e85d6SHong Zhang 200814e85d6SHong Zhang PetscFunctionBegin; 201814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 202*c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 203*c9ad9de0SHong Zhang PetscValidHeaderSpecific(Q,VEC_CLASSID,4); 204814e85d6SHong Zhang 205*c9ad9de0SHong Zhang ierr = PetscLogEventBegin(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr); 206814e85d6SHong Zhang if (ts->costintegrand) { 207814e85d6SHong Zhang PetscStackPush("TS user integrand in the cost function"); 208*c9ad9de0SHong Zhang ierr = (*ts->costintegrand)(ts,t,U,Q,ts->costintegrandctx);CHKERRQ(ierr); 209814e85d6SHong Zhang PetscStackPop; 210814e85d6SHong Zhang } else { 211*c9ad9de0SHong Zhang ierr = VecZeroEntries(Q);CHKERRQ(ierr); 212814e85d6SHong Zhang } 213814e85d6SHong Zhang 214*c9ad9de0SHong Zhang ierr = PetscLogEventEnd(TS_FunctionEval,ts,U,Q,0);CHKERRQ(ierr); 215814e85d6SHong Zhang PetscFunctionReturn(0); 216814e85d6SHong Zhang } 217814e85d6SHong Zhang 218814e85d6SHong Zhang /*@ 219*c9ad9de0SHong Zhang TSComputeDRDUFunction - Runs the user-defined DRDU function. 220814e85d6SHong Zhang 221814e85d6SHong Zhang Collective on TS 222814e85d6SHong Zhang 223814e85d6SHong Zhang Input Parameters: 224*c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate() 225*c9ad9de0SHong Zhang . t - current time 226*c9ad9de0SHong Zhang - U - stata vector 227*c9ad9de0SHong Zhang 228*c9ad9de0SHong Zhang Output Parameters: 229*c9ad9de0SHong Zhang . DRDU - vecotr array to hold the outputs 230814e85d6SHong Zhang 231814e85d6SHong Zhang Notes: 232*c9ad9de0SHong Zhang TSComputeDRDUFunction() is typically used for sensitivity implementation, 233814e85d6SHong Zhang so most users would not generally call this routine themselves. 234814e85d6SHong Zhang 235814e85d6SHong Zhang Level: developer 236814e85d6SHong Zhang 237814e85d6SHong Zhang .keywords: TS, sensitivity 238814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 239814e85d6SHong Zhang @*/ 240*c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDUFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 241814e85d6SHong Zhang { 242814e85d6SHong Zhang PetscErrorCode ierr; 243814e85d6SHong Zhang 244814e85d6SHong Zhang PetscFunctionBegin; 245814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 246*c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 247814e85d6SHong Zhang 248*c9ad9de0SHong Zhang PetscStackPush("TS user DRDU function for sensitivity analysis"); 249*c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 250814e85d6SHong Zhang PetscStackPop; 251814e85d6SHong Zhang PetscFunctionReturn(0); 252814e85d6SHong Zhang } 253814e85d6SHong Zhang 254814e85d6SHong Zhang /*@ 255814e85d6SHong Zhang TSComputeDRDPFunction - Runs the user-defined DRDP function. 256814e85d6SHong Zhang 257814e85d6SHong Zhang Collective on TS 258814e85d6SHong Zhang 259814e85d6SHong Zhang Input Parameters: 260*c9ad9de0SHong Zhang + ts - the TS context obtained from TSCreate() 261*c9ad9de0SHong Zhang . t - current time 262*c9ad9de0SHong Zhang - U - stata vector 263*c9ad9de0SHong Zhang 264*c9ad9de0SHong Zhang Output Parameters: 265*c9ad9de0SHong Zhang . DRDP - vecotr array to hold the outputs 266814e85d6SHong Zhang 267814e85d6SHong Zhang Notes: 268*c9ad9de0SHong Zhang TSComputeDRDPFunction() is typically used for sensitivity implementation, 269814e85d6SHong Zhang so most users would not generally call this routine themselves. 270814e85d6SHong Zhang 271814e85d6SHong Zhang Level: developer 272814e85d6SHong Zhang 273814e85d6SHong Zhang .keywords: TS, sensitivity 274814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 275814e85d6SHong Zhang @*/ 276*c9ad9de0SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 277814e85d6SHong Zhang { 278814e85d6SHong Zhang PetscErrorCode ierr; 279814e85d6SHong Zhang 280814e85d6SHong Zhang PetscFunctionBegin; 281814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 282*c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 283814e85d6SHong Zhang 284814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 285*c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 286814e85d6SHong Zhang PetscStackPop; 287814e85d6SHong Zhang PetscFunctionReturn(0); 288814e85d6SHong Zhang } 289814e85d6SHong Zhang 290814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 291814e85d6SHong Zhang 292814e85d6SHong Zhang /*@ 293814e85d6SHong 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 294814e85d6SHong Zhang for use by the TSAdjoint routines. 295814e85d6SHong Zhang 296814e85d6SHong Zhang Logically Collective on TS and Vec 297814e85d6SHong Zhang 298814e85d6SHong Zhang Input Parameters: 299814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 300814e85d6SHong 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 301814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 302814e85d6SHong Zhang 303814e85d6SHong Zhang Level: beginner 304814e85d6SHong Zhang 305814e85d6SHong Zhang Notes: 306814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 307814e85d6SHong Zhang 308814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 309814e85d6SHong Zhang 310814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values 311814e85d6SHong Zhang @*/ 312814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu) 313814e85d6SHong Zhang { 314814e85d6SHong Zhang PetscFunctionBegin; 315814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 316814e85d6SHong Zhang PetscValidPointer(lambda,2); 317814e85d6SHong Zhang ts->vecs_sensi = lambda; 318814e85d6SHong Zhang ts->vecs_sensip = mu; 319814e85d6SHong 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"); 320814e85d6SHong Zhang ts->numcost = numcost; 321814e85d6SHong Zhang PetscFunctionReturn(0); 322814e85d6SHong Zhang } 323814e85d6SHong Zhang 324814e85d6SHong Zhang /*@ 325814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 326814e85d6SHong Zhang 327814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 328814e85d6SHong Zhang 329814e85d6SHong Zhang Input Parameter: 330814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 331814e85d6SHong Zhang 332814e85d6SHong Zhang Output Parameter: 333814e85d6SHong Zhang + lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 334814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 335814e85d6SHong Zhang 336814e85d6SHong Zhang Level: intermediate 337814e85d6SHong Zhang 338814e85d6SHong Zhang .seealso: TSGetTimeStep() 339814e85d6SHong Zhang 340814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity 341814e85d6SHong Zhang @*/ 342814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu) 343814e85d6SHong Zhang { 344814e85d6SHong Zhang PetscFunctionBegin; 345814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 346814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 347814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 348814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 349814e85d6SHong Zhang PetscFunctionReturn(0); 350814e85d6SHong Zhang } 351814e85d6SHong Zhang 352814e85d6SHong Zhang /*@ 353814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 354814e85d6SHong Zhang of an adjoint solver 355814e85d6SHong Zhang 356814e85d6SHong Zhang Collective on TS 357814e85d6SHong Zhang 358814e85d6SHong Zhang Input Parameter: 359814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 360814e85d6SHong Zhang 361814e85d6SHong Zhang Level: advanced 362814e85d6SHong Zhang 363814e85d6SHong Zhang .keywords: TS, timestep, setup 364814e85d6SHong Zhang 365814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients() 366814e85d6SHong Zhang @*/ 367814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts) 368814e85d6SHong Zhang { 369814e85d6SHong Zhang PetscErrorCode ierr; 370814e85d6SHong Zhang 371814e85d6SHong Zhang PetscFunctionBegin; 372814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 373814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 374814e85d6SHong Zhang if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 375814e85d6SHong Zhang if (ts->vecs_sensip && !ts->Jacp) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSAdjointSetRHSJacobian() first"); 376814e85d6SHong Zhang 377814e85d6SHong Zhang if (ts->vec_costintegral) { /* if there is integral in the cost function */ 378*c9ad9de0SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr); 379814e85d6SHong Zhang if (ts->vecs_sensip){ 380814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr); 381814e85d6SHong Zhang } 382814e85d6SHong Zhang } 383814e85d6SHong Zhang 384814e85d6SHong Zhang if (ts->ops->adjointsetup) { 385814e85d6SHong Zhang ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr); 386814e85d6SHong Zhang } 387814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 388814e85d6SHong Zhang PetscFunctionReturn(0); 389814e85d6SHong Zhang } 390814e85d6SHong Zhang 391814e85d6SHong Zhang /*@ 392814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 393814e85d6SHong Zhang 394814e85d6SHong Zhang Logically Collective on TS 395814e85d6SHong Zhang 396814e85d6SHong Zhang Input Parameters: 397814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 398814e85d6SHong Zhang . steps - number of steps to use 399814e85d6SHong Zhang 400814e85d6SHong Zhang Level: intermediate 401814e85d6SHong Zhang 402814e85d6SHong Zhang Notes: 403814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 404814e85d6SHong Zhang so as to integrate back to less than the original timestep 405814e85d6SHong Zhang 406814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations 407814e85d6SHong Zhang 408814e85d6SHong Zhang .seealso: TSSetExactFinalTime() 409814e85d6SHong Zhang @*/ 410814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 411814e85d6SHong Zhang { 412814e85d6SHong Zhang PetscFunctionBegin; 413814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 414814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts,steps,2); 415814e85d6SHong Zhang if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 416814e85d6SHong Zhang if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 417814e85d6SHong Zhang ts->adjoint_max_steps = steps; 418814e85d6SHong Zhang PetscFunctionReturn(0); 419814e85d6SHong Zhang } 420814e85d6SHong Zhang 421814e85d6SHong Zhang /*@C 422814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 423814e85d6SHong Zhang 424814e85d6SHong Zhang Level: deprecated 425814e85d6SHong Zhang 426814e85d6SHong Zhang @*/ 427814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 428814e85d6SHong Zhang { 429814e85d6SHong Zhang PetscErrorCode ierr; 430814e85d6SHong Zhang 431814e85d6SHong Zhang PetscFunctionBegin; 432814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 433814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 434814e85d6SHong Zhang 435814e85d6SHong Zhang ts->rhsjacobianp = func; 436814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 437814e85d6SHong Zhang if(Amat) { 438814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 439814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 440814e85d6SHong Zhang ts->Jacp = Amat; 441814e85d6SHong Zhang } 442814e85d6SHong Zhang PetscFunctionReturn(0); 443814e85d6SHong Zhang } 444814e85d6SHong Zhang 445814e85d6SHong Zhang /*@C 446814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 447814e85d6SHong Zhang 448814e85d6SHong Zhang Level: deprecated 449814e85d6SHong Zhang 450814e85d6SHong Zhang @*/ 451*c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec U,Mat Amat) 452814e85d6SHong Zhang { 453814e85d6SHong Zhang PetscErrorCode ierr; 454814e85d6SHong Zhang 455814e85d6SHong Zhang PetscFunctionBegin; 456814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 457*c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 458814e85d6SHong Zhang PetscValidPointer(Amat,4); 459814e85d6SHong Zhang 460814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 461*c9ad9de0SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,U,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 462814e85d6SHong Zhang PetscStackPop; 463814e85d6SHong Zhang PetscFunctionReturn(0); 464814e85d6SHong Zhang } 465814e85d6SHong Zhang 466814e85d6SHong Zhang /*@ 467*c9ad9de0SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDUFunction() 468814e85d6SHong Zhang 469814e85d6SHong Zhang Level: deprecated 470814e85d6SHong Zhang 471814e85d6SHong Zhang @*/ 472*c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec U,Vec *DRDU) 473814e85d6SHong Zhang { 474814e85d6SHong Zhang PetscErrorCode ierr; 475814e85d6SHong Zhang 476814e85d6SHong Zhang PetscFunctionBegin; 477814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 478*c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 479814e85d6SHong Zhang 480814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 481*c9ad9de0SHong Zhang ierr = (*ts->drdufunction)(ts,t,U,DRDU,ts->costintegrandctx); CHKERRQ(ierr); 482814e85d6SHong Zhang PetscStackPop; 483814e85d6SHong Zhang PetscFunctionReturn(0); 484814e85d6SHong Zhang } 485814e85d6SHong Zhang 486814e85d6SHong Zhang /*@ 487814e85d6SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction() 488814e85d6SHong Zhang 489814e85d6SHong Zhang Level: deprecated 490814e85d6SHong Zhang 491814e85d6SHong Zhang @*/ 492*c9ad9de0SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec U,Vec *DRDP) 493814e85d6SHong Zhang { 494814e85d6SHong Zhang PetscErrorCode ierr; 495814e85d6SHong Zhang 496814e85d6SHong Zhang PetscFunctionBegin; 497814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 498*c9ad9de0SHong Zhang PetscValidHeaderSpecific(U,VEC_CLASSID,3); 499814e85d6SHong Zhang 500814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 501*c9ad9de0SHong Zhang ierr = (*ts->drdpfunction)(ts,t,U,DRDP,ts->costintegrandctx); CHKERRQ(ierr); 502814e85d6SHong Zhang PetscStackPop; 503814e85d6SHong Zhang PetscFunctionReturn(0); 504814e85d6SHong Zhang } 505814e85d6SHong Zhang 506814e85d6SHong Zhang /*@C 507814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 508814e85d6SHong Zhang 509814e85d6SHong Zhang Level: intermediate 510814e85d6SHong Zhang 511814e85d6SHong Zhang .keywords: TS, set, monitor 512814e85d6SHong Zhang 513814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 514814e85d6SHong Zhang @*/ 515814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 516814e85d6SHong Zhang { 517814e85d6SHong Zhang PetscErrorCode ierr; 518814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 519814e85d6SHong Zhang 520814e85d6SHong Zhang PetscFunctionBegin; 521814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 522814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 523814e85d6SHong Zhang ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 524814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 525814e85d6SHong Zhang PetscFunctionReturn(0); 526814e85d6SHong Zhang } 527814e85d6SHong Zhang 528814e85d6SHong Zhang /*@C 529814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 530814e85d6SHong Zhang 531814e85d6SHong Zhang Collective on TS 532814e85d6SHong Zhang 533814e85d6SHong Zhang Input Parameters: 534814e85d6SHong Zhang + ts - TS object you wish to monitor 535814e85d6SHong Zhang . name - the monitor type one is seeking 536814e85d6SHong Zhang . help - message indicating what monitoring is done 537814e85d6SHong Zhang . manual - manual page for the monitor 538814e85d6SHong Zhang . monitor - the monitor function 539814e85d6SHong 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 540814e85d6SHong Zhang 541814e85d6SHong Zhang Level: developer 542814e85d6SHong Zhang 543814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 544814e85d6SHong Zhang PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 545814e85d6SHong Zhang PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 546814e85d6SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 547814e85d6SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 548814e85d6SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 549814e85d6SHong Zhang PetscOptionsFList(), PetscOptionsEList() 550814e85d6SHong Zhang @*/ 551814e85d6SHong 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*)) 552814e85d6SHong Zhang { 553814e85d6SHong Zhang PetscErrorCode ierr; 554814e85d6SHong Zhang PetscViewer viewer; 555814e85d6SHong Zhang PetscViewerFormat format; 556814e85d6SHong Zhang PetscBool flg; 557814e85d6SHong Zhang 558814e85d6SHong Zhang PetscFunctionBegin; 55916413a6aSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 560814e85d6SHong Zhang if (flg) { 561814e85d6SHong Zhang PetscViewerAndFormat *vf; 562814e85d6SHong Zhang ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 563814e85d6SHong Zhang ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 564814e85d6SHong Zhang if (monitorsetup) { 565814e85d6SHong Zhang ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 566814e85d6SHong Zhang } 567814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 568814e85d6SHong Zhang } 569814e85d6SHong Zhang PetscFunctionReturn(0); 570814e85d6SHong Zhang } 571814e85d6SHong Zhang 572814e85d6SHong Zhang /*@C 573814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 574814e85d6SHong Zhang timestep to display the iteration's progress. 575814e85d6SHong Zhang 576814e85d6SHong Zhang Logically Collective on TS 577814e85d6SHong Zhang 578814e85d6SHong Zhang Input Parameters: 579814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 580814e85d6SHong Zhang . adjointmonitor - monitoring routine 581814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 582814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 583814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 584814e85d6SHong Zhang (may be NULL) 585814e85d6SHong Zhang 586814e85d6SHong Zhang Calling sequence of monitor: 587814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 588814e85d6SHong Zhang 589814e85d6SHong Zhang + ts - the TS context 590814e85d6SHong 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 591814e85d6SHong Zhang been interpolated to) 592814e85d6SHong Zhang . time - current time 593814e85d6SHong Zhang . u - current iterate 594814e85d6SHong Zhang . numcost - number of cost functionos 595814e85d6SHong Zhang . lambda - sensitivities to initial conditions 596814e85d6SHong Zhang . mu - sensitivities to parameters 597814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 598814e85d6SHong Zhang 599814e85d6SHong Zhang Notes: 600814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 601814e85d6SHong Zhang already has been loaded. 602814e85d6SHong Zhang 603814e85d6SHong Zhang Fortran Notes: 604814e85d6SHong Zhang Only a single monitor function can be set for each TS object 605814e85d6SHong Zhang 606814e85d6SHong Zhang Level: intermediate 607814e85d6SHong Zhang 608814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 609814e85d6SHong Zhang 610814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel() 611814e85d6SHong Zhang @*/ 612814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 613814e85d6SHong Zhang { 614814e85d6SHong Zhang PetscErrorCode ierr; 615814e85d6SHong Zhang PetscInt i; 616814e85d6SHong Zhang PetscBool identical; 617814e85d6SHong Zhang 618814e85d6SHong Zhang PetscFunctionBegin; 619814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 620814e85d6SHong Zhang for (i=0; i<ts->numbermonitors;i++) { 621814e85d6SHong Zhang ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 622814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 623814e85d6SHong Zhang } 624814e85d6SHong Zhang if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 625814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 626814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 627814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 628814e85d6SHong Zhang PetscFunctionReturn(0); 629814e85d6SHong Zhang } 630814e85d6SHong Zhang 631814e85d6SHong Zhang /*@C 632814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 633814e85d6SHong Zhang 634814e85d6SHong Zhang Logically Collective on TS 635814e85d6SHong Zhang 636814e85d6SHong Zhang Input Parameters: 637814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 638814e85d6SHong Zhang 639814e85d6SHong Zhang Notes: 640814e85d6SHong Zhang There is no way to remove a single, specific monitor. 641814e85d6SHong Zhang 642814e85d6SHong Zhang Level: intermediate 643814e85d6SHong Zhang 644814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 645814e85d6SHong Zhang 646814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 647814e85d6SHong Zhang @*/ 648814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts) 649814e85d6SHong Zhang { 650814e85d6SHong Zhang PetscErrorCode ierr; 651814e85d6SHong Zhang PetscInt i; 652814e85d6SHong Zhang 653814e85d6SHong Zhang PetscFunctionBegin; 654814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 655814e85d6SHong Zhang for (i=0; i<ts->numberadjointmonitors; i++) { 656814e85d6SHong Zhang if (ts->adjointmonitordestroy[i]) { 657814e85d6SHong Zhang ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 658814e85d6SHong Zhang } 659814e85d6SHong Zhang } 660814e85d6SHong Zhang ts->numberadjointmonitors = 0; 661814e85d6SHong Zhang PetscFunctionReturn(0); 662814e85d6SHong Zhang } 663814e85d6SHong Zhang 664814e85d6SHong Zhang /*@C 665814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 666814e85d6SHong Zhang 667814e85d6SHong Zhang Level: intermediate 668814e85d6SHong Zhang 669814e85d6SHong Zhang .keywords: TS, set, monitor 670814e85d6SHong Zhang 671814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 672814e85d6SHong Zhang @*/ 673814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 674814e85d6SHong Zhang { 675814e85d6SHong Zhang PetscErrorCode ierr; 676814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 677814e85d6SHong Zhang 678814e85d6SHong Zhang PetscFunctionBegin; 679814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 680814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 681814e85d6SHong Zhang ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 682814e85d6SHong 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); 683814e85d6SHong Zhang ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 684814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 685814e85d6SHong Zhang PetscFunctionReturn(0); 686814e85d6SHong Zhang } 687814e85d6SHong Zhang 688814e85d6SHong Zhang /*@C 689814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 690814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 691814e85d6SHong Zhang 692814e85d6SHong Zhang Collective on TS 693814e85d6SHong Zhang 694814e85d6SHong Zhang Input Parameters: 695814e85d6SHong Zhang + ts - the TS context 696814e85d6SHong Zhang . step - current time-step 697814e85d6SHong Zhang . ptime - current time 698814e85d6SHong Zhang . u - current state 699814e85d6SHong Zhang . numcost - number of cost functions 700814e85d6SHong Zhang . lambda - sensitivities to initial conditions 701814e85d6SHong Zhang . mu - sensitivities to parameters 702814e85d6SHong Zhang - dummy - either a viewer or NULL 703814e85d6SHong Zhang 704814e85d6SHong Zhang Level: intermediate 705814e85d6SHong Zhang 706814e85d6SHong Zhang .keywords: TS, vector, adjoint, monitor, view 707814e85d6SHong Zhang 708814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 709814e85d6SHong Zhang @*/ 710814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 711814e85d6SHong Zhang { 712814e85d6SHong Zhang PetscErrorCode ierr; 713814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 714814e85d6SHong Zhang PetscDraw draw; 715814e85d6SHong Zhang PetscReal xl,yl,xr,yr,h; 716814e85d6SHong Zhang char time[32]; 717814e85d6SHong Zhang 718814e85d6SHong Zhang PetscFunctionBegin; 719814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 720814e85d6SHong Zhang 721814e85d6SHong Zhang ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 722814e85d6SHong Zhang ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 723814e85d6SHong Zhang ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 724814e85d6SHong Zhang ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 725814e85d6SHong Zhang h = yl + .95*(yr - yl); 726814e85d6SHong Zhang ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 727814e85d6SHong Zhang ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 728814e85d6SHong Zhang PetscFunctionReturn(0); 729814e85d6SHong Zhang } 730814e85d6SHong Zhang 731814e85d6SHong Zhang /* 732814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 733814e85d6SHong Zhang 734814e85d6SHong Zhang Collective on TSAdjoint 735814e85d6SHong Zhang 736814e85d6SHong Zhang Input Parameter: 737814e85d6SHong Zhang . ts - the TS context 738814e85d6SHong Zhang 739814e85d6SHong Zhang Options Database Keys: 740814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 741814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 742814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 743814e85d6SHong Zhang 744814e85d6SHong Zhang Level: developer 745814e85d6SHong Zhang 746814e85d6SHong Zhang Notes: 747814e85d6SHong Zhang This is not normally called directly by users 748814e85d6SHong Zhang 749814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database 750814e85d6SHong Zhang 751814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 752814e85d6SHong Zhang */ 753814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 754814e85d6SHong Zhang { 755814e85d6SHong Zhang PetscBool tflg,opt; 756814e85d6SHong Zhang PetscErrorCode ierr; 757814e85d6SHong Zhang 758814e85d6SHong Zhang PetscFunctionBegin; 759814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,2); 760814e85d6SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 761814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 762814e85d6SHong Zhang ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 763814e85d6SHong Zhang if (opt) { 764814e85d6SHong Zhang ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 765814e85d6SHong Zhang ts->adjoint_solve = tflg; 766814e85d6SHong Zhang } 767814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 768814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 769814e85d6SHong Zhang opt = PETSC_FALSE; 770814e85d6SHong Zhang ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 771814e85d6SHong Zhang if (opt) { 772814e85d6SHong Zhang TSMonitorDrawCtx ctx; 773814e85d6SHong Zhang PetscInt howoften = 1; 774814e85d6SHong Zhang 775814e85d6SHong Zhang ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 776814e85d6SHong Zhang ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 777814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 778814e85d6SHong Zhang } 779814e85d6SHong Zhang PetscFunctionReturn(0); 780814e85d6SHong Zhang } 781814e85d6SHong Zhang 782814e85d6SHong Zhang /*@ 783814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 784814e85d6SHong Zhang 785814e85d6SHong Zhang Collective on TS 786814e85d6SHong Zhang 787814e85d6SHong Zhang Input Parameter: 788814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 789814e85d6SHong Zhang 790814e85d6SHong Zhang Level: intermediate 791814e85d6SHong Zhang 792814e85d6SHong Zhang .keywords: TS, adjoint, step 793814e85d6SHong Zhang 794814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve() 795814e85d6SHong Zhang @*/ 796814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts) 797814e85d6SHong Zhang { 798814e85d6SHong Zhang DM dm; 799814e85d6SHong Zhang PetscErrorCode ierr; 800814e85d6SHong Zhang 801814e85d6SHong Zhang PetscFunctionBegin; 802814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 803814e85d6SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 804814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 805814e85d6SHong Zhang 806814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vec_sol,(PetscObject)ts,"-ts_view_solution");CHKERRQ(ierr); 807814e85d6SHong Zhang 808814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 809814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 810814e85d6SHong 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); 811814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 812814e85d6SHong Zhang ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 813814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 814814e85d6SHong Zhang ts->adjoint_steps++; ts->steps--; 815814e85d6SHong Zhang 816814e85d6SHong Zhang if (ts->reason < 0) { 817814e85d6SHong Zhang if (ts->errorifstepfailed) { 818814e85d6SHong Zhang if (ts->reason == TS_DIVERGED_NONLINEAR_SOLVE) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_snes_failures or make negative to attempt recovery",TSConvergedReasons[ts->reason]); 819814e85d6SHong Zhang else if (ts->reason == TS_DIVERGED_STEP_REJECTED) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s, increase -ts_max_reject or make negative to attempt recovery",TSConvergedReasons[ts->reason]); 820814e85d6SHong Zhang else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]); 821814e85d6SHong Zhang } 822814e85d6SHong Zhang } else if (!ts->reason) { 823814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 824814e85d6SHong Zhang } 825814e85d6SHong Zhang PetscFunctionReturn(0); 826814e85d6SHong Zhang } 827814e85d6SHong Zhang 828814e85d6SHong Zhang /*@ 829814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 830814e85d6SHong Zhang 831814e85d6SHong Zhang Collective on TS 832814e85d6SHong Zhang 833814e85d6SHong Zhang Input Parameter: 834814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 835814e85d6SHong Zhang 836814e85d6SHong Zhang Options Database: 837814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 838814e85d6SHong Zhang 839814e85d6SHong Zhang Level: intermediate 840814e85d6SHong Zhang 841814e85d6SHong Zhang Notes: 842814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 843814e85d6SHong Zhang 844814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 845814e85d6SHong Zhang 846814e85d6SHong Zhang .keywords: TS, timestep, solve 847814e85d6SHong Zhang 848814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 849814e85d6SHong Zhang @*/ 850814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts) 851814e85d6SHong Zhang { 852814e85d6SHong Zhang PetscErrorCode ierr; 853814e85d6SHong Zhang 854814e85d6SHong Zhang PetscFunctionBegin; 855814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 856814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 857814e85d6SHong Zhang 858814e85d6SHong Zhang /* reset time step and iteration counters */ 859814e85d6SHong Zhang ts->adjoint_steps = 0; 860814e85d6SHong Zhang ts->ksp_its = 0; 861814e85d6SHong Zhang ts->snes_its = 0; 862814e85d6SHong Zhang ts->num_snes_failures = 0; 863814e85d6SHong Zhang ts->reject = 0; 864814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 865814e85d6SHong Zhang 866814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 867814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 868814e85d6SHong Zhang 869814e85d6SHong Zhang while (!ts->reason) { 870814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 871814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 872814e85d6SHong Zhang ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 873814e85d6SHong Zhang ierr = TSAdjointStep(ts);CHKERRQ(ierr); 874814e85d6SHong Zhang if (ts->vec_costintegral && !ts->costintegralfwd) { 875814e85d6SHong Zhang ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 876814e85d6SHong Zhang } 877814e85d6SHong Zhang } 878814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 879814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 880814e85d6SHong Zhang ts->solvetime = ts->ptime; 881814e85d6SHong Zhang ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 882814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 883814e85d6SHong Zhang ts->adjoint_max_steps = 0; 884814e85d6SHong Zhang PetscFunctionReturn(0); 885814e85d6SHong Zhang } 886814e85d6SHong Zhang 887814e85d6SHong Zhang /*@C 888814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 889814e85d6SHong Zhang 890814e85d6SHong Zhang Collective on TS 891814e85d6SHong Zhang 892814e85d6SHong Zhang Input Parameters: 893814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 894814e85d6SHong Zhang . step - step number that has just completed 895814e85d6SHong Zhang . ptime - model time of the state 896814e85d6SHong Zhang . u - state at the current model time 897814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 898814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 899814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 900814e85d6SHong Zhang 901814e85d6SHong Zhang Notes: 902814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 903814e85d6SHong Zhang Users would almost never call this routine directly. 904814e85d6SHong Zhang 905814e85d6SHong Zhang Level: developer 906814e85d6SHong Zhang 907814e85d6SHong Zhang .keywords: TS, timestep 908814e85d6SHong Zhang @*/ 909814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 910814e85d6SHong Zhang { 911814e85d6SHong Zhang PetscErrorCode ierr; 912814e85d6SHong Zhang PetscInt i,n = ts->numberadjointmonitors; 913814e85d6SHong Zhang 914814e85d6SHong Zhang PetscFunctionBegin; 915814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 916814e85d6SHong Zhang PetscValidHeaderSpecific(u,VEC_CLASSID,4); 9178860a134SJunchao Zhang ierr = VecLockReadPush(u);CHKERRQ(ierr); 918814e85d6SHong Zhang for (i=0; i<n; i++) { 919814e85d6SHong Zhang ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 920814e85d6SHong Zhang } 9218860a134SJunchao Zhang ierr = VecLockReadPop(u);CHKERRQ(ierr); 922814e85d6SHong Zhang PetscFunctionReturn(0); 923814e85d6SHong Zhang } 924814e85d6SHong Zhang 925814e85d6SHong Zhang /*@ 926814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 927814e85d6SHong Zhang 928814e85d6SHong Zhang Collective on TS 929814e85d6SHong Zhang 930814e85d6SHong Zhang Input Arguments: 931814e85d6SHong Zhang . ts - time stepping context 932814e85d6SHong Zhang 933814e85d6SHong Zhang Level: advanced 934814e85d6SHong Zhang 935814e85d6SHong Zhang Notes: 936814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 937814e85d6SHong Zhang 938814e85d6SHong Zhang .seealso: TSAdjointSolve(), TSAdjointStep 939814e85d6SHong Zhang @*/ 940814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts) 941814e85d6SHong Zhang { 942814e85d6SHong Zhang PetscErrorCode ierr; 943814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 944814e85d6SHong 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); 945814e85d6SHong Zhang ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 946814e85d6SHong Zhang PetscFunctionReturn(0); 947814e85d6SHong Zhang } 948814e85d6SHong Zhang 949814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 950814e85d6SHong Zhang 951814e85d6SHong Zhang /*@ 952814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 953814e85d6SHong Zhang of forward sensitivity analysis 954814e85d6SHong Zhang 955814e85d6SHong Zhang Collective on TS 956814e85d6SHong Zhang 957814e85d6SHong Zhang Input Parameter: 958814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 959814e85d6SHong Zhang 960814e85d6SHong Zhang Level: advanced 961814e85d6SHong Zhang 962814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup 963814e85d6SHong Zhang 964814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp() 965814e85d6SHong Zhang @*/ 966814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts) 967814e85d6SHong Zhang { 968814e85d6SHong Zhang PetscErrorCode ierr; 969814e85d6SHong Zhang 970814e85d6SHong Zhang PetscFunctionBegin; 971814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 972814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 973814e85d6SHong Zhang if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()"); 974814e85d6SHong Zhang if (ts->vecs_integral_sensip) { 975*c9ad9de0SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,ts->numcost,&ts->vecs_drdu);CHKERRQ(ierr); 976814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr); 977814e85d6SHong Zhang } 978814e85d6SHong Zhang 979814e85d6SHong Zhang if (ts->ops->forwardsetup) { 980814e85d6SHong Zhang ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 981814e85d6SHong Zhang } 982814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 983814e85d6SHong Zhang PetscFunctionReturn(0); 984814e85d6SHong Zhang } 985814e85d6SHong Zhang 986814e85d6SHong Zhang /*@ 987814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 988814e85d6SHong Zhang 989814e85d6SHong Zhang Input Parameter: 990814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 991814e85d6SHong Zhang . numfwdint- number of integrals 992814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 993814e85d6SHong Zhang 994814e85d6SHong Zhang Level: intermediate 995814e85d6SHong Zhang 996814e85d6SHong Zhang .keywords: TS, forward sensitivity 997814e85d6SHong Zhang 998814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 999814e85d6SHong Zhang @*/ 1000814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 1001814e85d6SHong Zhang { 1002814e85d6SHong Zhang PetscFunctionBegin; 1003814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1004814e85d6SHong 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()"); 1005814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 1006814e85d6SHong Zhang 1007814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 1008814e85d6SHong Zhang PetscFunctionReturn(0); 1009814e85d6SHong Zhang } 1010814e85d6SHong Zhang 1011814e85d6SHong Zhang /*@ 1012814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1013814e85d6SHong Zhang 1014814e85d6SHong Zhang Input Parameter: 1015814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1016814e85d6SHong Zhang 1017814e85d6SHong Zhang Output Parameter: 1018814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1019814e85d6SHong Zhang 1020814e85d6SHong Zhang Level: intermediate 1021814e85d6SHong Zhang 1022814e85d6SHong Zhang .keywords: TS, forward sensitivity 1023814e85d6SHong Zhang 1024814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1025814e85d6SHong Zhang @*/ 1026814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1027814e85d6SHong Zhang { 1028814e85d6SHong Zhang PetscFunctionBegin; 1029814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1030814e85d6SHong Zhang PetscValidPointer(vp,3); 1031814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1032814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1033814e85d6SHong Zhang PetscFunctionReturn(0); 1034814e85d6SHong Zhang } 1035814e85d6SHong Zhang 1036814e85d6SHong Zhang /*@ 1037814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1038814e85d6SHong Zhang 1039814e85d6SHong Zhang Collective on TS 1040814e85d6SHong Zhang 1041814e85d6SHong Zhang Input Arguments: 1042814e85d6SHong Zhang . ts - time stepping context 1043814e85d6SHong Zhang 1044814e85d6SHong Zhang Level: advanced 1045814e85d6SHong Zhang 1046814e85d6SHong Zhang Notes: 1047814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1048814e85d6SHong Zhang 1049814e85d6SHong Zhang .keywords: TS, forward sensitivity 1050814e85d6SHong Zhang 1051814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1052814e85d6SHong Zhang @*/ 1053814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts) 1054814e85d6SHong Zhang { 1055814e85d6SHong Zhang PetscErrorCode ierr; 1056814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1057814e85d6SHong Zhang if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1058814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1059814e85d6SHong Zhang ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1060814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1061814e85d6SHong Zhang PetscFunctionReturn(0); 1062814e85d6SHong Zhang } 1063814e85d6SHong Zhang 1064814e85d6SHong Zhang /*@ 1065814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1066814e85d6SHong Zhang 1067814e85d6SHong Zhang Logically Collective on TS and Vec 1068814e85d6SHong Zhang 1069814e85d6SHong Zhang Input Parameters: 1070814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1071814e85d6SHong Zhang . nump - number of parameters 1072814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1073814e85d6SHong Zhang 1074814e85d6SHong Zhang Level: beginner 1075814e85d6SHong Zhang 1076814e85d6SHong Zhang Notes: 1077814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1078814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1079814e85d6SHong Zhang You must call this function before TSSolve(). 1080814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1081814e85d6SHong Zhang 1082814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values 1083814e85d6SHong Zhang 1084814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1085814e85d6SHong Zhang @*/ 1086814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1087814e85d6SHong Zhang { 1088814e85d6SHong Zhang PetscErrorCode ierr; 1089814e85d6SHong Zhang 1090814e85d6SHong Zhang PetscFunctionBegin; 1091814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1092814e85d6SHong Zhang PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1093814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1094814e85d6SHong Zhang ts->num_parameters = nump; 1095814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1096814e85d6SHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1097814e85d6SHong Zhang ts->mat_sensip = Smat; 1098814e85d6SHong Zhang PetscFunctionReturn(0); 1099814e85d6SHong Zhang } 1100814e85d6SHong Zhang 1101814e85d6SHong Zhang /*@ 1102814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1103814e85d6SHong Zhang 1104814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1105814e85d6SHong Zhang 1106814e85d6SHong Zhang Output Parameter: 1107814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1108814e85d6SHong Zhang . nump - number of parameters 1109814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1110814e85d6SHong Zhang 1111814e85d6SHong Zhang Level: intermediate 1112814e85d6SHong Zhang 1113814e85d6SHong Zhang .keywords: TS, forward sensitivity 1114814e85d6SHong Zhang 1115814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1116814e85d6SHong Zhang @*/ 1117814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1118814e85d6SHong Zhang { 1119814e85d6SHong Zhang PetscFunctionBegin; 1120814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1121814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1122814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1123814e85d6SHong Zhang PetscFunctionReturn(0); 1124814e85d6SHong Zhang } 1125814e85d6SHong Zhang 1126814e85d6SHong Zhang /*@ 1127814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1128814e85d6SHong Zhang 1129814e85d6SHong Zhang Collective on TS 1130814e85d6SHong Zhang 1131814e85d6SHong Zhang Input Arguments: 1132814e85d6SHong Zhang . ts - time stepping context 1133814e85d6SHong Zhang 1134814e85d6SHong Zhang Level: advanced 1135814e85d6SHong Zhang 1136814e85d6SHong Zhang Notes: 1137814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1138814e85d6SHong Zhang 1139814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral() 1140814e85d6SHong Zhang @*/ 1141814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts) 1142814e85d6SHong Zhang { 1143814e85d6SHong Zhang PetscErrorCode ierr; 1144814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1145814e85d6SHong 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); 1146814e85d6SHong Zhang ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1147814e85d6SHong Zhang PetscFunctionReturn(0); 1148814e85d6SHong Zhang } 1149