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

Detailed Description

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

Definition at line 30 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)
 
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.
 
- 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
 
std::size_t const _element_id
 
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
 

Member Typedef Documentation

◆ GlobalDimNodalMatrixType

template<typename ShapeFunction >
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::GlobalDimNodalMatrixType
Initial value:
MatrixType< GlobalDim, ShapeFunction::NPOINTS > GlobalDimNodalMatrixType

Definition at line 41 of file HeatTransportBHELocalAssemblerSoil.h.

◆ NodalMatrixType

template<typename ShapeFunction >
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 36 of file HeatTransportBHELocalAssemblerSoil.h.

◆ NodalRowVectorType

template<typename ShapeFunction >
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::NodalRowVectorType = typename ShapeMatricesType::NodalRowVectorType

Definition at line 38 of file HeatTransportBHELocalAssemblerSoil.h.

◆ NodalVectorType

template<typename ShapeFunction >
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 37 of file HeatTransportBHELocalAssemblerSoil.h.

◆ ShapeMatrices

template<typename ShapeFunction >
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 40 of file HeatTransportBHELocalAssemblerSoil.h.

◆ ShapeMatricesType

template<typename ShapeFunction >
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::ShapeMatricesType

Constructor & Destructor Documentation

◆ HeatTransportBHELocalAssemblerSoil() [1/3]

template<typename ShapeFunction >
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::HeatTransportBHELocalAssemblerSoil ( HeatTransportBHELocalAssemblerSoil< ShapeFunction > const & )
delete

◆ HeatTransportBHELocalAssemblerSoil() [2/3]

template<typename ShapeFunction >
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::HeatTransportBHELocalAssemblerSoil ( HeatTransportBHELocalAssemblerSoil< ShapeFunction > && )
delete

◆ 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 )

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

36 : _process_data(process_data),
37 _integration_method(integration_method),
38 _element_id(e.getID())
39{
40 unsigned const n_integration_points =
42
43 _ip_data.reserve(n_integration_points);
44 _secondary_data.N.resize(n_integration_points);
45
48 3 /* GlobalDim */>(e, is_axially_symmetric,
50
52 x_position.setElementID(_element_id);
53
54 // ip data initialization
55 for (unsigned ip = 0; ip < n_integration_points; ip++)
56 {
57 x_position.setIntegrationPoint(ip);
58
59 // create the class IntegrationPointDataBHE in place
60 auto const& sm = _shape_matrices[ip];
61 double const w = _integration_method.getWeightedPoint(ip).getWeight() *
62 sm.integralMeasure * sm.detJ;
63 _ip_data.push_back({sm.N, sm.dNdx, w});
64
65 _secondary_data.N[ip] = sm.N;
66 }
67}
double getWeight() const
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
void setElementID(std::size_t element_id)
std::vector< ShapeMatrices, Eigen::aligned_allocator< ShapeMatrices > > _shape_matrices
std::vector< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType >, Eigen::aligned_allocator< IntegrationPointDataSoil< NodalRowVectorType, GlobalDimNodalMatrixType > > > _ip_data
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)
std::vector< ShapeMatrixType, Eigen::aligned_allocator< ShapeMatrixType > > N

References ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_element_id, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_integration_method, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_ip_data, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_secondary_data, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_shape_matrices, NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N, ParameterLib::SpatialPosition::setElementID(), and ParameterLib::SpatialPosition::setIntegrationPoint().

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 70 of file HeatTransportBHELocalAssemblerSoil-impl.h.

