OGS
ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim > Class Template Reference

Detailed Description

template<typename ShapeFunction, int GlobalDim>
class ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >

Definition at line 67 of file LiquidFlowLocalAssembler.h.

#include <LiquidFlowLocalAssembler.h>

Inheritance diagram for ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >:
[legend]
Collaboration diagram for ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >:
[legend]

Classes

struct  AnisotropicCalculator
 
struct  IsotropicCalculator
 

Public Member Functions

 LiquidFlowLocalAssembler (MeshLib::Element const &element, std::size_t const, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, LiquidFlowData 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 > &local_b_data) 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 > &velocity_cache) const override
 
- Public Member Functions inherited from ProcessLib::LocalAssemblerInterface
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_M_data, std::vector< double > &local_K_data, 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_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_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.
 
- Public Member Functions inherited from NumLib::ExtrapolatableElement
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 NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
 
using GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
 
using GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
 
using GlobalDimNodalMatrixType
 

Private Member Functions

template<typename LaplacianGravityVelocityCalculator >
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)
 
template<typename LaplacianGravityVelocityCalculator , typename VelocityCacheType >
void computeProjectedDarcyVelocity (const double t, const double dt, std::vector< double > const &local_x, ParameterLib::SpatialPosition const &pos, VelocityCacheType &darcy_velocity_at_ips) const
 
template<typename VelocityCacheType >
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
 

Private Attributes

MeshLib::Element const & _element
 
NumLib::GenericIntegrationMethod const & _integration_method
 
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
 
const LiquidFlowData_process_data
 

Member Typedef Documentation

◆ GlobalDimMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimMatrixType = typename ShapeMatricesType::GlobalDimMatrixType
private

Definition at line 79 of file LiquidFlowLocalAssembler.h.

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimNodalMatrixType
private
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 80 of file LiquidFlowLocalAssembler.h.

◆ GlobalDimVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::GlobalDimVectorType = typename ShapeMatricesType::GlobalDimVectorType
private

Definition at line 78 of file LiquidFlowLocalAssembler.h.

◆ LocalAssemblerTraits

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::LocalAssemblerTraits
private
Initial value:

Definition at line 72 of file LiquidFlowLocalAssembler.h.

◆ NodalMatrixType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::NodalMatrixType = typename LocalAssemblerTraits::LocalMatrix
private

Definition at line 75 of file LiquidFlowLocalAssembler.h.

◆ NodalRowVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
private

Definition at line 77 of file LiquidFlowLocalAssembler.h.

◆ NodalVectorType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::NodalVectorType = typename LocalAssemblerTraits::LocalVector
private

Definition at line 76 of file LiquidFlowLocalAssembler.h.

◆ ShapeMatrices

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
private

Definition at line 70 of file LiquidFlowLocalAssembler.h.

◆ ShapeMatricesType

template<typename ShapeFunction , int GlobalDim>
using ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::ShapeMatricesType = ShapeMatrixPolicyType<ShapeFunction, GlobalDim>
private

Definition at line 69 of file LiquidFlowLocalAssembler.h.

Constructor & Destructor Documentation

◆ LiquidFlowLocalAssembler()

template<typename ShapeFunction , int GlobalDim>
ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::LiquidFlowLocalAssembler ( MeshLib::Element const & element,
std::size_t const ,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
LiquidFlowData const & process_data )
inline

Definition at line 84 of file LiquidFlowLocalAssembler.h.

90 : _element(element),
91 _integration_method(integration_method),
92 _process_data(process_data)
93 {
94 unsigned const n_integration_points =
96 _ip_data.reserve(n_integration_points);
97
98 auto const& shape_matrices =
100 GlobalDim>(element, is_axially_symmetric,
102
105
106 double const aperture_size =
107 (_element.getDimension() == 3u)
108 ? 1.0
109 : _process_data.aperture_size(0.0, pos)[0];
110
111 for (unsigned ip = 0; ip < n_integration_points; ip++)
112 {
113 _ip_data.emplace_back(
114 shape_matrices[ip].dNdx,
116 shape_matrices[ip].integralMeasure *
117 shape_matrices[ip].detJ * aperture_size);
118 }
119 }
double getWeight() const
virtual constexpr unsigned getDimension() const =0
Get dimension of the mesh element.
std::size_t getID() const
Returns the ID of the element.
Definition Element.h:89
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
std::vector< IntegrationPointData< GlobalDimNodalMatrixType > > _ip_data
NumLib::GenericIntegrationMethod const & _integration_method
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)

References ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_element, ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_integration_method, ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_ip_data, ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_process_data, ProcessLib::LiquidFlow::LiquidFlowData::aperture_size, MeshLib::Element::getDimension(), MeshLib::Element::getID(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), and ParameterLib::SpatialPosition::setElementID().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction , int GlobalDim>
void ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::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 )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 26 of file LiquidFlowLocalAssembler-impl.h.

31{
34
35 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
37 vars.temperature =
39 .template value<double>(vars, pos, t, dt);
40 vars.liquid_phase_pressure = std::numeric_limits<double>::quiet_NaN();
42 MaterialPropertyLib::formEigenTensor<GlobalDim>(
44 vars, pos, t, dt));
45 // Note: For Inclined 1D in 2D/3D or 2D element in 3D, the first item in
46 // the assert must be changed to permeability.rows() ==
47 // _element->getDimension()
48 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
49
50 if (permeability.size() == 1)
51 { // isotropic or 1D problem.
52 assembleMatrixAndVector<IsotropicCalculator>(
53 t, dt, local_x, local_M_data, local_K_data, local_b_data);
54 }
55 else
56 {
57 assembleMatrixAndVector<AnisotropicCalculator>(
58 t, dt, local_x, local_M_data, local_K_data, local_b_data);
59 }
60}
typename ShapeMatricesType::GlobalDimMatrixType GlobalDimMatrixType
MaterialPropertyLib::MaterialSpatialDistributionMap media_map

References MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::permeability, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), and MaterialPropertyLib::VariableArray::temperature.

◆ assembleMatrixAndVector()

template<typename ShapeFunction , int GlobalDim>
template<typename LaplacianGravityVelocityCalculator >
void ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::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 )
private

Definition at line 110 of file LiquidFlowLocalAssembler-impl.h.

116{
117 auto const local_matrix_size = local_x.size();
118 assert(local_matrix_size == ShapeFunction::NPOINTS);
119
120 auto local_M = MathLib::createZeroedMatrix<NodalMatrixType>(
121 local_M_data, local_matrix_size, local_matrix_size);
122 auto local_K = MathLib::createZeroedMatrix<NodalMatrixType>(
123 local_K_data, local_matrix_size, local_matrix_size);
124 auto local_b = MathLib::createZeroedVector<NodalVectorType>(
125 local_b_data, local_matrix_size);
126
127 unsigned const n_integration_points =
129
132
133 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
134 auto const& liquid_phase = medium.phase("AqueousLiquid");
135
137 vars.temperature =
139 .template value<double>(vars, pos, t, dt);
140
141 GlobalDimVectorType const projected_body_force_vector =
145
146 auto const& Ns = _process_data.shape_matrix_cache
147 .NsHigherOrder<typename ShapeFunction::MeshElement>();
148
149 for (unsigned ip = 0; ip < n_integration_points; ip++)
150 {
151 auto const& ip_data = _ip_data[ip];
152 auto const& N = Ns[ip];
153
154 double p = 0.;
157
158 // Compute density:
159 auto const fluid_density =
161 .template value<double>(vars, pos, t, dt);
162 assert(fluid_density > 0.);
163 vars.density = fluid_density;
164
165 auto const ddensity_dpressure =
167 .template dValue<double>(
169 pos, t, dt);
170
171 auto const porosity =
173 .template value<double>(vars, pos, t, dt);
175 .template value<double>(vars, pos, t, dt);
176
177 // Assemble mass matrix, M
178 local_M.noalias() +=
179 (porosity * ddensity_dpressure / fluid_density + storage) *
180 N.transpose() * N * ip_data.integration_weight;
181
182 // Compute viscosity:
183 auto const viscosity =
185 .template value<double>(vars, pos, t, dt);
186
187 pos.setIntegrationPoint(ip);
189 MaterialPropertyLib::formEigenTensor<GlobalDim>(
191 vars, pos, t, dt));
192
193 // Assemble Laplacian, K, and RHS by the gravitational term
194 LaplacianGravityVelocityCalculator::calculateLaplacianAndGravityTerm(
195 local_K, local_b, ip_data, permeability, viscosity, fluid_density,
196 projected_body_force_vector, _process_data.has_gravity);
197 }
198}
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
auto const & NsHigherOrder() const
void setIntegrationPoint(unsigned integration_point)
typename ShapeMatricesType::GlobalDimVectorType GlobalDimVectorType
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)
double fluid_density(const double p, const double T, const double x)
NumLib::ShapeMatrixCache shape_matrix_cache
caches for each mesh element type the shape matrix
std::vector< Eigen::MatrixXd > const element_rotation_matrices
A vector of the rotation matrices for all elements.
Eigen::VectorXd const specific_body_force

References MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::liquid_phase_pressure, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::permeability, MaterialPropertyLib::porosity, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::storage, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

◆ computeDarcyVelocity()

template<typename ShapeFunction , int GlobalDim>
template<typename VelocityCacheType >
void ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::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
private

Definition at line 202 of file LiquidFlowLocalAssembler-impl.h.

207{
208 if (is_scalar_permeability)
209 { // isotropic or 1D problem.
210 computeProjectedDarcyVelocity<IsotropicCalculator>(
211 t, dt, local_x, pos, darcy_velocity_at_ips);
212 }
213 else
214 {
215 computeProjectedDarcyVelocity<AnisotropicCalculator>(
216 t, dt, local_x, pos, darcy_velocity_at_ips);
217 }
218}

◆ computeProjectedDarcyVelocity()

template<typename ShapeFunction , int GlobalDim>
template<typename LaplacianGravityVelocityCalculator , typename VelocityCacheType >
void ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::computeProjectedDarcyVelocity ( const double t,
const double dt,
std::vector< double > const & local_x,
ParameterLib::SpatialPosition const & pos,
VelocityCacheType & darcy_velocity_at_ips ) const
private

Definition at line 271 of file LiquidFlowLocalAssembler-impl.h.

276{
277 auto const local_matrix_size = local_x.size();
278 assert(local_matrix_size == ShapeFunction::NPOINTS);
279
280 const auto local_p_vec =
281 MathLib::toVector<NodalVectorType>(local_x, local_matrix_size);
282
283 unsigned const n_integration_points =
285
286 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
287 auto const& liquid_phase = medium.phase("AqueousLiquid");
288
290 vars.temperature =
292 .template value<double>(vars, pos, t, dt);
293
294 GlobalDimVectorType const projected_body_force_vector =
298
299 auto const& Ns = _process_data.shape_matrix_cache
300 .NsHigherOrder<typename ShapeFunction::MeshElement>();
301
302 for (unsigned ip = 0; ip < n_integration_points; ip++)
303 {
304 auto const& ip_data = _ip_data[ip];
305 auto const& N = Ns[ip];
306 double p = 0.;
309
310 // Compute density:
311 auto const fluid_density =
313 .template value<double>(vars, pos, t, dt);
314 vars.density = fluid_density;
315
316 // Compute viscosity:
317 auto const viscosity =
319 .template value<double>(vars, pos, t, dt);
320
322 MaterialPropertyLib::formEigenTensor<GlobalDim>(
324 vars, pos, t, dt));
325
326 darcy_velocity_at_ips.col(ip) =
327 LaplacianGravityVelocityCalculator::calculateVelocity(
328 local_p_vec, ip_data, permeability, viscosity, fluid_density,
329 projected_body_force_vector, _process_data.has_gravity);
330 }
331}

References MaterialPropertyLib::density, MaterialPropertyLib::VariableArray::density, MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::permeability, MaterialPropertyLib::reference_temperature, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::viscosity.

◆ getFlux()

template<typename ShapeFunction , int GlobalDim>
Eigen::Vector3d ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::getFlux ( MathLib::Point3d const & p_local_coords,
double const t,
std::vector< double > const & local_x ) const
overridevirtual

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 63 of file LiquidFlowLocalAssembler-impl.h.

66{
67 // TODO (tf) Temporary value not used by current material models. Need
68 // extension of getFlux interface
69 double const dt = std::numeric_limits<double>::quiet_NaN();
70
71 // Note: Axial symmetry is set to false here, because we only need dNdx
72 // here, which is not affected by axial symmetry.
73 auto const shape_matrices =
75 GlobalDim>(_element,
76 false /*is_axially_symmetric*/,
77 std::array{p_local_coords})[0];
78
79 // create pos object to access the correct media property
82
83 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
84 auto const& liquid_phase = medium.phase("AqueousLiquid");
85
87
88 double pressure = 0.0;
89 NumLib::shapeFunctionInterpolate(local_x, shape_matrices.N, pressure);
90 vars.liquid_phase_pressure = pressure;
91
92 GlobalDimMatrixType const intrinsic_permeability =
93 MaterialPropertyLib::formEigenTensor<GlobalDim>(
95 vars, pos, t, dt));
96 auto const viscosity =
98 .template value<double>(vars, pos, t, dt);
99
100 Eigen::Vector3d flux(0.0, 0.0, 0.0);
101 flux.head<GlobalDim>() =
102 -intrinsic_permeability / viscosity * shape_matrices.dNdx *
103 Eigen::Map<const NodalVectorType>(local_x.data(), local_x.size());
104
105 return flux;
106}
std::vector< typename ShapeMatricesType::ShapeMatrices, Eigen::aligned_allocator< typename ShapeMatricesType::ShapeMatrices > > computeShapeMatrices(MeshLib::Element const &e, bool const is_axially_symmetric, PointContainer const &points)

