1*814e85d6SHong Zhang #include <petsc/private/tsimpl.h> /*I "petscts.h" I*/ 2*814e85d6SHong Zhang #include <petscdraw.h> 3*814e85d6SHong Zhang 4*814e85d6SHong Zhang PetscLogEvent TS_AdjointStep, TS_ForwardStep; 5*814e85d6SHong Zhang 6*814e85d6SHong Zhang /* ------------------------ Sensitivity Context ---------------------------*/ 7*814e85d6SHong Zhang 8*814e85d6SHong Zhang /*@C 9*814e85d6SHong 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. 10*814e85d6SHong Zhang 11*814e85d6SHong Zhang Logically Collective on TS 12*814e85d6SHong Zhang 13*814e85d6SHong Zhang Input Parameters: 14*814e85d6SHong Zhang + ts - The TS context obtained from TSCreate() 15*814e85d6SHong Zhang - func - The function 16*814e85d6SHong Zhang 17*814e85d6SHong Zhang Calling sequence of func: 18*814e85d6SHong Zhang $ func (TS ts,PetscReal t,Vec y,Mat A,void *ctx); 19*814e85d6SHong Zhang + t - current timestep 20*814e85d6SHong Zhang . y - input vector (current ODE solution) 21*814e85d6SHong Zhang . A - output matrix 22*814e85d6SHong Zhang - ctx - [optional] user-defined function context 23*814e85d6SHong Zhang 24*814e85d6SHong Zhang Level: intermediate 25*814e85d6SHong Zhang 26*814e85d6SHong Zhang Notes: 27*814e85d6SHong 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 28*814e85d6SHong Zhang 29*814e85d6SHong Zhang .keywords: TS, sensitivity 30*814e85d6SHong Zhang .seealso: 31*814e85d6SHong Zhang @*/ 32*814e85d6SHong Zhang PetscErrorCode TSSetRHSJacobianP(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 33*814e85d6SHong Zhang { 34*814e85d6SHong Zhang PetscErrorCode ierr; 35*814e85d6SHong Zhang 36*814e85d6SHong Zhang PetscFunctionBegin; 37*814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 38*814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 39*814e85d6SHong Zhang 40*814e85d6SHong Zhang ts->rhsjacobianp = func; 41*814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 42*814e85d6SHong Zhang if(Amat) { 43*814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 44*814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 45*814e85d6SHong Zhang ts->Jacp = Amat; 46*814e85d6SHong Zhang } 47*814e85d6SHong Zhang PetscFunctionReturn(0); 48*814e85d6SHong Zhang } 49*814e85d6SHong Zhang 50*814e85d6SHong Zhang /*@C 51*814e85d6SHong Zhang TSComputeRHSJacobianP - Runs the user-defined JacobianP function. 52*814e85d6SHong Zhang 53*814e85d6SHong Zhang Collective on TS 54*814e85d6SHong Zhang 55*814e85d6SHong Zhang Input Parameters: 56*814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 57*814e85d6SHong Zhang 58*814e85d6SHong Zhang Level: developer 59*814e85d6SHong Zhang 60*814e85d6SHong Zhang .keywords: TS, sensitivity 61*814e85d6SHong Zhang .seealso: TSSetRHSJacobianP() 62*814e85d6SHong Zhang @*/ 63*814e85d6SHong Zhang PetscErrorCode TSComputeRHSJacobianP(TS ts,PetscReal t,Vec X,Mat Amat) 64*814e85d6SHong Zhang { 65*814e85d6SHong Zhang PetscErrorCode ierr; 66*814e85d6SHong Zhang 67*814e85d6SHong Zhang PetscFunctionBegin; 68*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 69*814e85d6SHong Zhang PetscValidHeaderSpecific(X,VEC_CLASSID,3); 70*814e85d6SHong Zhang PetscValidPointer(Amat,4); 71*814e85d6SHong Zhang 72*814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 73*814e85d6SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,X,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 74*814e85d6SHong Zhang PetscStackPop; 75*814e85d6SHong Zhang PetscFunctionReturn(0); 76*814e85d6SHong Zhang } 77*814e85d6SHong Zhang 78*814e85d6SHong Zhang /*@C 79*814e85d6SHong Zhang TSSetCostIntegrand - Sets the routine for evaluating the integral term in one or more cost functions 80*814e85d6SHong Zhang 81*814e85d6SHong Zhang Logically Collective on TS 82*814e85d6SHong Zhang 83*814e85d6SHong Zhang Input Parameters: 84*814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 85*814e85d6SHong Zhang . numcost - number of gradients to be computed, this is the number of cost functions 86*814e85d6SHong Zhang . costintegral - vector that stores the integral values 87*814e85d6SHong Zhang . rf - routine for evaluating the integrand function 88*814e85d6SHong Zhang . drdyf - function that computes the gradients of the r's with respect to y 89*814e85d6SHong 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) 90*814e85d6SHong Zhang . fwd - flag indicating whether to evaluate cost integral in the forward run or the adjoint run 91*814e85d6SHong Zhang - ctx - [optional] user-defined context for private data for the function evaluation routine (may be NULL) 92*814e85d6SHong Zhang 93*814e85d6SHong Zhang Calling sequence of rf: 94*814e85d6SHong Zhang $ PetscErrorCode rf(TS ts,PetscReal t,Vec y,Vec f,void *ctx); 95*814e85d6SHong Zhang 96*814e85d6SHong Zhang Calling sequence of drdyf: 97*814e85d6SHong Zhang $ PetscErroCode drdyf(TS ts,PetscReal t,Vec y,Vec *drdy,void *ctx); 98*814e85d6SHong Zhang 99*814e85d6SHong Zhang Calling sequence of drdpf: 100*814e85d6SHong Zhang $ PetscErroCode drdpf(TS ts,PetscReal t,Vec y,Vec *drdp,void *ctx); 101*814e85d6SHong Zhang 102*814e85d6SHong Zhang Level: intermediate 103*814e85d6SHong Zhang 104*814e85d6SHong Zhang Notes: 105*814e85d6SHong Zhang For optimization there is usually a single cost function (numcost = 1). For sensitivities there may be multiple cost functions 106*814e85d6SHong Zhang 107*814e85d6SHong Zhang .keywords: TS, sensitivity analysis, timestep, set, quadrature, function 108*814e85d6SHong Zhang 109*814e85d6SHong Zhang .seealso: TSSetRHSJacobianP(), TSGetCostGradients(), TSSetCostGradients() 110*814e85d6SHong Zhang @*/ 111*814e85d6SHong Zhang PetscErrorCode TSSetCostIntegrand(TS ts,PetscInt numcost,Vec costintegral,PetscErrorCode (*rf)(TS,PetscReal,Vec,Vec,void*), 112*814e85d6SHong Zhang PetscErrorCode (*drdyf)(TS,PetscReal,Vec,Vec*,void*), 113*814e85d6SHong Zhang PetscErrorCode (*drdpf)(TS,PetscReal,Vec,Vec*,void*), 114*814e85d6SHong Zhang PetscBool fwd,void *ctx) 115*814e85d6SHong Zhang { 116*814e85d6SHong Zhang PetscErrorCode ierr; 117*814e85d6SHong Zhang 118*814e85d6SHong Zhang PetscFunctionBegin; 119*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 120*814e85d6SHong Zhang if (costintegral) PetscValidHeaderSpecific(costintegral,VEC_CLASSID,3); 121*814e85d6SHong 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()"); 122*814e85d6SHong Zhang if (!ts->numcost) ts->numcost=numcost; 123*814e85d6SHong Zhang 124*814e85d6SHong Zhang if (costintegral) { 125*814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)costintegral);CHKERRQ(ierr); 126*814e85d6SHong Zhang ierr = VecDestroy(&ts->vec_costintegral);CHKERRQ(ierr); 127*814e85d6SHong Zhang ts->vec_costintegral = costintegral; 128*814e85d6SHong Zhang } else { 129*814e85d6SHong Zhang if (!ts->vec_costintegral) { /* Create a seq vec if user does not provide one */ 130*814e85d6SHong Zhang ierr = VecCreateSeq(PETSC_COMM_SELF,numcost,&ts->vec_costintegral);CHKERRQ(ierr); 131*814e85d6SHong Zhang } else { 132*814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegral,0.0);CHKERRQ(ierr); 133*814e85d6SHong Zhang } 134*814e85d6SHong Zhang } 135*814e85d6SHong Zhang if (!ts->vec_costintegrand) { 136*814e85d6SHong Zhang ierr = VecDuplicate(ts->vec_costintegral,&ts->vec_costintegrand);CHKERRQ(ierr); 137*814e85d6SHong Zhang } else { 138*814e85d6SHong Zhang ierr = VecSet(ts->vec_costintegrand,0.0);CHKERRQ(ierr); 139*814e85d6SHong Zhang } 140*814e85d6SHong Zhang ts->costintegralfwd = fwd; /* Evaluate the cost integral in forward run if fwd is true */ 141*814e85d6SHong Zhang ts->costintegrand = rf; 142*814e85d6SHong Zhang ts->costintegrandctx = ctx; 143*814e85d6SHong Zhang ts->drdyfunction = drdyf; 144*814e85d6SHong Zhang ts->drdpfunction = drdpf; 145*814e85d6SHong Zhang PetscFunctionReturn(0); 146*814e85d6SHong Zhang } 147*814e85d6SHong Zhang 148*814e85d6SHong Zhang /*@ 149*814e85d6SHong Zhang TSGetCostIntegral - Returns the values of the integral term in the cost functions. 150*814e85d6SHong Zhang It is valid to call the routine after a backward run. 151*814e85d6SHong Zhang 152*814e85d6SHong Zhang Not Collective 153*814e85d6SHong Zhang 154*814e85d6SHong Zhang Input Parameter: 155*814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 156*814e85d6SHong Zhang 157*814e85d6SHong Zhang Output Parameter: 158*814e85d6SHong Zhang . v - the vector containing the integrals for each cost function 159*814e85d6SHong Zhang 160*814e85d6SHong Zhang Level: intermediate 161*814e85d6SHong Zhang 162*814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 163*814e85d6SHong Zhang 164*814e85d6SHong Zhang .keywords: TS, sensitivity analysis 165*814e85d6SHong Zhang @*/ 166*814e85d6SHong Zhang PetscErrorCode TSGetCostIntegral(TS ts,Vec *v) 167*814e85d6SHong Zhang { 168*814e85d6SHong Zhang PetscFunctionBegin; 169*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 170*814e85d6SHong Zhang PetscValidPointer(v,2); 171*814e85d6SHong Zhang *v = ts->vec_costintegral; 172*814e85d6SHong Zhang PetscFunctionReturn(0); 173*814e85d6SHong Zhang } 174*814e85d6SHong Zhang 175*814e85d6SHong Zhang /*@ 176*814e85d6SHong Zhang TSComputeCostIntegrand - Evaluates the integral function in the cost functions. 177*814e85d6SHong Zhang 178*814e85d6SHong Zhang Input Parameters: 179*814e85d6SHong Zhang + ts - the TS context 180*814e85d6SHong Zhang . t - current time 181*814e85d6SHong Zhang - y - state vector, i.e. current solution 182*814e85d6SHong Zhang 183*814e85d6SHong Zhang Output Parameter: 184*814e85d6SHong Zhang . q - vector of size numcost to hold the outputs 185*814e85d6SHong Zhang 186*814e85d6SHong Zhang Note: 187*814e85d6SHong Zhang Most users should not need to explicitly call this routine, as it 188*814e85d6SHong Zhang is used internally within the sensitivity analysis context. 189*814e85d6SHong Zhang 190*814e85d6SHong Zhang Level: developer 191*814e85d6SHong Zhang 192*814e85d6SHong Zhang .keywords: TS, compute 193*814e85d6SHong Zhang 194*814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 195*814e85d6SHong Zhang @*/ 196*814e85d6SHong Zhang PetscErrorCode TSComputeCostIntegrand(TS ts,PetscReal t,Vec y,Vec q) 197*814e85d6SHong Zhang { 198*814e85d6SHong Zhang PetscErrorCode ierr; 199*814e85d6SHong Zhang 200*814e85d6SHong Zhang PetscFunctionBegin; 201*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 202*814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 203*814e85d6SHong Zhang PetscValidHeaderSpecific(q,VEC_CLASSID,4); 204*814e85d6SHong Zhang 205*814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_FunctionEval,ts,y,q,0);CHKERRQ(ierr); 206*814e85d6SHong Zhang if (ts->costintegrand) { 207*814e85d6SHong Zhang PetscStackPush("TS user integrand in the cost function"); 208*814e85d6SHong Zhang ierr = (*ts->costintegrand)(ts,t,y,q,ts->costintegrandctx);CHKERRQ(ierr); 209*814e85d6SHong Zhang PetscStackPop; 210*814e85d6SHong Zhang } else { 211*814e85d6SHong Zhang ierr = VecZeroEntries(q);CHKERRQ(ierr); 212*814e85d6SHong Zhang } 213*814e85d6SHong Zhang 214*814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_FunctionEval,ts,y,q,0);CHKERRQ(ierr); 215*814e85d6SHong Zhang PetscFunctionReturn(0); 216*814e85d6SHong Zhang } 217*814e85d6SHong Zhang 218*814e85d6SHong Zhang /*@ 219*814e85d6SHong Zhang TSComputeDRDYFunction - Runs the user-defined DRDY function. 220*814e85d6SHong Zhang 221*814e85d6SHong Zhang Collective on TS 222*814e85d6SHong Zhang 223*814e85d6SHong Zhang Input Parameters: 224*814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 225*814e85d6SHong Zhang 226*814e85d6SHong Zhang Notes: 227*814e85d6SHong Zhang TSAdjointComputeDRDYFunction() is typically used for sensitivity implementation, 228*814e85d6SHong Zhang so most users would not generally call this routine themselves. 229*814e85d6SHong Zhang 230*814e85d6SHong Zhang Level: developer 231*814e85d6SHong Zhang 232*814e85d6SHong Zhang .keywords: TS, sensitivity 233*814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 234*814e85d6SHong Zhang @*/ 235*814e85d6SHong Zhang PetscErrorCode TSComputeDRDYFunction(TS ts,PetscReal t,Vec y,Vec *drdy) 236*814e85d6SHong Zhang { 237*814e85d6SHong Zhang PetscErrorCode ierr; 238*814e85d6SHong Zhang 239*814e85d6SHong Zhang PetscFunctionBegin; 240*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 241*814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 242*814e85d6SHong Zhang 243*814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 244*814e85d6SHong Zhang ierr = (*ts->drdyfunction)(ts,t,y,drdy,ts->costintegrandctx); CHKERRQ(ierr); 245*814e85d6SHong Zhang PetscStackPop; 246*814e85d6SHong Zhang PetscFunctionReturn(0); 247*814e85d6SHong Zhang } 248*814e85d6SHong Zhang 249*814e85d6SHong Zhang /*@ 250*814e85d6SHong Zhang TSComputeDRDPFunction - Runs the user-defined DRDP function. 251*814e85d6SHong Zhang 252*814e85d6SHong Zhang Collective on TS 253*814e85d6SHong Zhang 254*814e85d6SHong Zhang Input Parameters: 255*814e85d6SHong Zhang . ts - The TS context obtained from TSCreate() 256*814e85d6SHong Zhang 257*814e85d6SHong Zhang Notes: 258*814e85d6SHong Zhang TSDRDPFunction() is typically used for sensitivity implementation, 259*814e85d6SHong Zhang so most users would not generally call this routine themselves. 260*814e85d6SHong Zhang 261*814e85d6SHong Zhang Level: developer 262*814e85d6SHong Zhang 263*814e85d6SHong Zhang .keywords: TS, sensitivity 264*814e85d6SHong Zhang .seealso: TSSetCostIntegrand() 265*814e85d6SHong Zhang @*/ 266*814e85d6SHong Zhang PetscErrorCode TSComputeDRDPFunction(TS ts,PetscReal t,Vec y,Vec *drdp) 267*814e85d6SHong Zhang { 268*814e85d6SHong Zhang PetscErrorCode ierr; 269*814e85d6SHong Zhang 270*814e85d6SHong Zhang PetscFunctionBegin; 271*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 272*814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 273*814e85d6SHong Zhang 274*814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 275*814e85d6SHong Zhang ierr = (*ts->drdpfunction)(ts,t,y,drdp,ts->costintegrandctx); CHKERRQ(ierr); 276*814e85d6SHong Zhang PetscStackPop; 277*814e85d6SHong Zhang PetscFunctionReturn(0); 278*814e85d6SHong Zhang } 279*814e85d6SHong Zhang 280*814e85d6SHong Zhang /* --------------------------- Adjoint sensitivity ---------------------------*/ 281*814e85d6SHong Zhang 282*814e85d6SHong Zhang /*@ 283*814e85d6SHong 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 284*814e85d6SHong Zhang for use by the TSAdjoint routines. 285*814e85d6SHong Zhang 286*814e85d6SHong Zhang Logically Collective on TS and Vec 287*814e85d6SHong Zhang 288*814e85d6SHong Zhang Input Parameters: 289*814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 290*814e85d6SHong 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 291*814e85d6SHong Zhang - mu - gradients with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 292*814e85d6SHong Zhang 293*814e85d6SHong Zhang Level: beginner 294*814e85d6SHong Zhang 295*814e85d6SHong Zhang Notes: 296*814e85d6SHong Zhang the entries in these vectors must be correctly initialized with the values lamda_i = df/dy|finaltime mu_i = df/dp|finaltime 297*814e85d6SHong Zhang 298*814e85d6SHong Zhang After TSAdjointSolve() is called the lamba and the mu contain the computed sensitivities 299*814e85d6SHong Zhang 300*814e85d6SHong Zhang .keywords: TS, timestep, set, sensitivity, initial values 301*814e85d6SHong Zhang @*/ 302*814e85d6SHong Zhang PetscErrorCode TSSetCostGradients(TS ts,PetscInt numcost,Vec *lambda,Vec *mu) 303*814e85d6SHong Zhang { 304*814e85d6SHong Zhang PetscFunctionBegin; 305*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 306*814e85d6SHong Zhang PetscValidPointer(lambda,2); 307*814e85d6SHong Zhang ts->vecs_sensi = lambda; 308*814e85d6SHong Zhang ts->vecs_sensip = mu; 309*814e85d6SHong 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"); 310*814e85d6SHong Zhang ts->numcost = numcost; 311*814e85d6SHong Zhang PetscFunctionReturn(0); 312*814e85d6SHong Zhang } 313*814e85d6SHong Zhang 314*814e85d6SHong Zhang /*@ 315*814e85d6SHong Zhang TSGetCostGradients - Returns the gradients from the TSAdjointSolve() 316*814e85d6SHong Zhang 317*814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 318*814e85d6SHong Zhang 319*814e85d6SHong Zhang Input Parameter: 320*814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 321*814e85d6SHong Zhang 322*814e85d6SHong Zhang Output Parameter: 323*814e85d6SHong Zhang + lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 324*814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 325*814e85d6SHong Zhang 326*814e85d6SHong Zhang Level: intermediate 327*814e85d6SHong Zhang 328*814e85d6SHong Zhang .seealso: TSGetTimeStep() 329*814e85d6SHong Zhang 330*814e85d6SHong Zhang .keywords: TS, timestep, get, sensitivity 331*814e85d6SHong Zhang @*/ 332*814e85d6SHong Zhang PetscErrorCode TSGetCostGradients(TS ts,PetscInt *numcost,Vec **lambda,Vec **mu) 333*814e85d6SHong Zhang { 334*814e85d6SHong Zhang PetscFunctionBegin; 335*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 336*814e85d6SHong Zhang if (numcost) *numcost = ts->numcost; 337*814e85d6SHong Zhang if (lambda) *lambda = ts->vecs_sensi; 338*814e85d6SHong Zhang if (mu) *mu = ts->vecs_sensip; 339*814e85d6SHong Zhang PetscFunctionReturn(0); 340*814e85d6SHong Zhang } 341*814e85d6SHong Zhang 342*814e85d6SHong Zhang /*@ 343*814e85d6SHong Zhang TSAdjointSetUp - Sets up the internal data structures for the later use 344*814e85d6SHong Zhang of an adjoint solver 345*814e85d6SHong Zhang 346*814e85d6SHong Zhang Collective on TS 347*814e85d6SHong Zhang 348*814e85d6SHong Zhang Input Parameter: 349*814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 350*814e85d6SHong Zhang 351*814e85d6SHong Zhang Level: advanced 352*814e85d6SHong Zhang 353*814e85d6SHong Zhang .keywords: TS, timestep, setup 354*814e85d6SHong Zhang 355*814e85d6SHong Zhang .seealso: TSCreate(), TSAdjointStep(), TSSetCostGradients() 356*814e85d6SHong Zhang @*/ 357*814e85d6SHong Zhang PetscErrorCode TSAdjointSetUp(TS ts) 358*814e85d6SHong Zhang { 359*814e85d6SHong Zhang PetscErrorCode ierr; 360*814e85d6SHong Zhang 361*814e85d6SHong Zhang PetscFunctionBegin; 362*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 363*814e85d6SHong Zhang if (ts->adjointsetupcalled) PetscFunctionReturn(0); 364*814e85d6SHong Zhang if (!ts->vecs_sensi) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSSetCostGradients() first"); 365*814e85d6SHong Zhang if (ts->vecs_sensip && !ts->Jacp) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_WRONGSTATE,"Must call TSAdjointSetRHSJacobian() first"); 366*814e85d6SHong Zhang 367*814e85d6SHong Zhang if (ts->vec_costintegral) { /* if there is integral in the cost function */ 368*814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensi[0],ts->numcost,&ts->vecs_drdy);CHKERRQ(ierr); 369*814e85d6SHong Zhang if (ts->vecs_sensip){ 370*814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr); 371*814e85d6SHong Zhang } 372*814e85d6SHong Zhang } 373*814e85d6SHong Zhang 374*814e85d6SHong Zhang if (ts->ops->adjointsetup) { 375*814e85d6SHong Zhang ierr = (*ts->ops->adjointsetup)(ts);CHKERRQ(ierr); 376*814e85d6SHong Zhang } 377*814e85d6SHong Zhang ts->adjointsetupcalled = PETSC_TRUE; 378*814e85d6SHong Zhang PetscFunctionReturn(0); 379*814e85d6SHong Zhang } 380*814e85d6SHong Zhang 381*814e85d6SHong Zhang /*@ 382*814e85d6SHong Zhang TSAdjointSetSteps - Sets the number of steps the adjoint solver should take backward in time 383*814e85d6SHong Zhang 384*814e85d6SHong Zhang Logically Collective on TS 385*814e85d6SHong Zhang 386*814e85d6SHong Zhang Input Parameters: 387*814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 388*814e85d6SHong Zhang . steps - number of steps to use 389*814e85d6SHong Zhang 390*814e85d6SHong Zhang Level: intermediate 391*814e85d6SHong Zhang 392*814e85d6SHong Zhang Notes: 393*814e85d6SHong Zhang Normally one does not call this and TSAdjointSolve() integrates back to the original timestep. One can call this 394*814e85d6SHong Zhang so as to integrate back to less than the original timestep 395*814e85d6SHong Zhang 396*814e85d6SHong Zhang .keywords: TS, timestep, set, maximum, iterations 397*814e85d6SHong Zhang 398*814e85d6SHong Zhang .seealso: TSSetExactFinalTime() 399*814e85d6SHong Zhang @*/ 400*814e85d6SHong Zhang PetscErrorCode TSAdjointSetSteps(TS ts,PetscInt steps) 401*814e85d6SHong Zhang { 402*814e85d6SHong Zhang PetscFunctionBegin; 403*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 404*814e85d6SHong Zhang PetscValidLogicalCollectiveInt(ts,steps,2); 405*814e85d6SHong Zhang if (steps < 0) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back a negative number of steps"); 406*814e85d6SHong Zhang if (steps > ts->steps) SETERRQ(PetscObjectComm((PetscObject)ts),PETSC_ERR_ARG_OUTOFRANGE,"Cannot step back more than the total number of forward steps"); 407*814e85d6SHong Zhang ts->adjoint_max_steps = steps; 408*814e85d6SHong Zhang PetscFunctionReturn(0); 409*814e85d6SHong Zhang } 410*814e85d6SHong Zhang 411*814e85d6SHong Zhang /*@C 412*814e85d6SHong Zhang TSAdjointSetRHSJacobian - Deprecated, use TSSetRHSJacobianP() 413*814e85d6SHong Zhang 414*814e85d6SHong Zhang Level: deprecated 415*814e85d6SHong Zhang 416*814e85d6SHong Zhang @*/ 417*814e85d6SHong Zhang PetscErrorCode TSAdjointSetRHSJacobian(TS ts,Mat Amat,PetscErrorCode (*func)(TS,PetscReal,Vec,Mat,void*),void *ctx) 418*814e85d6SHong Zhang { 419*814e85d6SHong Zhang PetscErrorCode ierr; 420*814e85d6SHong Zhang 421*814e85d6SHong Zhang PetscFunctionBegin; 422*814e85d6SHong Zhang PetscValidHeaderSpecific(ts, TS_CLASSID,1); 423*814e85d6SHong Zhang PetscValidHeaderSpecific(Amat,MAT_CLASSID,2); 424*814e85d6SHong Zhang 425*814e85d6SHong Zhang ts->rhsjacobianp = func; 426*814e85d6SHong Zhang ts->rhsjacobianpctx = ctx; 427*814e85d6SHong Zhang if(Amat) { 428*814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Amat);CHKERRQ(ierr); 429*814e85d6SHong Zhang ierr = MatDestroy(&ts->Jacp);CHKERRQ(ierr); 430*814e85d6SHong Zhang ts->Jacp = Amat; 431*814e85d6SHong Zhang } 432*814e85d6SHong Zhang PetscFunctionReturn(0); 433*814e85d6SHong Zhang } 434*814e85d6SHong Zhang 435*814e85d6SHong Zhang /*@C 436*814e85d6SHong Zhang TSAdjointComputeRHSJacobian - Deprecated, use TSComputeRHSJacobianP() 437*814e85d6SHong Zhang 438*814e85d6SHong Zhang Level: deprecated 439*814e85d6SHong Zhang 440*814e85d6SHong Zhang @*/ 441*814e85d6SHong Zhang PetscErrorCode TSAdjointComputeRHSJacobian(TS ts,PetscReal t,Vec X,Mat Amat) 442*814e85d6SHong Zhang { 443*814e85d6SHong Zhang PetscErrorCode ierr; 444*814e85d6SHong Zhang 445*814e85d6SHong Zhang PetscFunctionBegin; 446*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 447*814e85d6SHong Zhang PetscValidHeaderSpecific(X,VEC_CLASSID,3); 448*814e85d6SHong Zhang PetscValidPointer(Amat,4); 449*814e85d6SHong Zhang 450*814e85d6SHong Zhang PetscStackPush("TS user JacobianP function for sensitivity analysis"); 451*814e85d6SHong Zhang ierr = (*ts->rhsjacobianp)(ts,t,X,Amat,ts->rhsjacobianpctx); CHKERRQ(ierr); 452*814e85d6SHong Zhang PetscStackPop; 453*814e85d6SHong Zhang PetscFunctionReturn(0); 454*814e85d6SHong Zhang } 455*814e85d6SHong Zhang 456*814e85d6SHong Zhang /*@ 457*814e85d6SHong Zhang TSAdjointComputeDRDYFunction - Deprecated, use TSComputeDRDYFunction() 458*814e85d6SHong Zhang 459*814e85d6SHong Zhang Level: deprecated 460*814e85d6SHong Zhang 461*814e85d6SHong Zhang @*/ 462*814e85d6SHong Zhang PetscErrorCode TSAdjointComputeDRDYFunction(TS ts,PetscReal t,Vec y,Vec *drdy) 463*814e85d6SHong Zhang { 464*814e85d6SHong Zhang PetscErrorCode ierr; 465*814e85d6SHong Zhang 466*814e85d6SHong Zhang PetscFunctionBegin; 467*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 468*814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 469*814e85d6SHong Zhang 470*814e85d6SHong Zhang PetscStackPush("TS user DRDY function for sensitivity analysis"); 471*814e85d6SHong Zhang ierr = (*ts->drdyfunction)(ts,t,y,drdy,ts->costintegrandctx); CHKERRQ(ierr); 472*814e85d6SHong Zhang PetscStackPop; 473*814e85d6SHong Zhang PetscFunctionReturn(0); 474*814e85d6SHong Zhang } 475*814e85d6SHong Zhang 476*814e85d6SHong Zhang /*@ 477*814e85d6SHong Zhang TSAdjointComputeDRDPFunction - Deprecated, use TSComputeDRDPFunction() 478*814e85d6SHong Zhang 479*814e85d6SHong Zhang Level: deprecated 480*814e85d6SHong Zhang 481*814e85d6SHong Zhang @*/ 482*814e85d6SHong Zhang PetscErrorCode TSAdjointComputeDRDPFunction(TS ts,PetscReal t,Vec y,Vec *drdp) 483*814e85d6SHong Zhang { 484*814e85d6SHong Zhang PetscErrorCode ierr; 485*814e85d6SHong Zhang 486*814e85d6SHong Zhang PetscFunctionBegin; 487*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 488*814e85d6SHong Zhang PetscValidHeaderSpecific(y,VEC_CLASSID,3); 489*814e85d6SHong Zhang 490*814e85d6SHong Zhang PetscStackPush("TS user DRDP function for sensitivity analysis"); 491*814e85d6SHong Zhang ierr = (*ts->drdpfunction)(ts,t,y,drdp,ts->costintegrandctx); CHKERRQ(ierr); 492*814e85d6SHong Zhang PetscStackPop; 493*814e85d6SHong Zhang PetscFunctionReturn(0); 494*814e85d6SHong Zhang } 495*814e85d6SHong Zhang 496*814e85d6SHong Zhang /*@C 497*814e85d6SHong Zhang TSAdjointMonitorSensi - monitors the first lambda sensitivity 498*814e85d6SHong Zhang 499*814e85d6SHong Zhang Level: intermediate 500*814e85d6SHong Zhang 501*814e85d6SHong Zhang .keywords: TS, set, monitor 502*814e85d6SHong Zhang 503*814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 504*814e85d6SHong Zhang @*/ 505*814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSensi(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 506*814e85d6SHong Zhang { 507*814e85d6SHong Zhang PetscErrorCode ierr; 508*814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 509*814e85d6SHong Zhang 510*814e85d6SHong Zhang PetscFunctionBegin; 511*814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 512*814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 513*814e85d6SHong Zhang ierr = VecView(lambda[0],viewer);CHKERRQ(ierr); 514*814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 515*814e85d6SHong Zhang PetscFunctionReturn(0); 516*814e85d6SHong Zhang } 517*814e85d6SHong Zhang 518*814e85d6SHong Zhang /*@C 519*814e85d6SHong Zhang TSAdjointMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type indicated by the user 520*814e85d6SHong Zhang 521*814e85d6SHong Zhang Collective on TS 522*814e85d6SHong Zhang 523*814e85d6SHong Zhang Input Parameters: 524*814e85d6SHong Zhang + ts - TS object you wish to monitor 525*814e85d6SHong Zhang . name - the monitor type one is seeking 526*814e85d6SHong Zhang . help - message indicating what monitoring is done 527*814e85d6SHong Zhang . manual - manual page for the monitor 528*814e85d6SHong Zhang . monitor - the monitor function 529*814e85d6SHong 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 530*814e85d6SHong Zhang 531*814e85d6SHong Zhang Level: developer 532*814e85d6SHong Zhang 533*814e85d6SHong Zhang .seealso: PetscOptionsGetViewer(), PetscOptionsGetReal(), PetscOptionsHasName(), PetscOptionsGetString(), 534*814e85d6SHong Zhang PetscOptionsGetIntArray(), PetscOptionsGetRealArray(), PetscOptionsBool() 535*814e85d6SHong Zhang PetscOptionsInt(), PetscOptionsString(), PetscOptionsReal(), PetscOptionsBool(), 536*814e85d6SHong Zhang PetscOptionsName(), PetscOptionsBegin(), PetscOptionsEnd(), PetscOptionsHead(), 537*814e85d6SHong Zhang PetscOptionsStringArray(),PetscOptionsRealArray(), PetscOptionsScalar(), 538*814e85d6SHong Zhang PetscOptionsBoolGroupBegin(), PetscOptionsBoolGroup(), PetscOptionsBoolGroupEnd(), 539*814e85d6SHong Zhang PetscOptionsFList(), PetscOptionsEList() 540*814e85d6SHong Zhang @*/ 541*814e85d6SHong 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*)) 542*814e85d6SHong Zhang { 543*814e85d6SHong Zhang PetscErrorCode ierr; 544*814e85d6SHong Zhang PetscViewer viewer; 545*814e85d6SHong Zhang PetscViewerFormat format; 546*814e85d6SHong Zhang PetscBool flg; 547*814e85d6SHong Zhang 548*814e85d6SHong Zhang PetscFunctionBegin; 549*814e85d6SHong Zhang ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)ts),((PetscObject)ts)->prefix,name,&viewer,&format,&flg);CHKERRQ(ierr); 550*814e85d6SHong Zhang if (flg) { 551*814e85d6SHong Zhang PetscViewerAndFormat *vf; 552*814e85d6SHong Zhang ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 553*814e85d6SHong Zhang ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 554*814e85d6SHong Zhang if (monitorsetup) { 555*814e85d6SHong Zhang ierr = (*monitorsetup)(ts,vf);CHKERRQ(ierr); 556*814e85d6SHong Zhang } 557*814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,(PetscErrorCode (*)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);CHKERRQ(ierr); 558*814e85d6SHong Zhang } 559*814e85d6SHong Zhang PetscFunctionReturn(0); 560*814e85d6SHong Zhang } 561*814e85d6SHong Zhang 562*814e85d6SHong Zhang /*@C 563*814e85d6SHong Zhang TSAdjointMonitorSet - Sets an ADDITIONAL function that is to be used at every 564*814e85d6SHong Zhang timestep to display the iteration's progress. 565*814e85d6SHong Zhang 566*814e85d6SHong Zhang Logically Collective on TS 567*814e85d6SHong Zhang 568*814e85d6SHong Zhang Input Parameters: 569*814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 570*814e85d6SHong Zhang . adjointmonitor - monitoring routine 571*814e85d6SHong Zhang . adjointmctx - [optional] user-defined context for private data for the 572*814e85d6SHong Zhang monitor routine (use NULL if no context is desired) 573*814e85d6SHong Zhang - adjointmonitordestroy - [optional] routine that frees monitor context 574*814e85d6SHong Zhang (may be NULL) 575*814e85d6SHong Zhang 576*814e85d6SHong Zhang Calling sequence of monitor: 577*814e85d6SHong Zhang $ int adjointmonitor(TS ts,PetscInt steps,PetscReal time,Vec u,PetscInt numcost,Vec *lambda, Vec *mu,void *adjointmctx) 578*814e85d6SHong Zhang 579*814e85d6SHong Zhang + ts - the TS context 580*814e85d6SHong 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 581*814e85d6SHong Zhang been interpolated to) 582*814e85d6SHong Zhang . time - current time 583*814e85d6SHong Zhang . u - current iterate 584*814e85d6SHong Zhang . numcost - number of cost functionos 585*814e85d6SHong Zhang . lambda - sensitivities to initial conditions 586*814e85d6SHong Zhang . mu - sensitivities to parameters 587*814e85d6SHong Zhang - adjointmctx - [optional] adjoint monitoring context 588*814e85d6SHong Zhang 589*814e85d6SHong Zhang Notes: 590*814e85d6SHong Zhang This routine adds an additional monitor to the list of monitors that 591*814e85d6SHong Zhang already has been loaded. 592*814e85d6SHong Zhang 593*814e85d6SHong Zhang Fortran Notes: 594*814e85d6SHong Zhang Only a single monitor function can be set for each TS object 595*814e85d6SHong Zhang 596*814e85d6SHong Zhang Level: intermediate 597*814e85d6SHong Zhang 598*814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 599*814e85d6SHong Zhang 600*814e85d6SHong Zhang .seealso: TSAdjointMonitorCancel() 601*814e85d6SHong Zhang @*/ 602*814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorSet(TS ts,PetscErrorCode (*adjointmonitor)(TS,PetscInt,PetscReal,Vec,PetscInt,Vec*,Vec*,void*),void *adjointmctx,PetscErrorCode (*adjointmdestroy)(void**)) 603*814e85d6SHong Zhang { 604*814e85d6SHong Zhang PetscErrorCode ierr; 605*814e85d6SHong Zhang PetscInt i; 606*814e85d6SHong Zhang PetscBool identical; 607*814e85d6SHong Zhang 608*814e85d6SHong Zhang PetscFunctionBegin; 609*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 610*814e85d6SHong Zhang for (i=0; i<ts->numbermonitors;i++) { 611*814e85d6SHong Zhang ierr = PetscMonitorCompare((PetscErrorCode (*)(void))adjointmonitor,adjointmctx,adjointmdestroy,(PetscErrorCode (*)(void))ts->adjointmonitor[i],ts->adjointmonitorcontext[i],ts->adjointmonitordestroy[i],&identical);CHKERRQ(ierr); 612*814e85d6SHong Zhang if (identical) PetscFunctionReturn(0); 613*814e85d6SHong Zhang } 614*814e85d6SHong Zhang if (ts->numberadjointmonitors >= MAXTSMONITORS) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many adjoint monitors set"); 615*814e85d6SHong Zhang ts->adjointmonitor[ts->numberadjointmonitors] = adjointmonitor; 616*814e85d6SHong Zhang ts->adjointmonitordestroy[ts->numberadjointmonitors] = adjointmdestroy; 617*814e85d6SHong Zhang ts->adjointmonitorcontext[ts->numberadjointmonitors++] = (void*)adjointmctx; 618*814e85d6SHong Zhang PetscFunctionReturn(0); 619*814e85d6SHong Zhang } 620*814e85d6SHong Zhang 621*814e85d6SHong Zhang /*@C 622*814e85d6SHong Zhang TSAdjointMonitorCancel - Clears all the adjoint monitors that have been set on a time-step object. 623*814e85d6SHong Zhang 624*814e85d6SHong Zhang Logically Collective on TS 625*814e85d6SHong Zhang 626*814e85d6SHong Zhang Input Parameters: 627*814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 628*814e85d6SHong Zhang 629*814e85d6SHong Zhang Notes: 630*814e85d6SHong Zhang There is no way to remove a single, specific monitor. 631*814e85d6SHong Zhang 632*814e85d6SHong Zhang Level: intermediate 633*814e85d6SHong Zhang 634*814e85d6SHong Zhang .keywords: TS, timestep, set, adjoint, monitor 635*814e85d6SHong Zhang 636*814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 637*814e85d6SHong Zhang @*/ 638*814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorCancel(TS ts) 639*814e85d6SHong Zhang { 640*814e85d6SHong Zhang PetscErrorCode ierr; 641*814e85d6SHong Zhang PetscInt i; 642*814e85d6SHong Zhang 643*814e85d6SHong Zhang PetscFunctionBegin; 644*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 645*814e85d6SHong Zhang for (i=0; i<ts->numberadjointmonitors; i++) { 646*814e85d6SHong Zhang if (ts->adjointmonitordestroy[i]) { 647*814e85d6SHong Zhang ierr = (*ts->adjointmonitordestroy[i])(&ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 648*814e85d6SHong Zhang } 649*814e85d6SHong Zhang } 650*814e85d6SHong Zhang ts->numberadjointmonitors = 0; 651*814e85d6SHong Zhang PetscFunctionReturn(0); 652*814e85d6SHong Zhang } 653*814e85d6SHong Zhang 654*814e85d6SHong Zhang /*@C 655*814e85d6SHong Zhang TSAdjointMonitorDefault - the default monitor of adjoint computations 656*814e85d6SHong Zhang 657*814e85d6SHong Zhang Level: intermediate 658*814e85d6SHong Zhang 659*814e85d6SHong Zhang .keywords: TS, set, monitor 660*814e85d6SHong Zhang 661*814e85d6SHong Zhang .seealso: TSAdjointMonitorSet() 662*814e85d6SHong Zhang @*/ 663*814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDefault(TS ts,PetscInt step,PetscReal ptime,Vec v,PetscInt numcost,Vec *lambda,Vec *mu,PetscViewerAndFormat *vf) 664*814e85d6SHong Zhang { 665*814e85d6SHong Zhang PetscErrorCode ierr; 666*814e85d6SHong Zhang PetscViewer viewer = vf->viewer; 667*814e85d6SHong Zhang 668*814e85d6SHong Zhang PetscFunctionBegin; 669*814e85d6SHong Zhang PetscValidHeaderSpecific(viewer,PETSC_VIEWER_CLASSID,4); 670*814e85d6SHong Zhang ierr = PetscViewerPushFormat(viewer,vf->format);CHKERRQ(ierr); 671*814e85d6SHong Zhang ierr = PetscViewerASCIIAddTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 672*814e85d6SHong 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); 673*814e85d6SHong Zhang ierr = PetscViewerASCIISubtractTab(viewer,((PetscObject)ts)->tablevel);CHKERRQ(ierr); 674*814e85d6SHong Zhang ierr = PetscViewerPopFormat(viewer);CHKERRQ(ierr); 675*814e85d6SHong Zhang PetscFunctionReturn(0); 676*814e85d6SHong Zhang } 677*814e85d6SHong Zhang 678*814e85d6SHong Zhang /*@C 679*814e85d6SHong Zhang TSAdjointMonitorDrawSensi - Monitors progress of the adjoint TS solvers by calling 680*814e85d6SHong Zhang VecView() for the sensitivities to initial states at each timestep 681*814e85d6SHong Zhang 682*814e85d6SHong Zhang Collective on TS 683*814e85d6SHong Zhang 684*814e85d6SHong Zhang Input Parameters: 685*814e85d6SHong Zhang + ts - the TS context 686*814e85d6SHong Zhang . step - current time-step 687*814e85d6SHong Zhang . ptime - current time 688*814e85d6SHong Zhang . u - current state 689*814e85d6SHong Zhang . numcost - number of cost functions 690*814e85d6SHong Zhang . lambda - sensitivities to initial conditions 691*814e85d6SHong Zhang . mu - sensitivities to parameters 692*814e85d6SHong Zhang - dummy - either a viewer or NULL 693*814e85d6SHong Zhang 694*814e85d6SHong Zhang Level: intermediate 695*814e85d6SHong Zhang 696*814e85d6SHong Zhang .keywords: TS, vector, adjoint, monitor, view 697*814e85d6SHong Zhang 698*814e85d6SHong Zhang .seealso: TSAdjointMonitorSet(), TSAdjointMonitorDefault(), VecView() 699*814e85d6SHong Zhang @*/ 700*814e85d6SHong Zhang PetscErrorCode TSAdjointMonitorDrawSensi(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda,Vec *mu,void *dummy) 701*814e85d6SHong Zhang { 702*814e85d6SHong Zhang PetscErrorCode ierr; 703*814e85d6SHong Zhang TSMonitorDrawCtx ictx = (TSMonitorDrawCtx)dummy; 704*814e85d6SHong Zhang PetscDraw draw; 705*814e85d6SHong Zhang PetscReal xl,yl,xr,yr,h; 706*814e85d6SHong Zhang char time[32]; 707*814e85d6SHong Zhang 708*814e85d6SHong Zhang PetscFunctionBegin; 709*814e85d6SHong Zhang if (!(((ictx->howoften > 0) && (!(step % ictx->howoften))) || ((ictx->howoften == -1) && ts->reason))) PetscFunctionReturn(0); 710*814e85d6SHong Zhang 711*814e85d6SHong Zhang ierr = VecView(lambda[0],ictx->viewer);CHKERRQ(ierr); 712*814e85d6SHong Zhang ierr = PetscViewerDrawGetDraw(ictx->viewer,0,&draw);CHKERRQ(ierr); 713*814e85d6SHong Zhang ierr = PetscSNPrintf(time,32,"Timestep %d Time %g",(int)step,(double)ptime);CHKERRQ(ierr); 714*814e85d6SHong Zhang ierr = PetscDrawGetCoordinates(draw,&xl,&yl,&xr,&yr);CHKERRQ(ierr); 715*814e85d6SHong Zhang h = yl + .95*(yr - yl); 716*814e85d6SHong Zhang ierr = PetscDrawStringCentered(draw,.5*(xl+xr),h,PETSC_DRAW_BLACK,time);CHKERRQ(ierr); 717*814e85d6SHong Zhang ierr = PetscDrawFlush(draw);CHKERRQ(ierr); 718*814e85d6SHong Zhang PetscFunctionReturn(0); 719*814e85d6SHong Zhang } 720*814e85d6SHong Zhang 721*814e85d6SHong Zhang /* 722*814e85d6SHong Zhang TSAdjointSetFromOptions - Sets various TSAdjoint parameters from user options. 723*814e85d6SHong Zhang 724*814e85d6SHong Zhang Collective on TSAdjoint 725*814e85d6SHong Zhang 726*814e85d6SHong Zhang Input Parameter: 727*814e85d6SHong Zhang . ts - the TS context 728*814e85d6SHong Zhang 729*814e85d6SHong Zhang Options Database Keys: 730*814e85d6SHong Zhang + -ts_adjoint_solve <yes,no> After solving the ODE/DAE solve the adjoint problem (requires -ts_save_trajectory) 731*814e85d6SHong Zhang . -ts_adjoint_monitor - print information at each adjoint time step 732*814e85d6SHong Zhang - -ts_adjoint_monitor_draw_sensi - monitor the sensitivity of the first cost function wrt initial conditions (lambda[0]) graphically 733*814e85d6SHong Zhang 734*814e85d6SHong Zhang Level: developer 735*814e85d6SHong Zhang 736*814e85d6SHong Zhang Notes: 737*814e85d6SHong Zhang This is not normally called directly by users 738*814e85d6SHong Zhang 739*814e85d6SHong Zhang .keywords: TS, trajectory, timestep, set, options, database 740*814e85d6SHong Zhang 741*814e85d6SHong Zhang .seealso: TSSetSaveTrajectory(), TSTrajectorySetUp() 742*814e85d6SHong Zhang */ 743*814e85d6SHong Zhang PetscErrorCode TSAdjointSetFromOptions(PetscOptionItems *PetscOptionsObject,TS ts) 744*814e85d6SHong Zhang { 745*814e85d6SHong Zhang PetscBool tflg,opt; 746*814e85d6SHong Zhang PetscErrorCode ierr; 747*814e85d6SHong Zhang 748*814e85d6SHong Zhang PetscFunctionBegin; 749*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,2); 750*814e85d6SHong Zhang ierr = PetscOptionsHead(PetscOptionsObject,"TS Adjoint options");CHKERRQ(ierr); 751*814e85d6SHong Zhang tflg = ts->adjoint_solve ? PETSC_TRUE : PETSC_FALSE; 752*814e85d6SHong Zhang ierr = PetscOptionsBool("-ts_adjoint_solve","Solve the adjoint problem immediately after solving the forward problem","",tflg,&tflg,&opt);CHKERRQ(ierr); 753*814e85d6SHong Zhang if (opt) { 754*814e85d6SHong Zhang ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr); 755*814e85d6SHong Zhang ts->adjoint_solve = tflg; 756*814e85d6SHong Zhang } 757*814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor","Monitor adjoint timestep size","TSAdjointMonitorDefault",TSAdjointMonitorDefault,NULL);CHKERRQ(ierr); 758*814e85d6SHong Zhang ierr = TSAdjointMonitorSetFromOptions(ts,"-ts_adjoint_monitor_sensi","Monitor sensitivity in the adjoint computation","TSAdjointMonitorSensi",TSAdjointMonitorSensi,NULL);CHKERRQ(ierr); 759*814e85d6SHong Zhang opt = PETSC_FALSE; 760*814e85d6SHong Zhang ierr = PetscOptionsName("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",&opt);CHKERRQ(ierr); 761*814e85d6SHong Zhang if (opt) { 762*814e85d6SHong Zhang TSMonitorDrawCtx ctx; 763*814e85d6SHong Zhang PetscInt howoften = 1; 764*814e85d6SHong Zhang 765*814e85d6SHong Zhang ierr = PetscOptionsInt("-ts_adjoint_monitor_draw_sensi","Monitor adjoint sensitivities (lambda only) graphically","TSAdjointMonitorDrawSensi",howoften,&howoften,NULL);CHKERRQ(ierr); 766*814e85d6SHong Zhang ierr = TSMonitorDrawCtxCreate(PetscObjectComm((PetscObject)ts),0,0,PETSC_DECIDE,PETSC_DECIDE,300,300,howoften,&ctx);CHKERRQ(ierr); 767*814e85d6SHong Zhang ierr = TSAdjointMonitorSet(ts,TSAdjointMonitorDrawSensi,ctx,(PetscErrorCode (*)(void**))TSMonitorDrawCtxDestroy);CHKERRQ(ierr); 768*814e85d6SHong Zhang } 769*814e85d6SHong Zhang PetscFunctionReturn(0); 770*814e85d6SHong Zhang } 771*814e85d6SHong Zhang 772*814e85d6SHong Zhang /*@ 773*814e85d6SHong Zhang TSAdjointStep - Steps one time step backward in the adjoint run 774*814e85d6SHong Zhang 775*814e85d6SHong Zhang Collective on TS 776*814e85d6SHong Zhang 777*814e85d6SHong Zhang Input Parameter: 778*814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 779*814e85d6SHong Zhang 780*814e85d6SHong Zhang Level: intermediate 781*814e85d6SHong Zhang 782*814e85d6SHong Zhang .keywords: TS, adjoint, step 783*814e85d6SHong Zhang 784*814e85d6SHong Zhang .seealso: TSAdjointSetUp(), TSAdjointSolve() 785*814e85d6SHong Zhang @*/ 786*814e85d6SHong Zhang PetscErrorCode TSAdjointStep(TS ts) 787*814e85d6SHong Zhang { 788*814e85d6SHong Zhang DM dm; 789*814e85d6SHong Zhang PetscErrorCode ierr; 790*814e85d6SHong Zhang 791*814e85d6SHong Zhang PetscFunctionBegin; 792*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 793*814e85d6SHong Zhang ierr = TSGetDM(ts,&dm);CHKERRQ(ierr); 794*814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 795*814e85d6SHong Zhang 796*814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vec_sol,(PetscObject)ts,"-ts_view_solution");CHKERRQ(ierr); 797*814e85d6SHong Zhang 798*814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 799*814e85d6SHong Zhang ts->ptime_prev = ts->ptime; 800*814e85d6SHong 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); 801*814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 802*814e85d6SHong Zhang ierr = (*ts->ops->adjointstep)(ts);CHKERRQ(ierr); 803*814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_AdjointStep,ts,0,0,0);CHKERRQ(ierr); 804*814e85d6SHong Zhang ts->adjoint_steps++; ts->steps--; 805*814e85d6SHong Zhang 806*814e85d6SHong Zhang if (ts->reason < 0) { 807*814e85d6SHong Zhang if (ts->errorifstepfailed) { 808*814e85d6SHong 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]); 809*814e85d6SHong 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]); 810*814e85d6SHong Zhang else SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_NOT_CONVERGED,"TSStep has failed due to %s",TSConvergedReasons[ts->reason]); 811*814e85d6SHong Zhang } 812*814e85d6SHong Zhang } else if (!ts->reason) { 813*814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 814*814e85d6SHong Zhang } 815*814e85d6SHong Zhang PetscFunctionReturn(0); 816*814e85d6SHong Zhang } 817*814e85d6SHong Zhang 818*814e85d6SHong Zhang /*@ 819*814e85d6SHong Zhang TSAdjointSolve - Solves the discrete ajoint problem for an ODE/DAE 820*814e85d6SHong Zhang 821*814e85d6SHong Zhang Collective on TS 822*814e85d6SHong Zhang 823*814e85d6SHong Zhang Input Parameter: 824*814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 825*814e85d6SHong Zhang 826*814e85d6SHong Zhang Options Database: 827*814e85d6SHong Zhang . -ts_adjoint_view_solution <viewerinfo> - views the first gradient with respect to the initial values 828*814e85d6SHong Zhang 829*814e85d6SHong Zhang Level: intermediate 830*814e85d6SHong Zhang 831*814e85d6SHong Zhang Notes: 832*814e85d6SHong Zhang This must be called after a call to TSSolve() that solves the forward problem 833*814e85d6SHong Zhang 834*814e85d6SHong Zhang By default this will integrate back to the initial time, one can use TSAdjointSetSteps() to step back to a later time 835*814e85d6SHong Zhang 836*814e85d6SHong Zhang .keywords: TS, timestep, solve 837*814e85d6SHong Zhang 838*814e85d6SHong Zhang .seealso: TSCreate(), TSSetCostGradients(), TSSetSolution(), TSAdjointStep() 839*814e85d6SHong Zhang @*/ 840*814e85d6SHong Zhang PetscErrorCode TSAdjointSolve(TS ts) 841*814e85d6SHong Zhang { 842*814e85d6SHong Zhang PetscErrorCode ierr; 843*814e85d6SHong Zhang 844*814e85d6SHong Zhang PetscFunctionBegin; 845*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 846*814e85d6SHong Zhang ierr = TSAdjointSetUp(ts);CHKERRQ(ierr); 847*814e85d6SHong Zhang 848*814e85d6SHong Zhang /* reset time step and iteration counters */ 849*814e85d6SHong Zhang ts->adjoint_steps = 0; 850*814e85d6SHong Zhang ts->ksp_its = 0; 851*814e85d6SHong Zhang ts->snes_its = 0; 852*814e85d6SHong Zhang ts->num_snes_failures = 0; 853*814e85d6SHong Zhang ts->reject = 0; 854*814e85d6SHong Zhang ts->reason = TS_CONVERGED_ITERATING; 855*814e85d6SHong Zhang 856*814e85d6SHong Zhang if (!ts->adjoint_max_steps) ts->adjoint_max_steps = ts->steps; 857*814e85d6SHong Zhang if (ts->adjoint_steps >= ts->adjoint_max_steps) ts->reason = TS_CONVERGED_ITS; 858*814e85d6SHong Zhang 859*814e85d6SHong Zhang while (!ts->reason) { 860*814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 861*814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 862*814e85d6SHong Zhang ierr = TSAdjointEventHandler(ts);CHKERRQ(ierr); 863*814e85d6SHong Zhang ierr = TSAdjointStep(ts);CHKERRQ(ierr); 864*814e85d6SHong Zhang if (ts->vec_costintegral && !ts->costintegralfwd) { 865*814e85d6SHong Zhang ierr = TSAdjointCostIntegral(ts);CHKERRQ(ierr); 866*814e85d6SHong Zhang } 867*814e85d6SHong Zhang } 868*814e85d6SHong Zhang ierr = TSTrajectoryGet(ts->trajectory,ts,ts->steps,&ts->ptime);CHKERRQ(ierr); 869*814e85d6SHong Zhang ierr = TSAdjointMonitor(ts,ts->steps,ts->ptime,ts->vec_sol,ts->numcost,ts->vecs_sensi,ts->vecs_sensip);CHKERRQ(ierr); 870*814e85d6SHong Zhang ts->solvetime = ts->ptime; 871*814e85d6SHong Zhang ierr = TSTrajectoryViewFromOptions(ts->trajectory,NULL,"-ts_trajectory_view");CHKERRQ(ierr); 872*814e85d6SHong Zhang ierr = VecViewFromOptions(ts->vecs_sensi[0],(PetscObject) ts, "-ts_adjoint_view_solution");CHKERRQ(ierr); 873*814e85d6SHong Zhang ts->adjoint_max_steps = 0; 874*814e85d6SHong Zhang PetscFunctionReturn(0); 875*814e85d6SHong Zhang } 876*814e85d6SHong Zhang 877*814e85d6SHong Zhang /*@C 878*814e85d6SHong Zhang TSAdjointMonitor - Runs all user-provided adjoint monitor routines set using TSAdjointMonitorSet() 879*814e85d6SHong Zhang 880*814e85d6SHong Zhang Collective on TS 881*814e85d6SHong Zhang 882*814e85d6SHong Zhang Input Parameters: 883*814e85d6SHong Zhang + ts - time stepping context obtained from TSCreate() 884*814e85d6SHong Zhang . step - step number that has just completed 885*814e85d6SHong Zhang . ptime - model time of the state 886*814e85d6SHong Zhang . u - state at the current model time 887*814e85d6SHong Zhang . numcost - number of cost functions (dimension of lambda or mu) 888*814e85d6SHong Zhang . lambda - vectors containing the gradients of the cost functions with respect to the ODE/DAE solution variables 889*814e85d6SHong Zhang - mu - vectors containing the gradients of the cost functions with respect to the problem parameters 890*814e85d6SHong Zhang 891*814e85d6SHong Zhang Notes: 892*814e85d6SHong Zhang TSAdjointMonitor() is typically used automatically within the time stepping implementations. 893*814e85d6SHong Zhang Users would almost never call this routine directly. 894*814e85d6SHong Zhang 895*814e85d6SHong Zhang Level: developer 896*814e85d6SHong Zhang 897*814e85d6SHong Zhang .keywords: TS, timestep 898*814e85d6SHong Zhang @*/ 899*814e85d6SHong Zhang PetscErrorCode TSAdjointMonitor(TS ts,PetscInt step,PetscReal ptime,Vec u,PetscInt numcost,Vec *lambda, Vec *mu) 900*814e85d6SHong Zhang { 901*814e85d6SHong Zhang PetscErrorCode ierr; 902*814e85d6SHong Zhang PetscInt i,n = ts->numberadjointmonitors; 903*814e85d6SHong Zhang 904*814e85d6SHong Zhang PetscFunctionBegin; 905*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 906*814e85d6SHong Zhang PetscValidHeaderSpecific(u,VEC_CLASSID,4); 907*814e85d6SHong Zhang ierr = VecLockPush(u);CHKERRQ(ierr); 908*814e85d6SHong Zhang for (i=0; i<n; i++) { 909*814e85d6SHong Zhang ierr = (*ts->adjointmonitor[i])(ts,step,ptime,u,numcost,lambda,mu,ts->adjointmonitorcontext[i]);CHKERRQ(ierr); 910*814e85d6SHong Zhang } 911*814e85d6SHong Zhang ierr = VecLockPop(u);CHKERRQ(ierr); 912*814e85d6SHong Zhang PetscFunctionReturn(0); 913*814e85d6SHong Zhang } 914*814e85d6SHong Zhang 915*814e85d6SHong Zhang /*@ 916*814e85d6SHong Zhang TSAdjointCostIntegral - Evaluate the cost integral in the adjoint run. 917*814e85d6SHong Zhang 918*814e85d6SHong Zhang Collective on TS 919*814e85d6SHong Zhang 920*814e85d6SHong Zhang Input Arguments: 921*814e85d6SHong Zhang . ts - time stepping context 922*814e85d6SHong Zhang 923*814e85d6SHong Zhang Level: advanced 924*814e85d6SHong Zhang 925*814e85d6SHong Zhang Notes: 926*814e85d6SHong Zhang This function cannot be called until TSAdjointStep() has been completed. 927*814e85d6SHong Zhang 928*814e85d6SHong Zhang .seealso: TSAdjointSolve(), TSAdjointStep 929*814e85d6SHong Zhang @*/ 930*814e85d6SHong Zhang PetscErrorCode TSAdjointCostIntegral(TS ts) 931*814e85d6SHong Zhang { 932*814e85d6SHong Zhang PetscErrorCode ierr; 933*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 934*814e85d6SHong 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); 935*814e85d6SHong Zhang ierr = (*ts->ops->adjointintegral)(ts);CHKERRQ(ierr); 936*814e85d6SHong Zhang PetscFunctionReturn(0); 937*814e85d6SHong Zhang } 938*814e85d6SHong Zhang 939*814e85d6SHong Zhang /* ------------------ Forward (tangent linear) sensitivity ------------------*/ 940*814e85d6SHong Zhang 941*814e85d6SHong Zhang /*@ 942*814e85d6SHong Zhang TSForwardSetUp - Sets up the internal data structures for the later use 943*814e85d6SHong Zhang of forward sensitivity analysis 944*814e85d6SHong Zhang 945*814e85d6SHong Zhang Collective on TS 946*814e85d6SHong Zhang 947*814e85d6SHong Zhang Input Parameter: 948*814e85d6SHong Zhang . ts - the TS context obtained from TSCreate() 949*814e85d6SHong Zhang 950*814e85d6SHong Zhang Level: advanced 951*814e85d6SHong Zhang 952*814e85d6SHong Zhang .keywords: TS, forward sensitivity, setup 953*814e85d6SHong Zhang 954*814e85d6SHong Zhang .seealso: TSCreate(), TSDestroy(), TSSetUp() 955*814e85d6SHong Zhang @*/ 956*814e85d6SHong Zhang PetscErrorCode TSForwardSetUp(TS ts) 957*814e85d6SHong Zhang { 958*814e85d6SHong Zhang PetscErrorCode ierr; 959*814e85d6SHong Zhang 960*814e85d6SHong Zhang PetscFunctionBegin; 961*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 962*814e85d6SHong Zhang if (ts->forwardsetupcalled) PetscFunctionReturn(0); 963*814e85d6SHong Zhang if (ts->vec_costintegral && !ts->vecs_integral_sensip ) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"Must call TSForwardSetIntegralGradients() before TSSetCostIntegrand()"); 964*814e85d6SHong Zhang if (ts->vecs_integral_sensip) { 965*814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vec_sol,ts->numcost,&ts->vecs_drdy);CHKERRQ(ierr); 966*814e85d6SHong Zhang ierr = VecDuplicateVecs(ts->vecs_integral_sensip[0],ts->numcost,&ts->vecs_drdp);CHKERRQ(ierr); 967*814e85d6SHong Zhang } 968*814e85d6SHong Zhang 969*814e85d6SHong Zhang if (ts->ops->forwardsetup) { 970*814e85d6SHong Zhang ierr = (*ts->ops->forwardsetup)(ts);CHKERRQ(ierr); 971*814e85d6SHong Zhang } 972*814e85d6SHong Zhang ts->forwardsetupcalled = PETSC_TRUE; 973*814e85d6SHong Zhang PetscFunctionReturn(0); 974*814e85d6SHong Zhang } 975*814e85d6SHong Zhang 976*814e85d6SHong Zhang /*@ 977*814e85d6SHong Zhang TSForwardSetIntegralGradients - Set the vectors holding forward sensitivities of the integral term. 978*814e85d6SHong Zhang 979*814e85d6SHong Zhang Input Parameter: 980*814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 981*814e85d6SHong Zhang . numfwdint- number of integrals 982*814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 983*814e85d6SHong Zhang 984*814e85d6SHong Zhang Level: intermediate 985*814e85d6SHong Zhang 986*814e85d6SHong Zhang .keywords: TS, forward sensitivity 987*814e85d6SHong Zhang 988*814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 989*814e85d6SHong Zhang @*/ 990*814e85d6SHong Zhang PetscErrorCode TSForwardSetIntegralGradients(TS ts,PetscInt numfwdint,Vec *vp) 991*814e85d6SHong Zhang { 992*814e85d6SHong Zhang PetscFunctionBegin; 993*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 994*814e85d6SHong 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()"); 995*814e85d6SHong Zhang if (!ts->numcost) ts->numcost = numfwdint; 996*814e85d6SHong Zhang 997*814e85d6SHong Zhang ts->vecs_integral_sensip = vp; 998*814e85d6SHong Zhang PetscFunctionReturn(0); 999*814e85d6SHong Zhang } 1000*814e85d6SHong Zhang 1001*814e85d6SHong Zhang /*@ 1002*814e85d6SHong Zhang TSForwardGetIntegralGradients - Returns the forward sensitivities ofthe integral term. 1003*814e85d6SHong Zhang 1004*814e85d6SHong Zhang Input Parameter: 1005*814e85d6SHong Zhang . ts- the TS context obtained from TSCreate() 1006*814e85d6SHong Zhang 1007*814e85d6SHong Zhang Output Parameter: 1008*814e85d6SHong Zhang . vp = the vectors containing the gradients for each integral w.r.t. parameters 1009*814e85d6SHong Zhang 1010*814e85d6SHong Zhang Level: intermediate 1011*814e85d6SHong Zhang 1012*814e85d6SHong Zhang .keywords: TS, forward sensitivity 1013*814e85d6SHong Zhang 1014*814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1015*814e85d6SHong Zhang @*/ 1016*814e85d6SHong Zhang PetscErrorCode TSForwardGetIntegralGradients(TS ts,PetscInt *numfwdint,Vec **vp) 1017*814e85d6SHong Zhang { 1018*814e85d6SHong Zhang PetscFunctionBegin; 1019*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1020*814e85d6SHong Zhang PetscValidPointer(vp,3); 1021*814e85d6SHong Zhang if (numfwdint) *numfwdint = ts->numcost; 1022*814e85d6SHong Zhang if (vp) *vp = ts->vecs_integral_sensip; 1023*814e85d6SHong Zhang PetscFunctionReturn(0); 1024*814e85d6SHong Zhang } 1025*814e85d6SHong Zhang 1026*814e85d6SHong Zhang /*@ 1027*814e85d6SHong Zhang TSForwardStep - Compute the forward sensitivity for one time step. 1028*814e85d6SHong Zhang 1029*814e85d6SHong Zhang Collective on TS 1030*814e85d6SHong Zhang 1031*814e85d6SHong Zhang Input Arguments: 1032*814e85d6SHong Zhang . ts - time stepping context 1033*814e85d6SHong Zhang 1034*814e85d6SHong Zhang Level: advanced 1035*814e85d6SHong Zhang 1036*814e85d6SHong Zhang Notes: 1037*814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1038*814e85d6SHong Zhang 1039*814e85d6SHong Zhang .keywords: TS, forward sensitivity 1040*814e85d6SHong Zhang 1041*814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardSetUp() 1042*814e85d6SHong Zhang @*/ 1043*814e85d6SHong Zhang PetscErrorCode TSForwardStep(TS ts) 1044*814e85d6SHong Zhang { 1045*814e85d6SHong Zhang PetscErrorCode ierr; 1046*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1047*814e85d6SHong Zhang if (!ts->ops->forwardstep) SETERRQ1(PetscObjectComm((PetscObject)ts),PETSC_ERR_SUP,"%s does not provide forward sensitivity analysis",((PetscObject)ts)->type_name); 1048*814e85d6SHong Zhang ierr = PetscLogEventBegin(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1049*814e85d6SHong Zhang ierr = (*ts->ops->forwardstep)(ts);CHKERRQ(ierr); 1050*814e85d6SHong Zhang ierr = PetscLogEventEnd(TS_ForwardStep,ts,0,0,0);CHKERRQ(ierr); 1051*814e85d6SHong Zhang PetscFunctionReturn(0); 1052*814e85d6SHong Zhang } 1053*814e85d6SHong Zhang 1054*814e85d6SHong Zhang /*@ 1055*814e85d6SHong Zhang TSForwardSetSensitivities - Sets the initial value of the trajectory sensitivities of solution w.r.t. the problem parameters and initial values. 1056*814e85d6SHong Zhang 1057*814e85d6SHong Zhang Logically Collective on TS and Vec 1058*814e85d6SHong Zhang 1059*814e85d6SHong Zhang Input Parameters: 1060*814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1061*814e85d6SHong Zhang . nump - number of parameters 1062*814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1063*814e85d6SHong Zhang 1064*814e85d6SHong Zhang Level: beginner 1065*814e85d6SHong Zhang 1066*814e85d6SHong Zhang Notes: 1067*814e85d6SHong Zhang Forward sensitivity is also called 'trajectory sensitivity' in some fields such as power systems. 1068*814e85d6SHong Zhang This function turns on a flag to trigger TSSolve() to compute forward sensitivities automatically. 1069*814e85d6SHong Zhang You must call this function before TSSolve(). 1070*814e85d6SHong Zhang The entries in the sensitivity matrix must be correctly initialized with the values S = dy/dp|startingtime. 1071*814e85d6SHong Zhang 1072*814e85d6SHong Zhang .keywords: TS, timestep, set, forward sensitivity, initial values 1073*814e85d6SHong Zhang 1074*814e85d6SHong Zhang .seealso: TSForwardGetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1075*814e85d6SHong Zhang @*/ 1076*814e85d6SHong Zhang PetscErrorCode TSForwardSetSensitivities(TS ts,PetscInt nump,Mat Smat) 1077*814e85d6SHong Zhang { 1078*814e85d6SHong Zhang PetscErrorCode ierr; 1079*814e85d6SHong Zhang 1080*814e85d6SHong Zhang PetscFunctionBegin; 1081*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1082*814e85d6SHong Zhang PetscValidHeaderSpecific(Smat,MAT_CLASSID,3); 1083*814e85d6SHong Zhang ts->forward_solve = PETSC_TRUE; 1084*814e85d6SHong Zhang ts->num_parameters = nump; 1085*814e85d6SHong Zhang ierr = PetscObjectReference((PetscObject)Smat);CHKERRQ(ierr); 1086*814e85d6SHong Zhang ierr = MatDestroy(&ts->mat_sensip);CHKERRQ(ierr); 1087*814e85d6SHong Zhang ts->mat_sensip = Smat; 1088*814e85d6SHong Zhang PetscFunctionReturn(0); 1089*814e85d6SHong Zhang } 1090*814e85d6SHong Zhang 1091*814e85d6SHong Zhang /*@ 1092*814e85d6SHong Zhang TSForwardGetSensitivities - Returns the trajectory sensitivities 1093*814e85d6SHong Zhang 1094*814e85d6SHong Zhang Not Collective, but Vec returned is parallel if TS is parallel 1095*814e85d6SHong Zhang 1096*814e85d6SHong Zhang Output Parameter: 1097*814e85d6SHong Zhang + ts - the TS context obtained from TSCreate() 1098*814e85d6SHong Zhang . nump - number of parameters 1099*814e85d6SHong Zhang - Smat - sensitivities with respect to the parameters, the number of entries in these vectors is the same as the number of parameters 1100*814e85d6SHong Zhang 1101*814e85d6SHong Zhang Level: intermediate 1102*814e85d6SHong Zhang 1103*814e85d6SHong Zhang .keywords: TS, forward sensitivity 1104*814e85d6SHong Zhang 1105*814e85d6SHong Zhang .seealso: TSForwardSetSensitivities(), TSForwardSetIntegralGradients(), TSForwardGetIntegralGradients(), TSForwardStep() 1106*814e85d6SHong Zhang @*/ 1107*814e85d6SHong Zhang PetscErrorCode TSForwardGetSensitivities(TS ts,PetscInt *nump,Mat *Smat) 1108*814e85d6SHong Zhang { 1109*814e85d6SHong Zhang PetscFunctionBegin; 1110*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1111*814e85d6SHong Zhang if (nump) *nump = ts->num_parameters; 1112*814e85d6SHong Zhang if (Smat) *Smat = ts->mat_sensip; 1113*814e85d6SHong Zhang PetscFunctionReturn(0); 1114*814e85d6SHong Zhang } 1115*814e85d6SHong Zhang 1116*814e85d6SHong Zhang /*@ 1117*814e85d6SHong Zhang TSForwardCostIntegral - Evaluate the cost integral in the forward run. 1118*814e85d6SHong Zhang 1119*814e85d6SHong Zhang Collective on TS 1120*814e85d6SHong Zhang 1121*814e85d6SHong Zhang Input Arguments: 1122*814e85d6SHong Zhang . ts - time stepping context 1123*814e85d6SHong Zhang 1124*814e85d6SHong Zhang Level: advanced 1125*814e85d6SHong Zhang 1126*814e85d6SHong Zhang Notes: 1127*814e85d6SHong Zhang This function cannot be called until TSStep() has been completed. 1128*814e85d6SHong Zhang 1129*814e85d6SHong Zhang .seealso: TSSolve(), TSAdjointCostIntegral() 1130*814e85d6SHong Zhang @*/ 1131*814e85d6SHong Zhang PetscErrorCode TSForwardCostIntegral(TS ts) 1132*814e85d6SHong Zhang { 1133*814e85d6SHong Zhang PetscErrorCode ierr; 1134*814e85d6SHong Zhang PetscValidHeaderSpecific(ts,TS_CLASSID,1); 1135*814e85d6SHong 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); 1136*814e85d6SHong Zhang ierr = (*ts->ops->forwardintegral)(ts);CHKERRQ(ierr); 1137*814e85d6SHong Zhang PetscFunctionReturn(0); 1138*814e85d6SHong Zhang } 1139