Actual source code: theta.c

  1: /*
  2:   Code for timestepping with implicit Theta method
  3: */
  4: #include <petsc/private/tsimpl.h>
  5: #include <petscsnes.h>
  6: #include <petscdm.h>
  7: #include <petscmat.h>

  9: typedef struct {
 10:   /* context for time stepping */
 11:   PetscReal    stage_time;
 12:   Vec          Stages[2];   /* Storage for stage solutions */
 13:   Vec          X0, X, Xdot; /* Storage for u^n, u^n + dt a_{11} k_1, and time derivative u^{n+1}_t */
 14:   Vec          affine;      /* Affine vector needed for residual at beginning of step in endpoint formulation */
 15:   PetscReal    Theta;
 16:   PetscReal    shift; /* Shift parameter for SNES Jacobian, used by forward, TLM and adjoint */
 17:   PetscInt     order;
 18:   PetscBool    endpoint;
 19:   PetscBool    extrapolate;
 20:   TSStepStatus status;
 21:   Vec          VecCostIntegral0; /* Backup for roll-backs due to events, used by cost integral */
 22:   PetscReal    ptime0;           /* Backup for ts->ptime, the start time of current time step, used by TLM and cost integral */
 23:   PetscReal    time_step0;       /* Backup for ts->timestep, the step size of current time step, used by TLM and cost integral*/

 25:   /* context for sensitivity analysis */
 26:   PetscInt num_tlm;               /* Total number of tangent linear equations */
 27:   Vec     *VecsDeltaLam;          /* Increment of the adjoint sensitivity w.r.t IC at stage */
 28:   Vec     *VecsDeltaMu;           /* Increment of the adjoint sensitivity w.r.t P at stage */
 29:   Vec     *VecsSensiTemp;         /* Vector to be multiplied with Jacobian transpose */
 30:   Mat      MatFwdStages[2];       /* TLM Stages */
 31:   Mat      MatDeltaFwdSensip;     /* Increment of the forward sensitivity at stage */
 32:   Vec      VecDeltaFwdSensipCol;  /* Working vector for holding one column of the sensitivity matrix */
 33:   Mat      MatFwdSensip0;         /* backup for roll-backs due to events */
 34:   Mat      MatIntegralSensipTemp; /* Working vector for forward integral sensitivity */
 35:   Mat      MatIntegralSensip0;    /* backup for roll-backs due to events */
 36:   Vec     *VecsDeltaLam2;         /* Increment of the 2nd-order adjoint sensitivity w.r.t IC at stage */
 37:   Vec     *VecsDeltaMu2;          /* Increment of the 2nd-order adjoint sensitivity w.r.t P at stage */
 38:   Vec     *VecsSensi2Temp;        /* Working vectors that holds the residual for the second-order adjoint */
 39:   Vec     *VecsAffine;            /* Working vectors to store residuals */
 40:   /* context for error estimation */
 41:   Vec vec_sol_prev;
 42:   Vec vec_lte_work;
 43: } TS_Theta;

 45: static PetscErrorCode TSThetaGetX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
 46: {
 47:   TS_Theta *th = (TS_Theta *)ts->data;

 49:   PetscFunctionBegin;
 50:   if (X0) {
 51:     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_X0", X0));
 52:     else *X0 = ts->vec_sol;
 53:   }
 54:   if (Xdot) {
 55:     if (dm && dm != ts->dm) PetscCall(DMGetNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
 56:     else *Xdot = th->Xdot;
 57:   }
 58:   PetscFunctionReturn(PETSC_SUCCESS);
 59: }

 61: static PetscErrorCode TSThetaRestoreX0AndXdot(TS ts, DM dm, Vec *X0, Vec *Xdot)
 62: {
 63:   PetscFunctionBegin;
 64:   if (X0) {
 65:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_X0", X0));
 66:   }
 67:   if (Xdot) {
 68:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSTheta_Xdot", Xdot));
 69:   }
 70:   PetscFunctionReturn(PETSC_SUCCESS);
 71: }

 73: static PetscErrorCode DMCoarsenHook_TSTheta(DM fine, DM coarse, PetscCtx ctx)
 74: {
 75:   PetscFunctionBegin;
 76:   PetscFunctionReturn(PETSC_SUCCESS);
 77: }

 79: static PetscErrorCode DMRestrictHook_TSTheta(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, PetscCtx ctx)
 80: {
 81:   TS  ts = (TS)ctx;
 82:   Vec X0, Xdot, X0_c, Xdot_c;

 84:   PetscFunctionBegin;
 85:   PetscCall(TSThetaGetX0AndXdot(ts, fine, &X0, &Xdot));
 86:   PetscCall(TSThetaGetX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
 87:   PetscCall(MatRestrict(restrct, X0, X0_c));
 88:   PetscCall(MatRestrict(restrct, Xdot, Xdot_c));
 89:   PetscCall(VecPointwiseMult(X0_c, rscale, X0_c));
 90:   PetscCall(VecPointwiseMult(Xdot_c, rscale, Xdot_c));
 91:   PetscCall(TSThetaRestoreX0AndXdot(ts, fine, &X0, &Xdot));
 92:   PetscCall(TSThetaRestoreX0AndXdot(ts, coarse, &X0_c, &Xdot_c));
 93:   PetscFunctionReturn(PETSC_SUCCESS);
 94: }

 96: static PetscErrorCode DMSubDomainHook_TSTheta(DM dm, DM subdm, PetscCtx ctx)
 97: {
 98:   PetscFunctionBegin;
 99:   PetscFunctionReturn(PETSC_SUCCESS);
100: }

102: static PetscErrorCode DMSubDomainRestrictHook_TSTheta(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, PetscCtx ctx)
103: {
104:   TS  ts = (TS)ctx;
105:   Vec X0, Xdot, X0_sub, Xdot_sub;

107:   PetscFunctionBegin;
108:   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
109:   PetscCall(TSThetaGetX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));

111:   PetscCall(VecScatterBegin(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));
112:   PetscCall(VecScatterEnd(gscat, X0, X0_sub, INSERT_VALUES, SCATTER_FORWARD));

114:   PetscCall(VecScatterBegin(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));
115:   PetscCall(VecScatterEnd(gscat, Xdot, Xdot_sub, INSERT_VALUES, SCATTER_FORWARD));

117:   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
118:   PetscCall(TSThetaRestoreX0AndXdot(ts, subdm, &X0_sub, &Xdot_sub));
119:   PetscFunctionReturn(PETSC_SUCCESS);
120: }

122: static PetscErrorCode TSThetaEvaluateCostIntegral(TS ts)
123: {
124:   TS_Theta *th     = (TS_Theta *)ts->data;
125:   TS        quadts = ts->quadraturets;

127:   PetscFunctionBegin;
128:   if (th->endpoint) {
129:     /* Evolve ts->vec_costintegral to compute integrals */
130:     if (th->Theta != 1.0) {
131:       PetscCall(TSComputeRHSFunction(quadts, th->ptime0, th->X0, ts->vec_costintegrand));
132:       PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * (1.0 - th->Theta), ts->vec_costintegrand));
133:     }
134:     PetscCall(TSComputeRHSFunction(quadts, ts->ptime, ts->vec_sol, ts->vec_costintegrand));
135:     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0 * th->Theta, ts->vec_costintegrand));
136:   } else {
137:     PetscCall(TSComputeRHSFunction(quadts, th->stage_time, th->X, ts->vec_costintegrand));
138:     PetscCall(VecAXPY(quadts->vec_sol, th->time_step0, ts->vec_costintegrand));
139:   }
140:   PetscFunctionReturn(PETSC_SUCCESS);
141: }

143: static PetscErrorCode TSForwardCostIntegral_Theta(TS ts)
144: {
145:   TS_Theta *th     = (TS_Theta *)ts->data;
146:   TS        quadts = ts->quadraturets;

148:   PetscFunctionBegin;
149:   /* backup cost integral */
150:   PetscCall(VecCopy(quadts->vec_sol, th->VecCostIntegral0));
151:   PetscCall(TSThetaEvaluateCostIntegral(ts));
152:   PetscFunctionReturn(PETSC_SUCCESS);
153: }

155: static PetscErrorCode TSAdjointCostIntegral_Theta(TS ts)
156: {
157:   TS_Theta *th = (TS_Theta *)ts->data;

159:   PetscFunctionBegin;
160:   /* Like TSForwardCostIntegral(), the adjoint cost integral evaluation relies on ptime0 and time_step0. */
161:   th->ptime0     = ts->ptime + ts->time_step;
162:   th->time_step0 = -ts->time_step;
163:   PetscCall(TSThetaEvaluateCostIntegral(ts));
164:   PetscFunctionReturn(PETSC_SUCCESS);
165: }

