OGS
|
Definition at line 58 of file HeatConductionFEM.h.
#include <HeatConductionFEM.h>
Public Member Functions | |
LocalAssemblerData (MeshLib::Element const &element, std::size_t const local_matrix_size, bool is_axially_symmetric, unsigned const integration_order, HeatConductionProcessData const &process_data) | |
void | assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &) override |
void | assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_xdot, const double, const double, std::vector< double > &, std::vector< double > &, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override |
void | computeSecondaryVariableConcrete (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &) override |
Eigen::Map< const Eigen::RowVectorXd > | getShapeMatrix (const unsigned integration_point) const override |
Provides the shape matrix at the given integration point. More... | |
std::vector< double > const & | getIntPtHeatFluxX (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override |
std::vector< double > const & | getIntPtHeatFluxY (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override |
std::vector< double > const & | getIntPtHeatFluxZ (const double, std::vector< GlobalVector * > const &, std::vector< NumLib::LocalToGlobalIndexMap const * > const &, std::vector< double > &) const override |
Public Member Functions inherited from ProcessLib::LocalAssemblerInterface | |
virtual | ~LocalAssemblerInterface ()=default |
void | setInitialConditions (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, bool const use_monolithic_scheme, int const process_id) |
virtual void | initialize (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table) |
virtual void | preAssemble (double const, double const, std::vector< double > const &) |
virtual void | assembleForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) |
virtual void | assembleWithJacobianForStaggeredScheme (double const t, double const dt, Eigen::VectorXd const &local_x, Eigen::VectorXd const &local_xdot, const double dxdot_dx, const double dx_dx, int const process_id, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data, std::vector< double > &local_Jac_data) |
virtual void | computeSecondaryVariable (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, double const t, double const dt, std::vector< GlobalVector * > const &x, GlobalVector const &x_dot, int const process_id) |
virtual void | preTimestep (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, double const t, double const delta_t) |
virtual void | postTimestep (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, double const dt) |
void | postNonLinearSolver (std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table, GlobalVector const &x, GlobalVector const &xdot, double const t, double const dt, bool const use_monolithic_scheme, int const process_id) |
virtual std::vector< double > | interpolateNodalValuesToIntegrationPoints (std::vector< double > const &) |
virtual Eigen::Vector3d | getFlux (MathLib::Point3d const &, double const, std::vector< double > const &) const |
virtual Eigen::Vector3d | getFlux (MathLib::Point3d const &, double const, std::vector< std::vector< double >> const &) const |
Fits to staggered scheme. More... | |
Public Member Functions inherited from NumLib::ExtrapolatableElement | |
virtual | ~ExtrapolatableElement ()=default |
Private Types | |
using | ShapeMatricesType = ShapeMatrixPolicyType< ShapeFunction, GlobalDim > |
using | ShapeMatrices = typename ShapeMatricesType::ShapeMatrices |
using | LocalAssemblerTraits = ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim > |
using | NodalMatrixType = typename LocalAssemblerTraits::LocalMatrix |
using | NodalVectorType = typename LocalAssemblerTraits::LocalVector |
using | GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType |
Private Attributes | |
MeshLib::Element const & | _element |
HeatConductionProcessData const & | _process_data |
IntegrationMethod const | _integration_method |
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > | _shape_matrices |
std::vector< std::vector< double > > | _heat_fluxes |
|
private |
Definition at line 68 of file HeatConductionFEM.h.
|
private |
Definition at line 63 of file HeatConductionFEM.h.
|
private |
Definition at line 66 of file HeatConductionFEM.h.
|
private |
Definition at line 67 of file HeatConductionFEM.h.
|
private |
Definition at line 61 of file HeatConductionFEM.h.
|
private |
Definition at line 60 of file HeatConductionFEM.h.
|
inline |
The thermal_conductivity factor is directly integrated into the local element matrix.
Definition at line 73 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::NUM_NODAL_DOF.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 95 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_shape_matrices, MaterialPropertyLib::density, MeshLib::Element::getID(), MaterialPropertyLib::heat_capacity, ProcessLib::HeatConduction::HeatConductionProcessData::mass_lumping, ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ProcessLib::HeatConduction::NUM_NODAL_DOF, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::temperature, and MaterialPropertyLib::thermal_conductivity.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 160 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_shape_matrices, MaterialPropertyLib::density, MeshLib::Element::getID(), MaterialPropertyLib::heat_capacity, ProcessLib::HeatConduction::HeatConductionProcessData::mass_lumping, ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ProcessLib::HeatConduction::NUM_NODAL_DOF, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::storage, MaterialPropertyLib::temperature, and MaterialPropertyLib::thermal_conductivity.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 242 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_heat_fluxes, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_shape_matrices, MeshLib::Element::getID(), ProcessLib::HeatConduction::HeatConductionProcessData::media_map, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), MaterialPropertyLib::temperature, and MaterialPropertyLib::thermal_conductivity.
|
inlineoverridevirtual |
Implements ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface.
Definition at line 289 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_heat_fluxes.
|
inlineoverridevirtual |
Implements ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface.
Definition at line 299 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_heat_fluxes.
|
inlineoverridevirtual |
Implements ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface.
Definition at line 309 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::_heat_fluxes.
|
inlineoverridevirtual |
Provides the shape matrix at the given integration point.
Implements NumLib::ExtrapolatableElement.
Definition at line 280 of file HeatConductionFEM.h.
|
private |
Definition at line 320 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleWithJacobian(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete().
|
private |
Definition at line 327 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getIntPtHeatFluxX(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getIntPtHeatFluxY(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getIntPtHeatFluxZ().
|
private |
Definition at line 323 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleWithJacobian(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete().
|
private |
Definition at line 321 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleWithJacobian(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete().
|
private |
Definition at line 325 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::assembleWithJacobian(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::computeSecondaryVariableConcrete(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, IntegrationMethod, GlobalDim >::getShapeMatrix().