OGS
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction > Class Template Reference

Detailed Description

template<typename ShapeFunction>
class ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >

Definition at line 24 of file HeatTransportBHELocalAssemblerSoil.h.

#include <HeatTransportBHELocalAssemblerSoil.h>

Inheritance diagram for ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >:
[legend]
Collaboration diagram for ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >:
[legend]

Public Types

using ShapeMatricesType
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using GlobalDimNodalMatrixType

Public Member Functions

 HeatTransportBHELocalAssemblerSoil (HeatTransportBHELocalAssemblerSoil const &)=delete
 HeatTransportBHELocalAssemblerSoil (HeatTransportBHELocalAssemblerSoil &&)=delete
 HeatTransportBHELocalAssemblerSoil (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, bool const is_axially_symmetric, HeatTransportBHEProcessData &process_data, BHEMeshData const &bhe_mesh_data)
void assemble (double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, 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.
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 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 int getNumberOfVectorElementsForDeformation () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Private Attributes

HeatTransportBHEProcessData_process_data
std::vector< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType >, Eigen::aligned_allocator< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType > > > _ip_data
NumLib::GenericIntegrationMethod const & _integration_method
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
MeshLib::Element const & _element
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data

Member Typedef Documentation

◆ GlobalDimNodalMatrixType

Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 35 of file HeatTransportBHELocalAssemblerSoil.h.

◆ NodalMatrixType

◆ NodalRowVectorType

◆ NodalVectorType

◆ ShapeMatrices

◆ ShapeMatricesType

Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 28 of file HeatTransportBHELocalAssemblerSoil.h.

Constructor & Destructor Documentation

◆ HeatTransportBHELocalAssemblerSoil() [1/3]

◆ HeatTransportBHELocalAssemblerSoil() [2/3]

◆ HeatTransportBHELocalAssemblerSoil() [3/3]

template<typename ShapeFunction>
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::HeatTransportBHELocalAssemblerSoil ( MeshLib::Element const & e,
NumLib::GenericIntegrationMethod const & integration_method,
bool const is_axially_symmetric,
HeatTransportBHEProcessData & process_data,
BHEMeshData const & bhe_mesh_data )

Definition at line 25 of file HeatTransportBHELocalAssemblerSoil-impl.h.

34 _element(e)
35{
36 unsigned const n_integration_points =
37 _integration_method.getNumberOfPoints();
38
41
44 3 /* GlobalDim */>(e, is_axially_symmetric,
46
47 // ip data initialization
48 for (unsigned ip = 0; ip < n_integration_points; ip++)
49 {
50 // create the class IntegrationPointDataBHE in place
51 auto const& sm = _shape_matrices[ip];
52 double const w = _integration_method.getWeightedPoint(ip).getWeight() *
53 sm.integralMeasure * sm.detJ;
54 _ip_data.push_back({sm.N, sm.dNdx, w});
55
56 _secondary_data.N[ip] = sm.N;
57 }
58}
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
std::vector< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType >, Eigen::aligned_allocator< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType > > > _ip_data

References _element, _integration_method, _ip_data, _process_data, _secondary_data, _shape_matrices, and NumLib::initShapeMatrices().

Member Function Documentation

◆ assemble()

template<typename ShapeFunction>
void ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::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 > &  )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 61 of file HeatTransportBHELocalAssemblerSoil-impl.h.

66{
68 (void)local_x; // Avoid unused arg warning.
69
74
75 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
76 auto const& solid_phase =
78 auto const& liquid_phase =
80
82
83 unsigned const n_integration_points =
84 _integration_method.getNumberOfPoints();
85
86 for (unsigned ip = 0; ip < n_integration_points; ip++)
87 {
88 auto& ip_data = _ip_data[ip];
89 auto const& N = ip_data.N;
90 auto const& dNdx = ip_data.dNdx;
91 auto const& w = ip_data.integration_weight;
92
94 std::nullopt, _element.getID(),
97 _element, N))};
98
99 double T_int_pt = 0.0;
101
102 vars.temperature = T_int_pt;
103
104 // for now only using the solid and liquid phase parameters
105 auto const density_s =
107 .template value<double>(vars, pos, t, dt);
108
109 auto const heat_capacity_s =
111 .property(
113 .template value<double>(vars, pos, t, dt);
114
115 auto const density_f =
117 .template value<double>(vars, pos, t, dt);
118
119 auto const heat_capacity_f =
121 .property(
123 .template value<double>(vars, pos, t, dt);
124
125 auto const porosity =
127 .template value<double>(vars, pos, t, dt);
128
129 auto const velocity =
132 .template value<Eigen::Vector3d>(vars, pos, t, dt);
133
134 // calculate the hydrodynamic thermodispersion tensor
135 auto const thermal_conductivity =
137 medium
138 .property(
140 .value(vars, pos, t, dt));
141
143
144 double const velocity_magnitude = velocity.norm();
145
147 {
149 medium
152 .template value<double>();
154 medium
157 .template value<double>();
158
159 auto const thermal_dispersivity =
165 velocity_magnitude * velocity * velocity.transpose());
167 }
168
169 // assemble Conductance matrix
170 local_K.noalias() +=
172 N.transpose() * velocity.transpose() * dNdx * density_f *
174 w;
175
176 // assemble Mass matrix
177 local_M.noalias() += N.transpose() * N * w *
180 }
181
182 if (_process_data._mass_lumping)
183 {
184 // only mass lumping at the BHE connected soil elements
185 if (_process_data.mass_lumping_soil_elements[_element.getID()])
186 {
187 local_M = local_M.colwise().sum().eval().asDiagonal();
188 }
189 }
190
191 // debugging
192 // std::string sep = "\n----------------------------------------\n";
193 // Eigen::IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
194 // std::cout << local_K.format(CleanFmt) << sep;
195 // std::cout << local_M.format(CleanFmt) << sep;
196}
template Eigen::Matrix< double, 3, 3 > formEigenTensor< 3 >(MaterialPropertyLib::PropertyDataType const &values)
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)
void shapeFunctionInterpolate(const NodalValues &, const ShapeMatrix &)

References _element, _integration_method, _ip_data, _process_data, MaterialPropertyLib::AqueousLiquid, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor< 3 >(), NumLib::interpolateCoordinates(), MaterialPropertyLib::phase_velocity, MaterialPropertyLib::porosity, NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::Solid, MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.

Referenced by assembleWithJacobian().

◆ assembleWithJacobian()

template<typename ShapeFunction>
void ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::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 )
overridevirtual

Reimplemented from ProcessLib::LocalAssemblerInterface.

Definition at line 199 of file HeatTransportBHELocalAssemblerSoil-impl.h.

203{
205 auto const local_matrix_size = local_x.size();
206 // initialize x and x_prev
207 auto x =
211 // initialize local_Jac and local_rhs
216
220 local_rhs_data /*not going to be used*/);
221
222 // convert to matrix
227
228 // Jac matrix and rhs vector operation
229 local_Jac.noalias() += local_K + local_M / dt;
230 local_rhs.noalias() -= local_K * x + local_M * (x - x_prev) / dt;
231
232 local_M.setZero();
233 local_K.setZero();
234}
void assemble(double const, double const, std::vector< double > const &, std::vector< double > const &, std::vector< double > &, std::vector< double > &, std::vector< double > &) override
Eigen::Map< Vector > createZeroedVector(std::vector< double > &data, Eigen::VectorXd::Index size)
Eigen::Map< const Matrix > toMatrix(std::vector< double > const &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)

References assemble(), MathLib::createZeroedMatrix(), MathLib::createZeroedVector(), and MathLib::toMatrix().

◆ getShapeMatrix()

template<typename ShapeFunction>
Eigen::Map< const Eigen::RowVectorXd > ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::getShapeMatrix ( const unsigned integration_point) const
inlineoverridevirtual

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 63 of file HeatTransportBHELocalAssemblerSoil.h.

65 {
66 auto const& N = _secondary_data.N[integration_point];
67
68 // assumes N is stored contiguously in memory
69 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
70 }

References _secondary_data.

Member Data Documentation

◆ _element

◆ _integration_method

◆ _ip_data

◆ _process_data

◆ _secondary_data

template<typename ShapeFunction>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_secondary_data
private

◆ _shape_matrices

template<typename ShapeFunction>
std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices> > ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_shape_matrices
private

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