167: static PetscErrorCode TSTheta_SNESSolve(TS ts, Vec b, Vec x)
168: {
169:   PetscInt nits, lits;

171:   PetscFunctionBegin;
172:   PetscCall(SNESSolve(ts->snes, b, x));
173:   PetscCall(SNESGetIterationNumber(ts->snes, &nits));
174:   PetscCall(SNESGetLinearSolveIterations(ts->snes, &lits));
175:   ts->snes_its += nits;
176:   ts->ksp_its += lits;
177:   PetscFunctionReturn(PETSC_SUCCESS);
178: }

180: static PetscErrorCode TSResizeRegister_Theta(TS ts, PetscBool reg)
181: {
182:   TS_Theta *th = (TS_Theta *)ts->data;

184:   PetscFunctionBegin;
185:   if (reg) {
186:     PetscCall(TSResizeRegisterVec(ts, "ts:theta:sol_prev", th->vec_sol_prev));
187:     PetscCall(TSResizeRegisterVec(ts, "ts:theta:X0", th->X0));
188:   } else {
189:     PetscCall(TSResizeRetrieveVec(ts, "ts:theta:sol_prev", &th->vec_sol_prev));
190:     PetscCall(PetscObjectReference((PetscObject)th->vec_sol_prev));
191:     PetscCall(TSResizeRetrieveVec(ts, "ts:theta:X0", &th->X0));
192:     PetscCall(PetscObjectReference((PetscObject)th->X0));
193:   }
194:   PetscFunctionReturn(PETSC_SUCCESS);
195: }

197: static PetscErrorCode TSStep_Theta(TS ts)
198: {
199:   TS_Theta *th         = (TS_Theta *)ts->data;
200:   PetscInt  rejections = 0;
201:   PetscBool stageok, accept = PETSC_TRUE;
202:   PetscReal next_time_step = ts->time_step;

204:   PetscFunctionBegin;
205:   if (!ts->steprollback) {
206:     if (th->vec_sol_prev) PetscCall(VecCopy(th->X0, th->vec_sol_prev));
207:     PetscCall(VecCopy(ts->vec_sol, th->X0));
208:   }

210:   th->status = TS_STEP_INCOMPLETE;
211:   while (!ts->reason && th->status != TS_STEP_COMPLETE) {
212:     th->shift      = 1 / (th->Theta * ts->time_step);
213:     th->stage_time = ts->ptime + (th->endpoint ? (PetscReal)1 : th->Theta) * ts->time_step;
214:     PetscCall(VecCopy(th->X0, th->X));
215:     if (th->extrapolate && !ts->steprestart) PetscCall(VecAXPY(th->X, 1 / th->shift, th->Xdot));
216:     if (th->endpoint) { /* This formulation assumes linear time-independent mass matrix */
217:       if (!th->affine) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));
218:       PetscCall(VecZeroEntries(th->Xdot));
219:       PetscCall(TSComputeIFunction(ts, ts->ptime, th->X0, th->Xdot, th->affine, PETSC_FALSE));
220:       PetscCall(VecScale(th->affine, (th->Theta - 1) / th->Theta));
221:     }
222:     PetscCall(TSPreStage(ts, th->stage_time));
223:     PetscCall(TSTheta_SNESSolve(ts, th->endpoint ? th->affine : NULL, th->X));
224:     PetscCall(TSPostStage(ts, th->stage_time, 0, &th->X));
225:     PetscCall(TSAdaptCheckStage(ts->adapt, ts, th->stage_time, th->X, &stageok));
226:     if (!stageok) goto reject_step;

228:     if (th->endpoint) PetscCall(VecCopy(th->X, ts->vec_sol));
229:     else {
230:       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X)); /* th->Xdot is needed by TSInterpolate_Theta */
231:       if (th->Theta == 1.0) PetscCall(VecCopy(th->X, ts->vec_sol));              /* BEULER, stage already checked */
232:       else {
233:         PetscCall(VecAXPY(ts->vec_sol, ts->time_step, th->Xdot));
234:         PetscCall(TSAdaptCheckStage(ts->adapt, ts, ts->ptime + ts->time_step, ts->vec_sol, &stageok));
235:         if (!stageok) {
236:           PetscCall(VecCopy(th->X0, ts->vec_sol));
237:           goto reject_step;
238:         }
239:       }
240:     }

242:     th->status = TS_STEP_PENDING;
243:     PetscCall(TSAdaptChoose(ts->adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
244:     th->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
245:     if (!accept) {
246:       PetscCall(VecCopy(th->X0, ts->vec_sol));
247:       ts->time_step = next_time_step;
248:       goto reject_step;
249:     }

251:     if (ts->forward_solve || ts->costintegralfwd) { /* Save the info for the later use in cost integral evaluation */
252:       th->ptime0     = ts->ptime;
253:       th->time_step0 = ts->time_step;
254:     }
255:     ts->ptime += ts->time_step;
256:     ts->time_step = next_time_step;
257:     break;

259:   reject_step:
260:     ts->reject++;
261:     accept = PETSC_FALSE;
262:     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
263:       ts->reason = TS_DIVERGED_STEP_REJECTED;
264:       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
265:     }
266:   }
267:   PetscFunctionReturn(PETSC_SUCCESS);
268: }

270: static PetscErrorCode TSAdjointStepBEuler_Private(TS ts)
271: {
272:   TS_Theta      *th           = (TS_Theta *)ts->data;
273:   TS             quadts       = ts->quadraturets;
274:   Vec           *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
275:   Vec           *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
276:   Mat            J, Jpre, quadJ = NULL, quadJp = NULL;
277:   KSP            ksp;
278:   PetscScalar   *xarr;
279:   TSEquationType eqtype;
280:   PetscBool      isexplicitode = PETSC_FALSE;
281:   PetscReal      adjoint_time_step;

283:   PetscFunctionBegin;
284:   PetscCall(TSGetEquationType(ts, &eqtype));
285:   if (eqtype == TS_EQ_ODE_EXPLICIT) {
286:     isexplicitode = PETSC_TRUE;
287:     VecsDeltaLam  = ts->vecs_sensi;
288:     VecsDeltaLam2 = ts->vecs_sensi2;
289:   }
290:   th->status = TS_STEP_INCOMPLETE;
291:   PetscCall(SNESGetKSP(ts->snes, &ksp));
292:   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
293:   if (quadts) {
294:     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
295:     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
296:   }

298:   th->stage_time    = ts->ptime;
299:   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */

301:   /* Build RHS for first-order adjoint lambda_{n+1}/h + r_u^T(n+1) */
302:   if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));

304:   for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
305:     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
306:     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / adjoint_time_step)); /* lambda_{n+1}/h */
307:     if (quadJ) {
308:       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
309:       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
310:       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
311:       PetscCall(VecResetArray(ts->vec_drdu_col));
312:       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
313:     }
314:   }

316:   /* Build LHS for first-order adjoint */
317:   th->shift = 1. / adjoint_time_step;
318:   PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
319:   PetscCall(KSPSetOperators(ksp, J, Jpre));

321:   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
322:   for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
323:     KSPConvergedReason kspreason;
324:     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
325:     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
326:     if (kspreason < 0) {
327:       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
328:       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj));
329:     }
330:   }

332:   if (ts->vecs_sensi2) { /* U_{n+1} */
333:     /* Get w1 at t_{n+1} from TLM matrix */
334:     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
335:     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
336:     /* lambda_s^T F_UU w_1 */
337:     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
338:     /* lambda_s^T F_UP w_2 */
339:     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
340:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
341:       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
342:       PetscCall(VecScale(VecsSensi2Temp[nadj], 1. / adjoint_time_step));
343:       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
344:       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
345:     }
346:     /* Solve stage equation LHS X = RHS for second-order adjoint */
347:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
348:       KSPConvergedReason kspreason;
349:       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
350:       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
351:       if (kspreason < 0) {
352:         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
353:         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj));
354:       }
355:     }
356:   }

358:   /* Update sensitivities, and evaluate integrals if there is any */
359:   if (!isexplicitode) {
360:     th->shift = 0.0;
361:     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
362:     PetscCall(KSPSetOperators(ksp, J, Jpre));
363:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
364:       /* Add f_U \lambda_s to the original RHS */
365:       PetscCall(VecScale(VecsSensiTemp[nadj], -1.));
366:       PetscCall(MatMultTransposeAdd(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj], VecsSensiTemp[nadj]));
367:       PetscCall(VecScale(VecsSensiTemp[nadj], -adjoint_time_step));
368:       PetscCall(VecCopy(VecsSensiTemp[nadj], ts->vecs_sensi[nadj]));
369:       if (ts->vecs_sensi2) {
370:         PetscCall(MatMultTransposeAdd(J, VecsDeltaLam2[nadj], VecsSensi2Temp[nadj], VecsSensi2Temp[nadj]));
371:         PetscCall(VecScale(VecsSensi2Temp[nadj], -adjoint_time_step));
372:         PetscCall(VecCopy(VecsSensi2Temp[nadj], ts->vecs_sensi2[nadj]));
373:       }
374:     }
375:   }
376:   if (ts->vecs_sensip) {
377:     PetscCall(VecAXPBYPCZ(th->Xdot, -1. / adjoint_time_step, 1.0 / adjoint_time_step, 0, th->X0, ts->vec_sol));
378:     PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, 1. / adjoint_time_step, ts->Jacp, PETSC_FALSE)); /* get -f_p */
379:     if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
380:     if (ts->vecs_sensi2p) {
381:       /* lambda_s^T F_PU w_1 */
382:       PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
383:       /* lambda_s^T F_PP w_2 */
384:       PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
385:     }

387:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
388:       PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
389:       PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
390:       if (quadJp) {
391:         PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
392:         PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
393:         PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
394:         PetscCall(VecResetArray(ts->vec_drdp_col));
395:         PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
396:       }
397:       if (ts->vecs_sensi2p) {
398:         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
399:         PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, VecsDeltaMu2[nadj]));
400:         if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpu[nadj]));
401:         if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step, ts->vecs_fpp[nadj]));
402:       }
403:     }
404:   }

406:   if (ts->vecs_sensi2) {
407:     PetscCall(VecResetArray(ts->vec_sensip_col));
408:     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
409:   }
410:   th->status = TS_STEP_COMPLETE;
411:   PetscFunctionReturn(PETSC_SUCCESS);
412: }

414: static PetscErrorCode TSAdjointStep_Theta(TS ts)
415: {
416:   TS_Theta    *th           = (TS_Theta *)ts->data;
417:   TS           quadts       = ts->quadraturets;
418:   Vec         *VecsDeltaLam = th->VecsDeltaLam, *VecsDeltaMu = th->VecsDeltaMu, *VecsSensiTemp = th->VecsSensiTemp;
419:   Vec         *VecsDeltaLam2 = th->VecsDeltaLam2, *VecsDeltaMu2 = th->VecsDeltaMu2, *VecsSensi2Temp = th->VecsSensi2Temp;
420:   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
421:   KSP          ksp;
422:   PetscScalar *xarr;
423:   PetscReal    adjoint_time_step;
424:   PetscReal    adjoint_ptime; /* end time of the adjoint time step (ts->ptime is the start time, usually ts->ptime is larger than adjoint_ptime) */

426:   PetscFunctionBegin;
427:   if (th->Theta == 1.) {
428:     PetscCall(TSAdjointStepBEuler_Private(ts));
429:     PetscFunctionReturn(PETSC_SUCCESS);
430:   }
431:   th->status = TS_STEP_INCOMPLETE;
432:   PetscCall(SNESGetKSP(ts->snes, &ksp));
433:   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
434:   if (quadts) {
435:     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
436:     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
437:   }
438:   /* If endpoint=1, th->ptime and th->X0 will be used; if endpoint=0, th->stage_time and th->X will be used. */
439:   th->stage_time    = th->endpoint ? ts->ptime : (ts->ptime + (1. - th->Theta) * ts->time_step);
440:   adjoint_ptime     = ts->ptime + ts->time_step;
441:   adjoint_time_step = -ts->time_step; /* always positive since time_step is negative */

443:   if (!th->endpoint) {
444:     /* recover th->X0 using vec_sol and the stage value th->X */
445:     PetscCall(VecAXPBYPCZ(th->X0, 1.0 / (1.0 - th->Theta), th->Theta / (th->Theta - 1.0), 0, th->X, ts->vec_sol));
446:   }

448:   /* Build RHS for first-order adjoint */
449:   /* Cost function has an integral term */
450:   if (quadts) {
451:     if (th->endpoint) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
452:     else PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
453:   }

455:   for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
456:     PetscCall(VecCopy(ts->vecs_sensi[nadj], VecsSensiTemp[nadj]));
457:     PetscCall(VecScale(VecsSensiTemp[nadj], 1. / (th->Theta * adjoint_time_step)));
458:     if (quadJ) {
459:       PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
460:       PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
461:       PetscCall(VecAXPY(VecsSensiTemp[nadj], 1., ts->vec_drdu_col));
462:       PetscCall(VecResetArray(ts->vec_drdu_col));
463:       PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
464:     }
465:   }

467:   /* Build LHS for first-order adjoint */
468:   th->shift = 1. / (th->Theta * adjoint_time_step);
469:   if (th->endpoint) {
470:     PetscCall(TSComputeSNESJacobian(ts, ts->vec_sol, J, Jpre));
471:   } else {
472:     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre));
473:   }
474:   PetscCall(KSPSetOperators(ksp, J, Jpre));

476:   /* Solve stage equation LHS*lambda_s = RHS for first-order adjoint */
477:   for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
478:     KSPConvergedReason kspreason;
479:     PetscCall(KSPSolveTranspose(ksp, VecsSensiTemp[nadj], VecsDeltaLam[nadj]));
480:     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
481:     if (kspreason < 0) {
482:       ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
483:       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 1st-order adjoint solve\n", ts->steps, nadj));
484:     }
485:   }

487:   /* Second-order adjoint */
488:   if (ts->vecs_sensi2) { /* U_{n+1} */
489:     PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Operation not implemented in TS_Theta");
490:     /* Get w1 at t_{n+1} from TLM matrix */
491:     PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
492:     PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
493:     /* lambda_s^T F_UU w_1 */
494:     PetscCall(TSComputeIHessianProductFunctionUU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
495:     PetscCall(VecResetArray(ts->vec_sensip_col));
496:     PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
497:     /* lambda_s^T F_UP w_2 */
498:     PetscCall(TSComputeIHessianProductFunctionUP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
499:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) { /* compute the residual */
500:       PetscCall(VecCopy(ts->vecs_sensi2[nadj], VecsSensi2Temp[nadj]));
501:       PetscCall(VecScale(VecsSensi2Temp[nadj], th->shift));
502:       PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fuu[nadj]));
503:       if (ts->vecs_fup) PetscCall(VecAXPY(VecsSensi2Temp[nadj], -1., ts->vecs_fup[nadj]));
504:     }
505:     /* Solve stage equation LHS X = RHS for second-order adjoint */
506:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
507:       KSPConvergedReason kspreason;
508:       PetscCall(KSPSolveTranspose(ksp, VecsSensi2Temp[nadj], VecsDeltaLam2[nadj]));
509:       PetscCall(KSPGetConvergedReason(ksp, &kspreason));
510:       if (kspreason < 0) {
511:         ts->reason = TSADJOINT_DIVERGED_LINEAR_SOLVE;
512:         PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th cost function, transposed linear solve fails, stopping 2nd-order adjoint solve\n", ts->steps, nadj));
513:       }
514:     }
515:   }

517:   /* Update sensitivities, and evaluate integrals if there is any */
518:   if (th->endpoint) { /* two-stage Theta methods with th->Theta!=1, th->Theta==1 leads to BEuler */
519:     th->shift      = 1. / ((th->Theta - 1.) * adjoint_time_step);
520:     th->stage_time = adjoint_ptime;
521:     PetscCall(TSComputeSNESJacobian(ts, th->X0, J, Jpre));
522:     PetscCall(KSPSetOperators(ksp, J, Jpre));
523:     /* R_U at t_n */
524:     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, adjoint_ptime, th->X0, quadJ, NULL));
525:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
526:       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], ts->vecs_sensi[nadj]));
527:       if (quadJ) {
528:         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
529:         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
530:         PetscCall(VecAXPY(ts->vecs_sensi[nadj], -1., ts->vec_drdu_col));
531:         PetscCall(VecResetArray(ts->vec_drdu_col));
532:         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
533:       }
534:       PetscCall(VecScale(ts->vecs_sensi[nadj], 1. / th->shift));
535:     }

537:     /* Second-order adjoint */
538:     if (ts->vecs_sensi2) { /* U_n */
539:       /* Get w1 at t_n from TLM matrix */
540:       PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
541:       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
542:       /* lambda_s^T F_UU w_1 */
543:       PetscCall(TSComputeIHessianProductFunctionUU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fuu));
544:       PetscCall(VecResetArray(ts->vec_sensip_col));
545:       PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
546:       /* lambda_s^T F_UU w_2 */
547:       PetscCall(TSComputeIHessianProductFunctionUP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fup));
548:       for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
549:         /* M^T Lambda_s + h(1-theta) F_U^T Lambda_s + h(1-theta) lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2  */
550:         PetscCall(MatMultTranspose(J, VecsDeltaLam2[nadj], ts->vecs_sensi2[nadj]));
551:         PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fuu[nadj]));
552:         if (ts->vecs_fup) PetscCall(VecAXPY(ts->vecs_sensi2[nadj], 1., ts->vecs_fup[nadj]));
553:         PetscCall(VecScale(ts->vecs_sensi2[nadj], 1. / th->shift));
554:       }
555:     }

557:     th->stage_time = ts->ptime; /* recover the old value */

