Loading [MathJax]/extensions/MathMenu.js
OGS
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType > Class Template Reference

Detailed Description

template<typename ShapeFunction, typename BHEType>
class ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >

Definition at line 27 of file HeatTransportBHELocalAssemblerBHE.h.

#include <HeatTransportBHELocalAssemblerBHE.h>

Inheritance diagram for ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >:
[legend]
Collaboration diagram for ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >:
[legend]

Public Types

using ShapeMatricesType
using ShapeMatrices = typename ShapeMatricesType::ShapeMatrices
using NodalMatrixType = typename ShapeMatricesType::NodalMatrixType
using NodalVectorType = typename ShapeMatricesType::NodalVectorType
using BheLocalMatrixType
using BheLocalVectorType

Public Member Functions

 HeatTransportBHELocalAssemblerBHE (HeatTransportBHELocalAssemblerBHE const &)=delete
 HeatTransportBHELocalAssemblerBHE (HeatTransportBHELocalAssemblerBHE &&)=delete
 HeatTransportBHELocalAssemblerBHE (MeshLib::Element const &e, NumLib::GenericIntegrationMethod const &integration_method, BHEType const &bhe, 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 std::optional< VectorSegmentgetVectorDeformationSegment () const
Public Member Functions inherited from NumLib::ExtrapolatableElement
virtual ~ExtrapolatableElement ()=default

Private Attributes

HeatTransportBHEProcessData_process_data
std::vector< IntegrationPointDataBHE< ShapeMatricesType >, Eigen::aligned_allocator< IntegrationPointDataBHE< ShapeMatricesType > > > _ip_data
NumLib::GenericIntegrationMethod const & _integration_method
BHEType const & _bhe
std::size_t const _element_id
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
Eigen::Vector3d _element_direction
ShapeMatricesType::template MatrixType< bhe_unknowns_size, bhe_unknowns_size_R_matrix
ShapeMatricesType::template MatrixType< soil_temperature_size, soil_temperature_size_R_s_matrix
ShapeMatricesType::template MatrixType< bhe_unknowns_size, soil_temperature_size_R_pi_s_matrix

Static Private Attributes

static constexpr int bhe_unknowns = BHEType::number_of_unknowns
static constexpr int single_bhe_unknowns_size = ShapeFunction::NPOINTS
static constexpr int soil_temperature_size = ShapeFunction::NPOINTS
static constexpr int soil_temperature_index = 0
static constexpr int bhe_unknowns_size
static constexpr int bhe_unknowns_index = ShapeFunction::NPOINTS
static constexpr int local_matrix_size

Member Typedef Documentation

◆ BheLocalMatrixType

template<typename ShapeFunction, typename BHEType>
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::BheLocalMatrixType
Initial value:
typename ShapeMatricesType::template MatrixType<local_matrix_size,

Definition at line 48 of file HeatTransportBHELocalAssemblerBHE.h.

◆ BheLocalVectorType

template<typename ShapeFunction, typename BHEType>
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::BheLocalVectorType
Initial value:
typename ShapeMatricesType::template VectorType<local_matrix_size>

Definition at line 51 of file HeatTransportBHELocalAssemblerBHE.h.

◆ NodalMatrixType

template<typename ShapeFunction, typename BHEType>
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::NodalMatrixType = typename ShapeMatricesType::NodalMatrixType

Definition at line 44 of file HeatTransportBHELocalAssemblerBHE.h.

◆ NodalVectorType

template<typename ShapeFunction, typename BHEType>
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::NodalVectorType = typename ShapeMatricesType::NodalVectorType

Definition at line 45 of file HeatTransportBHELocalAssemblerBHE.h.

◆ ShapeMatrices

template<typename ShapeFunction, typename BHEType>
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::ShapeMatrices = typename ShapeMatricesType::ShapeMatrices

Definition at line 43 of file HeatTransportBHELocalAssemblerBHE.h.

◆ ShapeMatricesType

template<typename ShapeFunction, typename BHEType>
using ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::ShapeMatricesType
Initial value:
EigenFixedShapeMatrixPolicy< ShapeFunction, GlobalDim > ShapeMatrixPolicyType

Definition at line 41 of file HeatTransportBHELocalAssemblerBHE.h.

Constructor & Destructor Documentation

◆ HeatTransportBHELocalAssemblerBHE() [1/3]

template<typename ShapeFunction, typename BHEType>
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::HeatTransportBHELocalAssemblerBHE ( HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType > const & )
delete

◆ HeatTransportBHELocalAssemblerBHE() [2/3]

template<typename ShapeFunction, typename BHEType>
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::HeatTransportBHELocalAssemblerBHE ( HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType > && )
delete

◆ HeatTransportBHELocalAssemblerBHE() [3/3]

template<typename ShapeFunction, typename BHEType>
ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::HeatTransportBHELocalAssemblerBHE ( MeshLib::Element const & e,
NumLib::GenericIntegrationMethod const & integration_method,
BHEType const & bhe,
bool const is_axially_symmetric,
HeatTransportBHEProcessData & process_data )

Definition at line 22 of file HeatTransportBHELocalAssemblerBHE-impl.h.

31 _bhe(bhe),
32 _element_id(e.getID())
33{
34 // need to make sure that the BHE elements are one-dimensional
35 assert(e.getDimension() == 1);
36
37 unsigned const n_integration_points =
38 _integration_method.getNumberOfPoints();
39
42
43 auto const shape_matrices =
45 3 /* GlobalDim */>(e, is_axially_symmetric,
47
48 // ip data initialization
49 for (unsigned ip = 0; ip < n_integration_points; ip++)
50 {
51 auto const& sm = shape_matrices[ip];
52 // create the class IntegrationPointDataBHE in place
53 _ip_data.push_back(
54 {sm.N, sm.dNdx,
55 _integration_method.getWeightedPoint(ip).getWeight() *
56 sm.integralMeasure * sm.detJ});
57
58 _secondary_data.N[ip] = sm.N;
59 }
60
61 // calculate the element direction vector
62 auto const& p0 = e.getNode(0)->asEigenVector3d();
63 auto const& p1 = e.getNode(1)->asEigenVector3d();
64
66
70 static constexpr int max_num_thermal_exchange_terms = 5;
71 // formulate the local BHE R matrix
74 {
81 // Loop over Gauss points
82 for (unsigned ip = 0; ip < n_integration_points; ip++)
83 {
84 auto const& N = _ip_data[ip].N;
85 auto const& w = _ip_data[ip].integration_weight;
86
87 auto const& R = _bhe.thermalResistance(idx_bhe_unknowns);
88 // calculate mass matrix for current unknown
89 matBHE_loc_R += N.transpose() * N * (1 / R) * w;
90 } // end of loop over integration point
91
92 // The following assembly action is according to Diersch (2013) FEFLOW
93 // book please refer to M.127 and M.128 on page 955 and 956
94 // The if check is absolutely necessary because
95 // (i) In the CXA and CXC case, there are 3 exchange terms,
96 // and it is the same as the number of unknowns;
97 // (ii) In the 1U case, there are 4 exchange terms,
98 // and it is again same as the number of unknowns;
99 // (iii) In the 2U case, there are 5 exchange terms,
100 // and it is less than the number of unknowns (8).
102 {
106 }
107 } // end of loop over BHE unknowns
108
109 // debugging
110 // std::string sep =
111 // "\n----------------------------------------\n";
112 // Eigen::IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
113 // std::cout << "_R_matrix: \n" << sep;
114 // std::cout << _R_matrix.format(CleanFmt) << sep;
115 // std::cout << "_R_s_matrix: \n" << sep;
116 // std::cout << _R_s_matrix.format(CleanFmt) << sep;
117 // std::cout << "_R_pi_s_matrix: \n" << sep;
118 // std::cout << _R_pi_s_matrix.format(CleanFmt) << sep;
119}
ShapeMatricesType::template MatrixType< bhe_unknowns_size, bhe_unknowns_size > _R_matrix
ShapeMatricesType::template MatrixType< soil_temperature_size, soil_temperature_size > _R_s_matrix
std::vector< IntegrationPointDataBHE< ShapeMatricesType >, Eigen::aligned_allocator< IntegrationPointDataBHE< ShapeMatricesType > > > _ip_data
ShapeMatricesType::template MatrixType< bhe_unknowns_size, soil_temperature_size > _R_pi_s_matrix

References _bhe, _element_direction, _element_id, _integration_method, _ip_data, _process_data, _R_matrix, _R_pi_s_matrix, _R_s_matrix, _secondary_data, MathLib::Point3d::asEigenVector3d(), bhe_unknowns, bhe_unknowns_size, MeshLib::Element::getDimension(), MeshLib::Element::getNode(), NumLib::initShapeMatrices(), single_bhe_unknowns_size, and soil_temperature_size.

Member Function Documentation

◆ assemble()

template<typename ShapeFunction, typename BHEType>
void ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::assemble ( double const ,
double const ,
std::vector< double > const & ,
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 122 of file HeatTransportBHELocalAssemblerBHE-impl.h.

128{
133
134 unsigned const n_integration_points =
135 _integration_method.getNumberOfPoints();
136
137 auto const& pipe_heat_capacities = _bhe.pipeHeatCapacities();
138 auto const& pipe_heat_conductions = _bhe.pipeHeatConductions();
139 auto const& pipe_advection_vectors =
140 _bhe.pipeAdvectionVectors(_element_direction);
141 auto const& cross_section_areas = _bhe.crossSectionAreas();
142
143 // the mass and conductance matrix terms
144 for (unsigned ip = 0; ip < n_integration_points; ip++)
145 {
146 auto& ip_data = _ip_data[ip];
147
148 auto const& w = ip_data.integration_weight;
149 auto const& N = ip_data.N;
150 auto const& dNdx = ip_data.dNdx;
151
152 // looping over all unknowns.
155 {
156 // get coefficient of mass from corresponding BHE.
159 auto const& advection_vector =
162
163 int const single_bhe_unknowns_index =
166 // local M
167 local_M
171 .noalias() += N.transpose() * N * mass_coeff * A * w;
172
173 // local K
174 // laplace part
175 local_K
179 .noalias() += dNdx.transpose() * dNdx * lambda * A * w;
180 // advection part
181 local_K
185 .noalias() +=
186 N.transpose() * advection_vector.transpose() * dNdx * A * w;
187 }
188 }
189
190 // add the R matrix to local_K
193
194 // add the R_pi_s matrix to local_K
195 local_K
198 .noalias() += _R_pi_s_matrix;
199 local_K
202 .noalias() += _R_pi_s_matrix.transpose();
203
204 // add the R_s matrix to local_K
205 local_K
208 .noalias() += _bhe.number_of_grout_zones * _R_s_matrix;
209
210 // debugging
211 // std::string sep = "\n----------------------------------------\n";
212 // Eigen::IOFormat CleanFmt(4, 0, ", ", "\n", "[", "]");
213 // std::cout << local_K.format(CleanFmt) << sep;
214 // std::cout << local_M.format(CleanFmt) << sep;
215}
Eigen::Map< Matrix > createZeroedMatrix(std::vector< double > &data, Eigen::MatrixXd::Index rows, Eigen::MatrixXd::Index cols)

References _bhe, _element_direction, _integration_method, _ip_data, _R_matrix, _R_pi_s_matrix, _R_s_matrix, bhe_unknowns, bhe_unknowns_index, MathLib::createZeroedMatrix(), local_matrix_size, single_bhe_unknowns_size, and soil_temperature_index.

Referenced by assembleWithJacobian().

◆ assembleWithJacobian()

template<typename ShapeFunction, typename BHEType>
void ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::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 218 of file HeatTransportBHELocalAssemblerBHE-impl.h.

224{
225 auto const local_matrix_size = local_x.size();
226 // initialize x and x_prev
227 auto x =
231 // initialize local_Jac and local_rhs
236
240 local_rhs_data /*not going to be used*/);
241
242 // convert to matrix
247
248 // Jac matrix and rhs vector operation
249 local_Jac.noalias() += local_K + local_M / dt;
250 local_rhs.noalias() -= local_K * x + local_M * (x - x_prev) / dt;
251
252 local_M.setZero();
253 local_K.setZero();
254}
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(), local_matrix_size, and MathLib::toMatrix().

◆ getShapeMatrix()

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

Provides the shape matrix at the given integration point.

Implements NumLib::ExtrapolatableElement.

Definition at line 79 of file HeatTransportBHELocalAssemblerBHE.h.

81 {
82 auto const& N = _secondary_data.N[integration_point];
83
84 // assumes N is stored contiguously in memory
85 return Eigen::Map<const Eigen::RowVectorXd>(N.data(), N.size());
86 }

References _secondary_data.

Member Data Documentation

◆ _bhe

template<typename ShapeFunction, typename BHEType>
BHEType const& ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_bhe
private

◆ _element_direction

template<typename ShapeFunction, typename BHEType>
Eigen::Vector3d ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_element_direction
private

◆ _element_id

template<typename ShapeFunction, typename BHEType>
std::size_t const ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_element_id
private

◆ _integration_method

template<typename ShapeFunction, typename BHEType>
NumLib::GenericIntegrationMethod const& ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_integration_method
private

◆ _ip_data

template<typename ShapeFunction, typename BHEType>
std::vector< IntegrationPointDataBHE<ShapeMatricesType>, Eigen::aligned_allocator<IntegrationPointDataBHE<ShapeMatricesType> > > ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_ip_data
private

◆ _process_data

template<typename ShapeFunction, typename BHEType>
HeatTransportBHEProcessData& ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_process_data
private

◆ _R_matrix

template<typename ShapeFunction, typename BHEType>
ShapeMatricesType::template MatrixType<bhe_unknowns_size, bhe_unknowns_size> ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_R_matrix
private

◆ _R_pi_s_matrix

template<typename ShapeFunction, typename BHEType>
ShapeMatricesType::template MatrixType<bhe_unknowns_size, soil_temperature_size> ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_R_pi_s_matrix
private

◆ _R_s_matrix

template<typename ShapeFunction, typename BHEType>
ShapeMatricesType::template MatrixType<soil_temperature_size, soil_temperature_size> ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_R_s_matrix
private

◆ _secondary_data

template<typename ShapeFunction, typename BHEType>
SecondaryData<typename ShapeMatrices::ShapeType> ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_secondary_data
private

◆ bhe_unknowns

template<typename ShapeFunction, typename BHEType>
int ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::bhe_unknowns = BHEType::number_of_unknowns
staticconstexprprivate

◆ bhe_unknowns_index

template<typename ShapeFunction, typename BHEType>
int ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::bhe_unknowns_index = ShapeFunction::NPOINTS
staticconstexprprivate

Definition at line 36 of file HeatTransportBHELocalAssemblerBHE.h.

Referenced by assemble().

◆ bhe_unknowns_size

template<typename ShapeFunction, typename BHEType>
int ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::bhe_unknowns_size
staticconstexprprivate

◆ local_matrix_size

template<typename ShapeFunction, typename BHEType>
int ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::local_matrix_size
staticconstexprprivate
Initial value:

Definition at line 37 of file HeatTransportBHELocalAssemblerBHE.h.

Referenced by assemble(), and assembleWithJacobian().

◆ single_bhe_unknowns_size

template<typename ShapeFunction, typename BHEType>
int ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::single_bhe_unknowns_size = ShapeFunction::NPOINTS
staticconstexprprivate

◆ soil_temperature_index

template<typename ShapeFunction, typename BHEType>
int ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::soil_temperature_index = 0
staticconstexprprivate

Definition at line 33 of file HeatTransportBHELocalAssemblerBHE.h.

Referenced by assemble().

◆ soil_temperature_size

template<typename ShapeFunction, typename BHEType>
int ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::soil_temperature_size = ShapeFunction::NPOINTS
staticconstexprprivate

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