References NumLib::computeShapeMatrices(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::permeability, ParameterLib::SpatialPosition::setElementID(), NumLib::detail::shapeFunctionInterpolate(), and MaterialPropertyLib::viscosity.

◆ getIntPtDarcyVelocity()

template<typename ShapeFunction , int GlobalDim>
std::vector< double > const & ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::getIntPtDarcyVelocity ( const double t,
std::vector< GlobalVector * > const & x,
std::vector< NumLib::LocalToGlobalIndexMap const * > const & dof_table,
std::vector< double > & velocity_cache ) const
overridevirtual

Implements ProcessLib::LiquidFlow::LiquidFlowLocalAssemblerInterface.

Definition at line 222 of file LiquidFlowLocalAssembler-impl.h.

227{
228 // TODO (tf) Temporary value not used by current material models. Need
229 // extension of secondary variable interface.
230 double const dt = std::numeric_limits<double>::quiet_NaN();
231
232 constexpr int process_id = 0;
233 auto const indices =
234 NumLib::getIndices(_element.getID(), *dof_tables[process_id]);
235 auto const local_x = x[process_id]->get(indices);
236 auto const n_integration_points = _integration_method.getNumberOfPoints();
237 velocity_cache.clear();
238
241
242 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
244 vars.temperature =
246 .template value<double>(vars, pos, t, dt);
247 vars.liquid_phase_pressure = std::numeric_limits<double>::quiet_NaN();
248
250 MaterialPropertyLib::formEigenTensor<GlobalDim>(
252 vars, pos, t, dt));
253
254 assert(permeability.rows() == GlobalDim || permeability.rows() == 1);
255
256 bool const is_scalar_permeability = (permeability.size() == 1);
257
258 auto velocity_cache_vectors = MathLib::createZeroedMatrix<
259 Eigen::Matrix<double, GlobalDim, Eigen::Dynamic, Eigen::RowMajor>>(
260 velocity_cache, GlobalDim, n_integration_points);
261
262 computeDarcyVelocity(is_scalar_permeability, t, dt, local_x, pos,
263 velocity_cache_vectors);
264
265 return velocity_cache;
266}
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
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
std::vector< GlobalIndexType > getIndices(std::size_t const mesh_item_id, NumLib::LocalToGlobalIndexMap const &dof_table)

References MathLib::createZeroedMatrix(), NumLib::getIndices(), MaterialPropertyLib::VariableArray::liquid_phase_pressure, MaterialPropertyLib::permeability, MaterialPropertyLib::reference_temperature, ParameterLib::SpatialPosition::setElementID(), and MaterialPropertyLib::VariableArray::temperature.

◆ getShapeMatrix()

template<typename ShapeFunction , int GlobalDim>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 134 of file LiquidFlowLocalAssembler.h.

136 {
137 auto const& N =
139 .NsHigherOrder<typename ShapeFunction::MeshElement>();
140
141 // assumes N is stored contiguously in memory
142 return Eigen::Map<const Eigen::RowVectorXd>(
143 N[integration_point].data(), N[integration_point].size());
144 }
constexpr int size(int const displacement_dim)
Vectorized tensor size for given displacement dimension.

References ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_process_data, NumLib::ShapeMatrixCache::NsHigherOrder(), and ProcessLib::LiquidFlow::LiquidFlowData::shape_matrix_cache.

Member Data Documentation

◆ _element

template<typename ShapeFunction , int GlobalDim>
MeshLib::Element const& ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_element
private

◆ _integration_method

template<typename ShapeFunction , int GlobalDim>
NumLib::GenericIntegrationMethod const& ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_integration_method
private

◆ _ip_data

template<typename ShapeFunction , int GlobalDim>
std::vector<IntegrationPointData<GlobalDimNodalMatrixType> > ProcessLib::LiquidFlow::LiquidFlowLocalAssembler< ShapeFunction, GlobalDim >::_ip_data
private

◆ _process_data


The documentation for this class was generated from the following files: