Added class PolymerProperties with a constructor which takes argument from an eclipse parser.

This commit is contained in:
Xavier Raynaud 2012-02-29 11:27:59 +01:00
parent 25308b750e
commit 123762caa9
2 changed files with 109 additions and 0 deletions

View File

@ -39,6 +39,7 @@
// #include <opm/core/linalg/LinearSolverIstl.hpp>
#include <opm/polymer/TransportModelPolymer.hpp>
#include <opm/polymer/PolymerProperties.hpp>
#include <boost/filesystem/convenience.hpp>
#include <boost/scoped_ptr.hpp>

View File

@ -0,0 +1,108 @@
/*
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 <vector>
#include <opm/core/utility/linearInterpolation.hpp>
#include <opm/core/eclipse/SpecialEclipseFields.hpp>
#ifndef OPM_POLYMERPROPERTIES_HEADER_INCLUDED
#define OPM_POLYMERPROPERTIES_HEADER_INCLUDED
namespace Opm
{
class PolymerProperties
{
public:
PolymerProperties(double c_max_limit_arg, double omega_arg, double rhor_arg, double dps_arg,
std::vector<double> c_vals_visc, std::vector<double> visc_mult_vals,
std::vector<double> c_vals_ads, std::vector<double> ads_vals)
{
c_max_limit = c_max_limit_arg;
omega = omega_arg;
rhor = rhor_arg;
dps = dps_arg;
c_vals_visc_ = c_vals_visc;
visc_mult_vals_ = visc_mult_vals;
c_vals_ads_ = c_vals_ads;
ads_vals_ = ads_vals;
}
PolymerProperties(const EclipseGridParser& gridparser)
{
// We assume NTMISC=1
const std::vector<double>& plymax = gridparser.getPLYMAX().plymax_;
c_max_limit = plymax[0];
const std::vector<double>& tlmixpar = gridparser.getTLMIXPAR().tlmixpar_;
omega = tlmixpar[0];
rhor = gridparser.getFloatingPointValue("ROCKDEN")[0];
// We assume NTSFUN=1
const std::vector<double>& plyrock = gridparser.getPLYROCK().plyrock_;
dps = plyrock[0];
// We assume NTPVT=1
const PLYVISC& plyvisc = gridparser.getPLYVISC();
c_vals_visc_ = plyvisc.concentration_;
visc_mult_vals_ = plyvisc.factor_;
// We assume NTSFUN=1
const PLYADS& plyads = gridparser.getPLYADS();
c_vals_ads_ = plyads.local_concentration_;
ads_vals_ = plyads.adsorbed_concentration_;
}
double c_max_limit;
double omega;
double rhor;
double dps;
double viscMult(double c) const
{
return Opm::linearInterpolation(c_vals_visc_, visc_mult_vals_, c);
}
double viscMultWithDer(double c, double* der) const
{
*der = Opm::linearInterpolationDerivative(c_vals_visc_, visc_mult_vals_, c);
return Opm::linearInterpolation(c_vals_visc_, visc_mult_vals_, c);
}
double adsorbtion(double c) const
{
return Opm::linearInterpolation(c_vals_ads_, ads_vals_, c);
}
double adsorbtionWithDer(double c, double* der) const
{
*der = Opm::linearInterpolationDerivative(c_vals_ads_, ads_vals_, c);
return Opm::linearInterpolation(c_vals_ads_, ads_vals_, c);
}
private:
std::vector<double> c_vals_visc_;
std::vector<double> visc_mult_vals_;
std::vector<double> c_vals_ads_;
std::vector<double> ads_vals_;
};
} // namespace Opm
#endif // OPM_POLYMERPROPERTIES_HEADER_INCLUDED