template<typename ShapeFunction, typename BHEType>
class ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >
Definition at line 27 of file HeatTransportBHELocalAssemblerBHE.h.
|
| 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.
|
|
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.
|
|
template<typename ShapeFunction , typename BHEType >
Definition at line 22 of file HeatTransportBHELocalAssemblerBHE-impl.h.
33{
34
35 assert(e.getDimension() == 1);
36
37 unsigned const n_integration_points =
39
40 _ip_data.reserve(n_integration_points);
42
43 auto const shape_matrices =
45 3 >(e, is_axially_symmetric,
47
48
49 for (unsigned ip = 0; ip < n_integration_points; ip++)
50 {
51 auto const& sm = shape_matrices[ip];
52
54 {sm.N, sm.dNdx,
56 sm.integralMeasure * sm.detJ});
57
59 }
60
61
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
72 for (
int idx_bhe_unknowns = 0; idx_bhe_unknowns <
bhe_unknowns;
73 idx_bhe_unknowns++)
74 {
75 typename ShapeMatricesType::template MatrixType<
77 matBHE_loc_R = ShapeMatricesType::template MatrixType<
81
82 for (unsigned ip = 0; ip < n_integration_points; ip++)
83 {
85 auto const& w =
_ip_data[ip].integration_weight;
86
87 auto const& R =
_bhe.thermalResistance(idx_bhe_unknowns);
88
89 matBHE_loc_R +=
N.transpose() *
N * (1 / R) * w;
90 }
91
92
93
94
95
96
97
98
99
100
101 if (idx_bhe_unknowns < max_num_thermal_exchange_terms)
102 {
103 _bhe.template assembleRMatrices<ShapeFunction::NPOINTS>(
106 }
107 }
108
109
110
111
112
113
114
115
116
117
118
119}
MathLib::WeightedPoint const & getWeightedPoint(unsigned const igp) const
unsigned getNumberOfPoints() const
ShapeMatricesType::template MatrixType< bhe_unknowns_size, bhe_unknowns_size > _R_matrix
ShapeMatricesType::template MatrixType< soil_temperature_size, soil_temperature_size > _R_s_matrix
Eigen::Vector3d _element_direction
SecondaryData< typename ShapeMatrices::ShapeType > _secondary_data
ShapeMatrixPolicyType< ShapeFunction, 3 > ShapeMatricesType
static constexpr int soil_temperature_size
std::vector< IntegrationPointDataBHE< ShapeMatricesType >, Eigen::aligned_allocator< IntegrationPointDataBHE< ShapeMatricesType > > > _ip_data
ShapeMatricesType::template MatrixType< bhe_unknowns_size, soil_temperature_size > _R_pi_s_matrix
static constexpr int bhe_unknowns_size
static constexpr int bhe_unknowns
NumLib::GenericIntegrationMethod const & _integration_method
static constexpr int single_bhe_unknowns_size
std::size_t const _element_id
HeatTransportBHEProcessData & _process_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::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_bhe, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_element_direction, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_integration_method, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_ip_data, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_R_matrix, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_R_pi_s_matrix, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_R_s_matrix, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::_secondary_data, MathLib::Point3d::asEigenVector3d(), ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::bhe_unknowns, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::bhe_unknowns_size, MeshLib::Element::getDimension(), MeshLib::Element::getNode(), NumLib::GenericIntegrationMethod::getNumberOfPoints(), MathLib::WeightedPoint::getWeight(), NumLib::GenericIntegrationMethod::getWeightedPoint(), NumLib::initShapeMatrices(), ProcessLib::HeatTransportBHE::SecondaryData< ShapeMatrixType >::N, ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::single_bhe_unknowns_size, and ProcessLib::HeatTransportBHE::HeatTransportBHELocalAssemblerBHE< ShapeFunction, BHEType >::soil_temperature_size.