Tempus Version of the Day
Time Integration
Loading...
Searching...
No Matches
Tempus_StepperRKModifierXBase.hpp
Go to the documentation of this file.
1//@HEADER
2// *****************************************************************************
3// Tempus: Time Integration and Sensitivity Analysis Package
4//
5// Copyright 2017 NTESS and the Tempus contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8//@HEADER
9
10#ifndef Tempus_StepperRKModifierXBase_hpp
11#define Tempus_StepperRKModifierXBase_hpp
12
13#include "Tempus_config.hpp"
14#include "Tempus_SolutionHistory.hpp"
16
17#include "Teuchos_SerialDenseVector.hpp"
18
19namespace Tempus {
20
39template <class Scalar>
41 : virtual public Tempus::StepperRKAppAction<Scalar> {
42 private:
43 /* \brief Adaptor execute function
44 *
45 * This is an adaptor function to bridge between the AppAction
46 * interface and the ModifierX interface. It is meant to be private
47 * and non-virtual as deriving from this class should only need to
48 * implement the modify function.
49 *
50 * For the ModifierX interface, this adaptor maps the
51 * StepperRKAppAction::ACTION_LOCATION to the
52 * StepperRKModifierX::MODIFIERX_TYPE, and only pass the solution
53 * (\f$x\f$ and/or \f$\dot{x}\f$ and other parameters to the modify
54 * function.
55 */
56 void execute(
57 Teuchos::RCP<SolutionHistory<Scalar> > sh,
58 Teuchos::RCP<StepperRKBase<Scalar> > stepper,
60 {
61 using Teuchos::RCP;
62
64 const int stageNumber = stepper->getStageNumber();
65 Teuchos::SerialDenseVector<int, Scalar> c = stepper->getTableau()->c();
66 RCP<SolutionState<Scalar> > workingState = sh->getWorkingState();
67 const Scalar dt = workingState->getTimeStep();
68 Scalar time = sh->getCurrentState()->getTime();
69 if (stageNumber >= 0) time += c(stageNumber) * dt;
70 RCP<Thyra::VectorBase<Scalar> > x = workingState->getX();
71
72 switch (actLoc) {
74 modType = X_BEGIN_STEP;
75 break;
76 }
78 modType = X_BEGIN_STAGE;
79 break;
80 }
82 modType = X_BEFORE_SOLVE;
83 break;
84 }
86 modType = X_AFTER_SOLVE;
87 break;
88 }
90 modType = X_BEFORE_EXPLICIT_EVAL;
91 break;
92 }
94 modType = X_END_STAGE;
95 break;
96 }
98 modType = X_END_STEP;
99 time = workingState->getTime();
100 break;
101 }
102 default:
103 TEUCHOS_TEST_FOR_EXCEPTION(true, std::logic_error,
104 "Error - unknown action location.\n");
105 }
106
107 this->modify(x, time, dt, stageNumber, modType);
108 }
109
110 public:
121
123 virtual void modify(Teuchos::RCP<Thyra::VectorBase<Scalar> > /* x */,
124 const Scalar /* time */, const Scalar /* dt */,
125 const int /* stageNumber */,
126 const MODIFIER_TYPE modType) = 0;
127};
128
129} // namespace Tempus
130
131#endif // Tempus_StepperRKModifierXBase_hpp
SolutionHistory is basically a container of SolutionStates. SolutionHistory maintains a collection of...
Application Action for StepperRKBase.
ACTION_LOCATION
Indicates the location of application action (see algorithm).
Base class for Runge-Kutta methods, ExplicitRK, DIRK and IMEX.
MODIFIER_TYPE
Indicates the location of application action (see algorithm).
@ X_END_STEP
Modify at the end of the step.
@ X_BEGIN_STEP
Modify at the beginning of the step.
@ X_AFTER_SOLVE
Modify after the implicit solve.
@ X_END_STAGE
Modify at the end of the stage.
@ X_BEFORE_SOLVE
Modify before the implicit solve.
@ X_BEFORE_EXPLICIT_EVAL
Modify before the explicit evaluation.
@ X_BEGIN_STAGE
Modify at the beginning of the stage.
virtual void modify(Teuchos::RCP< Thyra::VectorBase< Scalar > >, const Scalar, const Scalar, const int, const MODIFIER_TYPE modType)=0
Modify solution based on the MODIFIER_TYPE.
void execute(Teuchos::RCP< SolutionHistory< Scalar > > sh, Teuchos::RCP< StepperRKBase< Scalar > > stepper, const typename StepperRKAppAction< Scalar >::ACTION_LOCATION actLoc)
Execute application action for RK Stepper.