20 #ifndef LIBMESH_RB_TEMPORAL_DISCRETIZATION_H 21 #define LIBMESH_RB_TEMPORAL_DISCRETIZATION_H 24 #include "libmesh/libmesh_common.h" 87 void set_control(
const std::vector<Real> & control);
135 #endif // LIBMESH_RB_TEMPORAL_DISCRETIZATION_H RBTemporalDiscretization & operator=(const RBTemporalDiscretization &)=default
std::vector< Real > _control
The RHS control (scalar function of time).
void set_delta_t(const Real delta_t_in)
void process_temporal_parameters_file(const std::string ¶meters_filename)
Read in and initialize parameters from parameters_filename.
Real _delta_t
The time-step size.
The libMesh namespace provides an interface to certain functionality in the library.
Real get_delta_t() const
Get/set delta_t, the time-step size.
unsigned int _current_time_step
The current time-step.
Real get_euler_theta() const
Get/set euler_theta, parameter that determines the temporal discretization.
void set_time_step(const unsigned int k)
void pull_temporal_discretization_data(RBTemporalDiscretization &other)
Pull the temporal discretization data from other.
virtual ~RBTemporalDiscretization()=default
unsigned int get_time_step() const
Get/set the current time-step.
RBTemporalDiscretization()
Constructor.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void set_euler_theta(const Real euler_theta_in)
unsigned int get_n_time_steps() const
Get/set the total number of time-steps.
Define a class that encapsulates the details of a "generalized Euler" temporal discretization to be u...
void set_n_time_steps(const unsigned int K)
unsigned int _n_time_steps
The number of time-steps.
Real get_control(const unsigned int k) const
Get/set the RHS control.
void set_control(const std::vector< Real > &control)
Real _euler_theta
The parameter that determines the generalized Euler scheme discretization that we employ...