![]() |
OGS
|
Definition at line 49 of file HeatConductionFEM.h.
#include <HeatConductionFEM.h>
Public Member Functions | |
LocalAssemblerData (MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, 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_x_prev, std::vector< double > &local_rhs_data, std::vector< double > &local_Jac_data) override |
Eigen::Map< const Eigen::RowVectorXd > | getShapeMatrix (const unsigned integration_point) const override |
Provides the shape matrix at the given integration point. | |
std::vector< double > const & | getIntPtHeatFlux (const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override |
![]() | |
virtual | ~LocalAssemblerInterface ()=default |
virtual void | setInitialConditions (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, double const t, 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_x_prev, 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_x_prev, int const process_id, 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_prev, 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, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) |
void | postNonLinearSolver (std::size_t const mesh_item_id, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< GlobalVector * > const &x, std::vector< GlobalVector * > const &x_prev, double const t, double const dt, int const process_id) |
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. | |
virtual std::optional< VectorSegment > | getVectorDeformationSegment () const |
![]() | |
virtual | ~ExtrapolatableElement ()=default |
Private Types | |
using | ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim> |
using | ShapeMatrices = typename ShapeMatricesType::ShapeMatrices |
using | LocalAssemblerTraits |
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 |
NumLib::GenericIntegrationMethod const & | _integration_method |
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > | _shape_matrices |
|
private |
Definition at line 59 of file HeatConductionFEM.h.
|
private |
Definition at line 54 of file HeatConductionFEM.h.
|
private |
Definition at line 57 of file HeatConductionFEM.h.
|
private |
Definition at line 58 of file HeatConductionFEM.h.
|
private |
Definition at line 52 of file HeatConductionFEM.h.
|
private |
Definition at line 51 of file HeatConductionFEM.h.
|
inline |
The thermal_conductivity factor is directly integrated into the local element matrix.
Definition at line 64 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::NUM_NODAL_DOF.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 84 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::interpolateCoordinates(), ProcessLib::HeatConduction::HeatConductionProcessData::mass_lumping, ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ProcessLib::HeatConduction::NUM_NODAL_DOF, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 152 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), ProcessLib::HeatConduction::HeatConductionProcessData::mass_lumping, ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ProcessLib::HeatConduction::NUM_NODAL_DOF, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.
|
inlineoverridevirtual |
Implements ProcessLib::HeatConduction::HeatConductionLocalAssemblerInterface.
Definition at line 238 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), NumLib::getIndices(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), ProcessLib::HeatConduction::HeatConductionProcessData::media_map, ParameterLib::SpatialPosition::setElementID(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.
|
inlineoverridevirtual |
Provides the shape matrix at the given integration point.
Implements NumLib::ExtrapolatableElement.
Definition at line 229 of file HeatConductionFEM.h.
References ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices.
|
private |
Definition at line 290 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleWithJacobian(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtHeatFlux().
|
private |
Definition at line 293 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleWithJacobian(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtHeatFlux().
|
private |
Definition at line 291 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleWithJacobian(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtHeatFlux().
|
private |
Definition at line 295 of file HeatConductionFEM.h.
Referenced by ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::assembleWithJacobian(), ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtHeatFlux(), and ProcessLib::HeatConduction::LocalAssemblerData< ShapeFunction, GlobalDim >::getShapeMatrix().