opm-simulators/opm/core/props/satfunc/SaturationPropsFromDeck.cpp
Andreas Lauser 750e9aedf4 use a .cpp instead of an _impl.hpp file for SaturationPropsFromDeck
this avoids having to include the "Evaluation.hpp" file as the first
thing in the morning.
2015-09-02 12:29:18 +02:00

255 lines
11 KiB
C++

/*
Copyright 2012 SINTEF ICT, Applied Mathematics.
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see <http://www.gnu.org/licenses/>.
*/
#include "config.h"
#include <opm/core/utility/UniformTableLinear.hpp>
#include <opm/core/utility/NonuniformTableLinear.hpp>
#include <opm/core/grid/GridHelpers.hpp>
#include <opm/core/simulator/ExplicitArraysFluidState.hpp>
#include <opm/core/simulator/ExplicitArraysSatDerivativesFluidState.hpp>
#include <opm/parser/eclipse/Utility/EndscaleWrapper.hpp>
#include <opm/parser/eclipse/Utility/ScalecrsWrapper.hpp>
#include <iostream>
#include <map>
#include "SaturationPropsFromDeck.hpp"
namespace Opm
{
// ----------- Methods of SaturationPropsFromDeck ---------
/// Default constructor.
SaturationPropsFromDeck::SaturationPropsFromDeck()
{
}
/// Initialize from deck.
void SaturationPropsFromDeck::init(Opm::DeckConstPtr deck,
Opm::EclipseStateConstPtr eclipseState,
std::shared_ptr<MaterialLawManager> materialLawManager,
const UnstructuredGrid& grid)
{
this->init(deck, eclipseState, materialLawManager, grid.number_of_cells,
grid.global_cell, grid.cell_centroids,
grid.dimensions);
}
/// Initialize from deck.
void SaturationPropsFromDeck::init(const PhaseUsage &phaseUsage,
std::shared_ptr<MaterialLawManager> materialLawManager)
{
phaseUsage_ = phaseUsage;
materialLawManager_ = materialLawManager;
}
/// \return P, the number of phases.
int SaturationPropsFromDeck::numPhases() const
{
return phaseUsage_.num_phases;
}
/// Relative permeability.
/// \param[in] n Number of data points.
/// \param[in] s Array of nP saturation values.
/// \param[in] cells Array of n cell indices to be associated with the s values.
/// \param[out] kr Array of nP relperm values, array must be valid before calling.
/// \param[out] dkrds If non-null: array of nP^2 relperm derivative values,
/// array must be valid before calling.
/// The P^2 derivative matrix is
/// m_{ij} = \frac{dkr_i}{ds^j},
/// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
void SaturationPropsFromDeck::relperm(const int n,
const double* s,
const int* cells,
double* kr,
double* dkrds) const
{
assert(cells != 0);
const int np = BlackoilPhases::MaxNumPhases;
if (dkrds) {
ExplicitArraysSatDerivativesFluidState fluidState(phaseUsage_);
fluidState.setSaturationArray(s);
typedef ExplicitArraysSatDerivativesFluidState::Evaluation Evaluation;
Evaluation relativePerms[BlackoilPhases::MaxNumPhases];
for (int i = 0; i < n; ++i) {
fluidState.setIndex(i);
const auto& params = materialLawManager_->materialLawParams(cells[i]);
MaterialLaw::relativePermeabilities(relativePerms, params, fluidState);
// copy the values calculated using opm-material to the target arrays
for (int krPhaseIdx = 0; krPhaseIdx < np; ++krPhaseIdx) {
kr[np*i + krPhaseIdx] = relativePerms[krPhaseIdx].value;
for (int satPhaseIdx = 0; satPhaseIdx < np; ++satPhaseIdx)
dkrds[np*np*i + satPhaseIdx*np + krPhaseIdx] = relativePerms[krPhaseIdx].derivatives[satPhaseIdx];
}
}
} else {
ExplicitArraysFluidState fluidState(phaseUsage_);
fluidState.setSaturationArray(s);
double relativePerms[BlackoilPhases::MaxNumPhases];
for (int i = 0; i < n; ++i) {
fluidState.setIndex(i);
const auto& params = materialLawManager_->materialLawParams(cells[i]);
MaterialLaw::relativePermeabilities(relativePerms, params, fluidState);
// copy the values calculated using opm-material to the target arrays
for (int krPhaseIdx = 0; krPhaseIdx < np; ++krPhaseIdx) {
kr[np*i + krPhaseIdx] = relativePerms[krPhaseIdx];
}
}
}
}
/// Capillary pressure.
/// \param[in] n Number of data points.
/// \param[in] s Array of nP saturation values.
/// \param[in] cells Array of n cell indices to be associated with the s values.
/// \param[out] pc Array of nP capillary pressure values, array must be valid before calling.
/// \param[out] dpcds If non-null: array of nP^2 derivative values,
/// array must be valid before calling.
/// The P^2 derivative matrix is
/// m_{ij} = \frac{dpc_i}{ds^j},
/// and is output in Fortran order (m_00 m_10 m_20 m01 ...)
void SaturationPropsFromDeck::capPress(const int n,
const double* s,
const int* cells,
double* pc,
double* dpcds) const
{
assert(cells != 0);
const int np = BlackoilPhases::MaxNumPhases;
if (dpcds) {
ExplicitArraysSatDerivativesFluidState fluidState(phaseUsage_);
typedef ExplicitArraysSatDerivativesFluidState::Evaluation Evaluation;
fluidState.setSaturationArray(s);
Evaluation capillaryPressures[BlackoilPhases::MaxNumPhases];
for (int i = 0; i < n; ++i) {
fluidState.setIndex(i);
const auto& params = materialLawManager_->materialLawParams(cells[i]);
MaterialLaw::capillaryPressures(capillaryPressures, params, fluidState);
// copy the values calculated using opm-material to the target arrays
for (int pcPhaseIdx = 0; pcPhaseIdx < np; ++pcPhaseIdx) {
double sign = (pcPhaseIdx == BlackoilPhases::Aqua)? -1.0 : 1.0;
pc[np*i + pcPhaseIdx] = sign*capillaryPressures[pcPhaseIdx].value;
for (int satPhaseIdx = 0; satPhaseIdx < np; ++satPhaseIdx)
dpcds[np*np*i + satPhaseIdx*np + pcPhaseIdx] = sign*capillaryPressures[pcPhaseIdx].derivatives[satPhaseIdx];
}
}
} else {
ExplicitArraysFluidState fluidState(phaseUsage_);
fluidState.setSaturationArray(s);
double capillaryPressures[BlackoilPhases::MaxNumPhases];
for (int i = 0; i < n; ++i) {
fluidState.setIndex(i);
const auto& params = materialLawManager_->materialLawParams(cells[i]);
MaterialLaw::capillaryPressures(capillaryPressures, params, fluidState);
// copy the values calculated using opm-material to the target arrays
for (int pcPhaseIdx = 0; pcPhaseIdx < np; ++pcPhaseIdx) {
double sign = (pcPhaseIdx == BlackoilPhases::Aqua)? -1.0 : 1.0;
pc[np*i + pcPhaseIdx] = sign*capillaryPressures[pcPhaseIdx];
}
}
}
}
/// Obtain the range of allowable saturation values.
/// \param[in] n Number of data points.
/// \param[in] cells Array of n cell indices.
/// \param[out] smin Array of nP minimum s values, array must be valid before calling.
/// \param[out] smax Array of nP maximum s values, array must be valid before calling.
void SaturationPropsFromDeck::satRange(const int n,
const int* cells,
double* smin,
double* smax) const
{
const int np = BlackoilPhases::MaxNumPhases;
for (int i = 0; i < n; ++i) {
const auto& scaledDrainageInfo =
materialLawManager_->oilWaterScaledEpsInfoDrainage(cells[i]);
smin[np*i + BlackoilPhases::Aqua] = scaledDrainageInfo.Swl;
smax[np*i + BlackoilPhases::Aqua] = scaledDrainageInfo.Swu;
smin[np*i + BlackoilPhases::Vapour] = scaledDrainageInfo.Sgl;
smax[np*i + BlackoilPhases::Vapour] = scaledDrainageInfo.Sgu;
smin[np*i + BlackoilPhases::Liquid] = 1.0;
smax[np*i + BlackoilPhases::Liquid] = 1.0;
smin[np*i + BlackoilPhases::Liquid] -= smax[np*i + BlackoilPhases::Aqua];
smax[np*i + BlackoilPhases::Liquid] -= smin[np*i + BlackoilPhases::Aqua];
smin[np*i + BlackoilPhases::Liquid] -= smax[np*i + BlackoilPhases::Vapour];
smax[np*i + BlackoilPhases::Liquid] -= smin[np*i + BlackoilPhases::Vapour];
smin[np*i + BlackoilPhases::Liquid] = std::max(0.0, smin[np*i + BlackoilPhases::Liquid]);
}
}
/// Update saturation state for the hysteresis tracking
/// \param[in] n Number of data points.
/// \param[in] s Array of nP saturation values.
void SaturationPropsFromDeck::updateSatHyst(const int n,
const int* cells,
const double* s)
{
assert(cells != 0);
if (materialLawManager_->enableHysteresis()) {
ExplicitArraysFluidState fluidState(phaseUsage_);
fluidState.setSaturationArray(s);
for (int i = 0; i < n; ++i) {
fluidState.setIndex(i);
materialLawManager_->updateHysteresis(fluidState, cells[i]);
}
}
}
/// Update capillary pressure scaling according to pressure diff. and initial water saturation.
/// \param[in] cell Cell index.
/// \param[in] pcow P_oil - P_water.
/// \param[in/out] swat Water saturation. / Possibly modified Water saturation.
void SaturationPropsFromDeck::swatInitScaling(const int cell,
const double pcow,
double& swat)
{
swat = materialLawManager_->applySwatinit(cell, pcow, swat);
}
} // namespace Opm