Skip to content
Snippets Groups Projects
Commit bb506c1e authored by HBShaoUFZ's avatar HBShaoUFZ Committed by Dmitri Naumov
Browse files

[PL] move commons into base class

parent 193a183d
No related branches found
No related tags found
No related merge requests found
Showing with 203 additions and 356 deletions
/**
* \file
*
* \copyright
* Copyright (c) 2012-2019, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#include "BHECommonCoaxial.h"
#include "FlowAndTemperatureControl.h"
#include "ThermoMechanicalFlowProperties.h"
namespace ProcessLib
{
namespace HeatTransportBHE
{
namespace BHE
{
constexpr std::pair<int, int>
BHECommonCoaxial::inflow_outflow_bc_component_ids[];
std::array<double, BHECommonCoaxial::number_of_unknowns>
BHECommonCoaxial::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 {{/*i*/ rho_r * specific_heat_capacity,
/*o*/ rho_r * specific_heat_capacity,
/*g*/ (1.0 - porosity_g) * rho_g * heat_cap_g}};
}
double BHECommonCoaxial::updateFlowRateAndTemperature(double const T_out,
double const current_time)
{
auto values = apply_visitor(
[&](auto const& control) { return control(T_out, current_time); },
flowAndTemperatureControl);
updateHeatTransferCoefficients(values.flow_rate);
return values.temperature;
}
std::array<double, BHECommonCoaxial::number_of_unknowns>
BHECommonCoaxial::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;
auto v = velocities();
// 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. 26-28, 23-25
return {{// pipe i, Eq. 26 and Eq. 23
(lambda_r + rho_r * Cp_r * alpha_L * std::abs(v[0])),
// pipe o, Eq. 27 and Eq. 24
(lambda_r + rho_r * Cp_r * alpha_L * std::abs(v[1])),
// pipe g, Eq. 28 and Eq. 25
(1.0 - porosity_g) * lambda_g}};
}
std::array<Eigen::Vector3d, BHECommonCoaxial::number_of_unknowns>
BHECommonCoaxial::pipeAdvectionVectors() const
{
double const& rho_r = refrigerant.density;
double const& Cp_r = refrigerant.specific_heat_capacity;
auto v = velocities();
return {{// pipe i, Eq. 26 and Eq. 23
{0, 0, -rho_r * Cp_r * std::abs(v[0])},
// pipe o, Eq. 27 and Eq. 24
{0, 0, rho_r * Cp_r * std::abs(v[1])},
// grout g, Eq. 28 and Eq. 25
{0, 0, 0}}};
}
} // namespace BHE
} // namespace HeatTransportBHE
} // namespace ProcessLib
......@@ -9,8 +9,11 @@
#pragma once
#include <Eigen/Eigen>
#include "BHECommon.h"
#include "FlowAndTemperatureControl.h"
#include "PipeConfigurationCoaxial.h"
#include "ThermoMechanicalFlowProperties.h"
namespace ProcessLib
{
......@@ -18,10 +21,10 @@ namespace HeatTransportBHE
{
namespace BHE
{
class BHECoaxialCommon : public BHECommon
class BHECommonCoaxial : public BHECommon
{
public:
BHECoaxialCommon(BoreholeGeometry const& borehole,
BHECommonCoaxial(BoreholeGeometry const& borehole,
RefrigerantProperties const& refrigerant,
GroutParameters const& grout,
FlowAndTemperatureControl const& flowAndTemperatureControl,
......@@ -30,12 +33,60 @@ public:
_pipes(pipes)
{
}
double thermalResistance(int const unknown_index) const
{
return _thermal_resistances[unknown_index];
};
double updateFlowRateAndTemperature(double T_out, double current_time);
static constexpr int number_of_unknowns = 3;
static constexpr int number_of_grout_zones = 1;
std::array<double, number_of_unknowns> pipeHeatCapacities() const;
static constexpr std::pair<int, int> inflow_outflow_bc_component_ids[] = {
{0, 1}};
std::array<double, number_of_unknowns> pipeHeatConductions() const;
std::array<Eigen::Vector3d, number_of_unknowns> pipeAdvectionVectors()
const;
protected:
virtual std::array<double, number_of_unknowns> calcThermalResistances(
double const Nu_inner_pipe, double const Nu_annulus_pipe) = 0;
void updateHeatTransferCoefficients(double const flow_rate)
{
auto const tm_flow_properties_annulus =
calculateThermoMechanicalFlowPropertiesAnnulus(
_pipes.inner_pipe,
_pipes.outer_pipe,
borehole_geometry.length,
refrigerant,
flow_rate);
_flow_velocity_annulus = tm_flow_properties_annulus.velocity;
auto const tm_flow_properties =
calculateThermoMechanicalFlowPropertiesPipe(
_pipes.inner_pipe,
borehole_geometry.length,
refrigerant,
flow_rate);
_flow_velocity_inner = tm_flow_properties.velocity;
_thermal_resistances =
calcThermalResistances(tm_flow_properties.nusselt_number,
tm_flow_properties_annulus.nusselt_number);
}
PipeConfigurationCoaxial const _pipes;
virtual std::array<double, 2> velocities() const = 0;
/// Here we store the thermal resistances needed for computation of the heat
/// exchange coefficients in the governing equations of BHE.
/// These governing equations can be found in
......@@ -44,8 +95,8 @@ protected:
std::array<double, number_of_unknowns> _thermal_resistances;
/// Flow velocity inside the pipes and annulus. Depends on the flow_rate.
double _flow_velocity, _flow_velocity_annulus;
double _flow_velocity_inner, _flow_velocity_annulus;
};
} // end of namespace BHE
} // end of namespace HeatTransportBHE
} // namespace BHE
} // namespace HeatTransportBHE
} // end of namespace ProcessLib
......@@ -8,10 +8,6 @@
*/
#include "BHE_CXA.h"
#include <boost/math/constants/constants.hpp>
#include "FlowAndTemperatureControl.h"
#include "Physics.h"
#include "ThermoMechanicalFlowProperties.h"
namespace ProcessLib
......@@ -25,7 +21,7 @@ BHE_CXA::BHE_CXA(BoreholeGeometry const& borehole,
GroutParameters const& grout,
FlowAndTemperatureControl const& flowAndTemperatureControl,
PipeConfigurationCoaxial const& pipes)
: BHECoaxialCommon{borehole, refrigerant, grout, flowAndTemperatureControl,
: BHECommonCoaxial{borehole, refrigerant, grout, flowAndTemperatureControl,
pipes}
{
// Initialize thermal resistances.
......@@ -38,91 +34,15 @@ BHE_CXA::BHE_CXA(BoreholeGeometry const& borehole,
updateHeatTransferCoefficients(values.flow_rate);
}
std::array<double, BHE_CXA::number_of_unknowns> BHE_CXA::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 {{/*i*/ rho_r * specific_heat_capacity,
/*o*/ rho_r * specific_heat_capacity,
/*g*/ (1.0 - porosity_g) * rho_g * heat_cap_g}};
}
std::array<double, BHE_CXA::number_of_unknowns> BHE_CXA::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_inner_pipe = std::abs(_flow_velocity);
double const velocity_norm_annulus = std::abs(_flow_velocity_annulus);
// 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. 23-25.
return {{// pipe i, Eq. 23
(lambda_r + rho_r * Cp_r * alpha_L * velocity_norm_annulus),
// pipe o, Eq. 24
(lambda_r + rho_r * Cp_r * alpha_L * velocity_norm_inner_pipe),
// pipe g, Eq. 25
(1.0 - porosity_g) * lambda_g}};
}
std::array<Eigen::Vector3d, BHE_CXA::number_of_unknowns>
BHE_CXA::pipeAdvectionVectors() const
{
double const& rho_r = refrigerant.density;
double const& Cp_r = refrigerant.specific_heat_capacity;
return {{// pipe i, Eq. 23
{0, 0, -rho_r * Cp_r * _flow_velocity_annulus},
// pipe o, Eq. 24
{0, 0, rho_r * Cp_r * _flow_velocity},
// grout g, Eq. 25
{0, 0, 0}}};
}
constexpr std::pair<int, int> BHE_CXA::inflow_outflow_bc_component_ids[];
void BHE_CXA::updateHeatTransferCoefficients(double const flow_rate)
{
auto const tm_flow_properties_annulus =
calculateThermoMechanicalFlowPropertiesAnnulus(_pipes.inner_pipe,
_pipes.outer_pipe,
borehole_geometry.length,
refrigerant,
flow_rate);
_flow_velocity_annulus = tm_flow_properties_annulus.velocity;
auto const tm_flow_properties = calculateThermoMechanicalFlowPropertiesPipe(
_pipes.inner_pipe, borehole_geometry.length, refrigerant, flow_rate);
_flow_velocity = tm_flow_properties.velocity;
_thermal_resistances =
calcThermalResistances(tm_flow_properties.nusselt_number,
tm_flow_properties_annulus.nusselt_number);
}
/// Nu_o is the Nusselt number of inner pipe, Nu_i is the Nusselt number of
/// annulus.
std::array<double, BHE_CXA::number_of_unknowns> BHE_CXA::calcThermalResistances(
double const Nu_inner_outflow_pipe, double const Nu_annulus)
double const Nu_inner_pipe, double const Nu_annulus_pipe)
{
// thermal resistances due to advective flow of refrigerant in the pipes
auto const R_advective = calculateAdvectiveThermalResistance(
_pipes.inner_pipe, _pipes.outer_pipe, refrigerant,
Nu_inner_outflow_pipe, Nu_annulus);
Nu_inner_pipe, Nu_annulus_pipe);
// thermal resistance due to thermal conductivity of the pipe wall material
auto const R_conductive = calculatePipeWallThermalResistance(
......@@ -142,25 +62,6 @@ std::array<double, BHE_CXA::number_of_unknowns> BHE_CXA::calcThermalResistances(
R_advective.b_annulus + R_conductive.annulus + R.conductive_b;
return {{R_fig, R_ff, R_gs}};
// keep the following lines------------------------------------------------
// when debuffing the code, printing the R and phi values are needed--------
// std::cout << "Rfig =" << R_fig << " Rff =" << R_ff << " Rgs =" << R_gs
// << "\n"; double phi_fig = 1.0 / (R_fig * S_i);
// double phi_ff = 1.0 / (R_ff * S_g1); double phi_gs = 1.0 / (R_gs * S_gs);
// std::cout << "phi_fig ="
// << phi_fig << " phi_ff =" << phi_ff << " phi_gs =" << phi_gs << "\n";
// -------------------------------------------------------------------------
}
double BHE_CXA::updateFlowRateAndTemperature(double const T_out,
double const current_time)
{
auto values = apply_visitor(
[&](auto const& control) { return control(T_out, current_time); },
flowAndTemperatureControl);
updateHeatTransferCoefficients(values.flow_rate);
return values.temperature;
}
} // namespace BHE
} // namespace HeatTransportBHE
......
......@@ -37,7 +37,7 @@ namespace BHE
* are calculated specifically during the initialization of the class.
*/
class BHE_CXA final : public BHECoaxialCommon
class BHE_CXA final : public BHECommonCoaxial
{
public:
BHE_CXA(BoreholeGeometry const& borehole,
......@@ -46,13 +46,6 @@ public:
FlowAndTemperatureControl const& flowAndTemperatureControl,
PipeConfigurationCoaxial const& pipes);
std::array<double, number_of_unknowns> pipeHeatCapacities() const;
std::array<double, number_of_unknowns> pipeHeatConductions() const;
std::array<Eigen::Vector3d, number_of_unknowns> pipeAdvectionVectors()
const;
template <int NPoints, typename SingleUnknownMatrixType,
typename RMatrixType, typename RPiSMatrixType,
typename RSMatrixType>
......@@ -105,17 +98,6 @@ public:
}
}
/// Return the inflow temperature for the boundary condition.
double updateFlowRateAndTemperature(double T_out, double current_time);
double thermalResistance(int const unknown_index) const
{
return _thermal_resistances[unknown_index];
}
static constexpr std::pair<int, int> inflow_outflow_bc_component_ids[] = {
{0, 1}};
public:
std::array<double, number_of_unknowns> const cross_section_areas = {
{_pipes.outer_pipe.area() - _pipes.inner_pipe.outsideArea(),
......@@ -123,10 +105,13 @@ public:
borehole_geometry.area() - _pipes.outer_pipe.outsideArea()}};
private:
void updateHeatTransferCoefficients(double const flow_rate);
std::array<double, number_of_unknowns> calcThermalResistances(
double const Nu_o, double const Nu_i);
double const Nu_inner_pipe, double const Nu_annulus_pipe) override;
std::array<double, 2> velocities() const override
{
return {_flow_velocity_annulus, _flow_velocity_inner};
}
};
} // namespace BHE
} // namespace HeatTransportBHE
......
......@@ -8,10 +8,6 @@
*/
#include "BHE_CXC.h"
#include <boost/math/constants/constants.hpp>
#include "FlowAndTemperatureControl.h"
#include "Physics.h"
#include "ThermoMechanicalFlowProperties.h"
namespace ProcessLib
......@@ -25,7 +21,7 @@ BHE_CXC::BHE_CXC(BoreholeGeometry const& borehole,
GroutParameters const& grout,
FlowAndTemperatureControl const& flowAndTemperatureControl,
PipeConfigurationCoaxial const& pipes)
: BHECoaxialCommon{borehole, refrigerant, grout, flowAndTemperatureControl,
: BHECommonCoaxial{borehole, refrigerant, grout, flowAndTemperatureControl,
pipes}
{
// Initialize thermal resistances.
......@@ -38,94 +34,18 @@ BHE_CXC::BHE_CXC(BoreholeGeometry const& borehole,
updateHeatTransferCoefficients(values.flow_rate);
}
std::array<double, BHE_CXC::number_of_unknowns> BHE_CXC::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 {{/*i*/ rho_r * specific_heat_capacity,
/*o*/ rho_r * specific_heat_capacity,
/*g*/ (1.0 - porosity_g) * rho_g * heat_cap_g}};
}
std::array<double, BHE_CXC::number_of_unknowns> BHE_CXC::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_inner_pipe = std::abs(_flow_velocity);
double const velocity_norm_annulus = std::abs(_flow_velocity_annulus);
// 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. 26-28.
return {{// pipe i, Eq. 26
(lambda_r + rho_r * Cp_r * alpha_L * velocity_norm_inner_pipe),
// pipe o, Eq. 27
(lambda_r + rho_r * Cp_r * alpha_L * velocity_norm_annulus),
// pipe g, Eq. 28
(1.0 - porosity_g) * lambda_g}};
}
std::array<Eigen::Vector3d, BHE_CXC::number_of_unknowns>
BHE_CXC::pipeAdvectionVectors() const
{
double const& rho_r = refrigerant.density;
double const& Cp_r = refrigerant.specific_heat_capacity;
return {{// pipe i, Eq. 26
{0, 0, -rho_r * Cp_r * _flow_velocity},
// pipe o, Eq. 27
{0, 0, rho_r * Cp_r * _flow_velocity_annulus},
// grout g, Eq. 28
{0, 0, 0}}};
}
constexpr std::pair<int, int> BHE_CXC::inflow_outflow_bc_component_ids[];
void BHE_CXC::updateHeatTransferCoefficients(double const flow_rate)
{
auto const tm_flow_properties_annulus =
calculateThermoMechanicalFlowPropertiesAnnulus(_pipes.inner_pipe,
_pipes.outer_pipe,
borehole_geometry.length,
refrigerant,
flow_rate);
_flow_velocity_annulus = tm_flow_properties_annulus.velocity;
auto const tm_flow_properties = calculateThermoMechanicalFlowPropertiesPipe(
_pipes.inner_pipe, borehole_geometry.length, refrigerant, flow_rate);
_flow_velocity = tm_flow_properties.velocity;
_thermal_resistances =
calcThermalResistances(tm_flow_properties_annulus.nusselt_number,
tm_flow_properties.nusselt_number);
}
/// Nu_o is the Nusselt number of inner pipe, Nu_i is the Nusselt number of
/// annulus.
std::array<double, BHE_CXC::number_of_unknowns> BHE_CXC::calcThermalResistances(
double const Nu_annulus, double const Nu_inner_inflow_pipe)
double const Nu_inner_pipe, double const Nu_annulus_pipe)
{
// thermal resistances due to advective flow of refrigerant in the pipes
auto const R_advective =
calculateAdvectiveThermalResistance(_pipes.inner_pipe,
_pipes.outer_pipe,
refrigerant,
Nu_inner_inflow_pipe,
Nu_annulus);
Nu_inner_pipe,
Nu_annulus_pipe);
// thermal resistance due to thermal conductivity of the pipe wall material
auto const R_conductive = calculatePipeWallThermalResistance(
......@@ -145,25 +65,6 @@ std::array<double, BHE_CXC::number_of_unknowns> BHE_CXC::calcThermalResistances(
R_advective.b_annulus + R_conductive.annulus + R.conductive_b;
return {{R_ff, R_fog, R_gs}};
// keep the following lines------------------------------------------------
// when debuffing the code, printing the R and phi values are needed--------
// std::cout << "Rfog =" << R_fog << " Rff =" << R_ff << " Rgs =" << R_gs
// << "\n"; double phi_fog = 1.0 / (R_fog * S_i);
// double phi_ff = 1.0 / (R_ff * S_g1); double phi_gs = 1.0 / (R_gs * S_gs);
// std::cout << "phi_fog ="
// << phi_fog << " phi_ff =" << phi_ff << " phi_gs =" << phi_gs << "\n";
// -------------------------------------------------------------------------
}
double BHE_CXC::updateFlowRateAndTemperature(double const T_out,
double const current_time)
{
auto values = apply_visitor(
[&](auto const& control) { return control(T_out, current_time); },
flowAndTemperatureControl);
updateHeatTransferCoefficients(values.flow_rate);
return values.temperature;
}
} // namespace BHE
} // namespace HeatTransportBHE
......
......@@ -36,7 +36,7 @@ namespace BHE
* surrounding soil is regulated through the thermal resistance values, which
* are calculated specifically during the initialization of the class.
*/
class BHE_CXC final : public BHECoaxialCommon
class BHE_CXC final : public BHECommonCoaxial
{
public:
BHE_CXC(BoreholeGeometry const& borehole,
......@@ -45,13 +45,6 @@ public:
FlowAndTemperatureControl const& flowAndTemperatureControl,
PipeConfigurationCoaxial const& pipes);
std::array<double, number_of_unknowns> pipeHeatCapacities() const;
std::array<double, number_of_unknowns> pipeHeatConductions() const;
std::array<Eigen::Vector3d, number_of_unknowns> pipeAdvectionVectors()
const;
template <int NPoints, typename SingleUnknownMatrixType,
typename RMatrixType, typename RPiSMatrixType,
typename RSMatrixType>
......@@ -104,17 +97,6 @@ public:
}
}
/// Return the inflow temperature for the boundary condition.
double updateFlowRateAndTemperature(double T_out, double current_time);
double thermalResistance(int const unknown_index) const
{
return _thermal_resistances[unknown_index];
}
static constexpr std::pair<int, int> inflow_outflow_bc_component_ids[] = {
{0, 1}};
public:
std::array<double, number_of_unknowns> const cross_section_areas = {
{_pipes.inner_pipe.area(),
......@@ -122,10 +104,13 @@ public:
borehole_geometry.area() - _pipes.outer_pipe.outsideArea()}};
private:
void updateHeatTransferCoefficients(double const flow_rate);
std::array<double, number_of_unknowns> calcThermalResistances(
double const Nu_o, double const Nu_i);
double const Nu_inner_pipe, double const Nu_annulus_pipe) override;
std::array<double, 2> velocities() const override
{
return {_flow_velocity_inner, _flow_velocity_annulus};
}
};
} // namespace BHE
} // namespace HeatTransportBHE
......
/**
* \file
*
* \copyright
* Copyright (c) 2012-2019, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#include "CreateBHECXC.h"
#include "BaseLib/ConfigTree.h"
#include "CreateFlowAndTemperatureControl.h"
namespace ProcessLib
{
namespace HeatTransportBHE
{
namespace BHE
{
BHE::BHE_CXC createBHECXC(
BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves)
{
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__borehole}
auto const borehole_geometry =
createBoreholeGeometry(config.getConfigSubtree("borehole"));
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__pipes}
auto const& pipes_config = config.getConfigSubtree("pipes");
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__pipes__outer}
Pipe const outer_pipe = createPipe(pipes_config.getConfigSubtree("outer"));
Pipe const inner_pipe =
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__pipes__inner}
createPipe(pipes_config.getConfigSubtree("inner"));
const double pipe_longitudinal_dispersion_length =
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__pipes__longitudinal_dispersion_length}
pipes_config.getConfigParameter<double>(
"longitudinal_dispersion_length");
PipeConfigurationCoaxial const pipes{inner_pipe, outer_pipe,
pipe_longitudinal_dispersion_length};
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__grout}
auto const grout = createGroutParameters(config.getConfigSubtree("grout"));
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__refrigerant}
auto const refrigerant =
createRefrigerantProperties(config.getConfigSubtree("refrigerant"));
auto const flowAndTemperatureControl = createFlowAndTemperatureControl(
//! \ogs_file_param{prj__processes__process__HEAT_TRANSPORT_BHE__borehole_heat_exchangers__borehole_heat_exchanger__flow_and_temperature_control}
config.getConfigSubtree("flow_and_temperature_control"),
curves,
refrigerant);
return {borehole_geometry, refrigerant, grout, flowAndTemperatureControl,
pipes};
}
} // namespace BHE
} // namespace HeatTransportBHE
} // namespace ProcessLib
/**
* \copyright
* Copyright (c) 2012-2019, OpenGeoSys Community (http://www.opengeosys.org)
* Distributed under a Modified BSD License.
* See accompanying file LICENSE.txt or
* http://www.opengeosys.org/project/license
*
*/
#pragma once
#include "BHE_CXC.h"
namespace BaseLib
{
class ConfigTree;
}
namespace ProcessLib
{
namespace HeatTransportBHE
{
namespace BHE
{
BHE::BHE_CXC createBHECXC(
BaseLib::ConfigTree const& bhe_conf,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves);
} // namespace BHE
} // namespace HeatTransportBHE
} // namespace ProcessLib
......@@ -9,7 +9,7 @@
*
*/
#include "CreateBHECXA.h"
#include "CreateBHECoaxial.h"
#include "BaseLib/ConfigTree.h"
#include "CreateFlowAndTemperatureControl.h"
namespace ProcessLib
......@@ -18,7 +18,12 @@ namespace HeatTransportBHE
{
namespace BHE
{
BHE::BHE_CXA createBHECXA(
static std::tuple<BoreholeGeometry,
RefrigerantProperties,
GroutParameters,
FlowAndTemperatureControl,
PipeConfigurationCoaxial>
parseBHECoaxialConfig(
BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
......@@ -58,6 +63,30 @@ BHE::BHE_CXA createBHECXA(
return {borehole_geometry, refrigerant, grout, flowAndTemperatureControl,
pipes};
}
template <typename T_BHE>
T_BHE createBHECoaxial(
BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves)
{
auto coaxial = parseBHECoaxialConfig(config, curves);
return {std::get<0>(coaxial), std::get<1>(coaxial), std::get<2>(coaxial),
std::get<3>(coaxial), std::get<4>(coaxial)};
}
template BHE_CXA createBHECoaxial<BHE_CXA>(
BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves);
template BHE_CXC createBHECoaxial<BHE_CXC>(
BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves);
} // namespace BHE
} // namespace HeatTransportBHE
} // namespace ProcessLib
......@@ -9,7 +9,9 @@
#pragma once
#include "BHECommonCoaxial.h"
#include "BHE_CXA.h"
#include "BHE_CXC.h"
namespace BaseLib
{
......@@ -21,8 +23,9 @@ namespace HeatTransportBHE
{
namespace BHE
{
BHE::BHE_CXA createBHECXA(
BaseLib::ConfigTree const& bhe_conf,
template <typename T_BHE>
T_BHE createBHECoaxial(
BaseLib::ConfigTree const& config,
std::map<std::string,
std::unique_ptr<MathLib::PiecewiseLinearInterpolation>> const&
curves);
......
......@@ -16,8 +16,7 @@
#include "BHE/BHETypes.h"
#include "BHE/CreateBHE1U.h"
#include "BHE/CreateBHECXA.h"
#include "BHE/CreateBHECXC.h"
#include "BHE/CreateBHECoaxial.h"
#include "HeatTransportBHEProcess.h"
#include "HeatTransportBHEProcessData.h"
......@@ -191,13 +190,13 @@ std::unique_ptr<Process> createHeatTransportBHEProcess(
if (bhe_type == "CXA")
{
bhes.push_back(BHE::createBHECXA(bhe_config, curves));
bhes.push_back(BHE::createBHECoaxial<BHE::BHE_CXA>(bhe_config, curves));
continue;
}
if (bhe_type == "CXC")
{
bhes.push_back(BHE::createBHECXC(bhe_config, curves));
bhes.push_back(BHE::createBHECoaxial<BHE::BHE_CXC>(bhe_config, curves));
continue;
}
OGS_FATAL("Unknown BHE type '%s'.", bhe_type.c_str());
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment