ug4
|
Time Discretization Interface. More...
#include <time_disc_interface.h>
Public Types | |
typedef TAlgebra | algebra_type |
Algebra type. More... | |
typedef IDomainDiscretization< TAlgebra > | domain_discretization_type |
Domain Discretization type. More... | |
typedef algebra_type::vector_type | vector_type |
Vector type. More... | |
Public Types inherited from ug::IAssemble< TAlgebra > | |
typedef TAlgebra | algebra_type |
Algebra type. More... | |
typedef TAlgebra::matrix_type | matrix_type |
Type of algebra matrix. More... | |
typedef TAlgebra::vector_type | vector_type |
Type of algebra vector. More... | |
Public Member Functions | |
virtual SmartPtr< IConstraint< TAlgebra > > | constraint (size_t i) |
returns the i'th constraint More... | |
SmartPtr< IDomainDiscretization< TAlgebra > > | domain_disc () const |
return underlying domain disc More... | |
virtual void | finish_step (SmartPtr< VectorTimeSeries< vector_type > > currSol)=0 |
virtual number | future_time () const =0 |
returns the future time point (i.e. the one that will be computed) More... | |
ITimeDiscretization (SmartPtr< IDomainDiscretization< TAlgebra > > spDD) | |
create and set domain discretization More... | |
virtual size_t | num_constraints () const |
returns the number of constraint More... | |
virtual size_t | num_prev_steps () const =0 |
returns number of previous time steps needed More... | |
virtual size_t | num_stages () const =0 |
returns the number of stages More... | |
virtual void | prepare_step (SmartPtr< VectorTimeSeries< vector_type > > prevSol, number dt)=0 |
prepares the assembling of Defect/Jacobian for a time step More... | |
virtual void | set_stage (size_t stage)=0 |
sets the stage More... | |
virtual void | prepare_step_elem (SmartPtr< VectorTimeSeries< vector_type > > prevSol, number dt, const GridLevel &gl)=0 |
prepares the assembling of Defect/Jacobian for a time step More... | |
void | prepare_step_elem (SmartPtr< VectorTimeSeries< vector_type > > prevSol, number dt) |
prepares the assembling of Defect/Jacobian for a time step More... | |
virtual void | finish_step_elem (SmartPtr< VectorTimeSeries< vector_type > > currSol, const GridLevel &gl)=0 |
void | finish_step_elem (SmartPtr< VectorTimeSeries< vector_type > > currSol) |
Public Member Functions inherited from ug::IAssemble< TAlgebra > | |
void | adjust_solution (vector_type &u) |
virtual void | adjust_solution (vector_type &u, const GridLevel &gl)=0 |
sets dirichlet values in solution vector More... | |
void | assemble_defect (vector_type &d, const vector_type &u) |
virtual void | assemble_defect (vector_type &d, const vector_type &u, const GridLevel &gl)=0 |
assembles Defect More... | |
void | assemble_jacobian (matrix_type &J, const vector_type &u) |
virtual void | assemble_jacobian (matrix_type &J, const vector_type &u, const GridLevel &gl)=0 |
assembles Jacobian (or Approximation of Jacobian) More... | |
void | assemble_linear (matrix_type &A, vector_type &b) |
virtual void | assemble_linear (matrix_type &A, vector_type &b, const GridLevel &gl)=0 |
Assembles Matrix and Right-Hand-Side for a linear problem. More... | |
void | assemble_mass_matrix (matrix_type &M, const vector_type &u) |
virtual void | assemble_mass_matrix (matrix_type &M, const vector_type &u, const GridLevel &gl) |
assembles mass matrix More... | |
void | assemble_rhs (vector_type &b) |
virtual void | assemble_rhs (vector_type &b, const GridLevel &gl)=0 |
Assembles Right-Hand-Side for a linear problem. More... | |
virtual void | assemble_rhs (vector_type &rhs, const vector_type &u) |
virtual void | assemble_rhs (vector_type &rhs, const vector_type &u, const GridLevel &gl)=0 |
assembles rhs More... | |
void | assemble_stiffness_matrix (matrix_type &A, const vector_type &u) |
virtual void | assemble_stiffness_matrix (matrix_type &A, const vector_type &u, const GridLevel &gl) |
assembles stiffness matrix More... | |
virtual | ~IAssemble () |
Virtual Destructor. More... | |
Protected Attributes | |
SmartPtr< IDomainDiscretization< TAlgebra > > | m_spDomDisc |
Domain Discretization. More... | |
Private Member Functions | |
SmartPtr< AssemblingTuner< TAlgebra > > | ass_tuner () |
ConstSmartPtr< AssemblingTuner< TAlgebra > > | ass_tuner () const |
Time Discretization Interface.
Defines the time discretization interface.
This class uses a ISpatialDiscretization in order to implement the IAssemble interface.
After the method prepare step has been called, Jacobian/Defect can be computed.
TAlgebra | Algebra Type |
typedef TAlgebra ug::ITimeDiscretization< TAlgebra >::algebra_type |
Algebra type.
typedef IDomainDiscretization<TAlgebra> ug::ITimeDiscretization< TAlgebra >::domain_discretization_type |
Domain Discretization type.
typedef algebra_type::vector_type ug::ITimeDiscretization< TAlgebra >::vector_type |
Vector type.
|
inline |
create and set domain discretization
[in] | dd | Domain Discretization |
|
inlineprivatevirtual |
returns the assemble adapter
Implements ug::IAssemble< TAlgebra >.
References ug::ITimeDiscretization< TAlgebra >::m_spDomDisc.
|
inlineprivatevirtual |
returns the assemble adapter
Implements ug::IAssemble< TAlgebra >.
References ug::ITimeDiscretization< TAlgebra >::m_spDomDisc.
|
inlinevirtual |
returns the i'th constraint
Implements ug::IAssemble< TAlgebra >.
References ug::ITimeDiscretization< TAlgebra >::m_spDomDisc.
|
inline |
return underlying domain disc
References ug::ITimeDiscretization< TAlgebra >::m_spDomDisc.
|
pure virtual |
finishes a time step and allows to adapt data depending on the current solution elemDisc-wise This function is called after the assembling routines at the end of a time step.
[in] | currSol | the solution at the previous time steps |
[in] | dt | size of time step |
Implemented in ug::MultiStepTimeDiscretization< TAlgebra >.
|
inline |
finishes a time step and allows to adapt data depending on the current solution element-wise This function is called after the assembling routines at the end of a time step. Within this function "fsh_timestep_elem" is called which allows modifying data depending on the current solution at element-level.
[in] | currSol | the current solution |
[in] | dd | DoF Distribution |
References ug::ITimeDiscretization< TAlgebra >::finish_step_elem().
|
pure virtual |
finishes a time step and allows to adapt data depending on the current solution element-wise This function is called after the assembling routines at the end of a time step. Within this function "fsh_timestep_elem" is called which allows modifying data depending on the current solution at element-level.
[in] | currSol | the current solution |
[in] | dd | DoF Distribution |
Implemented in ug::MultiStepTimeDiscretization< TAlgebra >.
Referenced by ug::ITimeDiscretization< TAlgebra >::finish_step_elem().
|
pure virtual |
returns the future time point (i.e. the one that will be computed)
Implemented in ug::MultiStepTimeDiscretization< TAlgebra >, and ug::LinearImplicitEuler< class >.
|
inlinevirtual |
returns the number of constraint
Implements ug::IAssemble< TAlgebra >.
References ug::ITimeDiscretization< TAlgebra >::m_spDomDisc.
|
pure virtual |
returns number of previous time steps needed
Implemented in ug::MultiStepTimeDiscretization< TAlgebra >, and ug::LinearImplicitEuler< class >.
|
pure virtual |
returns the number of stages
Implemented in ug::SDIRK< TAlgebra >, ug::BDF< TAlgebra >, ug::ThetaTimeStep< TAlgebra >, and ug::LinearImplicitEuler< class >.
|
pure virtual |
prepares the assembling of Defect/Jacobian for a time step
This function supplies the TimeDiscretization with previous time steps and step size before the assembling routines can be called.
[in] | prevSol | the solution at the previous time steps |
[in] | dt | size of time step |
Implemented in ug::SDIRK< TAlgebra >, and ug::MultiStepTimeDiscretization< TAlgebra >.
|
inline |
prepares the assembling of Defect/Jacobian for a time step
This function supplies the TimeDiscretization with previous time steps and step size before the assembling routines can be called. A sub-routine at element-level ("prep_timestep_elem") is called within this function.
[in] | prevSol | the solution at the previous time steps |
[in] | dt | size of time step |
[in] | dd | DoF Distribution |
References ug::ITimeDiscretization< TAlgebra >::prepare_step_elem().
|
pure virtual |
prepares the assembling of Defect/Jacobian for a time step
This function supplies the TimeDiscretization with previous time steps and step size before the assembling routines can be called. A sub-routine at element-level ("prep_timestep_elem") is called within this function.
[in] | prevSol | the solution at the previous time steps |
[in] | dt | size of time step |
[in] | dd | DoF Distribution |
Implemented in ug::MultiStepTimeDiscretization< TAlgebra >.
Referenced by ug::ITimeDiscretization< TAlgebra >::prepare_step_elem().
|
pure virtual |
sets the stage
Implemented in ug::SDIRK< TAlgebra >, ug::BDF< TAlgebra >, ug::ThetaTimeStep< TAlgebra >, and ug::LinearImplicitEuler< class >.
|
protected |
Domain Discretization.
Referenced by ug::ITimeDiscretization< TAlgebra >::ass_tuner(), ug::ITimeDiscretization< TAlgebra >::constraint(), ug::ITimeDiscretization< TAlgebra >::domain_disc(), ug::MultiStepTimeDiscretization< TAlgebra >::invalidate_error(), ug::MultiStepTimeDiscretization< TAlgebra >::is_error_valid(), and ug::ITimeDiscretization< TAlgebra >::num_constraints().