2012-08-27 10:56:01 -05:00
|
|
|
/*
|
|
|
|
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/>.
|
|
|
|
*/
|
2015-08-31 03:59:03 -05:00
|
|
|
#include "config.h"
|
2012-08-27 10:56:01 -05:00
|
|
|
|
2015-09-14 07:58:05 -05:00
|
|
|
#include <opm/core/props/satfunc/SaturationPropsFromDeck.hpp>
|
|
|
|
|
|
|
|
#include <opm/material/fluidmatrixinteractions/EclMaterialLawManager.hpp>
|
|
|
|
|
2012-08-27 10:56:01 -05:00
|
|
|
#include <opm/core/utility/UniformTableLinear.hpp>
|
|
|
|
#include <opm/core/utility/NonuniformTableLinear.hpp>
|
2014-02-25 08:45:50 -06:00
|
|
|
#include <opm/core/grid/GridHelpers.hpp>
|
2015-08-31 03:59:03 -05:00
|
|
|
#include <opm/core/simulator/ExplicitArraysFluidState.hpp>
|
|
|
|
#include <opm/core/simulator/ExplicitArraysSatDerivativesFluidState.hpp>
|
2012-08-27 10:56:01 -05:00
|
|
|
|
2013-08-28 07:19:12 -05:00
|
|
|
#include <iostream>
|
2014-09-01 12:09:03 -05:00
|
|
|
#include <map>
|
2013-08-28 07:19:12 -05:00
|
|
|
|
2012-08-27 10:56:01 -05:00
|
|
|
namespace Opm
|
|
|
|
{
|
|
|
|
|
2015-09-14 07:58:05 -05:00
|
|
|
typedef SaturationPropsFromDeck::MaterialLawManager::MaterialLaw MaterialLaw;
|
2012-08-27 10:56:01 -05:00
|
|
|
|
|
|
|
// ----------- Methods of SaturationPropsFromDeck ---------
|
|
|
|
|
|
|
|
|
|
|
|
/// Default constructor.
|
2015-06-26 06:23:37 -05:00
|
|
|
SaturationPropsFromDeck::SaturationPropsFromDeck()
|
2012-08-27 10:56:01 -05:00
|
|
|
{
|
|
|
|
}
|
|
|
|
|
2015-07-28 10:28:51 -05:00
|
|
|
/// Initialize from deck.
|
|
|
|
void SaturationPropsFromDeck::init(const PhaseUsage &phaseUsage,
|
|
|
|
std::shared_ptr<MaterialLawManager> materialLawManager)
|
|
|
|
{
|
|
|
|
phaseUsage_ = phaseUsage;
|
|
|
|
materialLawManager_ = materialLawManager;
|
|
|
|
}
|
|
|
|
|
2012-08-27 10:56:01 -05:00
|
|
|
/// \return P, the number of phases.
|
2015-06-26 06:23:37 -05:00
|
|
|
int SaturationPropsFromDeck::numPhases() const
|
2012-08-27 10:56:01 -05:00
|
|
|
{
|
2015-07-28 10:28:51 -05:00
|
|
|
return phaseUsage_.num_phases;
|
2012-08-27 10:56:01 -05:00
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
/// 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 ...)
|
2015-06-26 06:23:37 -05:00
|
|
|
void SaturationPropsFromDeck::relperm(const int n,
|
2012-08-27 10:56:01 -05:00
|
|
|
const double* s,
|
|
|
|
const int* cells,
|
|
|
|
double* kr,
|
|
|
|
double* dkrds) const
|
|
|
|
{
|
2013-08-28 07:00:35 -05:00
|
|
|
assert(cells != 0);
|
2012-08-27 10:56:01 -05:00
|
|
|
|
2015-10-21 08:27:52 -05:00
|
|
|
const int np = numPhases();
|
2012-08-27 10:56:01 -05:00
|
|
|
if (dkrds) {
|
2015-07-28 10:28:51 -05:00
|
|
|
ExplicitArraysSatDerivativesFluidState fluidState(phaseUsage_);
|
|
|
|
fluidState.setSaturationArray(s);
|
|
|
|
|
|
|
|
typedef ExplicitArraysSatDerivativesFluidState::Evaluation Evaluation;
|
|
|
|
Evaluation relativePerms[BlackoilPhases::MaxNumPhases];
|
2012-08-27 10:56:01 -05:00
|
|
|
for (int i = 0; i < n; ++i) {
|
2015-06-26 07:07:05 -05:00
|
|
|
fluidState.setIndex(i);
|
2015-07-28 10:28:51 -05:00
|
|
|
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];
|
2012-11-29 09:36:13 -06:00
|
|
|
}
|
2012-08-27 10:56:01 -05:00
|
|
|
}
|
|
|
|
} else {
|
2015-07-28 10:28:51 -05:00
|
|
|
ExplicitArraysFluidState fluidState(phaseUsage_);
|
|
|
|
fluidState.setSaturationArray(s);
|
|
|
|
|
|
|
|
double relativePerms[BlackoilPhases::MaxNumPhases];
|
2012-08-27 10:56:01 -05:00
|
|
|
for (int i = 0; i < n; ++i) {
|
2015-08-12 06:38:10 -05:00
|
|
|
fluidState.setIndex(i);
|
2015-07-28 10:28:51 -05:00
|
|
|
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];
|
2012-11-29 09:36:13 -06:00
|
|
|
}
|
2012-08-27 10:56:01 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
|
/// 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 ...)
|
2015-06-26 06:23:37 -05:00
|
|
|
void SaturationPropsFromDeck::capPress(const int n,
|
2012-08-27 10:56:01 -05:00
|
|
|
const double* s,
|
|
|
|
const int* cells,
|
|
|
|
double* pc,
|
|
|
|
double* dpcds) const
|
|
|
|
{
|
2013-08-28 07:00:35 -05:00
|
|
|
assert(cells != 0);
|
2012-08-27 10:56:01 -05:00
|
|
|
|
2015-10-21 08:27:52 -05:00
|
|
|
const int np = numPhases();
|
2012-08-27 10:56:01 -05:00
|
|
|
if (dpcds) {
|
2015-07-28 10:28:51 -05:00
|
|
|
ExplicitArraysSatDerivativesFluidState fluidState(phaseUsage_);
|
|
|
|
typedef ExplicitArraysSatDerivativesFluidState::Evaluation Evaluation;
|
|
|
|
fluidState.setSaturationArray(s);
|
|
|
|
|
|
|
|
Evaluation capillaryPressures[BlackoilPhases::MaxNumPhases];
|
2012-08-27 10:56:01 -05:00
|
|
|
for (int i = 0; i < n; ++i) {
|
2015-06-26 07:07:05 -05:00
|
|
|
fluidState.setIndex(i);
|
2015-07-28 10:28:51 -05:00
|
|
|
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];
|
2014-01-28 09:36:55 -06:00
|
|
|
}
|
2012-08-27 10:56:01 -05:00
|
|
|
}
|
|
|
|
} else {
|
2015-07-28 10:28:51 -05:00
|
|
|
ExplicitArraysFluidState fluidState(phaseUsage_);
|
|
|
|
fluidState.setSaturationArray(s);
|
|
|
|
|
|
|
|
double capillaryPressures[BlackoilPhases::MaxNumPhases];
|
2014-01-28 09:36:55 -06:00
|
|
|
for (int i = 0; i < n; ++i) {
|
2015-06-26 07:07:05 -05:00
|
|
|
fluidState.setIndex(i);
|
2015-07-28 10:28:51 -05:00
|
|
|
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];
|
2014-01-28 09:36:55 -06:00
|
|
|
}
|
2012-08-27 10:56:01 -05:00
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
|
|
|
|
|
|
|
/// 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.
|
2015-06-26 06:23:37 -05:00
|
|
|
void SaturationPropsFromDeck::satRange(const int n,
|
2012-08-27 10:56:01 -05:00
|
|
|
const int* cells,
|
2012-09-05 04:28:54 -05:00
|
|
|
double* smin,
|
|
|
|
double* smax) const
|
2015-06-26 07:07:05 -05:00
|
|
|
{
|
2015-09-01 06:18:44 -05:00
|
|
|
int wpos = phaseUsage_.phase_pos[BlackoilPhases::Aqua];
|
|
|
|
int gpos = phaseUsage_.phase_pos[BlackoilPhases::Vapour];
|
|
|
|
int opos = phaseUsage_.phase_pos[BlackoilPhases::Liquid];
|
|
|
|
|
2015-10-21 08:27:52 -05:00
|
|
|
const int np = numPhases();
|
2015-07-28 10:28:51 -05:00
|
|
|
for (int i = 0; i < n; ++i) {
|
|
|
|
const auto& scaledDrainageInfo =
|
|
|
|
materialLawManager_->oilWaterScaledEpsInfoDrainage(cells[i]);
|
|
|
|
|
2015-09-01 06:18:44 -05:00
|
|
|
if (phaseUsage_.phase_used[BlackoilPhases::Aqua]) {
|
|
|
|
smin[np*i + wpos] = scaledDrainageInfo.Swl;
|
|
|
|
smax[np*i + wpos] = scaledDrainageInfo.Swu;
|
|
|
|
}
|
|
|
|
|
|
|
|
if (phaseUsage_.phase_used[BlackoilPhases::Vapour]) {
|
|
|
|
smin[np*i + gpos] = scaledDrainageInfo.Sgl;
|
|
|
|
smax[np*i + gpos] = scaledDrainageInfo.Sgu;
|
|
|
|
}
|
|
|
|
|
|
|
|
if (phaseUsage_.phase_used[BlackoilPhases::Liquid]) {
|
|
|
|
smin[np*i + opos] = 1.0;
|
|
|
|
smax[np*i + opos] = 1.0;
|
|
|
|
if (phaseUsage_.phase_used[BlackoilPhases::Aqua]) {
|
|
|
|
smin[np*i + opos] -= smax[np*i + wpos];
|
|
|
|
smax[np*i + opos] -= smin[np*i + wpos];
|
|
|
|
}
|
|
|
|
if (phaseUsage_.phase_used[BlackoilPhases::Vapour]) {
|
|
|
|
smin[np*i + opos] -= smax[np*i + gpos];
|
|
|
|
smax[np*i + opos] -= smin[np*i + gpos];
|
|
|
|
}
|
|
|
|
smin[np*i + opos] = std::max(0.0, smin[np*i + opos]);
|
|
|
|
}
|
2012-09-05 04:28:54 -05:00
|
|
|
}
|
2012-08-27 10:56:01 -05:00
|
|
|
}
|
|
|
|
|
2014-01-28 09:36:55 -06:00
|
|
|
/// Update saturation state for the hysteresis tracking
|
|
|
|
/// \param[in] n Number of data points.
|
|
|
|
/// \param[in] s Array of nP saturation values.
|
2015-06-26 06:23:37 -05:00
|
|
|
void SaturationPropsFromDeck::updateSatHyst(const int n,
|
2014-01-28 09:36:55 -06:00
|
|
|
const int* cells,
|
|
|
|
const double* s)
|
|
|
|
{
|
|
|
|
assert(cells != 0);
|
|
|
|
|
2015-07-28 10:28:51 -05:00
|
|
|
if (materialLawManager_->enableHysteresis()) {
|
|
|
|
ExplicitArraysFluidState fluidState(phaseUsage_);
|
|
|
|
fluidState.setSaturationArray(s);
|
|
|
|
|
2014-01-28 09:36:55 -06:00
|
|
|
for (int i = 0; i < n; ++i) {
|
2015-07-28 10:28:51 -05:00
|
|
|
fluidState.setIndex(i);
|
|
|
|
materialLawManager_->updateHysteresis(fluidState, cells[i]);
|
2014-01-28 09:36:55 -06:00
|
|
|
}
|
2015-07-28 10:28:51 -05:00
|
|
|
}
|
2014-01-28 09:36:55 -06:00
|
|
|
}
|
2012-08-27 10:56:01 -05:00
|
|
|
|
2014-06-10 07:02:22 -05:00
|
|
|
|
|
|
|
/// Update capillary pressure scaling according to pressure diff. and initial water saturation.
|
2014-08-04 12:18:43 -05:00
|
|
|
/// \param[in] cell Cell index.
|
2014-06-10 07:02:22 -05:00
|
|
|
/// \param[in] pcow P_oil - P_water.
|
2014-08-04 12:18:43 -05:00
|
|
|
/// \param[in/out] swat Water saturation. / Possibly modified Water saturation.
|
2015-06-26 06:23:37 -05:00
|
|
|
void SaturationPropsFromDeck::swatInitScaling(const int cell,
|
2014-08-04 12:18:43 -05:00
|
|
|
const double pcow,
|
2014-08-04 12:23:03 -05:00
|
|
|
double& swat)
|
2014-06-10 07:02:22 -05:00
|
|
|
{
|
2015-07-28 10:28:51 -05:00
|
|
|
swat = materialLawManager_->applySwatinit(cell, pcow, swat);
|
2014-06-10 07:02:22 -05:00
|
|
|
}
|
2012-08-27 10:56:01 -05:00
|
|
|
} // namespace Opm
|