559:     if (ts->vecs_sensip) { /* sensitivities wrt parameters */
560:       /* U_{n+1} */
561:       th->shift = 1.0 / (adjoint_time_step * th->Theta);
562:       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
563:       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, -1. / (th->Theta * adjoint_time_step), ts->Jacp, PETSC_FALSE));
564:       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
565:       for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
566:         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
567:         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu[nadj]));
568:         if (quadJp) {
569:           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
570:           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
571:           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * th->Theta, ts->vec_drdp_col));
572:           PetscCall(VecResetArray(ts->vec_drdp_col));
573:           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
574:         }
575:       }
576:       if (ts->vecs_sensi2p) { /* second-order */
577:         /* Get w1 at t_{n+1} from TLM matrix */
578:         PetscCall(MatDenseGetColumn(ts->mat_sensip, 0, &xarr));
579:         PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
580:         /* lambda_s^T F_PU w_1 */
581:         PetscCall(TSComputeIHessianProductFunctionPU(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
582:         PetscCall(VecResetArray(ts->vec_sensip_col));
583:         PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));

585:         /* lambda_s^T F_PP w_2 */
586:         PetscCall(TSComputeIHessianProductFunctionPP(ts, th->stage_time, ts->vec_sol, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
587:         for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
588:           /* Mu2 <- Mu2 + h theta F_P^T Lambda_s + h theta (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2)  */
589:           PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
590:           PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, VecsDeltaMu2[nadj]));
591:           if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpu[nadj]));
592:           if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * th->Theta, ts->vecs_fpp[nadj]));
593:         }
594:       }

596:       /* U_s */
597:       PetscCall(VecZeroEntries(th->Xdot));
598:       PetscCall(TSComputeIJacobianP(ts, adjoint_ptime, th->X0, th->Xdot, 1. / ((th->Theta - 1.0) * adjoint_time_step), ts->Jacp, PETSC_FALSE));
599:       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, adjoint_ptime, th->X0, quadJp));
600:       for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
601:         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
602:         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu[nadj]));
603:         if (quadJp) {
604:           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
605:           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
606:           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step * (1.0 - th->Theta), ts->vec_drdp_col));
607:           PetscCall(VecResetArray(ts->vec_drdp_col));
608:           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
609:         }
610:         if (ts->vecs_sensi2p) { /* second-order */
611:           /* Get w1 at t_n from TLM matrix */
612:           PetscCall(MatDenseGetColumn(th->MatFwdSensip0, 0, &xarr));
613:           PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
614:           /* lambda_s^T F_PU w_1 */
615:           PetscCall(TSComputeIHessianProductFunctionPU(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_sensip_col, ts->vecs_fpu));
616:           PetscCall(VecResetArray(ts->vec_sensip_col));
617:           PetscCall(MatDenseRestoreColumn(th->MatFwdSensip0, &xarr));
618:           /* lambda_s^T F_PP w_2 */
619:           PetscCall(TSComputeIHessianProductFunctionPP(ts, adjoint_ptime, th->X0, VecsDeltaLam, ts->vec_dir, ts->vecs_fpp));
620:           for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
621:             /* Mu2 <- Mu2 + h(1-theta) F_P^T Lambda_s + h(1-theta) (lambda_s^T F_UU w_1 + lambda_s^T F_UP w_2) */
622:             PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam2[nadj], VecsDeltaMu2[nadj]));
623:             PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), VecsDeltaMu2[nadj]));
624:             if (ts->vecs_fpu) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpu[nadj]));
625:             if (ts->vecs_fpp) PetscCall(VecAXPY(ts->vecs_sensi2p[nadj], -adjoint_time_step * (1.0 - th->Theta), ts->vecs_fpp[nadj]));
626:           }
627:         }
628:       }
629:     }
630:   } else { /* one-stage case */
631:     th->shift = 0.0;
632:     PetscCall(TSComputeSNESJacobian(ts, th->X, J, Jpre)); /* get -f_y */
633:     PetscCall(KSPSetOperators(ksp, J, Jpre));
634:     if (quadts) PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
635:     for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
636:       PetscCall(MatMultTranspose(J, VecsDeltaLam[nadj], VecsSensiTemp[nadj]));
637:       PetscCall(VecAXPY(ts->vecs_sensi[nadj], -adjoint_time_step, VecsSensiTemp[nadj]));
638:       if (quadJ) {
639:         PetscCall(MatDenseGetColumn(quadJ, nadj, &xarr));
640:         PetscCall(VecPlaceArray(ts->vec_drdu_col, xarr));
641:         PetscCall(VecAXPY(ts->vecs_sensi[nadj], adjoint_time_step, ts->vec_drdu_col));
642:         PetscCall(VecResetArray(ts->vec_drdu_col));
643:         PetscCall(MatDenseRestoreColumn(quadJ, &xarr));
644:       }
645:     }
646:     if (ts->vecs_sensip) {
647:       th->shift = 1.0 / (adjoint_time_step * th->Theta);
648:       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
649:       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
650:       if (quadts) PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
651:       for (PetscInt nadj = 0; nadj < ts->numcost; nadj++) {
652:         PetscCall(MatMultTranspose(ts->Jacp, VecsDeltaLam[nadj], VecsDeltaMu[nadj]));
653:         PetscCall(VecAXPY(ts->vecs_sensip[nadj], -adjoint_time_step, VecsDeltaMu[nadj]));
654:         if (quadJp) {
655:           PetscCall(MatDenseGetColumn(quadJp, nadj, &xarr));
656:           PetscCall(VecPlaceArray(ts->vec_drdp_col, xarr));
657:           PetscCall(VecAXPY(ts->vecs_sensip[nadj], adjoint_time_step, ts->vec_drdp_col));
658:           PetscCall(VecResetArray(ts->vec_drdp_col));
659:           PetscCall(MatDenseRestoreColumn(quadJp, &xarr));
660:         }
661:       }
662:     }
663:   }

665:   th->status = TS_STEP_COMPLETE;
666:   PetscFunctionReturn(PETSC_SUCCESS);
667: }

669: static PetscErrorCode TSInterpolate_Theta(TS ts, PetscReal t, Vec X)
670: {
671:   TS_Theta *th = (TS_Theta *)ts->data;
672:   PetscReal dt = t - ts->ptime;

674:   PetscFunctionBegin;
675:   PetscCall(VecCopy(ts->vec_sol, th->X));
676:   if (th->endpoint) dt *= th->Theta;
677:   PetscCall(VecWAXPY(X, dt, th->Xdot, th->X));
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: static PetscErrorCode TSEvaluateWLTE_Theta(TS ts, NormType wnormtype, PetscInt *order, PetscReal *wlte)
682: {
683:   TS_Theta *th = (TS_Theta *)ts->data;
684:   Vec       X  = ts->vec_sol;      /* X = solution */
685:   Vec       Y  = th->vec_lte_work; /* Y = X + LTE  */
686:   PetscReal wltea, wlter;

688:   PetscFunctionBegin;
689:   if (!th->vec_sol_prev) {
690:     *wlte = -1;
691:     PetscFunctionReturn(PETSC_SUCCESS);
692:   }
693:   /* Cannot compute LTE in first step or in restart after event */
694:   if (ts->steprestart) {
695:     *wlte = -1;
696:     PetscFunctionReturn(PETSC_SUCCESS);
697:   }
698:   /* Compute LTE using backward differences with non-constant time step */
699:   {
700:     PetscReal   h = ts->time_step, h_prev = ts->ptime - ts->ptime_prev;
701:     PetscReal   a = 1 + h_prev / h;
702:     PetscScalar scal[3];
703:     Vec         vecs[3];

705:     scal[0] = -1 / a;
706:     scal[1] = +1 / (a - 1);
707:     scal[2] = -1 / (a * (a - 1));
708:     vecs[0] = X;
709:     vecs[1] = th->X0;
710:     vecs[2] = th->vec_sol_prev;
711:     PetscCall(VecCopy(X, Y));
712:     PetscCall(VecMAXPY(Y, 3, scal, vecs));
713:     PetscCall(TSErrorWeightedNorm(ts, X, Y, wnormtype, wlte, &wltea, &wlter));
714:   }
715:   if (order) *order = 2;
716:   PetscFunctionReturn(PETSC_SUCCESS);
717: }

719: static PetscErrorCode TSRollBack_Theta(TS ts)
720: {
721:   TS_Theta *th     = (TS_Theta *)ts->data;
722:   TS        quadts = ts->quadraturets;

724:   PetscFunctionBegin;
725:   if (quadts && ts->costintegralfwd) PetscCall(VecCopy(th->VecCostIntegral0, quadts->vec_sol));
726:   th->status = TS_STEP_INCOMPLETE;
727:   if (ts->mat_sensip) PetscCall(MatCopy(th->MatFwdSensip0, ts->mat_sensip, SAME_NONZERO_PATTERN));
728:   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(th->MatIntegralSensip0, quadts->mat_sensip, SAME_NONZERO_PATTERN));
729:   PetscFunctionReturn(PETSC_SUCCESS);
730: }

