31namespace RichardsComponentTransport
33template <
typename NodalRowVectorType,
typename GlobalDimNodalMatrixType,
34 typename NodalMatrixType>
38 GlobalDimNodalMatrixType
const& dNdx_,
39 double const& integration_weight_,
40 NodalMatrixType
const mass_operator_)
48 NodalRowVectorType
const N;
49 GlobalDimNodalMatrixType
const dNdx;
65 std::vector<GlobalVector*>
const& x,
66 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
67 std::vector<double>& cache)
const = 0;
71 std::vector<GlobalVector*>
const& x,
72 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
73 std::vector<double>& cache)
const = 0;
76template <
typename ShapeFunction,
int GlobalDim>
87 typename ShapeMatricesType::template VectorType<
NUM_NODAL_DOF *
88 ShapeFunction::NPOINTS>;
102 std::size_t
const local_matrix_size,
104 bool is_axially_symmetric,
108 void assemble(
double const t,
double const dt,
109 std::vector<double>
const& local_x,
110 std::vector<double>
const& local_x_prev,
111 std::vector<double>& local_M_data,
112 std::vector<double>& local_K_data,
113 std::vector<double>& local_b_data)
override;
117 std::vector<GlobalVector*>
const& x,
118 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
119 std::vector<double>& cache)
const override;
122 const unsigned integration_point)
const override;
126 std::vector<GlobalVector*>
const& x,
127 std::vector<NumLib::LocalToGlobalIndexMap const*>
const& dof_table,
128 std::vector<double>& cache)
const override;
ProcessVariable const & _transport_process_variable
std::vector< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType >, Eigen::aligned_allocator< IntegrationPointData< NodalRowVectorType, GlobalDimNodalMatrixType, NodalMatrixType > > > _ip_data
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
typename ShapeMatricesType::GlobalDimNodalMatrixType GlobalDimNodalMatrixType
typename ShapeMatricesType::NodalMatrixType NodalMatrixType
typename ShapeMatricesType::template VectorType< NUM_NODAL_DOF * ShapeFunction::NPOINTS > LocalVectorType
typename ShapeMatricesType::NodalVectorType NodalVectorType
static const int pressure_index
static const int concentration_index
LocalAssemblerData(MeshLib::Element const &element, std::size_t const local_matrix_size, NumLib::GenericIntegrationMethod const &integration_method, bool is_axially_symmetric, RichardsComponentTransportProcessData const &process_data, ProcessVariable const &transport_process_variable)
std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const override
unsigned const _element_id
typename ShapeMatricesType::NodalRowVectorType NodalRowVectorType
typename ShapeMatricesType::template MatrixType< NUM_NODAL_DOF *ShapeFunction::NPOINTS, NUM_NODAL_DOF *ShapeFunction::NPOINTS > LocalMatrixType
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
static const int concentration_size
static const int pressure_size
NumLib::GenericIntegrationMethod const & _integration_method
typename ShapeMatricesType::ShapeMatrices ShapeMatrices
void assemble(double const t, double const dt, std::vector< double > const &local_x, std::vector< double > const &local_x_prev, std::vector< double > &local_M_data, std::vector< double > &local_K_data, std::vector< double > &local_b_data) override
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
RichardsComponentTransportProcessData const & _process_data
Eigen::Map< const Eigen::RowVectorXd > getShapeMatrix(const unsigned integration_point) const override
Provides the shape matrix at the given integration point.
virtual 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 =0
virtual std::vector< double > const & getIntPtSaturation(const double t, std::vector< GlobalVector * > const &x, std::vector< NumLib::LocalToGlobalIndexMap const * > const &dof_table, std::vector< double > &cache) const =0
const unsigned NUM_NODAL_DOF
NumLib::ShapeMatrices< NodalRowVectorType, DimNodalMatrixType, DimMatrixType, GlobalDimNodalMatrixType > ShapeMatrices
MatrixType< ShapeFunction::NPOINTS, ShapeFunction::NPOINTS > NodalMatrixType
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType
MatrixType< GlobalDim, GlobalDim > GlobalDimMatrixType
VectorType< GlobalDim > GlobalDimVectorType
VectorType< ShapeFunction::NPOINTS > NodalVectorType
RowVectorType< ShapeFunction::NPOINTS > NodalRowVectorType
NodalRowVectorType const N
double const integration_weight
IntegrationPointData(NodalRowVectorType const &N_, GlobalDimNodalMatrixType const &dNdx_, double const &integration_weight_, NodalMatrixType const mass_operator_)
GlobalDimNodalMatrixType const dNdx
NodalMatrixType const mass_operator
EIGEN_MAKE_ALIGNED_OPERATOR_NEW