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

Detailed Description

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

Definition at line 23 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.
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 34 of file HeatTransportBHELocalAssemblerSoil.h.

◆ NodalMatrixType

◆ NodalRowVectorType

◆ NodalVectorType

◆ ShapeMatrices

◆ ShapeMatricesType

Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 27 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 )

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

32 _element(e)
33{
34 unsigned const n_integration_points =
35 _integration_method.getNumberOfPoints();
36
39
42 3 /* GlobalDim */>(e, is_axially_symmetric,
44
45 // ip data initialization
46 for (unsigned ip = 0; ip < n_integration_points; ip++)
47 {
48 // create the class IntegrationPointDataBHE in place
49 auto const& sm = _shape_matrices[ip];
50 double const w = _integration_method.getWeightedPoint(ip).getWeight() *
51 sm.integralMeasure * sm.detJ;
52 _ip_data.push_back({sm.N, sm.dNdx, w});
53
54 _secondary_data.N[ip] = sm.N;
55 }
56}
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 59 of file HeatTransportBHELocalAssemblerSoil-impl.h.

64{
66 (void)local_x; // Avoid unused arg warning.
67
72
73 auto const& medium = *_process_data.media_map.getMedium(_element.getID());
74 auto const& solid_phase = medium.phase("Solid");
75 auto const& liquid_phase = medium.phase("AqueousLiquid");
76
78
79 unsigned const n_integration_points =
80 _integration_method.getNumberOfPoints();
81
82 for (unsigned ip = 0; ip < n_integration_points; ip++)
83 {
84 auto& ip_data = _ip_data[ip];
85 auto const& N = ip_data.N;
86 auto const& dNdx = ip_data.dNdx;
87 auto const& w = ip_data.integration_weight;
88
90 std::nullopt, _element.getID(),
93 _element, N))};
94
95 double T_int_pt = 0.0;
97
98 vars.temperature = T_int_pt;
99
100 // for now only using the solid and liquid phase parameters
101 auto const density_s =
103 .template value<double>(vars, pos, t, dt);
104
105 auto const heat_capacity_s =
107 .property(
109 .template value<double>(vars, pos, t, dt);
110
111 auto const density_f =
113 .template value<double>(vars, pos, t, dt);
114
115 auto const heat_capacity_f =
117 .property(
119 .template value<double>(vars, pos, t, dt);
120
121 auto const porosity =
123 .template value<double>(vars, pos, t, dt);
124
125 auto const velocity =
128 .template value<Eigen::Vector3d>(vars, pos, t, dt);
129
130 // calculate the hydrodynamic thermodispersion tensor
131 auto const thermal_conductivity =
133 medium
134 .property(
136 .value(vars, pos, t, dt));
137
139
140 double const velocity_magnitude = velocity.norm();
141
143 {
145 medium
148 .template value<double>();
150 medium
153 .template value<double>();
154
155 auto const thermal_dispersivity =
161 velocity_magnitude * velocity * velocity.transpose());
163 }
164
165 // assemble Conductance matrix
166 local_K.noalias() +=
168 N.transpose() * velocity.transpose() * dNdx * density_f *
170 w;
171
172 // assemble Mass matrix
173 local_M.noalias() += N.transpose() * N * w *
176 }
177
178 if (_process_data._mass_lumping)
179 {
180 // only mass lumping at the BHE connected soil elements
181 if (_process_data.mass_lumping_soil_elements[_element.getID()])
182 {
183 local_M = local_M.colwise().sum().eval().asDiagonal();
184 }
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}
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, MathLib::createZeroedMatrix(), MaterialPropertyLib::density, MaterialPropertyLib::formEigenTensor< 3 >(), NumLib::interpolateCoordinates(), MaterialPropertyLib::phase_velocity, MaterialPropertyLib::porosity, NumLib::detail::shapeFunctionInterpolate(), 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 195 of file HeatTransportBHELocalAssemblerSoil-impl.h.

199{
201 auto const local_matrix_size = local_x.size();
202 // initialize x and x_prev
203 auto x =
207 // initialize local_Jac and local_rhs
212
216 local_rhs_data /*not going to be used*/);
217
218 // convert to matrix
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 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 61 of file HeatTransportBHELocalAssemblerSoil.h.

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

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: