From 08fff44b70d34d218e1a455716356cf526a5fb07 Mon Sep 17 00:00:00 2001 From: Haibing Shao <haibing.shao@ufz.de> Date: Mon, 5 Nov 2018 10:50:39 +0100 Subject: [PATCH] [PL] Code cleanup and typos. 1) move function bodies from .h to .cpp file. 2) add documentation of the BHE_1U class. 3) add the default case for the switch. --- .../BHEInflowDirichletBoundaryCondition.h | 2 +- ProcessLib/HeatTransportBHE/BHE/BHE_1U.cpp | 88 ++++++++++++++- ProcessLib/HeatTransportBHE/BHE/BHE_1U.h | 102 +++++------------- ProcessLib/HeatTransportBHE/BHE/MeshUtils.cpp | 2 +- .../BHE/PipeConfiguration1U.h | 2 +- .../HeatTransportBHEProcess.cpp | 7 +- .../HeatTransportBHELocalAssemblerBHE-impl.h | 2 +- .../HeatTransportBHELocalAssemblerBHE.h | 2 +- .../HeatTransportBHELocalAssemblerSoil-impl.h | 6 +- .../HeatTransportBHELocalAssemblerSoil.h | 2 +- .../LocalAssemblers/SecondaryData.h | 2 +- 11 files changed, 126 insertions(+), 91 deletions(-) diff --git a/ProcessLib/BoundaryCondition/BHEInflowDirichletBoundaryCondition.h b/ProcessLib/BoundaryCondition/BHEInflowDirichletBoundaryCondition.h index b955493cda3..87d22fc3471 100644 --- a/ProcessLib/BoundaryCondition/BHEInflowDirichletBoundaryCondition.h +++ b/ProcessLib/BoundaryCondition/BHEInflowDirichletBoundaryCondition.h @@ -36,7 +36,7 @@ public: bc_values.ids[0] = _in_out_global_indices.first; // here call the corresponding BHE functions auto const T_out = x[_in_out_global_indices.second]; - bc_values.values[0] = _bhe.getTinByTout(T_out, t); + bc_values.values[0] = _bhe.updateFlowRateAndTemperature(T_out, t); } private: diff --git a/ProcessLib/HeatTransportBHE/BHE/BHE_1U.cpp b/ProcessLib/HeatTransportBHE/BHE/BHE_1U.cpp index 1fd83dc95d1..76603585b82 100644 --- a/ProcessLib/HeatTransportBHE/BHE/BHE_1U.cpp +++ b/ProcessLib/HeatTransportBHE/BHE/BHE_1U.cpp @@ -13,10 +13,87 @@ #include "FlowAndTemperatureControl.h" #include "Physics.h" -using namespace ProcessLib::HeatTransportBHE::BHE; +namespace ProcessLib +{ +namespace HeatTransportBHE +{ +namespace BHE +{ +BHE_1U::BHE_1U(BoreholeGeometry const& borehole, + RefrigerantProperties const& refrigerant, + GroutParameters const& grout, + FlowAndTemperatureControl const& flowAndTemperatureControl, + PipeConfiguration1U const& pipes) + : BHECommon{borehole, refrigerant, grout, flowAndTemperatureControl}, + _pipes(pipes) +{ + // Initialize thermal resistances. + auto values = apply_visitor( + [&](auto const& control) { + return control(refrigerant.reference_temperature, + 0. /* initial time */); + }, + flowAndTemperatureControl); + updateHeatTransferCoefficients(values.flow_rate); +} + +std::array<double, BHE_1U::number_of_unknowns> BHE_1U::pipeHeatCapacities() + const +{ + double const& rho_r = refrigerant.density; + double const& specific_heat_capacity = refrigerant.specific_heat_capacity; + double const& rho_g = grout.rho_g; + double const& porosity_g = grout.porosity_g; + double const& heat_cap_g = grout.heat_cap_g; + + return {{/*i1*/ rho_r * specific_heat_capacity, + /*o1*/ rho_r * specific_heat_capacity, + /*g1*/ (1.0 - porosity_g) * rho_g * heat_cap_g, + /*g2*/ (1.0 - porosity_g) * rho_g * heat_cap_g}}; +} -namespace +std::array<double, BHE_1U::number_of_unknowns> BHE_1U::pipeHeatConductions() + const { + double const& lambda_r = refrigerant.thermal_conductivity; + double const& rho_r = refrigerant.density; + double const& Cp_r = refrigerant.specific_heat_capacity; + double const& alpha_L = _pipes.longitudinal_dispersion_length; + double const& porosity_g = grout.porosity_g; + double const& lambda_g = grout.lambda_g; + + double const velocity_norm = std::abs(_flow_velocity) * std::sqrt(2); + + // Here we calculate the laplace coefficients in the governing + // equations of BHE. These governing equations can be found in + // 1) Diersch (2013) FEFLOW book on page 952, M.120-122, or + // 2) Diersch (2011) Comp & Geosci 37:1122-1135, Eq. 19-22. + return {{// pipe i1, Eq. 19 + (lambda_r + rho_r * Cp_r * alpha_L * velocity_norm), + // pipe o1, Eq. 20 + (lambda_r + rho_r * Cp_r * alpha_L * velocity_norm), + // pipe g1, Eq. 21 + (1.0 - porosity_g) * lambda_g, + // pipe g2, Eq. 22 + (1.0 - porosity_g) * lambda_g}}; +} + +std::array<Eigen::Vector3d, BHE_1U::number_of_unknowns> +BHE_1U::pipeAdvectionVectors() const +{ + double const& rho_r = refrigerant.density; + double const& Cp_r = refrigerant.specific_heat_capacity; + + return {{// pipe i1, Eq. 19 + {0, 0, -rho_r * Cp_r * _flow_velocity}, + // pipe o1, Eq. 20 + {0, 0, rho_r * Cp_r * _flow_velocity}, + // grout g1, Eq. 21 + {0, 0, 0}, + // grout g2, Eq. 22 + {0, 0, 0}}}; +} + double compute_R_gs(double const chi, double const R_g) { return (1 - chi) * R_g; @@ -84,7 +161,6 @@ std::pair<double, double> thermalResistancesGroutSoil(double chi, } return {R_gg, R_gs}; } -} // namespace constexpr std::pair<int, int> BHE_1U::inflow_outflow_bc_component_ids[]; @@ -168,7 +244,8 @@ std::array<double, BHE_1U::number_of_unknowns> BHE_1U::calcThermalResistances( // ------------------------------------------------------------------------- } -double BHE_1U::getTinByTout(double const T_out, double const current_time) +double BHE_1U::updateFlowRateAndTemperature(double const T_out, + double const current_time) { auto values = apply_visitor( [&](auto const& control) { return control(T_out, current_time); }, @@ -176,3 +253,6 @@ double BHE_1U::getTinByTout(double const T_out, double const current_time) updateHeatTransferCoefficients(values.flow_rate); return values.temperature; } +} // namespace BHE +} // namespace HeatTransportBHE +} // namespace ProcessLib diff --git a/ProcessLib/HeatTransportBHE/BHE/BHE_1U.h b/ProcessLib/HeatTransportBHE/BHE/BHE_1U.h index 1aaa0c91eb3..0063864dce7 100644 --- a/ProcessLib/HeatTransportBHE/BHE/BHE_1U.h +++ b/ProcessLib/HeatTransportBHE/BHE/BHE_1U.h @@ -23,6 +23,19 @@ namespace HeatTransportBHE { namespace BHE { +/** + * The BHE_1U class is the realization of 1U type of Borehole Heate Exchanger. + * In this class, the pipe heat capacity, pipe heat conductiion, pie advection + * vectors are intialized according to the geometry of 1U type of BHE. + * For 1U type of BHE, 4 primary unknowns are assigned on the 1D BHE elements. + * They are the temperature in inflow pipe T_in, temperature in outflow pipe + * T_out temperature of the two grout zones sorrounding the inflow and outflow + * pipe T_g1 and T_g2. These primary varaibles are solved according to heat + * convection and conduction equations on the pipes and also in the grout zones. + * The interaction of the 1U type of BHE and the sorrounding soil is regulated + * through the thermal resistance values, which are calculated specifically + * during the initialization of the class. + */ class BHE_1U final : public BHECommon { public: @@ -30,76 +43,16 @@ public: RefrigerantProperties const& refrigerant, GroutParameters const& grout, FlowAndTemperatureControl const& flowAndTemperatureControl, - PipeConfiguration1U const& pipes) - : BHECommon{borehole, refrigerant, grout, flowAndTemperatureControl}, - _pipes(pipes) - { - // Initialize thermal resistances. - auto values = apply_visitor( - [&](auto const& control) { - return control(refrigerant.reference_temperature, - 0. /* initial time */); - }, - flowAndTemperatureControl); - updateHeatTransferCoefficients(values.flow_rate); - } + PipeConfiguration1U const& pipes); static constexpr int number_of_unknowns = 4; - std::array<double, number_of_unknowns> pipeHeatCapacities() const - { - double const& rho_r = refrigerant.density; - double const& specific_heat_capacity = - refrigerant.specific_heat_capacity; - double const& rho_g = grout.rho_g; - double const& porosity_g = grout.porosity_g; - double const& heat_cap_g = grout.heat_cap_g; - - return {{/*i1*/ rho_r * specific_heat_capacity, - /*o1*/ rho_r * specific_heat_capacity, - /*g1*/ (1.0 - porosity_g) * rho_g * heat_cap_g, - /*g2*/ (1.0 - porosity_g) * rho_g * heat_cap_g}}; - } + std::array<double, number_of_unknowns> pipeHeatCapacities() const; - std::array<double, number_of_unknowns> pipeHeatConductions() const - { - double const& lambda_r = refrigerant.thermal_conductivity; - double const& rho_r = refrigerant.density; - double const& Cp_r = refrigerant.specific_heat_capacity; - double const& alpha_L = _pipes.longitudinal_despersion_length; - double const& porosity_g = grout.porosity_g; - double const& lambda_g = grout.lambda_g; - - double const velocity_norm = std::abs(_flow_velocity) * std::sqrt(2); - - // Here we calculate the laplace coefficients in the governing - // equations of BHE. These governing equations can be found in - // 1) Diersch (2013) FEFLOW book on page 952, M.120-122, or - // 2) Diersch (2011) Comp & Geosci 37:1122-1135, Eq. 19-22. - return {{// pipe i1, Eq. 19 - (lambda_r + rho_r * Cp_r * alpha_L * velocity_norm), - // pipe o1, Eq. 20 - (lambda_r + rho_r * Cp_r * alpha_L * velocity_norm), - // pipe g1, Eq. 21 - (1.0 - porosity_g) * lambda_g, - // pipe g2, Eq. 22 - (1.0 - porosity_g) * lambda_g}}; - } + std::array<double, number_of_unknowns> pipeHeatConductions() const; - std::array<Eigen::Vector3d, number_of_unknowns> pipeAdvectionVectors() const - { - double const& rho_r = refrigerant.density; - double const& Cp_r = refrigerant.specific_heat_capacity; - - return {{// pipe i1, Eq. 19 - {0, 0, -rho_r * Cp_r * _flow_velocity}, - // pipe o1, Eq. 20 - {0, 0, rho_r * Cp_r * _flow_velocity}, - // grout g1, Eq. 21 - {0, 0, 0}, - // grout g2, Eq. 22 - {0, 0, 0}}}; - } + std::array<Eigen::Vector3d, number_of_unknowns> pipeAdvectionVectors() + const; template <int NPoints, typename SingleUnknownMatrixType, @@ -125,7 +78,7 @@ public: 1.0 * matBHE_loc_R; // K_i1 R_matrix.block(2 * NPoints, 2 * NPoints, NPoints, NPoints) += 1.0 * matBHE_loc_R; // K_ig - break; + return; case 1: // PHI_fog R_matrix.block(NPoints, 3 * NPoints, NPoints, NPoints) += -1.0 * matBHE_loc_R; @@ -136,7 +89,7 @@ public: 1.0 * matBHE_loc_R; // K_o1 R_matrix.block(3 * NPoints, 3 * NPoints, NPoints, NPoints) += 1.0 * matBHE_loc_R; // K_og - break; + return; case 2: // PHI_gg R_matrix.block(2 * NPoints, 3 * NPoints, NPoints, NPoints) += -1.0 * matBHE_loc_R; @@ -149,7 +102,7 @@ public: R_matrix.block(3 * NPoints, 3 * NPoints, NPoints, NPoints) += 1.0 * matBHE_loc_R; // K_og // see Diersch 2013 FEFLOW // book page 954 Table M.2 - break; + return; case 3: // PHI_gs R_s_matrix.template block<NPoints, NPoints>(0, 0).noalias() += 1.0 * matBHE_loc_R; @@ -162,15 +115,16 @@ public: 1.0 * matBHE_loc_R; // K_ig R_matrix.block(3 * NPoints, 3 * NPoints, NPoints, NPoints) += 1.0 * matBHE_loc_R; // K_og - break; + return; + default: + OGS_FATAL( + "Error!!! In the function BHE_1U::assembleRMatrices, " + "the index of bhe unknowns is out of range! "); } } - /** - * return the inflow temperature based on outflow temperature and fixed - * power. - */ - double getTinByTout(double T_out, double current_time); + /// Return the inflow temperature for the boundary condition. + double updateFlowRateAndTemperature(double T_out, double current_time); double thermalResistance(int const unknown_index) const { diff --git a/ProcessLib/HeatTransportBHE/BHE/MeshUtils.cpp b/ProcessLib/HeatTransportBHE/BHE/MeshUtils.cpp index 54059b5f75b..1243f4bed72 100644 --- a/ProcessLib/HeatTransportBHE/BHE/MeshUtils.cpp +++ b/ProcessLib/HeatTransportBHE/BHE/MeshUtils.cpp @@ -106,4 +106,4 @@ BHEMeshData getBHEDataInMesh(MeshLib::Mesh const& mesh) return {bhe_material_ids, bhe_elements, bhe_nodes}; } } // end of namespace HeatTransportBHE -} // namespace ProcessLib \ No newline at end of file +} // namespace ProcessLib diff --git a/ProcessLib/HeatTransportBHE/BHE/PipeConfiguration1U.h b/ProcessLib/HeatTransportBHE/BHE/PipeConfiguration1U.h index 72f4d8b4383..b71c3115f72 100644 --- a/ProcessLib/HeatTransportBHE/BHE/PipeConfiguration1U.h +++ b/ProcessLib/HeatTransportBHE/BHE/PipeConfiguration1U.h @@ -30,7 +30,7 @@ struct PipeConfiguration1U /// Distance between pipes. double const distance; - double const longitudinal_despersion_length; + double const longitudinal_dispersion_length; }; } // namespace BHE } // namespace HeatTransportBHE diff --git a/ProcessLib/HeatTransportBHE/HeatTransportBHEProcess.cpp b/ProcessLib/HeatTransportBHE/HeatTransportBHEProcess.cpp index 2d609b204ad..e940f8f73b1 100644 --- a/ProcessLib/HeatTransportBHE/HeatTransportBHEProcess.cpp +++ b/ProcessLib/HeatTransportBHE/HeatTransportBHEProcess.cpp @@ -44,7 +44,7 @@ HeatTransportBHEProcess::HeatTransportBHEProcess( { OGS_FATAL( "The number of the given BHE properties (%d) are not consistent " - "with the number of BHE groups in a mesh (%d).", + "with the number of BHE groups in the mesh (%d).", _process_data._vec_BHE_property.size(), _bheMeshData.BHE_mat_IDs.size()); } @@ -54,6 +54,9 @@ HeatTransportBHEProcess::HeatTransportBHEProcess( { OGS_FATAL("Not able to get material IDs! "); } + + _process_data._mesh_prop_materialIDs = material_ids; + // create a map from a material ID to a BHE ID for (int i = 0; i < static_cast<int>(_bheMeshData.BHE_mat_IDs.size()); i++) { @@ -61,8 +64,6 @@ HeatTransportBHEProcess::HeatTransportBHEProcess( _process_data._map_materialID_to_BHE_ID[_bheMeshData.BHE_mat_IDs[i]] = i; } - - _process_data._mesh_prop_materialIDs = material_ids; } void HeatTransportBHEProcess::constructDofTable() diff --git a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE-impl.h b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE-impl.h index 3fb23405916..dadaaef4337 100644 --- a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE-impl.h +++ b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE-impl.h @@ -28,7 +28,7 @@ HeatTransportBHELocalAssemblerBHE<ShapeFunction, IntegrationMethod, BHEType>:: : _process_data(process_data), _integration_method(integration_order), _bhe(bhe), - element_id(e.getID()) + _element_id(e.getID()) { // need to make sure that the BHE elements are one-dimensional assert(e.getDimension() == 1); diff --git a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE.h b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE.h index ed74e2601ef..0942665922d 100644 --- a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE.h +++ b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerBHE.h @@ -84,7 +84,7 @@ private: BHEType const& _bhe; - std::size_t const element_id; + std::size_t const _element_id; SecondaryData<typename ShapeMatrices::ShapeType> _secondary_data; diff --git a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil-impl.h b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil-impl.h index 058044c0dac..5184ff8db39 100644 --- a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil-impl.h +++ b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil-impl.h @@ -33,7 +33,7 @@ HeatTransportBHELocalAssemblerSoil<ShapeFunction, IntegrationMethod>:: HeatTransportBHEProcessData& process_data) : _process_data(process_data), _integration_method(integration_order), - element_id(e.getID()) + _element_id(e.getID()) { unsigned const n_integration_points = _integration_method.getNumberOfPoints(); @@ -48,7 +48,7 @@ HeatTransportBHELocalAssemblerSoil<ShapeFunction, IntegrationMethod>:: e, is_axially_symmetric, _integration_method); SpatialPosition x_position; - x_position.setElementID(element_id); + x_position.setElementID(_element_id); // ip data initialization for (unsigned ip = 0; ip < n_integration_points; ip++) @@ -86,7 +86,7 @@ void HeatTransportBHELocalAssemblerSoil<ShapeFunction, IntegrationMethod>:: _integration_method.getNumberOfPoints(); SpatialPosition pos; - pos.setElementID(element_id); + pos.setElementID(_element_id); for (unsigned ip = 0; ip < n_integration_points; ip++) { diff --git a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil.h b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil.h index 51088a08481..74f2d7fe668 100644 --- a/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil.h +++ b/ProcessLib/HeatTransportBHE/LocalAssemblers/HeatTransportBHELocalAssemblerSoil.h @@ -75,7 +75,7 @@ private: std::vector<ShapeMatrices, Eigen::aligned_allocator<ShapeMatrices>> _shape_matrices; - std::size_t const element_id; + std::size_t const _element_id; SecondaryData<typename ShapeMatrices::ShapeType> _secondary_data; }; diff --git a/ProcessLib/HeatTransportBHE/LocalAssemblers/SecondaryData.h b/ProcessLib/HeatTransportBHE/LocalAssemblers/SecondaryData.h index 4a7484928e3..19f87a2c59b 100644 --- a/ProcessLib/HeatTransportBHE/LocalAssemblers/SecondaryData.h +++ b/ProcessLib/HeatTransportBHE/LocalAssemblers/SecondaryData.h @@ -19,7 +19,7 @@ namespace ProcessLib { namespace HeatTransportBHE { -/// Used by for extrapolation of the integration point values. It is ordered +/// Used for extrapolation of the integration point values. It is ordered /// (and stored) by integration points. template <typename ShapeMatrixType> struct SecondaryData -- GitLab