/* 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 . */ #ifndef OPM_POLYMERPROPERTIES_HEADER_INCLUDED #define OPM_POLYMERPROPERTIES_HEADER_INCLUDED #include #include #include namespace Opm { class PolymerProperties { public: PolymerProperties() { } enum AdsorptionBehaviour { Desorption = 1, NoDesorption = 2 }; PolymerProperties(double c_max, double mix_param, double rock_density, double dead_pore_vol, double res_factor, double c_max_ads, AdsorptionBehaviour ads_index, const std::vector& c_vals_visc, const std::vector& visc_mult_vals, const std::vector& c_vals_ads, const std::vector& ads_vals ) : c_max_(c_max), mix_param_(mix_param), rock_density_(rock_density), dead_pore_vol_(dead_pore_vol), res_factor_(res_factor), c_max_ads_(c_max_ads), ads_index_(ads_index), 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) { readFromDeck(gridparser); } void set(double c_max, double mix_param, double rock_density, double dead_pore_vol, double res_factor, double c_max_ads, AdsorptionBehaviour ads_index, const std::vector& c_vals_visc, const std::vector& visc_mult_vals, const std::vector& c_vals_ads, const std::vector& ads_vals ) { c_max_ = c_max; mix_param_ = mix_param; rock_density_ = rock_density; dead_pore_vol_ = dead_pore_vol; res_factor_ = res_factor; c_max_ads_ = c_max_ads; c_vals_visc_ = c_vals_visc; visc_mult_vals_ = visc_mult_vals; c_vals_ads_ = c_vals_ads; ads_vals_ = ads_vals; ads_index_ = ads_index; } void readFromDeck(const EclipseGridParser& gridparser) { // We assume NTMISC=1 const std::vector& plymax = gridparser.getPLYMAX().plymax_; c_max_ = plymax[0]; const std::vector& tlmixpar = gridparser.getTLMIXPAR().tlmixpar_; mix_param_ = tlmixpar[0]; // We assume NTSFUN=1 const std::vector& plyrock = gridparser.getPLYROCK().plyrock_; ASSERT(plyrock.size() == 5); dead_pore_vol_ = plyrock[0]; res_factor_ = plyrock[1]; rock_density_ = plyrock[2]; ads_index_ = static_cast(plyrock[3]); c_max_ads_ = plyrock[4]; // 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 cMax() const; double mixParam() const; double rockDensity() const; double deadPoreVol() const; double resFactor() const; double cMaxAds() const; int adsIndex() const; double viscMult(double c) const; double viscMultWithDer(double c, double* der) const; void simpleAdsorption(double c, double& c_ads) const; void simpleAdsorptionWithDer(double c, double& c_ads, double& dc_ads_dc) const; void adsorption(double c, double cmax, double& c_ads) const; void adsorptionWithDer(double c, double cmax, double& c_ads, double& dc_ads_dc) const; void effectiveVisc(const double c, const double* visc, double* visc_eff) const; void effectiveInvVisc(const double c, const double* visc, double* inv_visc_eff) const; void effectiveViscWithDer(const double c, const double* visc, double* visc_eff, double* dvisc_eff_dc) const; void effectiveRelperm(const double c, const double cmax, const double* relperm, double& eff_relperm_wat) const; void effectiveRelpermWithDer (const double c, const double cmax, const double* relperm, const double* drelperm_ds, double& eff_relperm_wat, double& deff_relperm_wat_ds, double& deff_relperm_wat_dc) const; void effectiveMobilities(const double c, const double cmax, const double* visc, const double* relperm, std::vector& mob) const; void effectiveMobilitiesWithDer(const double c, const double cmax, const double* visc, const double* relperm, const double* drelpermds, std::vector& mob, std::vector& dmobds, double& dmobwatdc) const; void effectiveMobilitiesBoth(const double c, const double cmax, const double* visc, const double* relperm, const double* drelperm_ds, std::vector& mob, std::vector& dmob_ds, double& dmobwat_dc, bool if_with_der) const; void computeMcWithDer(const double& c, double& mc) const; void computeMcWithDer(const double& c, double& mc, double& dmc_dc) const; void computeMcBoth(const double& c, double& mc, double& dmc_dc, bool if_with_der) const; private: double c_max_; double mix_param_; double rock_density_; double dead_pore_vol_; double res_factor_; double c_max_ads_; AdsorptionBehaviour ads_index_; std::vector c_vals_visc_; std::vector visc_mult_vals_; std::vector c_vals_ads_; std::vector ads_vals_; void simpleAdsorptionBoth(double c, double& c_ads, double& dc_ads_dc, bool if_with_der) const; void adsorptionBoth(double c, double cmax, double& c_ads, double& dc_ads_dc, bool if_with_der) const; void effectiveViscBoth(const double c, const double* visc, double* visc_eff, double* dvisc_eff_dc, bool if_with_der) const; void effectiveRelpermBoth(const double c, const double cmax, const double* relperm, const double* drelperm_ds, double& eff_relperm_wat, double& deff_relperm_wat_ds, double& deff_relperm_wat_dc, bool if_with_der) const; }; } // namespace Opm #endif // OPM_POLYMERPROPERTIES_HEADER_INCLUDED