75{
76 assert(local_x.size() == ShapeFunction::NPOINTS);
77 (void)local_x; // Avoid unused arg warning.
78
80 local_M_data, ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
82 local_K_data, ShapeFunction::NPOINTS, ShapeFunction::NPOINTS);
83
86
87 auto const& medium = *_process_data.media_map.getMedium(_element_id);
88 auto const& solid_phase = medium.phase("Solid");
89 auto const& liquid_phase = medium.phase("AqueousLiquid");
90
92
93 unsigned const n_integration_points =
95
96 for (unsigned ip = 0; ip < n_integration_points; ip++)
97 {
98 pos.setIntegrationPoint(ip);
99 auto& ip_data = _ip_data[ip];
100 auto const& N = ip_data.N;
101 auto const& dNdx = ip_data.dNdx;
102 auto const& w = ip_data.integration_weight;
103
104 double T_int_pt = 0.0;
105 NumLib::shapeFunctionInterpolate(local_x, N, T_int_pt);
106
107 vars.temperature = T_int_pt;
108
109 // for now only using the solid and liquid phase parameters
110 auto const density_s =
112 .template value<double>(vars, pos, t, dt);
113
114 auto const heat_capacity_s =
115 solid_phase
116 .property(
118 .template value<double>(vars, pos, t, dt);
119
120 auto const density_f =
121 liquid_phase.property(MaterialPropertyLib::PropertyType::density)
122 .template value<double>(vars, pos, t, dt);
123
124 auto const heat_capacity_f =
125 liquid_phase
126 .property(
128 .template value<double>(vars, pos, t, dt);
129
130 auto const porosity =
132 .template value<double>(vars, pos, t, dt);
133
134 auto const velocity =
135 liquid_phase
137 .template value<Eigen::Vector3d>(vars, pos, t, dt);
138
139 // calculate the hydrodynamic thermodispersion tensor
140 auto const thermal_conductivity =
142 medium
143 .property(
145 .value(vars, pos, t, dt));
146
147 auto thermal_conductivity_dispersivity = thermal_conductivity;
148
149 double const velocity_magnitude = velocity.norm();
150
151 if (velocity_magnitude >= std::numeric_limits<double>::epsilon())
152 {
153 auto const thermal_dispersivity_longitudinal =
154 medium
156 thermal_longitudinal_dispersivity)
157 .template value<double>();
158 auto const thermal_dispersivity_transversal =
159 medium
161 thermal_transversal_dispersivity)
162 .template value<double>();
163
164 auto const thermal_dispersivity =
165 density_f * heat_capacity_f *
166 (thermal_dispersivity_transversal * velocity_magnitude *
167 Eigen::Matrix3d::Identity() +
168 (thermal_dispersivity_longitudinal -
169 thermal_dispersivity_transversal) /
170 velocity_magnitude * velocity * velocity.transpose());
171 thermal_conductivity_dispersivity += thermal_dispersivity;
172 }
173
174 // assemble Conductance matrix
175 local_K.noalias() +=
176 (dNdx.transpose() * thermal_conductivity_dispersivity * dNdx +
177 N.transpose() * velocity.transpose() * dNdx * density_f *
178 heat_capacity_f) *
179 w;
180
181 // assemble Mass matrix
182 local_M.noalias() += N.transpose() * N * w *
183 (density_s * heat_capacity_s * (1 - porosity) +
184 density_f * heat_capacity_f * porosity);
185 }
186
187 // debugging
188 // std::string sep = "\n----------------------------------------\n";
189 // Eigen::IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
190 // std::cout << local_K.format(CleanFmt) << sep;
191 // std::cout << local_M.format(CleanFmt) << sep;
192}
Phase const & phase(std::size_t index) const
Definition Medium.cpp:33
void setIntegrationPoint(unsigned integration_point)
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 &)
MaterialPropertyLib::MaterialSpatialDistributionMap media_map

References MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor< 3 >(), MaterialPropertyLib::phase_velocity, MaterialPropertyLib::porosity, ParameterLib::SpatialPosition::setElementID(), ParameterLib::SpatialPosition::setIntegrationPoint(), NumLib::detail::shapeFunctionInterpolate(), MaterialPropertyLib::specific_heat_capacity, MaterialPropertyLib::VariableArray::temperature, and MaterialPropertyLib::thermal_conductivity.

◆ 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 195 of file HeatTransportBHELocalAssemblerSoil-impl.h.

199{
200 assert(local_x.size() == ShapeFunction::NPOINTS);
201 auto const local_matrix_size = local_x.size();
202 // initialize x and x_prev
203 auto x =
204 Eigen::Map<NodalVectorType const>(local_x.data(), local_matrix_size);
205 auto x_prev = Eigen::Map<NodalVectorType const>(local_x_prev.data(),
206 local_matrix_size);
207 // initialize local_Jac and local_rhs
209 local_Jac_data, local_matrix_size, local_matrix_size);
211 local_rhs_data, local_matrix_size);
212
213 std::vector<double> local_M_data;
214 std::vector<double> local_K_data;
215 assemble(t, dt, local_x, local_x_prev, local_M_data, local_K_data,
216 local_rhs_data /*not going to be used*/);
217
218 // convert to matrix
220 local_M_data, local_matrix_size, local_matrix_size);
222 local_K_data, local_matrix_size, local_matrix_size);
223
224 // Jac matrix and rhs vector operation
225 local_Jac.noalias() += local_K + local_M / dt;
226 local_rhs.noalias() -= local_K * x + local_M * (x - x_prev) / dt;
227
228 local_M.setZero();
229 local_K.setZero();
230}
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 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 68 of file HeatTransportBHELocalAssemblerSoil.h.

70 {
71 auto const& N = _secondary_data.N[integration_point];
72
73 // assumes N is stored contiguously in memory
74 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
75 }

References ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_secondary_data, and ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N.

Member Data Documentation

◆ _element_id

◆ _integration_method

◆ _ip_data

◆ _process_data

template<typename ShapeFunction >
HeatTransportBHEProcessData& ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerSoil< ShapeFunction >::_process_data
private

Definition at line 78 of file HeatTransportBHELocalAssemblerSoil.h.

◆ _secondary_data

◆ _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: