28template <
typename GlobalDimNodalMatrixType>
32 double const& integration_weight_)
37 GlobalDimNodalMatrixType
const dNdx;
52 std::vector<GlobalVector*>
const& x,
53 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_tables,
54 std::vector<double>& cache)
const = 0;
57template <
typename ShapeFunction,
int GlobalDim>
79 bool const is_axially_symmetric,
85 unsigned const n_integration_points =
87 _ip_data.reserve(n_integration_points);
89 auto const& shape_matrices =
91 GlobalDim>(element, is_axially_symmetric,
97 double const aperture_size =
102 for (
unsigned ip = 0; ip < n_integration_points; ip++)
105 shape_matrices[ip].dNdx,
107 shape_matrices[ip].integralMeasure *
108 shape_matrices[ip].detJ * aperture_size);
112 void assemble(
double const t,
double const dt,
113 std::vector<double>
const& local_x,
114 std::vector<double>
const& ,
115 std::vector<double>& local_M_data,
116 std::vector<double>& local_K_data,
117 std::vector<double>& local_b_data)
override;
123 std::vector<double>
const& local_x)
const override;
126 const unsigned integration_point)
const override
130 .NsHigherOrder<
typename ShapeFunction::MeshElement>();
133 return Eigen::Map<const Eigen::RowVectorXd>(
134 N[integration_point].data(), N[integration_point].size());
139 std::vector<GlobalVector*>
const& x,
140 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
141 std::vector<double>& velocity_cache)
const override;
147 std::vector<IntegrationPointData<GlobalDimNodalMatrixType>>
_ip_data;
156 Eigen::Map<NodalMatrixType>& local_K,
157 Eigen::Map<NodalVectorType>& local_b,
160 double const mu,
double const rho_L,
162 bool const has_gravity);
165 Eigen::Map<const NodalVectorType>
const& local_p,
169 bool const has_gravity);
179 Eigen::Map<NodalMatrixType>& local_K,
180 Eigen::Map<NodalVectorType>& local_b,
183 double const mu,
double const rho_L,
185 bool const has_gravity);
188 Eigen::Map<const NodalVectorType>
const& local_p,
192 bool const has_gravity);
195 template <
typename LaplacianGravityVelocityCalculator>
197 std::vector<double>
const& local_x,
198 std::vector<double>& local_M_data,
199 std::vector<double>& local_K_data,
200 std::vector<double>& local_b_data);
202 template <
typename LaplacianGravityVelocityCalculator,
203 typename VelocityCacheType>
205 const double t,
const double dt, std::vector<double>
const& local_x,
207 VelocityCacheType& darcy_velocity_at_ips)
const;
209 template <
typename VelocityCacheType>
212 std::vector<double>
const& local_x,
214 VelocityCacheType& darcy_velocity_at_ips)
const;
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType
void setElementID(std::size_t element_id)
virtual std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_tables, std::vector< double > &cache) const =0
Eigen::Vector3d getFlux(MathLib::Point3d const &p_local_coords, double const t, std::vector< double > const &local_x) const override
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
typename LocalAssemblerTraits::LocalMatrix NodalMatrixType
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
LiquidFlowLocalAssembler(MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LiquidFlowData const &process_data)
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
const LiquidFlowData & _process_data
ShapeMatrixPolicyType< ShapeFunction, GlobalDim > ShapeMatricesType
void computeProjectedDarcyVelocity(const double t, const double dt, std::vector< double > const &local_x, ParameterLib::SpatialPosition const &pos, VelocityCacheType &darcy_velocity_at_ips) const
typename LocalAssemblerTraits::LocalVector NodalVectorType
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 > &local_b_data) override
MeshLib::Element const & _element
NumLib::GenericIntegrationMethod const & _integration_method
ProcessLib::LocalAssemblerTraits< ShapeMatricesType, ShapeFunction::NPOINTS, NUM_NODAL_DOF, GlobalDim > LocalAssemblerTraits
void assembleMatrixAndVector(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data)
void computeDarcyVelocity(bool const is_scalar_permeability, const double t, const double dt, std::vector< double > const &local_x, ParameterLib::SpatialPosition const &pos, VelocityCacheType &darcy_velocity_at_ips) const
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
std::vector< double > const & getIntPtDarcyVelocity(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &velocity_cache) const override
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > initShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, IntegrationMethod const &integration_method)
const unsigned NUM_NODAL_DOF
detail::LocalAssemblerTraitsFixed< ShpPol, NNodes, NodalDOF, Dim > LocalAssemblerTraits
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
GlobalDimNodalMatrixType const dNdx
EIGEN_MAKE_ALIGNED_OPERATOR_NEW
IntegrationPointData(GlobalDimNodalMatrixType const &dNdx_, double const &integration_weight_)
double const integration_weight
static void calculateLaplacianAndGravityTerm(Eigen::Map< NodalMatrixType > &local_K, Eigen::Map< NodalVectorType > &local_b, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability_with_density_factor, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static Eigen::Matrix< double, GlobalDim, 1 > calculateVelocity(Eigen::Map< const NodalVectorType > const &local_p, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static void calculateLaplacianAndGravityTerm(Eigen::Map< NodalMatrixType > &local_K, Eigen::Map< NodalVectorType > &local_b, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability_with_density_factor, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
static Eigen::Matrix< double, GlobalDim, 1 > calculateVelocity(Eigen::Map< const NodalVectorType > const &local_p, IntegrationPointData< GlobalDimNodalMatrixType > const &ip_data, GlobalDimMatrixType const &permeability, double const mu, double const rho_L, GlobalDimVectorType const &specific_body_force, bool const has_gravity)
Vector< NNodes *NodalDOF > LocalVector
Matrix< NNodes *NodalDOF, NNodes *NodalDOF > LocalMatrix