732: static PetscErrorCode TSForwardStep_Theta(TS ts)
733: {
734:   TS_Theta    *th                   = (TS_Theta *)ts->data;
735:   TS           quadts               = ts->quadraturets;
736:   Mat          MatDeltaFwdSensip    = th->MatDeltaFwdSensip;
737:   Vec          VecDeltaFwdSensipCol = th->VecDeltaFwdSensipCol;
738:   KSP          ksp;
739:   Mat          J, Jpre, quadJ = NULL, quadJp = NULL;
740:   PetscScalar *barr, *xarr;
741:   PetscReal    previous_shift;

743:   PetscFunctionBegin;
744:   previous_shift = th->shift;
745:   PetscCall(MatCopy(ts->mat_sensip, th->MatFwdSensip0, SAME_NONZERO_PATTERN));

747:   if (quadts && quadts->mat_sensip) PetscCall(MatCopy(quadts->mat_sensip, th->MatIntegralSensip0, SAME_NONZERO_PATTERN));
748:   PetscCall(SNESGetKSP(ts->snes, &ksp));
749:   PetscCall(TSGetIJacobian(ts, &J, &Jpre, NULL, NULL));
750:   if (quadts) {
751:     PetscCall(TSGetRHSJacobian(quadts, &quadJ, NULL, NULL, NULL));
752:     PetscCall(TSGetRHSJacobianP(quadts, &quadJp, NULL, NULL));
753:   }

755:   /* Build RHS */
756:   if (th->endpoint) { /* 2-stage method*/
757:     th->shift = 1. / ((th->Theta - 1.) * th->time_step0);
758:     PetscCall(TSComputeIJacobian(ts, th->ptime0, th->X0, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
759:     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DETERMINE, &MatDeltaFwdSensip));
760:     PetscCall(MatScale(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta));

762:     /* Add the f_p forcing terms */
763:     if (ts->Jacp) {
764:       PetscCall(VecZeroEntries(th->Xdot));
765:       PetscCall(TSComputeIJacobianP(ts, th->ptime0, th->X0, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
766:       PetscCall(MatAXPY(MatDeltaFwdSensip, (th->Theta - 1.) / th->Theta, ts->Jacp, SUBSET_NONZERO_PATTERN));
767:       th->shift = previous_shift;
768:       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, ts->vec_sol));
769:       PetscCall(TSComputeIJacobianP(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
770:       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
771:     }
772:   } else { /* 1-stage method */
773:     th->shift = 0.0;
774:     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
775:     PetscCall(MatMatMult(J, ts->mat_sensip, MAT_REUSE_MATRIX, PETSC_DETERMINE, &MatDeltaFwdSensip));
776:     PetscCall(MatScale(MatDeltaFwdSensip, -1.));

778:     /* Add the f_p forcing terms */
779:     if (ts->Jacp) {
780:       th->shift = previous_shift;
781:       PetscCall(VecAXPBYPCZ(th->Xdot, -th->shift, th->shift, 0, th->X0, th->X));
782:       PetscCall(TSComputeIJacobianP(ts, th->stage_time, th->X, th->Xdot, th->shift, ts->Jacp, PETSC_FALSE));
783:       PetscCall(MatAXPY(MatDeltaFwdSensip, -1., ts->Jacp, SUBSET_NONZERO_PATTERN));
784:     }
785:   }

787:   /* Build LHS */
788:   th->shift = previous_shift; /* recover the previous shift used in TSStep_Theta() */
789:   if (th->endpoint) {
790:     PetscCall(TSComputeIJacobian(ts, th->stage_time, ts->vec_sol, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
791:   } else {
792:     PetscCall(TSComputeIJacobian(ts, th->stage_time, th->X, th->Xdot, th->shift, J, Jpre, PETSC_FALSE));
793:   }
794:   PetscCall(KSPSetOperators(ksp, J, Jpre));

796:   /*
797:     Evaluate the first stage of integral gradients with the 2-stage method:
798:     drdu|t_n*S(t_n) + drdp|t_n
799:     This is done before the linear solve because the sensitivity variable S(t_n) will be propagated to S(t_{n+1})
800:   */
801:   if (th->endpoint) { /* 2-stage method only */
802:     if (quadts && quadts->mat_sensip) {
803:       PetscCall(TSComputeRHSJacobian(quadts, th->ptime0, th->X0, quadJ, NULL));
804:       PetscCall(TSComputeRHSJacobianP(quadts, th->ptime0, th->X0, quadJp));
805:       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp));
806:       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
807:       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * (1. - th->Theta), th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
808:     }
809:   }

811:   /* Solve the tangent linear equation for forward sensitivities to parameters */
812:   for (PetscInt ntlm = 0; ntlm < th->num_tlm; ntlm++) {
813:     KSPConvergedReason kspreason;
814:     PetscCall(MatDenseGetColumn(MatDeltaFwdSensip, ntlm, &barr));
815:     PetscCall(VecPlaceArray(VecDeltaFwdSensipCol, barr));
816:     if (th->endpoint) {
817:       PetscCall(MatDenseGetColumn(ts->mat_sensip, ntlm, &xarr));
818:       PetscCall(VecPlaceArray(ts->vec_sensip_col, xarr));
819:       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, ts->vec_sensip_col));
820:       PetscCall(VecResetArray(ts->vec_sensip_col));
821:       PetscCall(MatDenseRestoreColumn(ts->mat_sensip, &xarr));
822:     } else {
823:       PetscCall(KSPSolve(ksp, VecDeltaFwdSensipCol, VecDeltaFwdSensipCol));
824:     }
825:     PetscCall(KSPGetConvergedReason(ksp, &kspreason));
826:     if (kspreason < 0) {
827:       ts->reason = TSFORWARD_DIVERGED_LINEAR_SOLVE;
828:       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", %" PetscInt_FMT "th tangent linear solve, linear solve fails, stopping tangent linear solve\n", ts->steps, ntlm));
829:     }
830:     PetscCall(VecResetArray(VecDeltaFwdSensipCol));
831:     PetscCall(MatDenseRestoreColumn(MatDeltaFwdSensip, &barr));
832:   }

834:   /*
835:     Evaluate the second stage of integral gradients with the 2-stage method:
836:     drdu|t_{n+1}*S(t_{n+1}) + drdp|t_{n+1}
837:   */
838:   if (quadts && quadts->mat_sensip) {
839:     if (!th->endpoint) {
840:       PetscCall(MatAXPY(ts->mat_sensip, 1, MatDeltaFwdSensip, SAME_NONZERO_PATTERN)); /* stage sensitivity */
841:       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, th->X, quadJ, NULL));
842:       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, th->X, quadJp));
843:       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp));
844:       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
845:       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
846:       PetscCall(MatAXPY(ts->mat_sensip, (1. - th->Theta) / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
847:     } else {
848:       PetscCall(TSComputeRHSJacobian(quadts, th->stage_time, ts->vec_sol, quadJ, NULL));
849:       PetscCall(TSComputeRHSJacobianP(quadts, th->stage_time, ts->vec_sol, quadJp));
850:       PetscCall(MatTransposeMatMult(ts->mat_sensip, quadJ, MAT_REUSE_MATRIX, PETSC_DETERMINE, &th->MatIntegralSensipTemp));
851:       PetscCall(MatAXPY(th->MatIntegralSensipTemp, 1, quadJp, SAME_NONZERO_PATTERN));
852:       PetscCall(MatAXPY(quadts->mat_sensip, th->time_step0 * th->Theta, th->MatIntegralSensipTemp, SAME_NONZERO_PATTERN));
853:     }
854:   } else {
855:     if (!th->endpoint) PetscCall(MatAXPY(ts->mat_sensip, 1. / th->Theta, MatDeltaFwdSensip, SAME_NONZERO_PATTERN));
856:   }
857:   PetscFunctionReturn(PETSC_SUCCESS);
858: }

860: static PetscErrorCode TSForwardGetStages_Theta(TS ts, PetscInt *ns, Mat *stagesensip[])
861: {
862:   TS_Theta *th = (TS_Theta *)ts->data;

864:   PetscFunctionBegin;
865:   if (ns) {
866:     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
867:     else *ns = 2;                                   /* endpoint form */
868:   }
869:   if (stagesensip) {
870:     if (!th->endpoint && th->Theta != 1.0) {
871:       th->MatFwdStages[0] = th->MatDeltaFwdSensip;
872:     } else {
873:       th->MatFwdStages[0] = th->MatFwdSensip0;
874:       th->MatFwdStages[1] = ts->mat_sensip; /* stiffly accurate */
875:     }
876:     *stagesensip = th->MatFwdStages;
877:   }
878:   PetscFunctionReturn(PETSC_SUCCESS);
879: }

881: static PetscErrorCode TSReset_Theta(TS ts)
882: {
883:   TS_Theta *th = (TS_Theta *)ts->data;

885:   PetscFunctionBegin;
886:   PetscCall(VecDestroy(&th->X));
887:   PetscCall(VecDestroy(&th->Xdot));
888:   PetscCall(VecDestroy(&th->X0));
889:   PetscCall(VecDestroy(&th->affine));

891:   PetscCall(VecDestroy(&th->vec_sol_prev));
892:   PetscCall(VecDestroy(&th->vec_lte_work));

894:   PetscCall(VecDestroy(&th->VecCostIntegral0));
895:   PetscFunctionReturn(PETSC_SUCCESS);
896: }

