Skip to content
Snippets Groups Projects
Commit 5b01d89a authored by Kilian Weishaupt's avatar Kilian Weishaupt
Browse files

[PNM] Add diffusive enthalpy transport

parent 9d46fcf1
No related branches found
No related tags found
1 merge request!2568[pnm] Various imrprovements
......@@ -26,13 +26,22 @@
#define DUMUX_FLUX_PNM_FOURIERS_LAW_HH
#include <dumux/common/math.hh>
#include <dumux/flux/referencesystemformulation.hh>
#include <type_traits>
namespace Dumux::PoreNetwork {
namespace Detail {
struct NoDiffusionType {};
} // end namespace Detail
/*!
* \ingroup PoreNetworkModels
* \brief Specialization of Fourier's Law for the pore-network model.
*/
template<class MolecularDiffusionType = Detail::NoDiffusionType>
struct PNMFouriersLaw
{
......@@ -68,7 +77,41 @@ struct PNMFouriersLaw
const Scalar gradT = deltaT/fluxVarsCache.throatLength();
heatflux += thermalConductivity*gradT*area;
if constexpr (!std::is_same_v<MolecularDiffusionType, Detail::NoDiffusionType>)
heatflux += componentEnthalpyFlux_(problem, element, fvGeometry, elemVolVars, scvf, elemFluxVarsCache, phaseIdx);
}
return heatflux;
}
private:
template<class Problem, class Element, class FVElementGeometry,
class ElementVolumeVariables, class ElementFluxVariablesCache>
static auto componentEnthalpyFlux_(const Problem& problem,
const Element& element,
const FVElementGeometry& fvGeometry,
const ElementVolumeVariables& elemVolVars,
const typename FVElementGeometry::SubControlVolumeFace& scvf,
const ElementFluxVariablesCache& elemFluxVarsCache,
const int phaseIdx)
{
using Scalar = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables::value_type;
Scalar heatflux = 0.0;
using FluidSystem = typename ElementVolumeVariables::VolumeVariables::FluidSystem;
const auto diffusiveFlux = MolecularDiffusionType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarsCache);
for (int compIdx = 0; compIdx < ElementVolumeVariables::VolumeVariables::numFluidComponents(); ++compIdx)
{
const bool insideIsUpstream = diffusiveFlux[compIdx] > 0.0;
const auto& upstreamVolVars = insideIsUpstream ? elemVolVars[scvf.insideScvIdx()] : elemVolVars[scvf.outsideScvIdx()];
const Scalar componentEnthalpy = FluidSystem::componentEnthalpy(upstreamVolVars.fluidState(), phaseIdx, compIdx);
if (MolecularDiffusionType::referenceSystemFormulation() == ReferenceSystemFormulation::massAveraged)
heatflux += diffusiveFlux[compIdx] * componentEnthalpy;
else
heatflux += diffusiveFlux[compIdx] * FluidSystem::molarMass(compIdx) * componentEnthalpy;
}
return heatflux;
}
};
......
......@@ -242,5 +242,12 @@ struct ThermalConductivityModel<TypeTag, TTag::PNMOnePNCNI>
using type = ThermalConductivityAverage<GetPropType<TypeTag, Properties::Scalar>>;
}; //!< Use the average for effective conductivities
// template<class TypeTag>
// struct HeatConductionType<TypeTag, TTag::PNMOnePNCNI>
// {
// TODO uncomment this as soon as there is a generalized approach for component enthalpies in all fluid systems
// using type = Dumux::PoreNetwork::PNMFouriersLaw<GetPropType<TypeTag, MolecularDiffusionType>>;
// }; //!< Use Fourier's law and also consider enthalpy transport by molecular diffusion
} // end namespace Dumux::Properties
#endif
......@@ -68,7 +68,7 @@ public:
};
template<class TypeTag>
struct HeatConductionType<TypeTag, TTag::PoreNetworkModel> { using type = Dumux::PoreNetwork::PNMFouriersLaw; };
struct HeatConductionType<TypeTag, TTag::PoreNetworkModel> { using type = Dumux::PoreNetwork::PNMFouriersLaw<>; };
//! The labels
template<class TypeTag>
......
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