25#ifndef DUMUX_FREEFLOW_NC_VOLUMEVARIABLES_HH
26#define DUMUX_FREEFLOW_NC_VOLUMEVARIABLES_HH
29#include <dune/common/exceptions.hh>
38template <
class Traits>
44 using Scalar =
typename Traits::PrimaryVariables::value_type;
46 static constexpr bool useMoles = Traits::ModelTraits::useMoles();
65 template<
class ElementSolution,
class Problem,
class Element,
class SubControlVolume>
66 void update(
const ElementSolution &elemSol,
67 const Problem &problem,
68 const Element &element,
69 const SubControlVolume& scv)
71 ParentType::update(elemSol, problem, element, scv);
75 typename FluidSystem::ParameterCache paramCache;
78 for (
unsigned int compIIdx = 0; compIIdx < ParentType::numFluidComponents(); ++compIIdx)
80 for (
unsigned int compJIdx = 0; compJIdx < ParentType::numFluidComponents(); ++compJIdx)
83 if(compIIdx != compJIdx)
86 = FluidSystem::binaryDiffusionCoefficient(
fluidState_,
99 template<
class ElementSolution,
class Problem,
class Element,
class SubControlVolume>
101 const Problem& problem,
102 const Element& element,
103 const SubControlVolume& scv,
107 fluidState.setPressure(0, elemSol[0][Indices::pressureIdx]);
114 Scalar sumFracMinorComp = 0.0;
116 for(
int compIdx = 1; compIdx < ParentType::numFluidComponents(); ++compIdx)
120 Scalar moleOrMassFraction = elemSol[0][Indices::conti0EqIdx+compIdx] + 1.0;
121 moleOrMassFraction = moleOrMassFraction - 1.0;
123 fluidState.setMoleFraction(0, compIdx, moleOrMassFraction);
125 fluidState.setMassFraction(0, compIdx, moleOrMassFraction);
126 sumFracMinorComp += moleOrMassFraction;
129 fluidState.setMoleFraction(0, 0, 1.0 - sumFracMinorComp);
131 fluidState.setMassFraction(0, 0, 1.0 - sumFracMinorComp);
133 typename FluidSystem::ParameterCache paramCache;
146 const Scalar h = ParentType::enthalpy(
fluidState, paramCache);
245 return FluidSystem::molarMass(compIdx);
264 if (compIIdx == compJIdx)
265 DUNE_THROW(Dune::InvalidStateException,
"Diffusion coefficient called for compIIdx = compJIdx");
289 std::array<std::array<Scalar, ParentType::numFluidComponents()>, ParentType::numFluidComponents()>
diffCoefficient_ = {};
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Property tag Indices
Definition: porousmediumflow/sequential/properties.hh:59
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
std::string viscosity(int phaseIdx) noexcept
I/O name of viscosity for multiphase systems.
Definition: name.hh:74
std::string molarDensity(int phaseIdx) noexcept
I/O name of molar density for multiphase systems.
Definition: name.hh:83
std::string density(int phaseIdx) noexcept
I/O name of density for multiphase systems.
Definition: name.hh:65
Volume variables for the single-phase, multi-component free-flow model.
Definition: freeflow/compositional/volumevariables.hh:40
Scalar moleFraction(int phaseIdx, int compIdx) const
Returns the mole fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:215
const FluidState & fluidState() const
Return the fluid state of the control volume.
Definition: freeflow/compositional/volumevariables.hh:284
Scalar diffusionCoefficient(int compIIdx, int compJIdx=0) const
Returns the diffusion coefficient .
Definition: freeflow/compositional/volumevariables.hh:262
Scalar massFraction(int phaseIdx, int compIdx) const
Returns the mass fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:194
std::array< std::array< Scalar, ParentType::numFluidComponents()>, ParentType::numFluidComponents()> diffCoefficient_
Definition: freeflow/compositional/volumevariables.hh:289
Scalar effectiveViscosity() const
Return the effective dynamic viscosity of the fluid within the control volume.
Definition: freeflow/compositional/volumevariables.hh:178
FluidState fluidState_
Definition: freeflow/compositional/volumevariables.hh:288
typename Traits::ModelTraits::Indices Indices
export the indices type
Definition: freeflow/compositional/volumevariables.hh:54
Scalar temperature() const
Return temperature inside the sub-control volume.
Definition: freeflow/compositional/volumevariables.hh:171
static void completeFluidState(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv, FluidState &fluidState)
Update the fluid state.
Definition: freeflow/compositional/volumevariables.hh:100
typename Traits::FluidState FluidState
export the fluid state type
Definition: freeflow/compositional/volumevariables.hh:52
void update(const ElementSolution &elemSol, const Problem &problem, const Element &element, const SubControlVolume &scv)
Update all quantities for a given control volume.
Definition: freeflow/compositional/volumevariables.hh:66
Scalar viscosity(int phaseIdx=0) const
Return the dynamic viscosity of the fluid within the control volume.
Definition: freeflow/compositional/volumevariables.hh:185
Scalar pressure(int phaseIdx=0) const
Return the effective pressure of a given phase within the control volume.
Definition: freeflow/compositional/volumevariables.hh:154
Scalar molarMass(int compIdx) const
Returns the molar mass of a given component.
Definition: freeflow/compositional/volumevariables.hh:243
Scalar moleFraction(int compIdx) const
Returns the mole fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:225
Scalar density(int phaseIdx=0) const
Return the mass density of a given phase within the control volume.
Definition: freeflow/compositional/volumevariables.hh:161
Scalar massFraction(int compIdx) const
Returns the mass fraction of a component in the phase .
Definition: freeflow/compositional/volumevariables.hh:204
Scalar effectiveDiffusivity(int compIIdx, int compJIdx=0) const
Returns the effective diffusion coefficient .
Definition: freeflow/compositional/volumevariables.hh:275
Scalar molarDensity(int phaseIdx=0) const
Returns the mass density of a given phase .
Definition: freeflow/compositional/volumevariables.hh:233
typename Traits::FluidSystem FluidSystem
export the underlying fluid system
Definition: freeflow/compositional/volumevariables.hh:50
Scalar averageMolarMass(const int phaseIdx=0) const
Returns the average molar mass of the fluid phase.
Definition: freeflow/compositional/volumevariables.hh:253
Definition: freeflow/volumevariables.hh:34