898: static PetscErrorCode TSAdjointReset_Theta(TS ts)
899: {
900:   TS_Theta *th = (TS_Theta *)ts->data;

902:   PetscFunctionBegin;
903:   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam));
904:   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu));
905:   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaLam2));
906:   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsDeltaMu2));
907:   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensiTemp));
908:   PetscCall(VecDestroyVecs(ts->numcost, &th->VecsSensi2Temp));
909:   PetscFunctionReturn(PETSC_SUCCESS);
910: }

912: static PetscErrorCode TSDestroy_Theta(TS ts)
913: {
914:   PetscFunctionBegin;
915:   PetscCall(TSReset_Theta(ts));
916:   if (ts->dm) {
917:     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
918:     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));
919:   }
920:   PetscCall(PetscFree(ts->data));
921:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", NULL));
922:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", NULL));
923:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", NULL));
924:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", NULL));
925:   PetscFunctionReturn(PETSC_SUCCESS);
926: }

928: /*
929:   This defines the nonlinear equation that is to be solved with SNES
930:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0

932:   Note that U here is the stage argument. This means that U = U_{n+1} only if endpoint = true,
933:   otherwise U = theta U_{n+1} + (1 - theta) U0, which for the case of implicit midpoint is
934:   U = (U_{n+1} + U0)/2
935: */
936: static PetscErrorCode SNESTSFormFunction_Theta(SNES snes, Vec x, Vec y, TS ts)
937: {
938:   TS_Theta *th = (TS_Theta *)ts->data;
939:   Vec       X0, Xdot;
940:   DM        dm, dmsave;
941:   PetscReal shift = th->shift;

943:   PetscFunctionBegin;
944:   PetscCall(SNESGetDM(snes, &dm));
945:   /* When using the endpoint variant, this is actually 1/Theta * Xdot */
946:   PetscCall(TSThetaGetX0AndXdot(ts, dm, &X0, &Xdot));
947:   if (x != X0) PetscCall(VecAXPBYPCZ(Xdot, -shift, shift, 0, X0, x));
948:   else PetscCall(VecZeroEntries(Xdot));
949:   /* DM monkey-business allows user code to call TSGetDM() inside of functions evaluated on levels of FAS */
950:   dmsave = ts->dm;
951:   ts->dm = dm;
952:   PetscCall(TSComputeIFunction(ts, th->stage_time, x, Xdot, y, PETSC_FALSE));
953:   ts->dm = dmsave;
954:   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, &X0, &Xdot));
955:   PetscFunctionReturn(PETSC_SUCCESS);
956: }

958: static PetscErrorCode SNESTSFormJacobian_Theta(SNES snes, Vec x, Mat A, Mat B, TS ts)
959: {
960:   TS_Theta *th = (TS_Theta *)ts->data;
961:   Vec       Xdot;
962:   DM        dm, dmsave;
963:   PetscReal shift = th->shift;

965:   PetscFunctionBegin;
966:   PetscCall(SNESGetDM(snes, &dm));
967:   /* Xdot has already been computed in SNESTSFormFunction_Theta (SNES guarantees this) */
968:   PetscCall(TSThetaGetX0AndXdot(ts, dm, NULL, &Xdot));

970:   dmsave = ts->dm;
971:   ts->dm = dm;
972:   PetscCall(TSComputeIJacobian(ts, th->stage_time, x, Xdot, shift, A, B, PETSC_FALSE));
973:   ts->dm = dmsave;
974:   PetscCall(TSThetaRestoreX0AndXdot(ts, dm, NULL, &Xdot));
975:   PetscFunctionReturn(PETSC_SUCCESS);
976: }

978: static PetscErrorCode TSForwardSetUp_Theta(TS ts)
979: {
980:   TS_Theta *th     = (TS_Theta *)ts->data;
981:   TS        quadts = ts->quadraturets;

983:   PetscFunctionBegin;
984:   /* combine sensitivities to parameters and sensitivities to initial values into one array */
985:   th->num_tlm = ts->num_parameters;
986:   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatDeltaFwdSensip));
987:   if (quadts && quadts->mat_sensip) {
988:     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensipTemp));
989:     PetscCall(MatDuplicate(quadts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatIntegralSensip0));
990:   }
991:   /* backup sensitivity results for roll-backs */
992:   PetscCall(MatDuplicate(ts->mat_sensip, MAT_DO_NOT_COPY_VALUES, &th->MatFwdSensip0));

994:   PetscCall(VecDuplicate(ts->vec_sol, &th->VecDeltaFwdSensipCol));
995:   PetscFunctionReturn(PETSC_SUCCESS);
996: }

998: static PetscErrorCode TSForwardReset_Theta(TS ts)
999: {
1000:   TS_Theta *th     = (TS_Theta *)ts->data;
1001:   TS        quadts = ts->quadraturets;

1003:   PetscFunctionBegin;
1004:   if (quadts && quadts->mat_sensip) {
1005:     PetscCall(MatDestroy(&th->MatIntegralSensipTemp));
1006:     PetscCall(MatDestroy(&th->MatIntegralSensip0));
1007:   }
1008:   PetscCall(VecDestroy(&th->VecDeltaFwdSensipCol));
1009:   PetscCall(MatDestroy(&th->MatDeltaFwdSensip));
1010:   PetscCall(MatDestroy(&th->MatFwdSensip0));
1011:   PetscFunctionReturn(PETSC_SUCCESS);
1012: }

1014: static PetscErrorCode TSSetUp_Theta(TS ts)
1015: {
1016:   TS_Theta *th     = (TS_Theta *)ts->data;
1017:   TS        quadts = ts->quadraturets;
1018:   PetscBool match;

1020:   PetscFunctionBegin;
1021:   if (!th->VecCostIntegral0 && quadts && ts->costintegralfwd) { /* back up cost integral */
1022:     PetscCall(VecDuplicate(quadts->vec_sol, &th->VecCostIntegral0));
1023:   }
1024:   if (!th->X) PetscCall(VecDuplicate(ts->vec_sol, &th->X));
1025:   if (!th->Xdot) PetscCall(VecDuplicate(ts->vec_sol, &th->Xdot));
1026:   if (!th->X0) PetscCall(VecDuplicate(ts->vec_sol, &th->X0));
1027:   if (th->endpoint) PetscCall(VecDuplicate(ts->vec_sol, &th->affine));

1029:   th->order = (th->Theta == 0.5) ? 2 : 1;
1030:   th->shift = 1 / (th->Theta * ts->time_step);

1032:   PetscCall(TSGetDM(ts, &ts->dm));
1033:   PetscCall(DMCoarsenHookAdd(ts->dm, DMCoarsenHook_TSTheta, DMRestrictHook_TSTheta, ts));
1034:   PetscCall(DMSubDomainHookAdd(ts->dm, DMSubDomainHook_TSTheta, DMSubDomainRestrictHook_TSTheta, ts));

1036:   PetscCall(TSGetAdapt(ts, &ts->adapt));
1037:   PetscCall(TSAdaptCandidatesClear(ts->adapt));
1038:   PetscCall(PetscObjectTypeCompare((PetscObject)ts->adapt, TSADAPTNONE, &match));
1039:   if (!match) {
1040:     if (!th->vec_sol_prev) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_sol_prev));
1041:     if (!th->vec_lte_work) PetscCall(VecDuplicate(ts->vec_sol, &th->vec_lte_work));
1042:   }
1043:   PetscCall(TSGetSNES(ts, &ts->snes));

1045:   ts->stifflyaccurate = (!th->endpoint && th->Theta != 1.0) ? PETSC_FALSE : PETSC_TRUE;
1046:   PetscFunctionReturn(PETSC_SUCCESS);
1047: }

1049: static PetscErrorCode TSAdjointSetUp_Theta(TS ts)
1050: {
1051:   TS_Theta *th = (TS_Theta *)ts->data;

1053:   PetscFunctionBegin;
1054:   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam));
1055:   PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsSensiTemp));
1056:   if (ts->vecs_sensip) PetscCall(VecDuplicateVecs(ts->vecs_sensip[0], ts->numcost, &th->VecsDeltaMu));
1057:   if (ts->vecs_sensi2) {
1058:     PetscCall(VecDuplicateVecs(ts->vecs_sensi[0], ts->numcost, &th->VecsDeltaLam2));
1059:     PetscCall(VecDuplicateVecs(ts->vecs_sensi2[0], ts->numcost, &th->VecsSensi2Temp));
1060:     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1061:     if (!ts->ihessianproduct_fuu) ts->vecs_fuu = ts->vecs_guu;
1062:     if (!ts->ihessianproduct_fup) ts->vecs_fup = ts->vecs_gup;
1063:   }
1064:   if (ts->vecs_sensi2p) {
1065:     PetscCall(VecDuplicateVecs(ts->vecs_sensi2p[0], ts->numcost, &th->VecsDeltaMu2));
1066:     /* hack ts to make implicit TS solver work when users provide only explicit versions of callbacks (RHSFunction,RHSJacobian,RHSHessian etc.) */
1067:     if (!ts->ihessianproduct_fpu) ts->vecs_fpu = ts->vecs_gpu;
1068:     if (!ts->ihessianproduct_fpp) ts->vecs_fpp = ts->vecs_gpp;
1069:   }
1070:   PetscFunctionReturn(PETSC_SUCCESS);
1071: }

