polymerUtilities.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2012 SINTEF ICT, Applied Mathematics.
3 
4  This file is part of the Open Porous Media project (OPM).
5 
6  OPM is free software: you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation, either version 3 of the License, or
9  (at your option) any later version.
10 
11  OPM is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with OPM. If not, see <http://www.gnu.org/licenses/>.
18 */
19 
20 #ifndef OPM_POLYMERUTILITIES_HEADER_INCLUDED
21 #define OPM_POLYMERUTILITIES_HEADER_INCLUDED
22 
23 
24 #include <opm/core/grid.h>
25 #include <opm/core/props/IncompPropertiesInterface.hpp>
26 #include <opm/core/props/BlackoilPropertiesInterface.hpp>
30 #include <opm/core/props/rock/RockCompressibility.hpp>
31 #include <opm/core/utility/SparseVector.hpp>
32 #include <vector>
33 
34 
35 namespace Opm
36 {
37 
45  void computeTotalMobility(const Opm::IncompPropertiesInterface& props,
46  const Opm::PolymerProperties& polyprops,
47  const std::vector<int>& cells,
48  const std::vector<double>& s,
49  const std::vector<double>& c,
50  const std::vector<double>& cmax,
51  std::vector<double>& totmob);
52 
63  void computeTotalMobilityOmega(const Opm::IncompPropertiesInterface& props,
64  const Opm::PolymerProperties& polyprops,
65  const std::vector<int>& cells,
66  const std::vector<double>& s,
67  const std::vector<double>& c,
68  const std::vector<double>& cmax,
69  std::vector<double>& totmob,
70  std::vector<double>& omega);
71 
80  void computeFractionalFlow(const Opm::IncompPropertiesInterface& props,
81  const Opm::PolymerProperties& polyprops,
82  const std::vector<int>& cells,
83  const std::vector<double>& s,
84  const std::vector<double>& c,
85  const std::vector<double>& cmax,
86  std::vector<double>& fractional_flows);
87 
99  void computeFractionalFlow(const Opm::BlackoilPropertiesInterface& props,
100  const Opm::PolymerProperties& polyprops,
101  const std::vector<int>& cells,
102  const std::vector<double>& p,
103  const std::vector<double>& T,
104  const std::vector<double>& z,
105  const std::vector<double>& s,
106  const std::vector<double>& c,
107  const std::vector<double>& cmax,
108  std::vector<double>& fractional_flows);
109 
128  void computeInjectedProduced(const IncompPropertiesInterface& props,
129  const Opm::PolymerProperties& polyprops,
130  const PolymerState& state,
131  const std::vector<double>& transport_src,
132  const std::vector<double>& inj_c,
133  const double dt,
134  double* injected,
135  double* produced,
136  double& polyinj,
137  double& polyprod);
138 
157  void computeInjectedProduced(const BlackoilPropertiesInterface& props,
158  const Opm::PolymerProperties& polyprops,
159  const PolymerBlackoilState& state,
160  const std::vector<double>& transport_src,
161  const std::vector<double>& inj_c,
162  const double dt,
163  double* injected,
164  double* produced,
165  double& polyinj,
166  double& polyprod);
167 
174  double computePolymerMass(const std::vector<double>& pv,
175  const std::vector<double>& s,
176  const std::vector<double>& c,
177  const double dps);
178 
185  double computePolymerAdsorbed(const IncompPropertiesInterface& props,
186  const Opm::PolymerProperties& polyprops,
187  const std::vector<double>& pv,
188  const std::vector<double>& cmax);
189 
198  double computePolymerAdsorbed(const UnstructuredGrid& grid,
199  const BlackoilPropertiesInterface& props,
200  const Opm::PolymerProperties& polyprops,
201  const PolymerBlackoilState& state,
202  const RockCompressibility* rock_comp);
203 
204 
205 
206 
207 } // namespace Opm
208 
209 
210 #endif // OPM_POLYMERUTILITIES_HEADER_INCLUDED
Definition: CompressibleTpfaPolymer.hpp:32
std::vector< double > & cmax
Definition: GravityColumnSolverPolymer_impl.hpp:78
double computePolymerMass(const std::vector< double > &pv, const std::vector< double > &s, const std::vector< double > &c, const double dps)
Computes total (free) polymer mass over all grid cells.
Definition: PolymerProperties.hpp:34
void computeTotalMobility(const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &totmob)
Computes total mobility for a set of s/c values.
double computePolymerAdsorbed(const IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< double > &pv, const std::vector< double > &cmax)
Computes total absorbed polymer mass over all grid cells.
void computeFractionalFlow(const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &fractional_flows)
void computeInjectedProduced(const IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const PolymerState &state, const std::vector< double > &transport_src, const std::vector< double > &inj_c, const double dt, double *injected, double *produced, double &polyinj, double &polyprod)
Computes injected and produced volumes of all phases, and injected and produced polymer mass...
void computeTotalMobilityOmega(const Opm::IncompPropertiesInterface &props, const Opm::PolymerProperties &polyprops, const std::vector< int > &cells, const std::vector< double > &s, const std::vector< double > &c, const std::vector< double > &cmax, std::vector< double > &totmob, std::vector< double > &omega)
Computes total mobility and omega for a set of s/c values.