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 9814e85d6SHong Zhang TSSetRHSJacobianP - Sets the function that computes the Jacobian of G w.r.t. the parameters p where y_t = G(y,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 20814e85d6SHong Zhang . y - 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 @*/ 63814e85d6SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec X,Mat Amat) 64814e85d6SHong Zhang { 65814e85d6SHong Zhang PetscErrorCode ierr; 66814e85d6SHong Zhang 67814e85d6SHong Zhang PetscFunctionBegin; 68814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 69814e85d6SHong Zhang PetscValidHeaderSpecific(X,VEC_CLASSID,3); 70814e85d6SHong Zhang PetscValidPointer(Amat,4); 71814e85d6SHong Zhang 72814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 73814e85d6SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,X,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 88814e85d6SHong Zhang . drdyf - function that computes the gradients of the r's with respect to y 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: 94814e85d6SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec y,Vec f,void *ctx); 95814e85d6SHong Zhang 96814e85d6SHong Zhang Calling sequence of drdyf: 97814e85d6SHong Zhang $ PetscErroCode drdyf(TS ts,PetscReal t,Vec y,Vec *drdy,void *ctx); 98814e85d6SHong Zhang 99814e85d6SHong Zhang Calling sequence of drdpf: 100814e85d6SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec y,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*), 112814e85d6SHong Zhang PetscErrorCode (*drdyf)(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; 143814e85d6SHong Zhang ts->drdyfunction = drdyf; 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 181814e85d6SHong Zhang - y - state vector, i.e. current solution 182814e85d6SHong Zhang 183814e85d6SHong Zhang Output Parameter: 184814e85d6SHong 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 @*/ 196814e85d6SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec y,Vec q) 197814e85d6SHong Zhang { 198814e85d6SHong Zhang PetscErrorCode ierr; 199814e85d6SHong Zhang 200814e85d6SHong Zhang PetscFunctionBegin; 201814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 202814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 203814e85d6SHong Zhang PetscValidHeaderSpecific(q,VEC_CLASSID,4); 204814e85d6SHong Zhang 205814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_FunctionEval,ts,y,q,0);CHKERRQ(ierr); 206814e85d6SHong Zhang if (ts->costintegrand) { 207814e85d6SHong Zhang PetscStackPush("TS user integrand in the cost function"); 208814e85d6SHong Zhang ierr = (*ts->costintegrand)(ts,t,y,q,ts->costintegrandctx);CHKERRQ(ierr); 209814e85d6SHong Zhang PetscStackPop; 210814e85d6SHong Zhang } else { 211814e85d6SHong Zhang ierr = VecZeroEntries(q);CHKERRQ(ierr); 212814e85d6SHong Zhang } 213814e85d6SHong Zhang 214814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_FunctionEval,ts,y,q,0);CHKERRQ(ierr); 215814e85d6SHong Zhang PetscFunctionReturn(0); 216814e85d6SHong Zhang } 217814e85d6SHong Zhang 218814e85d6SHong Zhang /*@ 219814e85d6SHong Zhang TSComputeDRDYFunction - Runs the user-defined DRDY function. 220814e85d6SHong Zhang 221814e85d6SHong Zhang Collective on TS 222814e85d6SHong Zhang 223814e85d6SHong Zhang Input Parameters: 224814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 225814e85d6SHong Zhang 226814e85d6SHong Zhang Notes: 227814e85d6SHong Zhang TSAdjointComputeDRDYFunction() is typically used for sensitivity implementation, 228814e85d6SHong Zhang so most users would not generally call this routine themselves. 229814e85d6SHong Zhang 230814e85d6SHong Zhang Level: developer 231814e85d6SHong Zhang 232814e85d6SHong Zhang .keywords: TS, sensitivity 233814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 234814e85d6SHong Zhang @*/ 235814e85d6SHong Zhang PetscErrorCode TSComputeDRDYFunction(TS ts,PetscReal t,Vec y,Vec *drdy) 236814e85d6SHong Zhang { 237814e85d6SHong Zhang PetscErrorCode ierr; 238814e85d6SHong Zhang 239814e85d6SHong Zhang PetscFunctionBegin; 240814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 241814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 242814e85d6SHong Zhang 243814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 244814e85d6SHong Zhang ierr = (*ts->drdyfunction)(ts,t,y,drdy,ts->costintegrandctx); CHKERRQ(ierr); 245814e85d6SHong Zhang PetscStackPop; 246814e85d6SHong Zhang PetscFunctionReturn(0); 247814e85d6SHong Zhang } 248814e85d6SHong Zhang 249814e85d6SHong Zhang /*@ 250814e85d6SHong Zhang TSComputeDRDPFunction - Runs the user-defined DRDP function. 251814e85d6SHong Zhang 252814e85d6SHong Zhang Collective on TS 253814e85d6SHong Zhang 254814e85d6SHong Zhang Input Parameters: 255814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 256814e85d6SHong Zhang 257814e85d6SHong Zhang Notes: 258814e85d6SHong Zhang TSDRDPFunction() is typically used for sensitivity implementation, 259814e85d6SHong Zhang so most users would not generally call this routine themselves. 260814e85d6SHong Zhang 261814e85d6SHong Zhang Level: developer 262814e85d6SHong Zhang 263814e85d6SHong Zhang .keywords: TS, sensitivity 264814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 265814e85d6SHong Zhang @*/ 266814e85d6SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec y,Vec *drdp) 267814e85d6SHong Zhang { 268814e85d6SHong Zhang PetscErrorCode ierr; 269814e85d6SHong Zhang 270814e85d6SHong Zhang PetscFunctionBegin; 271814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 272814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 273814e85d6SHong Zhang 274814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 275814e85d6SHong Zhang ierr = (*ts->drdpfunction)(ts,t,y,drdp,ts->costintegrandctx); CHKERRQ(ierr); 276814e85d6SHong Zhang PetscStackPop; 277814e85d6SHong Zhang PetscFunctionReturn(0); 278814e85d6SHong Zhang } 279814e85d6SHong Zhang 280814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 281814e85d6SHong Zhang 282814e85d6SHong Zhang /*@ 283814e85d6SHong 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 284814e85d6SHong Zhang for use by the TSAdjoint routines. 285814e85d6SHong Zhang 286814e85d6SHong Zhang Logically Collective on TS and Vec 287814e85d6SHong Zhang 288814e85d6SHong Zhang Input Parameters: 289814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 290814e85d6SHong 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 291814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 292814e85d6SHong Zhang 293814e85d6SHong Zhang Level: beginner 294814e85d6SHong Zhang 295814e85d6SHong Zhang Notes: 296814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 297814e85d6SHong Zhang 298814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 299814e85d6SHong Zhang 300814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values 301814e85d6SHong Zhang @*/ 302814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu) 303814e85d6SHong Zhang { 304814e85d6SHong Zhang PetscFunctionBegin; 305814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 306814e85d6SHong Zhang PetscValidPointer(lambda,2); 307814e85d6SHong Zhang ts->vecs_sensi = lambda; 308814e85d6SHong Zhang ts->vecs_sensip = mu; 309814e85d6SHong 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"); 310814e85d6SHong Zhang ts->numcost = numcost; 311814e85d6SHong Zhang PetscFunctionReturn(0); 312814e85d6SHong Zhang } 313814e85d6SHong Zhang 314814e85d6SHong Zhang /*@ 315814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 316814e85d6SHong Zhang 317814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 318814e85d6SHong Zhang 319814e85d6SHong Zhang Input Parameter: 320814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 321814e85d6SHong Zhang 322814e85d6SHong Zhang Output Parameter: 323814e85d6SHong Zhang + lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 324814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 325814e85d6SHong Zhang 326814e85d6SHong Zhang Level: intermediate 327814e85d6SHong Zhang 328814e85d6SHong Zhang .seealso: TSGetTimeStep() 329814e85d6SHong Zhang 330814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity 331814e85d6SHong Zhang @*/ 332814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu) 333814e85d6SHong Zhang { 334814e85d6SHong Zhang PetscFunctionBegin; 335814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 336814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 337814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 338814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 339814e85d6SHong Zhang PetscFunctionReturn(0); 340814e85d6SHong Zhang } 341814e85d6SHong Zhang 342814e85d6SHong Zhang /*@ 343814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 344814e85d6SHong Zhang of an adjoint solver 345814e85d6SHong Zhang 346814e85d6SHong Zhang Collective on TS 347814e85d6SHong Zhang 348814e85d6SHong Zhang Input Parameter: 349814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 350814e85d6SHong Zhang 351814e85d6SHong Zhang Level: advanced 352814e85d6SHong Zhang 353814e85d6SHong Zhang .keywords: TS, timestep, setup 354814e85d6SHong Zhang 355814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients() 356814e85d6SHong Zhang @*/ 357814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts) 358814e85d6SHong Zhang { 359814e85d6SHong Zhang PetscErrorCode ierr; 360814e85d6SHong Zhang 361814e85d6SHong Zhang PetscFunctionBegin; 362814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 363814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 364814e85d6SHong Zhang if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 365814e85d6SHong Zhang if (ts->vecs_sensip && !ts->Jacp) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSAdjointSetRHSJacobian() first"); 366814e85d6SHong Zhang 367814e85d6SHong Zhang if (ts->vec_costintegral) { /* if there is integral in the cost function */ 368814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&ts->vecs_drdy);CHKERRQ(ierr); 369814e85d6SHong Zhang if (ts->vecs_sensip){ 370814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr); 371814e85d6SHong Zhang } 372814e85d6SHong Zhang } 373814e85d6SHong Zhang 374814e85d6SHong Zhang if (ts->ops->adjointsetup) { 375814e85d6SHong Zhang ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr); 376814e85d6SHong Zhang } 377814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 378814e85d6SHong Zhang PetscFunctionReturn(0); 379814e85d6SHong Zhang } 380814e85d6SHong Zhang 381814e85d6SHong Zhang /*@ 382814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 383814e85d6SHong Zhang 384814e85d6SHong Zhang Logically Collective on TS 385814e85d6SHong Zhang 386814e85d6SHong Zhang Input Parameters: 387814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 388*a2b725a8SWilliam Gropp - steps - number of steps to use 389814e85d6SHong Zhang 390814e85d6SHong Zhang Level: intermediate 391814e85d6SHong Zhang 392814e85d6SHong Zhang Notes: 393814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 394814e85d6SHong Zhang so as to integrate back to less than the original timestep 395814e85d6SHong Zhang 396814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations 397814e85d6SHong Zhang 398814e85d6SHong Zhang .seealso: TSSetExactFinalTime() 399814e85d6SHong Zhang @*/ 400814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 401814e85d6SHong Zhang { 402814e85d6SHong Zhang PetscFunctionBegin; 403814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 404814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts,steps,2); 405814e85d6SHong Zhang if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 406814e85d6SHong Zhang if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 407814e85d6SHong Zhang ts->adjoint_max_steps = steps; 408814e85d6SHong Zhang PetscFunctionReturn(0); 409814e85d6SHong Zhang } 410814e85d6SHong Zhang 411814e85d6SHong Zhang /*@C 412814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 413814e85d6SHong Zhang 414814e85d6SHong Zhang Level: deprecated 415814e85d6SHong Zhang 416814e85d6SHong Zhang @*/ 417814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 418814e85d6SHong Zhang { 419814e85d6SHong Zhang PetscErrorCode ierr; 420814e85d6SHong Zhang 421814e85d6SHong Zhang PetscFunctionBegin; 422814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 423814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 424814e85d6SHong Zhang 425814e85d6SHong Zhang ts->rhsjacobianp = func; 426814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 427814e85d6SHong Zhang if(Amat) { 428814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 429814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 430814e85d6SHong Zhang ts->Jacp = Amat; 431814e85d6SHong Zhang } 432814e85d6SHong Zhang PetscFunctionReturn(0); 433814e85d6SHong Zhang } 434814e85d6SHong Zhang 435814e85d6SHong Zhang /*@C 436814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 437814e85d6SHong Zhang 438814e85d6SHong Zhang Level: deprecated 439814e85d6SHong Zhang 440814e85d6SHong Zhang @*/ 441814e85d6SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat) 442814e85d6SHong Zhang { 443814e85d6SHong Zhang PetscErrorCode ierr; 444814e85d6SHong Zhang 445814e85d6SHong Zhang PetscFunctionBegin; 446814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 447814e85d6SHong Zhang PetscValidHeaderSpecific(X,VEC_CLASSID,3); 448814e85d6SHong Zhang PetscValidPointer(Amat,4); 449814e85d6SHong Zhang 450814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 451814e85d6SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,X,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 452814e85d6SHong Zhang PetscStackPop; 453814e85d6SHong Zhang PetscFunctionReturn(0); 454814e85d6SHong Zhang } 455814e85d6SHong Zhang 456814e85d6SHong Zhang /*@ 457814e85d6SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDYFunction() 458814e85d6SHong Zhang 459814e85d6SHong Zhang Level: deprecated 460814e85d6SHong Zhang 461814e85d6SHong Zhang @*/ 462814e85d6SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec y,Vec *drdy) 463814e85d6SHong Zhang { 464814e85d6SHong Zhang PetscErrorCode ierr; 465814e85d6SHong Zhang 466814e85d6SHong Zhang PetscFunctionBegin; 467814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 468814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 469814e85d6SHong Zhang 470814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 471814e85d6SHong Zhang ierr = (*ts->drdyfunction)(ts,t,y,drdy,ts->costintegrandctx); CHKERRQ(ierr); 472814e85d6SHong Zhang PetscStackPop; 473814e85d6SHong Zhang PetscFunctionReturn(0); 474814e85d6SHong Zhang } 475814e85d6SHong Zhang 476814e85d6SHong Zhang /*@ 477814e85d6SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction() 478814e85d6SHong Zhang 479814e85d6SHong Zhang Level: deprecated 480814e85d6SHong Zhang 481814e85d6SHong Zhang @*/ 482814e85d6SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec y,Vec *drdp) 483814e85d6SHong Zhang { 484814e85d6SHong Zhang PetscErrorCode ierr; 485814e85d6SHong Zhang 486814e85d6SHong Zhang PetscFunctionBegin; 487814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 488814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 489814e85d6SHong Zhang 490814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 491814e85d6SHong Zhang ierr = (*ts->drdpfunction)(ts,t,y,drdp,ts->costintegrandctx); CHKERRQ(ierr); 492814e85d6SHong Zhang PetscStackPop; 493814e85d6SHong Zhang PetscFunctionReturn(0); 494814e85d6SHong Zhang } 495814e85d6SHong Zhang 496814e85d6SHong Zhang /*@C 497814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 498814e85d6SHong Zhang 499814e85d6SHong Zhang Level: intermediate 500814e85d6SHong Zhang 501814e85d6SHong Zhang .keywords: TS, set, monitor 502814e85d6SHong Zhang 503814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 504814e85d6SHong Zhang @*/ 505814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 506814e85d6SHong Zhang { 507814e85d6SHong Zhang PetscErrorCode ierr; 508814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 509814e85d6SHong Zhang 510814e85d6SHong Zhang PetscFunctionBegin; 511814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 512814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 513814e85d6SHong Zhang ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 514814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 515814e85d6SHong Zhang PetscFunctionReturn(0); 516814e85d6SHong Zhang } 517814e85d6SHong Zhang 518814e85d6SHong Zhang /*@C 519814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 520814e85d6SHong Zhang 521814e85d6SHong Zhang Collective on TS 522814e85d6SHong Zhang 523814e85d6SHong Zhang Input Parameters: 524814e85d6SHong Zhang + ts - TS object you wish to monitor 525814e85d6SHong Zhang . name - the monitor type one is seeking 526814e85d6SHong Zhang . help - message indicating what monitoring is done 527814e85d6SHong Zhang . manual - manual page for the monitor 528814e85d6SHong Zhang . monitor - the monitor function 529814e85d6SHong 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 530814e85d6SHong Zhang 531814e85d6SHong Zhang Level: developer 532814e85d6SHong Zhang 533814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 534814e85d6SHong Zhang PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 535814e85d6SHong Zhang PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 536814e85d6SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 537814e85d6SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 538814e85d6SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 539814e85d6SHong Zhang PetscOptionsFList(), PetscOptionsEList() 540814e85d6SHong Zhang @*/ 541814e85d6SHong 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*)) 542814e85d6SHong Zhang { 543814e85d6SHong Zhang PetscErrorCode ierr; 544814e85d6SHong Zhang PetscViewer viewer; 545814e85d6SHong Zhang PetscViewerFormat format; 546814e85d6SHong Zhang PetscBool flg; 547814e85d6SHong Zhang 548814e85d6SHong Zhang PetscFunctionBegin; 54916413a6aSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject) ts)->options,((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 550814e85d6SHong Zhang if (flg) { 551814e85d6SHong Zhang PetscViewerAndFormat *vf; 552814e85d6SHong Zhang ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 553814e85d6SHong Zhang ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 554814e85d6SHong Zhang if (monitorsetup) { 555814e85d6SHong Zhang ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 556814e85d6SHong Zhang } 557814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 558814e85d6SHong Zhang } 559814e85d6SHong Zhang PetscFunctionReturn(0); 560814e85d6SHong Zhang } 561814e85d6SHong Zhang 562814e85d6SHong Zhang /*@C 563814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 564814e85d6SHong Zhang timestep to display the iteration's progress. 565814e85d6SHong Zhang 566814e85d6SHong Zhang Logically Collective on TS 567814e85d6SHong Zhang 568814e85d6SHong Zhang Input Parameters: 569814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 570814e85d6SHong Zhang . adjointmonitor - monitoring routine 571814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 572814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 573814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 574814e85d6SHong Zhang (may be NULL) 575814e85d6SHong Zhang 576814e85d6SHong Zhang Calling sequence of monitor: 577814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 578814e85d6SHong Zhang 579814e85d6SHong Zhang + ts - the TS context 580814e85d6SHong 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 581814e85d6SHong Zhang been interpolated to) 582814e85d6SHong Zhang . time - current time 583814e85d6SHong Zhang . u - current iterate 584814e85d6SHong Zhang . numcost - number of cost functionos 585814e85d6SHong Zhang . lambda - sensitivities to initial conditions 586814e85d6SHong Zhang . mu - sensitivities to parameters 587814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 588814e85d6SHong Zhang 589814e85d6SHong Zhang Notes: 590814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 591814e85d6SHong Zhang already has been loaded. 592814e85d6SHong Zhang 593814e85d6SHong Zhang Fortran Notes: 594814e85d6SHong Zhang Only a single monitor function can be set for each TS object 595814e85d6SHong Zhang 596814e85d6SHong Zhang Level: intermediate 597814e85d6SHong Zhang 598814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 599814e85d6SHong Zhang 600814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel() 601814e85d6SHong Zhang @*/ 602814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 603814e85d6SHong Zhang { 604814e85d6SHong Zhang PetscErrorCode ierr; 605814e85d6SHong Zhang PetscInt i; 606814e85d6SHong Zhang PetscBool identical; 607814e85d6SHong Zhang 608814e85d6SHong Zhang PetscFunctionBegin; 609814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 610814e85d6SHong Zhang for (i=0; i<ts->numbermonitors;i++) { 611814e85d6SHong Zhang ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 612814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 613814e85d6SHong Zhang } 614814e85d6SHong Zhang if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 615814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 616814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 617814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 618814e85d6SHong Zhang PetscFunctionReturn(0); 619814e85d6SHong Zhang } 620814e85d6SHong Zhang 621814e85d6SHong Zhang /*@C 622814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 623814e85d6SHong Zhang 624814e85d6SHong Zhang Logically Collective on TS 625814e85d6SHong Zhang 626814e85d6SHong Zhang Input Parameters: 627814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 628814e85d6SHong Zhang 629814e85d6SHong Zhang Notes: 630814e85d6SHong Zhang There is no way to remove a single, specific monitor. 631814e85d6SHong Zhang 632814e85d6SHong Zhang Level: intermediate 633814e85d6SHong Zhang 634814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 635814e85d6SHong Zhang 636814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 637814e85d6SHong Zhang @*/ 638814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts) 639814e85d6SHong Zhang { 640814e85d6SHong Zhang PetscErrorCode ierr; 641814e85d6SHong Zhang PetscInt i; 642814e85d6SHong Zhang 643814e85d6SHong Zhang PetscFunctionBegin; 644814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 645814e85d6SHong Zhang for (i=0; i<ts->numberadjointmonitors; i++) { 646814e85d6SHong Zhang if (ts->adjointmonitordestroy[i]) { 647814e85d6SHong Zhang ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 648814e85d6SHong Zhang } 649814e85d6SHong Zhang } 650814e85d6SHong Zhang ts->numberadjointmonitors = 0; 651814e85d6SHong Zhang PetscFunctionReturn(0); 652814e85d6SHong Zhang } 653814e85d6SHong Zhang 654814e85d6SHong Zhang /*@C 655814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 656814e85d6SHong Zhang 657814e85d6SHong Zhang Level: intermediate 658814e85d6SHong Zhang 659814e85d6SHong Zhang .keywords: TS, set, monitor 660814e85d6SHong Zhang 661814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 662814e85d6SHong Zhang @*/ 663814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 664814e85d6SHong Zhang { 665814e85d6SHong Zhang PetscErrorCode ierr; 666814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 667814e85d6SHong Zhang 668814e85d6SHong Zhang PetscFunctionBegin; 669814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 670814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 671814e85d6SHong Zhang ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 672814e85d6SHong 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); 673814e85d6SHong Zhang ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 674814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 675814e85d6SHong Zhang PetscFunctionReturn(0); 676814e85d6SHong Zhang } 677814e85d6SHong Zhang 678814e85d6SHong Zhang /*@C 679814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 680814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 681814e85d6SHong Zhang 682814e85d6SHong Zhang Collective on TS 683814e85d6SHong Zhang 684814e85d6SHong Zhang Input Parameters: 685814e85d6SHong Zhang + ts - the TS context 686814e85d6SHong Zhang . step - current time-step 687814e85d6SHong Zhang . ptime - current time 688814e85d6SHong Zhang . u - current state 689814e85d6SHong Zhang . numcost - number of cost functions 690814e85d6SHong Zhang . lambda - sensitivities to initial conditions 691814e85d6SHong Zhang . mu - sensitivities to parameters 692814e85d6SHong Zhang - dummy - either a viewer or NULL 693814e85d6SHong Zhang 694814e85d6SHong Zhang Level: intermediate 695814e85d6SHong Zhang 696814e85d6SHong Zhang .keywords: TS, vector, adjoint, monitor, view 697814e85d6SHong Zhang 698814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 699814e85d6SHong Zhang @*/ 700814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 701814e85d6SHong Zhang { 702814e85d6SHong Zhang PetscErrorCode ierr; 703814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 704814e85d6SHong Zhang PetscDraw draw; 705814e85d6SHong Zhang PetscReal xl,yl,xr,yr,h; 706814e85d6SHong Zhang char time[32]; 707814e85d6SHong Zhang 708814e85d6SHong Zhang PetscFunctionBegin; 709814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 710814e85d6SHong Zhang 711814e85d6SHong Zhang ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 712814e85d6SHong Zhang ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 713814e85d6SHong Zhang ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 714814e85d6SHong Zhang ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 715814e85d6SHong Zhang h = yl + .95*(yr - yl); 716814e85d6SHong Zhang ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 717814e85d6SHong Zhang ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 718814e85d6SHong Zhang PetscFunctionReturn(0); 719814e85d6SHong Zhang } 720814e85d6SHong Zhang 721814e85d6SHong Zhang /* 722814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 723814e85d6SHong Zhang 724814e85d6SHong Zhang Collective on TSAdjoint 725814e85d6SHong Zhang 726814e85d6SHong Zhang Input Parameter: 727814e85d6SHong Zhang . ts - the TS context 728814e85d6SHong Zhang 729814e85d6SHong Zhang Options Database Keys: 730814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 731814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 732814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 733814e85d6SHong Zhang 734814e85d6SHong Zhang Level: developer 735814e85d6SHong Zhang 736814e85d6SHong Zhang Notes: 737814e85d6SHong Zhang This is not normally called directly by users 738814e85d6SHong Zhang 739814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database 740814e85d6SHong Zhang 741814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 742814e85d6SHong Zhang */ 743814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 744814e85d6SHong Zhang { 745814e85d6SHong Zhang PetscBool tflg,opt; 746814e85d6SHong Zhang PetscErrorCode ierr; 747814e85d6SHong Zhang 748814e85d6SHong Zhang PetscFunctionBegin; 749814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,2); 750814e85d6SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 751814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 752814e85d6SHong Zhang ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 753814e85d6SHong Zhang if (opt) { 754814e85d6SHong Zhang ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 755814e85d6SHong Zhang ts->adjoint_solve = tflg; 756814e85d6SHong Zhang } 757814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 758814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 759814e85d6SHong Zhang opt = PETSC_FALSE; 760814e85d6SHong Zhang ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 761814e85d6SHong Zhang if (opt) { 762814e85d6SHong Zhang TSMonitorDrawCtx ctx; 763814e85d6SHong Zhang PetscInt howoften = 1; 764814e85d6SHong Zhang 765814e85d6SHong Zhang ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 766814e85d6SHong Zhang ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 767814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 768814e85d6SHong Zhang } 769814e85d6SHong Zhang PetscFunctionReturn(0); 770814e85d6SHong Zhang } 771814e85d6SHong Zhang 772814e85d6SHong Zhang /*@ 773814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 774814e85d6SHong Zhang 775814e85d6SHong Zhang Collective on TS 776814e85d6SHong Zhang 777814e85d6SHong Zhang Input Parameter: 778814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 779814e85d6SHong Zhang 780814e85d6SHong Zhang Level: intermediate 781814e85d6SHong Zhang 782814e85d6SHong Zhang .keywords: TS, adjoint, step 783814e85d6SHong Zhang 784814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve() 785814e85d6SHong Zhang @*/ 786814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts) 787814e85d6SHong Zhang { 788814e85d6SHong Zhang DM dm; 789814e85d6SHong Zhang PetscErrorCode ierr; 790814e85d6SHong Zhang 791814e85d6SHong Zhang PetscFunctionBegin; 792814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 793814e85d6SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 794814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 795814e85d6SHong Zhang 796814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vec_sol,(PetscObject)ts,"-ts_view_solution");CHKERRQ(ierr); 797814e85d6SHong Zhang 798814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 799814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 800814e85d6SHong 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); 801814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 802814e85d6SHong Zhang ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 803814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 804814e85d6SHong Zhang ts->adjoint_steps++; ts->steps--; 805814e85d6SHong Zhang 806814e85d6SHong Zhang if (ts->reason < 0) { 807814e85d6SHong Zhang if (ts->errorifstepfailed) { 808814e85d6SHong 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]); 809814e85d6SHong 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]); 810814e85d6SHong Zhang else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]); 811814e85d6SHong Zhang } 812814e85d6SHong Zhang } else if (!ts->reason) { 813814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 814814e85d6SHong Zhang } 815814e85d6SHong Zhang PetscFunctionReturn(0); 816814e85d6SHong Zhang } 817814e85d6SHong Zhang 818814e85d6SHong Zhang /*@ 819814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 820814e85d6SHong Zhang 821814e85d6SHong Zhang Collective on TS 822814e85d6SHong Zhang 823814e85d6SHong Zhang Input Parameter: 824814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 825814e85d6SHong Zhang 826814e85d6SHong Zhang Options Database: 827814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 828814e85d6SHong Zhang 829814e85d6SHong Zhang Level: intermediate 830814e85d6SHong Zhang 831814e85d6SHong Zhang Notes: 832814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 833814e85d6SHong Zhang 834814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 835814e85d6SHong Zhang 836814e85d6SHong Zhang .keywords: TS, timestep, solve 837814e85d6SHong Zhang 838814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 839814e85d6SHong Zhang @*/ 840814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts) 841814e85d6SHong Zhang { 842814e85d6SHong Zhang PetscErrorCode ierr; 843814e85d6SHong Zhang 844814e85d6SHong Zhang PetscFunctionBegin; 845814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 846814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 847814e85d6SHong Zhang 848814e85d6SHong Zhang /* reset time step and iteration counters */ 849814e85d6SHong Zhang ts->adjoint_steps = 0; 850814e85d6SHong Zhang ts->ksp_its = 0; 851814e85d6SHong Zhang ts->snes_its = 0; 852814e85d6SHong Zhang ts->num_snes_failures = 0; 853814e85d6SHong Zhang ts->reject = 0; 854814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 855814e85d6SHong Zhang 856814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 857814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 858814e85d6SHong Zhang 859814e85d6SHong Zhang while (!ts->reason) { 860814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 861814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 862814e85d6SHong Zhang ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 863814e85d6SHong Zhang ierr = TSAdjointStep(ts);CHKERRQ(ierr); 864814e85d6SHong Zhang if (ts->vec_costintegral && !ts->costintegralfwd) { 865814e85d6SHong Zhang ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 866814e85d6SHong Zhang } 867814e85d6SHong Zhang } 868814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 869814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 870814e85d6SHong Zhang ts->solvetime = ts->ptime; 871814e85d6SHong Zhang ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 872814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 873814e85d6SHong Zhang ts->adjoint_max_steps = 0; 874814e85d6SHong Zhang PetscFunctionReturn(0); 875814e85d6SHong Zhang } 876814e85d6SHong Zhang 877814e85d6SHong Zhang /*@C 878814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 879814e85d6SHong Zhang 880814e85d6SHong Zhang Collective on TS 881814e85d6SHong Zhang 882814e85d6SHong Zhang Input Parameters: 883814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 884814e85d6SHong Zhang . step - step number that has just completed 885814e85d6SHong Zhang . ptime - model time of the state 886814e85d6SHong Zhang . u - state at the current model time 887814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 888814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 889814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 890814e85d6SHong Zhang 891814e85d6SHong Zhang Notes: 892814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 893814e85d6SHong Zhang Users would almost never call this routine directly. 894814e85d6SHong Zhang 895814e85d6SHong Zhang Level: developer 896814e85d6SHong Zhang 897814e85d6SHong Zhang .keywords: TS, timestep 898814e85d6SHong Zhang @*/ 899814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 900814e85d6SHong Zhang { 901814e85d6SHong Zhang PetscErrorCode ierr; 902814e85d6SHong Zhang PetscInt i,n = ts->numberadjointmonitors; 903814e85d6SHong Zhang 904814e85d6SHong Zhang PetscFunctionBegin; 905814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 906814e85d6SHong Zhang PetscValidHeaderSpecific(u,VEC_CLASSID,4); 9078860a134SJunchao Zhang ierr = VecLockReadPush(u);CHKERRQ(ierr); 908814e85d6SHong Zhang for (i=0; i<n; i++) { 909814e85d6SHong Zhang ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 910814e85d6SHong Zhang } 9118860a134SJunchao Zhang ierr = VecLockReadPop(u);CHKERRQ(ierr); 912814e85d6SHong Zhang PetscFunctionReturn(0); 913814e85d6SHong Zhang } 914814e85d6SHong Zhang 915814e85d6SHong Zhang /*@ 916814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 917814e85d6SHong Zhang 918814e85d6SHong Zhang Collective on TS 919814e85d6SHong Zhang 920814e85d6SHong Zhang Input Arguments: 921814e85d6SHong Zhang . ts - time stepping context 922814e85d6SHong Zhang 923814e85d6SHong Zhang Level: advanced 924814e85d6SHong Zhang 925814e85d6SHong Zhang Notes: 926814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 927814e85d6SHong Zhang 928814e85d6SHong Zhang .seealso: TSAdjointSolve(), TSAdjointStep 929814e85d6SHong Zhang @*/ 930814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts) 931814e85d6SHong Zhang { 932814e85d6SHong Zhang PetscErrorCode ierr; 933814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 934814e85d6SHong 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); 935814e85d6SHong Zhang ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 936814e85d6SHong Zhang PetscFunctionReturn(0); 937814e85d6SHong Zhang } 938814e85d6SHong Zhang 939814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 940814e85d6SHong Zhang 941814e85d6SHong Zhang /*@ 942814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 943814e85d6SHong Zhang of forward sensitivity analysis 944814e85d6SHong Zhang 945814e85d6SHong Zhang Collective on TS 946814e85d6SHong Zhang 947814e85d6SHong Zhang Input Parameter: 948814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 949814e85d6SHong Zhang 950814e85d6SHong Zhang Level: advanced 951814e85d6SHong Zhang 952814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup 953814e85d6SHong Zhang 954814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp() 955814e85d6SHong Zhang @*/ 956814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts) 957814e85d6SHong Zhang { 958814e85d6SHong Zhang PetscErrorCode ierr; 959814e85d6SHong Zhang 960814e85d6SHong Zhang PetscFunctionBegin; 961814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 962814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 963814e85d6SHong Zhang if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()"); 964814e85d6SHong Zhang if (ts->vecs_integral_sensip) { 965814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,ts->numcost,&ts->vecs_drdy);CHKERRQ(ierr); 966814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr); 967814e85d6SHong Zhang } 968814e85d6SHong Zhang 969814e85d6SHong Zhang if (ts->ops->forwardsetup) { 970814e85d6SHong Zhang ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 971814e85d6SHong Zhang } 972814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 973814e85d6SHong Zhang PetscFunctionReturn(0); 974814e85d6SHong Zhang } 975814e85d6SHong Zhang 976814e85d6SHong Zhang /*@ 977814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 978814e85d6SHong Zhang 979814e85d6SHong Zhang Input Parameter: 980*a2b725a8SWilliam Gropp + ts- the TS context obtained from TSCreate() 981814e85d6SHong Zhang . numfwdint- number of integrals 982*a2b725a8SWilliam Gropp - vp = the vectors containing the gradients for each integral w.r.t. parameters 983814e85d6SHong Zhang 984814e85d6SHong Zhang Level: intermediate 985814e85d6SHong Zhang 986814e85d6SHong Zhang .keywords: TS, forward sensitivity 987814e85d6SHong Zhang 988814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 989814e85d6SHong Zhang @*/ 990814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 991814e85d6SHong Zhang { 992814e85d6SHong Zhang PetscFunctionBegin; 993814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 994814e85d6SHong 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()"); 995814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 996814e85d6SHong Zhang 997814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 998814e85d6SHong Zhang PetscFunctionReturn(0); 999814e85d6SHong Zhang } 1000814e85d6SHong Zhang 1001814e85d6SHong Zhang /*@ 1002814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1003814e85d6SHong Zhang 1004814e85d6SHong Zhang Input Parameter: 1005814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1006814e85d6SHong Zhang 1007814e85d6SHong Zhang Output Parameter: 1008814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1009814e85d6SHong Zhang 1010814e85d6SHong Zhang Level: intermediate 1011814e85d6SHong Zhang 1012814e85d6SHong Zhang .keywords: TS, forward sensitivity 1013814e85d6SHong Zhang 1014814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1015814e85d6SHong Zhang @*/ 1016814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1017814e85d6SHong Zhang { 1018814e85d6SHong Zhang PetscFunctionBegin; 1019814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1020814e85d6SHong Zhang PetscValidPointer(vp,3); 1021814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1022814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1023814e85d6SHong Zhang PetscFunctionReturn(0); 1024814e85d6SHong Zhang } 1025814e85d6SHong Zhang 1026814e85d6SHong Zhang /*@ 1027814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1028814e85d6SHong Zhang 1029814e85d6SHong Zhang Collective on TS 1030814e85d6SHong Zhang 1031814e85d6SHong Zhang Input Arguments: 1032814e85d6SHong Zhang . ts - time stepping context 1033814e85d6SHong Zhang 1034814e85d6SHong Zhang Level: advanced 1035814e85d6SHong Zhang 1036814e85d6SHong Zhang Notes: 1037814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1038814e85d6SHong Zhang 1039814e85d6SHong Zhang .keywords: TS, forward sensitivity 1040814e85d6SHong Zhang 1041814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1042814e85d6SHong Zhang @*/ 1043814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts) 1044814e85d6SHong Zhang { 1045814e85d6SHong Zhang PetscErrorCode ierr; 1046814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1047814e85d6SHong Zhang if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1048814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1049814e85d6SHong Zhang ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1050814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1051814e85d6SHong Zhang PetscFunctionReturn(0); 1052814e85d6SHong Zhang } 1053814e85d6SHong Zhang 1054814e85d6SHong Zhang /*@ 1055814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1056814e85d6SHong Zhang 1057814e85d6SHong Zhang Logically Collective on TS and Vec 1058814e85d6SHong Zhang 1059814e85d6SHong Zhang Input Parameters: 1060814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1061814e85d6SHong Zhang . nump - number of parameters 1062814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1063814e85d6SHong Zhang 1064814e85d6SHong Zhang Level: beginner 1065814e85d6SHong Zhang 1066814e85d6SHong Zhang Notes: 1067814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1068814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1069814e85d6SHong Zhang You must call this function before TSSolve(). 1070814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1071814e85d6SHong Zhang 1072814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values 1073814e85d6SHong Zhang 1074814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1075814e85d6SHong Zhang @*/ 1076814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1077814e85d6SHong Zhang { 1078814e85d6SHong Zhang PetscErrorCode ierr; 1079814e85d6SHong Zhang 1080814e85d6SHong Zhang PetscFunctionBegin; 1081814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1082814e85d6SHong Zhang PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1083814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1084814e85d6SHong Zhang ts->num_parameters = nump; 1085814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1086814e85d6SHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1087814e85d6SHong Zhang ts->mat_sensip = Smat; 1088814e85d6SHong Zhang PetscFunctionReturn(0); 1089814e85d6SHong Zhang } 1090814e85d6SHong Zhang 1091814e85d6SHong Zhang /*@ 1092814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1093814e85d6SHong Zhang 1094814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1095814e85d6SHong Zhang 1096814e85d6SHong Zhang Output Parameter: 1097814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1098814e85d6SHong Zhang . nump - number of parameters 1099814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1100814e85d6SHong Zhang 1101814e85d6SHong Zhang Level: intermediate 1102814e85d6SHong Zhang 1103814e85d6SHong Zhang .keywords: TS, forward sensitivity 1104814e85d6SHong Zhang 1105814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1106814e85d6SHong Zhang @*/ 1107814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1108814e85d6SHong Zhang { 1109814e85d6SHong Zhang PetscFunctionBegin; 1110814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1111814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1112814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1113814e85d6SHong Zhang PetscFunctionReturn(0); 1114814e85d6SHong Zhang } 1115814e85d6SHong Zhang 1116814e85d6SHong Zhang /*@ 1117814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1118814e85d6SHong Zhang 1119814e85d6SHong Zhang Collective on TS 1120814e85d6SHong Zhang 1121814e85d6SHong Zhang Input Arguments: 1122814e85d6SHong Zhang . ts - time stepping context 1123814e85d6SHong Zhang 1124814e85d6SHong Zhang Level: advanced 1125814e85d6SHong Zhang 1126814e85d6SHong Zhang Notes: 1127814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1128814e85d6SHong Zhang 1129814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral() 1130814e85d6SHong Zhang @*/ 1131814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts) 1132814e85d6SHong Zhang { 1133814e85d6SHong Zhang PetscErrorCode ierr; 1134814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1135814e85d6SHong 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); 1136814e85d6SHong Zhang ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1137814e85d6SHong Zhang PetscFunctionReturn(0); 1138814e85d6SHong Zhang } 1139