1073: static PetscErrorCode TSSetFromOptions_Theta(TS ts, PetscOptionItems PetscOptionsObject)
1074: {
1075:   TS_Theta *th = (TS_Theta *)ts->data;

1077:   PetscFunctionBegin;
1078:   PetscOptionsHeadBegin(PetscOptionsObject, "Theta ODE solver options");
1079:   {
1080:     PetscCall(PetscOptionsReal("-ts_theta_theta", "Location of stage (0<Theta<=1)", "TSThetaSetTheta", th->Theta, &th->Theta, NULL));
1081:     PetscCall(PetscOptionsBool("-ts_theta_endpoint", "Use the endpoint instead of midpoint form of the Theta method", "TSThetaSetEndpoint", th->endpoint, &th->endpoint, NULL));
1082:     PetscCall(PetscOptionsBool("-ts_theta_initial_guess_extrapolate", "Extrapolate stage initial guess from previous solution (sometimes unstable)", "TSThetaSetExtrapolate", th->extrapolate, &th->extrapolate, NULL));
1083:   }
1084:   PetscOptionsHeadEnd();
1085:   PetscFunctionReturn(PETSC_SUCCESS);
1086: }

1088: static PetscErrorCode TSView_Theta(TS ts, PetscViewer viewer)
1089: {
1090:   TS_Theta *th = (TS_Theta *)ts->data;
1091:   PetscBool isascii;

1093:   PetscFunctionBegin;
1094:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
1095:   if (isascii) {
1096:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Theta=%g\n", (double)th->Theta));
1097:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Extrapolation=%s\n", th->extrapolate ? "yes" : "no"));
1098:   }
1099:   PetscFunctionReturn(PETSC_SUCCESS);
1100: }

1102: static PetscErrorCode TSThetaGetTheta_Theta(TS ts, PetscReal *theta)
1103: {
1104:   TS_Theta *th = (TS_Theta *)ts->data;

1106:   PetscFunctionBegin;
1107:   *theta = th->Theta;
1108:   PetscFunctionReturn(PETSC_SUCCESS);
1109: }

1111: static PetscErrorCode TSThetaSetTheta_Theta(TS ts, PetscReal theta)
1112: {
1113:   TS_Theta *th = (TS_Theta *)ts->data;

1115:   PetscFunctionBegin;
1116:   PetscCheck(theta > 0 && theta <= 1, PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_OUTOFRANGE, "Theta %g not in range (0,1]", (double)theta);
1117:   th->Theta = theta;
1118:   th->order = (th->Theta == 0.5) ? 2 : 1;
1119:   PetscFunctionReturn(PETSC_SUCCESS);
1120: }

1122: static PetscErrorCode TSThetaGetEndpoint_Theta(TS ts, PetscBool *endpoint)
1123: {
1124:   TS_Theta *th = (TS_Theta *)ts->data;

1126:   PetscFunctionBegin;
1127:   *endpoint = th->endpoint;
1128:   PetscFunctionReturn(PETSC_SUCCESS);
1129: }

1131: static PetscErrorCode TSThetaSetEndpoint_Theta(TS ts, PetscBool flg)
1132: {
1133:   TS_Theta *th = (TS_Theta *)ts->data;

1135:   PetscFunctionBegin;
1136:   th->endpoint = flg;
1137:   PetscFunctionReturn(PETSC_SUCCESS);
1138: }

1140: #if defined(PETSC_HAVE_COMPLEX)
1141: static PetscErrorCode TSComputeLinearStability_Theta(TS ts, PetscReal xr, PetscReal xi, PetscReal *yr, PetscReal *yi)
1142: {
1143:   PetscComplex z  = xr + xi * PETSC_i, f;
1144:   TS_Theta    *th = (TS_Theta *)ts->data;

1146:   PetscFunctionBegin;
1147:   f   = (1.0 + (1.0 - th->Theta) * z) / (1.0 - th->Theta * z);
1148:   *yr = PetscRealPartComplex(f);
1149:   *yi = PetscImaginaryPartComplex(f);
1150:   PetscFunctionReturn(PETSC_SUCCESS);
1151: }
1152: #endif

1154: static PetscErrorCode TSGetStages_Theta(TS ts, PetscInt *ns, Vec *Y[])
1155: {
1156:   TS_Theta *th = (TS_Theta *)ts->data;

1158:   PetscFunctionBegin;
1159:   if (ns) {
1160:     if (!th->endpoint && th->Theta != 1.0) *ns = 1; /* midpoint form */
1161:     else *ns = 2;                                   /* endpoint form */
1162:   }
1163:   if (Y) {
1164:     if (!th->endpoint && th->Theta != 1.0) {
1165:       th->Stages[0] = th->X;
1166:     } else {
1167:       th->Stages[0] = th->X0;
1168:       th->Stages[1] = ts->vec_sol; /* stiffly accurate */
1169:     }
1170:     *Y = th->Stages;
1171:   }
1172:   PetscFunctionReturn(PETSC_SUCCESS);
1173: }

1175: /*MC
1176:       TSTHETA - DAE solver using the implicit Theta method

1178:    Level: beginner

1180:    Options Database Keys:
1181: +  -ts_theta_theta Theta                   - Location of stage (0<Theta<=1)
1182: .  -ts_theta_endpoint flag                 - Use the endpoint (like Crank-Nicholson) instead of midpoint form of the Theta method
1183: -  -ts_theta_initial_guess_extrapolate flg - Extrapolate stage initial guess from previous solution (sometimes unstable)

1185:    Notes:
1186: .vb
1187:   -ts_type theta -ts_theta_theta 1.0 corresponds to backward Euler (TSBEULER)
1188:   -ts_type theta -ts_theta_theta 0.5 corresponds to the implicit midpoint rule
1189:   -ts_type theta -ts_theta_theta 0.5 -ts_theta_endpoint corresponds to Crank-Nicholson (TSCN)
1190: .ve

1192:    The endpoint variant of the Theta method and backward Euler can be applied to DAE. The midpoint variant is not suitable for DAEs because it is not stiffly accurate.

1194:    The midpoint variant is cast as a 1-stage implicit Runge-Kutta method.

1196: .vb
1197:   Theta | Theta
1198:   -------------
1199:         |  1
1200: .ve

1202:    For the default Theta=0.5, this is also known as the implicit midpoint rule.

1204:    When the endpoint variant is chosen, the method becomes a 2-stage method with first stage explicit:

1206: .vb
1207:   0 | 0         0
1208:   1 | 1-Theta   Theta
1209:   -------------------
1210:     | 1-Theta   Theta
1211: .ve

1213:    For the default Theta=0.5, this is the trapezoid rule (also known as Crank-Nicolson, see TSCN).

1215:    To apply a diagonally implicit RK method to DAE, the stage formula
1216: .vb
1217:   Y_i = X + h sum_j a_ij Y'_j
1218: .ve
1219:    is interpreted as a formula for Y'_i in terms of Y_i and known values (Y'_j, j<i)

1221: .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSCN`, `TSBEULER`, `TSThetaSetTheta()`, `TSThetaSetEndpoint()`
1222: M*/
1223: PETSC_EXTERN PetscErrorCode TSCreate_Theta(TS ts)
1224: {
1225:   TS_Theta *th;

1227:   PetscFunctionBegin;
1228:   ts->ops->reset          = TSReset_Theta;
1229:   ts->ops->adjointreset   = TSAdjointReset_Theta;
1230:   ts->ops->destroy        = TSDestroy_Theta;
1231:   ts->ops->view           = TSView_Theta;
1232:   ts->ops->setup          = TSSetUp_Theta;
1233:   ts->ops->adjointsetup   = TSAdjointSetUp_Theta;
1234:   ts->ops->adjointreset   = TSAdjointReset_Theta;
1235:   ts->ops->step           = TSStep_Theta;
1236:   ts->ops->interpolate    = TSInterpolate_Theta;
1237:   ts->ops->evaluatewlte   = TSEvaluateWLTE_Theta;
1238:   ts->ops->rollback       = TSRollBack_Theta;
1239:   ts->ops->resizeregister = TSResizeRegister_Theta;
1240:   ts->ops->setfromoptions = TSSetFromOptions_Theta;
1241:   ts->ops->snesfunction   = SNESTSFormFunction_Theta;
1242:   ts->ops->snesjacobian   = SNESTSFormJacobian_Theta;
1243: #if defined(PETSC_HAVE_COMPLEX)
1244:   ts->ops->linearstability = TSComputeLinearStability_Theta;
1245: #endif
1246:   ts->ops->getstages       = TSGetStages_Theta;
1247:   ts->ops->adjointstep     = TSAdjointStep_Theta;
1248:   ts->ops->adjointintegral = TSAdjointCostIntegral_Theta;
1249:   ts->ops->forwardintegral = TSForwardCostIntegral_Theta;
1250:   ts->default_adapt_type   = TSADAPTNONE;

1252:   ts->ops->forwardsetup     = TSForwardSetUp_Theta;
1253:   ts->ops->forwardreset     = TSForwardReset_Theta;
1254:   ts->ops->forwardstep      = TSForwardStep_Theta;
1255:   ts->ops->forwardgetstages = TSForwardGetStages_Theta;

1257:   ts->usessnes = PETSC_TRUE;

1259:   PetscCall(PetscNew(&th));
1260:   ts->data = (void *)th;

1262:   th->VecsDeltaLam   = NULL;
1263:   th->VecsDeltaMu    = NULL;
1264:   th->VecsSensiTemp  = NULL;
1265:   th->VecsSensi2Temp = NULL;

1267:   th->extrapolate = PETSC_FALSE;
1268:   th->Theta       = 0.5;
1269:   th->order       = 2;
1270:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetTheta_C", TSThetaGetTheta_Theta));
1271:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetTheta_C", TSThetaSetTheta_Theta));
1272:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaGetEndpoint_C", TSThetaGetEndpoint_Theta));
1273:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSThetaSetEndpoint_C", TSThetaSetEndpoint_Theta));
1274:   PetscFunctionReturn(PETSC_SUCCESS);
1275: }

