SolventPvt: move some more code to compile unit

This commit is contained in:
Arne Morten Kvarving
2023-01-05 09:16:39 +01:00
parent 68336940e6
commit 8e78bec581
3 changed files with 56 additions and 40 deletions

View File

@@ -55,6 +55,7 @@ list (APPEND MAIN_SOURCE_FILES
src/opm/material/fluidsystems/blackoilpvt/DryGasPvt.cpp
src/opm/material/fluidsystems/blackoilpvt/DryHumidGasPvt.cpp
src/opm/material/fluidsystems/blackoilpvt/LiveOilPvt.cpp
src/opm/material/fluidsystems/blackoilpvt/SolventPvt.cpp
)
if(ENABLE_ECL_INPUT)
list(APPEND MAIN_SOURCE_FILES
@@ -261,7 +262,6 @@ if(ENABLE_ECL_INPUT)
src/opm/material/fluidsystems/blackoilpvt/GasPvtThermal.cpp
src/opm/material/fluidsystems/blackoilpvt/OilPvtMultiplexer.cpp
src/opm/material/fluidsystems/blackoilpvt/OilPvtThermal.cpp
src/opm/material/fluidsystems/blackoilpvt/SolventPvt.cpp
src/opm/material/fluidsystems/blackoilpvt/WaterPvtMultiplexer.cpp
src/opm/material/fluidsystems/blackoilpvt/WaterPvtThermal.cpp
src/opm/material/fluidsystems/blackoilpvt/WetGasPvt.cpp

View File

@@ -59,13 +59,7 @@ public:
void initFromState(const EclipseState& eclState, const Schedule&);
#endif
void setNumRegions(size_t numRegions)
{
solventReferenceDensity_.resize(numRegions);
inverseSolventB_.resize(numRegions);
inverseSolventBMu_.resize(numRegions);
solventMu_.resize(numRegions);
}
void setNumRegions(size_t numRegions);
/*!
* \brief Initialize the reference density of the solvent gas for a given PVT region
@@ -86,42 +80,13 @@ public:
*
* \param samplePoints A container of \f$(p_g, B_s)\f$ values
*/
void setSolventFormationVolumeFactor(unsigned regionIdx, const SamplingPoints& samplePoints)
{
SamplingPoints tmp(samplePoints);
auto it = tmp.begin();
const auto& endIt = tmp.end();
for (; it != endIt; ++ it)
std::get<1>(*it) = 1.0/std::get<1>(*it);
inverseSolventB_[regionIdx].setContainerOfTuples(tmp);
assert(inverseSolventB_[regionIdx].monotonic());
}
void setSolventFormationVolumeFactor(unsigned regionIdx,
const SamplingPoints& samplePoints);
/*!
* \brief Finish initializing the oil phase PVT properties.
*/
void initEnd()
{
// calculate the final 2D functions which are used for interpolation.
size_t numRegions = solventMu_.size();
for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
// calculate the table which stores the inverse of the product of the solvent
// formation volume factor and its viscosity
const auto& solventMu = solventMu_[regionIdx];
const auto& invSolventB = inverseSolventB_[regionIdx];
assert(solventMu.numSamples() == invSolventB.numSamples());
std::vector<Scalar> pressureValues(solventMu.numSamples());
std::vector<Scalar> invSolventBMuValues(solventMu.numSamples());
for (unsigned pIdx = 0; pIdx < solventMu.numSamples(); ++pIdx) {
pressureValues[pIdx] = invSolventB.xAt(pIdx);
invSolventBMuValues[pIdx] = invSolventB.valueAt(pIdx) * (1.0/solventMu.valueAt(pIdx));
}
inverseSolventBMu_[regionIdx].setXYContainers(pressureValues, invSolventBMuValues);
}
}
void initEnd();
/*!
* \brief Return the number of PVT regions which are considered by this PVT-object.

View File

@@ -26,14 +26,17 @@
#include <opm/common/ErrorMacros.hpp>
#if HAVE_ECL_INPUT
#include <opm/input/eclipse/EclipseState/EclipseState.hpp>
#include <opm/input/eclipse/EclipseState/Tables/TableManager.hpp>
#include <opm/input/eclipse/EclipseState/Tables/PvdsTable.hpp>
#endif
#include <fmt/format.h>
namespace Opm {
#if HAVE_ECL_INPUT
template<class Scalar>
void SolventPvt<Scalar>::
initFromState(const EclipseState& eclState, const Schedule&)
@@ -72,6 +75,54 @@ initFromState(const EclipseState& eclState, const Schedule&)
initEnd();
}
#endif
template<class Scalar>
void SolventPvt<Scalar>::setNumRegions(size_t numRegions)
{
solventReferenceDensity_.resize(numRegions);
inverseSolventB_.resize(numRegions);
inverseSolventBMu_.resize(numRegions);
solventMu_.resize(numRegions);
}
template<class Scalar>
void SolventPvt<Scalar>::
setSolventFormationVolumeFactor(unsigned regionIdx,
const SamplingPoints& samplePoints)
{
SamplingPoints tmp(samplePoints);
auto it = tmp.begin();
const auto& endIt = tmp.end();
for (; it != endIt; ++ it)
std::get<1>(*it) = 1.0/std::get<1>(*it);
inverseSolventB_[regionIdx].setContainerOfTuples(tmp);
assert(inverseSolventB_[regionIdx].monotonic());
}
template<class Scalar>
void SolventPvt<Scalar>::initEnd()
{
// calculate the final 2D functions which are used for interpolation.
size_t numRegions = solventMu_.size();
for (unsigned regionIdx = 0; regionIdx < numRegions; ++ regionIdx) {
// calculate the table which stores the inverse of the product of the solvent
// formation volume factor and its viscosity
const auto& solventMu = solventMu_[regionIdx];
const auto& invSolventB = inverseSolventB_[regionIdx];
assert(solventMu.numSamples() == invSolventB.numSamples());
std::vector<Scalar> pressureValues(solventMu.numSamples());
std::vector<Scalar> invSolventBMuValues(solventMu.numSamples());
for (unsigned pIdx = 0; pIdx < solventMu.numSamples(); ++pIdx) {
pressureValues[pIdx] = invSolventB.xAt(pIdx);
invSolventBMuValues[pIdx] = invSolventB.valueAt(pIdx) * (1.0/solventMu.valueAt(pIdx));
}
inverseSolventBMu_[regionIdx].setXYContainers(pressureValues, invSolventBMuValues);
}
}
template class SolventPvt<double>;
template class SolventPvt<float>;