Skip to content
Snippets Groups Projects
Commit d8e07bac authored by Katharina Heck's avatar Katharina Heck Committed by Kilian Weishaupt
Browse files

[1pnc] make compatible with maxwellstefan

parent 4c41c116
No related branches found
No related tags found
2 merge requests!1337WIP Fix/dirichlet caching v2,!1020Stokes-darcy coupling with maxwell stefan difusion
...@@ -423,6 +423,11 @@ public: ...@@ -423,6 +423,11 @@ public:
{ {
if (phaseIdx == liquidPhaseIdx) if (phaseIdx == liquidPhaseIdx)
{ {
if (compIIdx > compJIdx)
{
using std::swap;
swap(compIIdx, compJIdx);
}
//! \todo TODO implement binary coefficients //! \todo TODO implement binary coefficients
// http://webserver.dmt.upm.es/~isidoro/dat1/Mass%20diffusivity%20data.htm // http://webserver.dmt.upm.es/~isidoro/dat1/Mass%20diffusivity%20data.htm
// The link above was given as a reference in brine_air fluid system. // The link above was given as a reference in brine_air fluid system.
......
...@@ -99,11 +99,21 @@ public: ...@@ -99,11 +99,21 @@ public:
// Could be avoided if diffusion coefficients also // Could be avoided if diffusion coefficients also
// became part of the fluid state. // became part of the fluid state.
typename FluidSystem::ParameterCache paramCache; typename FluidSystem::ParameterCache paramCache;
paramCache.updatePhase(fluidState_, 0); paramCache.updatePhase(fluidState_, 0);
diffCoeff_[0] = 0.0; // the main component with itself doesn't have a binary diffusion coefficient for (unsigned int compIIdx = 0; compIIdx < numFluidComps; ++compIIdx)
for (unsigned int compJIdx = 1; compJIdx < numFluidComps; ++compJIdx) {
diffCoeff_[compJIdx] = FluidSystem::binaryDiffusionCoefficient(fluidState_, paramCache, 0, 0, compJIdx); for (unsigned int compJIdx = 0; compJIdx < numFluidComps; ++compJIdx)
{
if(compIIdx != compJIdx)
diffCoeff_[compIIdx][compJIdx] = FluidSystem::binaryDiffusionCoefficient(fluidState_,
paramCache,
0,
compIIdx,
compJIdx);
}
}
} }
/*! /*!
...@@ -294,7 +304,7 @@ public: ...@@ -294,7 +304,7 @@ public:
Scalar diffusionCoefficient(int phaseIdx, int compIdx) const Scalar diffusionCoefficient(int phaseIdx, int compIdx) const
{ {
assert(compIdx < numFluidComps); assert(compIdx < numFluidComps);
return diffCoeff_[compIdx]; return diffCoeff_[phaseIdx][compIdx];
} }
/*! /*!
...@@ -330,8 +340,8 @@ protected: ...@@ -330,8 +340,8 @@ protected:
SolidState solidState_; SolidState solidState_;
private: private:
PermeabilityType permeability_; //!< Effective permeability within the control volume PermeabilityType permeability_;
Dune::FieldVector<Scalar, numFluidComps> diffCoeff_; //!< Binary diffusion coefficients std::array<std::array<Scalar, numFluidComps>, numFluidComps> diffCoeff_;
}; };
} // end namespace Dumux } // end namespace Dumux
......
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