From c2d568514f1088b394e1592d8ec51293f67754b3 Mon Sep 17 00:00:00 2001 From: Katherina Baber <katherina.baber@gmail.com> Date: Tue, 15 May 2012 11:45:23 +0000 Subject: [PATCH] further name changes of the indices transportCompIdx (former x1Idx) is now massOrMoleFracIdx comp0Idx is now phaseCompIdx comp1Idx is now transportCompIdx reviewed by Klaus git-svn-id: svn://svn.iws.uni-stuttgart.de/DUMUX/dumux/trunk@8313 2fb0f335-1f38-0410-981e-8018bf24f1b0 --- dumux/boxmodels/1p2c/1p2cfluxvariables.hh | 12 ++++----- dumux/boxmodels/1p2c/1p2cindices.hh | 19 +++++++------- dumux/boxmodels/1p2c/1p2clocalresidual.hh | 24 +++++++++--------- dumux/boxmodels/1p2c/1p2cproperties.hh | 1 + dumux/boxmodels/1p2c/1p2cpropertydefaults.hh | 5 +++- dumux/boxmodels/1p2c/1p2cvolumevariables.hh | 26 ++++++++++---------- test/boxmodels/1p2c/1p2coutflowproblem.hh | 6 ++--- 7 files changed, 49 insertions(+), 44 deletions(-) diff --git a/dumux/boxmodels/1p2c/1p2cfluxvariables.hh b/dumux/boxmodels/1p2c/1p2cfluxvariables.hh index 21c3dd7b26..376b88bfcd 100644 --- a/dumux/boxmodels/1p2c/1p2cfluxvariables.hh +++ b/dumux/boxmodels/1p2c/1p2cfluxvariables.hh @@ -62,7 +62,7 @@ class OnePTwoCFluxVariables typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices; enum { - comp1Idx = Indices::comp1Idx + transportCompIdx = Indices::transportCompIdx }; typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView; @@ -361,17 +361,17 @@ protected: // the concentration gradient [mol/m^3/m] tmp = feGrad; - tmp *= elemVolVars[volVarsIdx].molarity(comp1Idx); + tmp *= elemVolVars[volVarsIdx].molarity(transportCompIdx); concentrationGrad_ += tmp; // the mole-fraction gradient tmp = feGrad; - tmp *= elemVolVars[volVarsIdx].moleFraction(comp1Idx); + tmp *= elemVolVars[volVarsIdx].moleFraction(transportCompIdx); moleFractionGrad_ += tmp; // the mass-fraction gradient tmp = feGrad; - tmp *= elemVolVars[volVarsIdx].massFraction(comp1Idx); + tmp *= elemVolVars[volVarsIdx].massFraction(transportCompIdx); massFractionGrad_ += tmp; // phase viscosity @@ -398,9 +398,9 @@ protected: potentialGrad_ = tmp; potentialGrad_ *= volVarsJ.pressure() - volVarsI.pressure(); concentrationGrad_ = tmp; - concentrationGrad_ *= volVarsJ.molarity(comp1Idx) - volVarsI.molarity(comp1Idx); + concentrationGrad_ *= volVarsJ.molarity(transportCompIdx) - volVarsI.molarity(transportCompIdx); moleFractionGrad_ = tmp; - moleFractionGrad_ *= volVarsJ.moleFraction(comp1Idx) - volVarsI.moleFraction(comp1Idx); + moleFractionGrad_ *= volVarsJ.moleFraction(transportCompIdx) - volVarsI.moleFraction(transportCompIdx); } /////////////// diff --git a/dumux/boxmodels/1p2c/1p2cindices.hh b/dumux/boxmodels/1p2c/1p2cindices.hh index 44e72b5db7..f525910cfc 100644 --- a/dumux/boxmodels/1p2c/1p2cindices.hh +++ b/dumux/boxmodels/1p2c/1p2cindices.hh @@ -40,17 +40,18 @@ namespace Dumux * \ingroup BoxIndices * \brief The indices for the isothermal single-phase, two-component model. */ -template <int PVOffset = 0> +template <class TypeTag, int PVOffset = 0> struct OnePTwoCIndices { //! Set the default phase used by the fluid system to the first one - static const int phaseIdx = 0; + static const int phaseIdx = GET_PROP_VALUE(TypeTag, PhaseIdx);; - //! Set the default for the first component the fluid system's first component - static const int comp0Idx = 0; - //! Set the default for the second component the fluid system's second component - static const int comp1Idx = 1; + //! Component indices + static const int phaseCompIdx = phaseIdx;//!< The index of the main component of the considered phase + static const int transportCompIdx = (unsigned int)(1-phaseIdx); //!< The index of the transported (minor) component; ASSUMES phase indices of 0 and 1 + static const int comp0Idx = phaseCompIdx;//!< \deprecated use phaseCompIdx + static const int comp1Idx = transportCompIdx;//!< \deprecated use phaseCompIdx // Equation indices static const int conti0EqIdx = PVOffset + 0; //!< continuity equation index @@ -61,9 +62,9 @@ struct OnePTwoCIndices // primary variable indices static const int pressureIdx = PVOffset + 0; //!< pressure - static const int transportCompIdx = PVOffset + 1; //!< mole fraction of the second component - static const int x1Idx = transportCompIdx; // \deprecated use transportCompIdx instead - static const int massOrMoleFractionIdx = transportCompIdx; // \deprecated use transportCompIdx instead + static const int massOrMoleFracIdx = PVOffset + 1; //!< mole fraction of the second component + static const int x1Idx = massOrMoleFracIdx; // \deprecated use massOrMoleFracIdx instead + }; // \} diff --git a/dumux/boxmodels/1p2c/1p2clocalresidual.hh b/dumux/boxmodels/1p2c/1p2clocalresidual.hh index a4f7a7b5b4..3a17d7f6d2 100644 --- a/dumux/boxmodels/1p2c/1p2clocalresidual.hh +++ b/dumux/boxmodels/1p2c/1p2clocalresidual.hh @@ -77,12 +77,12 @@ protected: //phase index phaseIdx = Indices::phaseIdx, - comp1Idx = Indices::comp1Idx, + transportCompIdx = Indices::transportCompIdx, }; // indices of the primary variables enum { pressuerIdx = Indices::pressureIdx, - transportCompIdx = Indices::transportCompIdx, + massOrMoleFracIdx = Indices::massOrMoleFracIdx, }; // indices of the equations enum { @@ -133,7 +133,7 @@ public: volVars.fluidState().density(phaseIdx)*volVars.porosity(); //storage term of the transport equation - massfractions storage[transportEqIdx] += - volVars.fluidState().density(phaseIdx) * volVars.fluidState().massFraction(phaseIdx, comp1Idx) * volVars.porosity(); + volVars.fluidState().density(phaseIdx) * volVars.fluidState().massFraction(phaseIdx, transportCompIdx) * volVars.porosity(); } else { @@ -142,7 +142,7 @@ public: storage[conti0EqIdx] += volVars.molarDensity()*volVars.porosity(); // storage term of the transport equation - molefractions storage[transportEqIdx] += - volVars.fluidState().molarDensity(phaseIdx)*volVars.fluidState().moleFraction(phaseIdx, comp1Idx) * + volVars.fluidState().molarDensity(phaseIdx)*volVars.fluidState().moleFraction(phaseIdx, transportCompIdx) * volVars.porosity(); } @@ -204,9 +204,9 @@ public: // advective flux of the second component - massfraction flux[transportEqIdx] += fluxVars.KmvpNormal() * - (( upwindWeight_)*up.fluidState().density(phaseIdx) * up.fluidState().massFraction(phaseIdx, comp1Idx)/up.viscosity() + (( upwindWeight_)*up.fluidState().density(phaseIdx) * up.fluidState().massFraction(phaseIdx, transportCompIdx)/up.viscosity() + - (1 - upwindWeight_)*dn.fluidState().density(phaseIdx)*dn.fluidState().massFraction(phaseIdx, comp1Idx)/dn.viscosity()); + (1 - upwindWeight_)*dn.fluidState().density(phaseIdx)*dn.fluidState().massFraction(phaseIdx, transportCompIdx)/dn.viscosity()); } else { @@ -221,9 +221,9 @@ public: // advective flux of the second component -molefraction flux[transportEqIdx] += fluxVars.KmvpNormal() * - (( upwindWeight_)*up.molarDensity() * up.fluidState().moleFraction(phaseIdx, comp1Idx)/up.viscosity() + (( upwindWeight_)*up.molarDensity() * up.fluidState().moleFraction(phaseIdx, transportCompIdx)/up.viscosity() + - (1 - upwindWeight_)*dn.molarDensity() * dn.fluidState().moleFraction(phaseIdx, comp1Idx)/dn.viscosity()); + (1 - upwindWeight_)*dn.molarDensity() * dn.fluidState().moleFraction(phaseIdx, transportCompIdx)/dn.viscosity()); } } @@ -243,22 +243,22 @@ public: if(!useMoles) { // diffusive flux of the second component - massfraction - tmp = -(fluxVars.moleFractionGrad(comp1Idx)*fluxVars.face().normal); + tmp = -(fluxVars.moleFractionGrad(transportCompIdx)*fluxVars.face().normal); tmp *= fluxVars.porousDiffCoeff() * fluxVars.molarDensity(); // convert it to a mass flux and add it - flux[transportEqIdx] += tmp * FluidSystem::molarMass(comp1Idx); + flux[transportEqIdx] += tmp * FluidSystem::molarMass(transportCompIdx); } else { // diffusive flux of the second component - molefraction - tmp = -(fluxVars.moleFractionGrad(comp1Idx)*fluxVars.face().normal); + tmp = -(fluxVars.moleFractionGrad(transportCompIdx)*fluxVars.face().normal); tmp *= fluxVars.porousDiffCoeff() * fluxVars.molarDensity(); // dispersive flux of second component - molefraction // Vector normalDisp; // fluxVars.dispersionTensor().mv(fluxVars.face().normal, normalDisp); // tmp -= fluxVars.molarDensity()* - // (normalDisp * fluxVars.moleFractionGrad(comp1Idx)); + // (normalDisp * fluxVars.moleFractionGrad(transportCompIdx)); flux[transportEqIdx] += tmp; } diff --git a/dumux/boxmodels/1p2c/1p2cproperties.hh b/dumux/boxmodels/1p2c/1p2cproperties.hh index 3ce9b8f980..7eae384522 100644 --- a/dumux/boxmodels/1p2c/1p2cproperties.hh +++ b/dumux/boxmodels/1p2c/1p2cproperties.hh @@ -55,6 +55,7 @@ NEW_TYPE_TAG(BoxOnePTwoC, INHERITS_FROM(BoxModel)); ////////////////////////////////////////////////////////////////// NEW_PROP_TAG(NumPhases); //!< Number of fluid phases in the system +NEW_PROP_TAG(PhaseIdx); //!< A phase index in to allow that a two-phase fluidsystem is used NEW_PROP_TAG(NumComponents); //!< Number of fluid components in the system NEW_PROP_TAG(OnePTwoCIndices); //!< DEPRECATED Enumerations for the 1p2c models NEW_PROP_TAG(Indices); //!< Enumerations for the model diff --git a/dumux/boxmodels/1p2c/1p2cpropertydefaults.hh b/dumux/boxmodels/1p2c/1p2cpropertydefaults.hh index 0b8157e9f7..df040c82ce 100644 --- a/dumux/boxmodels/1p2c/1p2cpropertydefaults.hh +++ b/dumux/boxmodels/1p2c/1p2cpropertydefaults.hh @@ -73,13 +73,16 @@ SET_TYPE_PROP(BoxOnePTwoC, FluxVariables, OnePTwoCFluxVariables<TypeTag>); SET_SCALAR_PROP(BoxOnePTwoC, UpwindWeight, 1.0); //! Set the indices used by the 1p2c model -SET_TYPE_PROP(BoxOnePTwoC, Indices, Dumux::OnePTwoCIndices<0>); +SET_TYPE_PROP(BoxOnePTwoC, Indices, Dumux::OnePTwoCIndices<TypeTag, 0>); //! DEPRECATED OnePTwoCIndices property SET_TYPE_PROP(BoxOnePTwoC, OnePTwoCIndices, typename GET_PROP_TYPE(TypeTag, Indices)); //! DEPRECATED SpatialParameters property SET_TYPE_PROP(BoxOnePTwoC, SpatialParameters, typename GET_PROP_TYPE(TypeTag, SpatialParams)); + +//! Set the phaseIndex per default to zero (important for two-phase fluidsystems). +SET_INT_PROP(BoxOnePTwoC, PhaseIdx, 0); } // \} } diff --git a/dumux/boxmodels/1p2c/1p2cvolumevariables.hh b/dumux/boxmodels/1p2c/1p2cvolumevariables.hh index fb67a8ff46..dd04927683 100644 --- a/dumux/boxmodels/1p2c/1p2cvolumevariables.hh +++ b/dumux/boxmodels/1p2c/1p2cvolumevariables.hh @@ -58,13 +58,13 @@ class OnePTwoCVolumeVariables : public BoxVolumeVariables<TypeTag> typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices; enum { phaseIdx = Indices::phaseIdx, - comp0Idx = Indices::comp0Idx, - comp1Idx = Indices::comp1Idx, + phaseCompIdx = Indices::phaseCompIdx, + transportCompIdx = Indices::transportCompIdx, }; //indices of primary variables enum{ pressureIdx = Indices::pressureIdx, - transportCompIdx = Indices::transportCompIdx + massOrMoleFracIdx = Indices::massOrMoleFracIdx }; typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView; @@ -113,8 +113,8 @@ public: diffCoeff_ = FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, phaseIdx, - comp0Idx, - comp1Idx); + phaseCompIdx, + transportCompIdx); Valgrind::CheckDefined(porosity_); Valgrind::CheckDefined(tortuosity_); @@ -141,19 +141,19 @@ public: fluidState.setPressure(phaseIdx, priVars[pressureIdx]); - Scalar x1 = priVars[transportCompIdx]; //mole or mass fraction of component 1 + Scalar x1 = priVars[massOrMoleFracIdx]; //mole or mass fraction of component 1 if(!useMoles) //mass-fraction formulation { // convert mass to mole fractions - Scalar M0 = FluidSystem::molarMass(comp0Idx); - Scalar M1 = FluidSystem::molarMass(comp1Idx); + Scalar M0 = FluidSystem::molarMass(phaseCompIdx); + Scalar M1 = FluidSystem::molarMass(transportCompIdx); //meanMolarMass if x1_ is a massfraction Scalar meanMolarMass = M0*M1/(M1 + x1*(M0 - M1)); x1 *= meanMolarMass/M1; } - fluidState.setMoleFraction(phaseIdx, comp0Idx, 1 - x1); - fluidState.setMoleFraction(phaseIdx, comp1Idx, x1); + fluidState.setMoleFraction(phaseIdx, phaseCompIdx, 1 - x1); + fluidState.setMoleFraction(phaseIdx, transportCompIdx, x1); typename FluidSystem::ParameterCache paramCache; paramCache.updatePhase(fluidState, phaseIdx); @@ -190,7 +190,7 @@ public: * \param compIdx The index of the component */ Scalar moleFraction(int compIdx) const - { return fluidState_.moleFraction(phaseIdx, (compIdx==0)?comp0Idx:comp1Idx); } + { return fluidState_.moleFraction(phaseIdx, (compIdx==0)?phaseCompIdx:transportCompIdx); } /*! * \brief Return mass fraction \f$\mathrm{[kg/kg]}\f$ of a component in the phase. @@ -198,7 +198,7 @@ public: * \param compIdx The index of the component */ Scalar massFraction(int compIdx) const - { return fluidState_.massFraction(phaseIdx, (compIdx==0)?comp0Idx:comp1Idx); } + { return fluidState_.massFraction(phaseIdx, (compIdx==0)?phaseCompIdx:transportCompIdx); } /*! * \brief Return concentration \f$\mathrm{[mol/m^3]}\f$ of a component in the phase. @@ -206,7 +206,7 @@ public: * \param compIdx The index of the component */ Scalar molarity(int compIdx) const - { return fluidState_.molarity(phaseIdx, (compIdx==0)?comp0Idx:comp1Idx); } + { return fluidState_.molarity(phaseIdx, (compIdx==0)?phaseCompIdx:transportCompIdx); } /*! * \brief Return the effective pressure \f$\mathrm{[Pa]}\f$ of a given phase within diff --git a/test/boxmodels/1p2c/1p2coutflowproblem.hh b/test/boxmodels/1p2c/1p2coutflowproblem.hh index f7f5346406..6fd1e2af28 100644 --- a/test/boxmodels/1p2c/1p2coutflowproblem.hh +++ b/test/boxmodels/1p2c/1p2coutflowproblem.hh @@ -136,7 +136,7 @@ class OnePTwoCOutflowProblem : public PorousMediaBoxProblem<TypeTag> enum { // indices of the primary variables pressureIdx = Indices::pressureIdx, - transportCompIdx = Indices::transportCompIdx, + massOrMoleFracIdx = Indices::massOrMoleFracIdx, }; enum { // indices of the equations @@ -220,7 +220,7 @@ public: initial_(priVars, globalPos); //condition for the N2 molefraction at left boundary if(globalPos[0] < eps_) - priVars[transportCompIdx] = 2.0e-5; + priVars[massOrMoleFracIdx] = 2.0e-5; } /*! @@ -289,7 +289,7 @@ private: const GlobalPosition &globalPos) const { priVars[pressureIdx] = 2e5 - 1e5*globalPos[0];//0.0; //initial condition for the pressure - priVars[transportCompIdx] = 0.0; //initial condition for the N2 molefraction + priVars[massOrMoleFracIdx] = 0.0; //initial condition for the N2 molefraction } const Scalar eps_; -- GitLab