1277: /*@
1278:   TSThetaGetTheta - Get the abscissa of the stage in (0,1] for `TSTHETA`

1280:   Not Collective

1282:   Input Parameter:
1283: . ts - timestepping context

1285:   Output Parameter:
1286: . theta - stage abscissa

1288:   Level: advanced

1290:   Note:
1291:   Use of this function is normally only required to hack `TSTHETA` to use a modified integration scheme.

1293: .seealso: [](ch_ts), `TSThetaSetTheta()`, `TSTHETA`
1294: @*/
1295: PetscErrorCode TSThetaGetTheta(TS ts, PetscReal *theta)
1296: {
1297:   PetscFunctionBegin;
1299:   PetscAssertPointer(theta, 2);
1300:   PetscUseMethod(ts, "TSThetaGetTheta_C", (TS, PetscReal *), (ts, theta));
1301:   PetscFunctionReturn(PETSC_SUCCESS);
1302: }

1304: /*@
1305:   TSThetaSetTheta - Set the abscissa of the stage in (0,1]  for `TSTHETA`

1307:   Not Collective

1309:   Input Parameters:
1310: + ts    - timestepping context
1311: - theta - stage abscissa

1313:   Options Database Key:
1314: . -ts_theta_theta theta - set theta

1316:   Level: intermediate

1318: .seealso: [](ch_ts), `TSThetaGetTheta()`, `TSTHETA`, `TSCN`
1319: @*/
1320: PetscErrorCode TSThetaSetTheta(TS ts, PetscReal theta)
1321: {
1322:   PetscFunctionBegin;
1324:   PetscTryMethod(ts, "TSThetaSetTheta_C", (TS, PetscReal), (ts, theta));
1325:   PetscFunctionReturn(PETSC_SUCCESS);
1326: }

1328: /*@
1329:   TSThetaGetEndpoint - Gets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`

1331:   Not Collective

1333:   Input Parameter:
1334: . ts - timestepping context

1336:   Output Parameter:
1337: . endpoint - `PETSC_TRUE` when using the endpoint variant

1339:   Level: advanced

1341: .seealso: [](ch_ts), `TSThetaSetEndpoint()`, `TSTHETA`, `TSCN`
1342: @*/
1343: PetscErrorCode TSThetaGetEndpoint(TS ts, PetscBool *endpoint)
1344: {
1345:   PetscFunctionBegin;
1347:   PetscAssertPointer(endpoint, 2);
1348:   PetscUseMethod(ts, "TSThetaGetEndpoint_C", (TS, PetscBool *), (ts, endpoint));
1349:   PetscFunctionReturn(PETSC_SUCCESS);
1350: }

1352: /*@
1353:   TSThetaSetEndpoint - Sets whether to use the endpoint variant of the method (e.g. trapezoid/Crank-Nicolson instead of midpoint rule) for `TSTHETA`

1355:   Not Collective

1357:   Input Parameters:
1358: + ts  - timestepping context
1359: - flg - `PETSC_TRUE` to use the endpoint variant

1361:   Options Database Key:
1362: . -ts_theta_endpoint flg - use the endpoint variant

1364:   Level: intermediate

1366: .seealso: [](ch_ts), `TSTHETA`, `TSCN`
1367: @*/
1368: PetscErrorCode TSThetaSetEndpoint(TS ts, PetscBool flg)
1369: {
1370:   PetscFunctionBegin;
1372:   PetscTryMethod(ts, "TSThetaSetEndpoint_C", (TS, PetscBool), (ts, flg));
1373:   PetscFunctionReturn(PETSC_SUCCESS);
1374: }

1376: /*
1377:  * TSBEULER and TSCN are straightforward specializations of TSTHETA.
1378:  * The creation functions for these specializations are below.
1379:  */

1381: static PetscErrorCode TSSetUp_BEuler(TS ts)
1382: {
1383:   TS_Theta *th = (TS_Theta *)ts->data;

1385:   PetscFunctionBegin;
1386:   PetscCheck(th->Theta == 1.0, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change the default value (1) of theta when using backward Euler");
1387:   PetscCheck(!th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change to the endpoint form of the Theta methods when using backward Euler");
1388:   PetscCall(TSSetUp_Theta(ts));
1389:   PetscFunctionReturn(PETSC_SUCCESS);
1390: }

1392: static PetscErrorCode TSView_BEuler(TS ts, PetscViewer viewer)
1393: {
1394:   PetscFunctionBegin;
1395:   PetscFunctionReturn(PETSC_SUCCESS);
1396: }

1398: /*MC
1399:       TSBEULER - ODE solver using the implicit backward Euler method

1401:   Level: beginner

1403:   Note:
1404:   `TSBEULER` is equivalent to `TSTHETA` with Theta=1.0 or `-ts_type theta -ts_theta_theta 1.0`

1406: .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSEULER`, `TSCN`, `TSTHETA`
1407: M*/
1408: PETSC_EXTERN PetscErrorCode TSCreate_BEuler(TS ts)
1409: {
1410:   PetscFunctionBegin;
1411:   PetscCall(TSCreate_Theta(ts));
1412:   PetscCall(TSThetaSetTheta(ts, 1.0));
1413:   PetscCall(TSThetaSetEndpoint(ts, PETSC_FALSE));
1414:   ts->ops->setup = TSSetUp_BEuler;
1415:   ts->ops->view  = TSView_BEuler;
1416:   PetscFunctionReturn(PETSC_SUCCESS);
1417: }

1419: static PetscErrorCode TSSetUp_CN(TS ts)
1420: {
1421:   TS_Theta *th = (TS_Theta *)ts->data;

1423:   PetscFunctionBegin;
1424:   PetscCheck(th->Theta == 0.5, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change the default value (0.5) of theta when using Crank-Nicolson");
1425:   PetscCheck(th->endpoint, PetscObjectComm((PetscObject)ts), PETSC_ERR_OPT_OVERWRITE, "Can not change to the midpoint form of the Theta methods when using Crank-Nicolson");
1426:   PetscCall(TSSetUp_Theta(ts));
1427:   PetscFunctionReturn(PETSC_SUCCESS);
1428: }

1430: static PetscErrorCode TSView_CN(TS ts, PetscViewer viewer)
1431: {
1432:   PetscFunctionBegin;
1433:   PetscFunctionReturn(PETSC_SUCCESS);
1434: }

1436: /*MC
1437:   TSCN - ODE solver using the implicit Crank-Nicolson method.

1439:   Level: beginner

1441:   Notes:
1442:   `TSCN` is equivalent to `TSTHETA` with Theta=0.5 and the "endpoint" option set. I.e.
1443: .vb
1444:   -ts_type theta
1445:   -ts_theta_theta 0.5
1446:   -ts_theta_endpoint
1447: .ve

1449: .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSBEULER`, `TSTHETA`, `TSType`
1450: M*/
1451: PETSC_EXTERN PetscErrorCode TSCreate_CN(TS ts)
1452: {
1453:   PetscFunctionBegin;
1454:   PetscCall(TSCreate_Theta(ts));
1455:   PetscCall(TSThetaSetTheta(ts, 0.5));
1456:   PetscCall(TSThetaSetEndpoint(ts, PETSC_TRUE));
1457:   ts->ops->setup = TSSetUp_CN;
1458:   ts->ops->view  = TSView_CN;
1459:   PetscFunctionReturn(PETSC_SUCCESS);
1460: }