![]() |
OGS
|
Definition at line 50 of file SteadyStateDiffusionFEM.h.
#include <SteadyStateDiffusionFEM.h>
Public Member Functions | |
LocalAssemblerData (MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, SteadyStateDiffusionData const &process_data) | |
void | assemble (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &, std::vector< double > &, std::vector< double > &local_K_data, std::vector< double > &) override |
Eigen::Vector3d | getFlux (MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const 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 & | getIntPtDarcyVelocity (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 | assembleWithJacobian (double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_b_data, std::vector< double > &local_Jac_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< 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 |
SteadyStateDiffusionData const & | _process_data |
NumLib::GenericIntegrationMethod const & | _integration_method |
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > | _shape_matrices |
|
private |
Definition at line 60 of file SteadyStateDiffusionFEM.h.
|
private |
Definition at line 55 of file SteadyStateDiffusionFEM.h.
|
private |
Definition at line 58 of file SteadyStateDiffusionFEM.h.
|
private |
Definition at line 59 of file SteadyStateDiffusionFEM.h.
|
private |
Definition at line 53 of file SteadyStateDiffusionFEM.h.
|
private |
Definition at line 52 of file SteadyStateDiffusionFEM.h.
|
inline |
The hydraulic_conductivity factor is directly integrated into the local element matrix.
Definition at line 65 of file SteadyStateDiffusionFEM.h.
|
inlineoverridevirtual |
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 81 of file SteadyStateDiffusionFEM.h.
References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::diffusion, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), NumLib::GenericIntegrationMethod::getWeightedPoint(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, ProcessLib::SteadyStateDiffusion::NUM_NODAL_DOF, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.
|
inlineoverridevirtual |
Computes the flux in the point p_local_coords
that is given in local coordinates using the values from local_x
.
Reimplemented from ProcessLib::LocalAssemblerInterface.
Definition at line 130 of file SteadyStateDiffusionFEM.h.
References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, NumLib::computeShapeMatrices(), MaterialPropertyLib::diffusion, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::VariableArray::temperature.
|
inlineoverridevirtual |
Implements ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionLocalAssemblerInterface.
Definition at line 184 of file SteadyStateDiffusionFEM.h.
References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_element, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_process_data, ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices, MathLib::createZeroedMatrix(), MaterialPropertyLib::diffusion, MaterialPropertyLib::formEigenTensor(), MeshLib::Element::getID(), NumLib::getIndices(), MaterialPropertyLib::MaterialSpatialDistributionMap::getMedium(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, ProcessLib::SteadyStateDiffusion::SteadyStateDiffusionData::media_map, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MathLib::toVector().
|
inlineoverridevirtual |
Provides the shape matrix at the given integration point.
Implements NumLib::ExtrapolatableElement.
Definition at line 175 of file SteadyStateDiffusionFEM.h.
References ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::_shape_matrices.
|
private |
Definition at line 242 of file SteadyStateDiffusionFEM.h.
Referenced by ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getFlux(), and ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity().
|
private |
|
private |
Definition at line 243 of file SteadyStateDiffusionFEM.h.
Referenced by ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getFlux(), and ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity().
|
private |
Definition at line 247 of file SteadyStateDiffusionFEM.h.
Referenced by ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::assemble(), ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity(), and ProcessLib::SteadyStateDiffusion::LocalAssemblerData< ShapeFunction, GlobalDim >::getShapeMatrix().