SimpleFluid2pWrappingProps_impl.hpp
Go to the documentation of this file.
1 /*===========================================================================
2 //
3 // File: SimpleFluid2pWrappingProps.cpp
4 //
5 // Author: hnil <hnil@sintef.no>
6 //
7 // Created: 15 Nov 2012
8 //==========================================================================*/
9 /*
10  Copyright 2011 SINTEF ICT, Applied Mathematics.
11  Copyright 2011 Statoil ASA.
12 
13  This file is part of the Open Porous Media Project (OPM).
14 
15  OPM is free software: you can redistribute it and/or modify
16  it under the terms of the GNU General Public License as published by
17  the Free Software Foundation, either version 3 of the License, or
18  (at your option) any later version.
19 
20  OPM is distributed in the hope that it will be useful,
21  but WITHOUT ANY WARRANTY; without even the implied warranty of
22  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
23  GNU General Public License for more details.
24 
25  You should have received a copy of the GNU General Public License
26  along with OPM. If not, see <http://www.gnu.org/licenses/>.
27 */
28 #ifndef SIMPLEFLUID2PWRAPPINGPROPS_IMPL_HPP
29 #define SIMPLEFLUID2PWRAPPINGPROPS_IMPL_HPP
30 
32 #include <cassert>
33 #include <opm/common/ErrorMacros.hpp>
34 namespace Opm{
35 
37  : props_(props),
38  smin_(props.numCells()*props.numPhases()),
39  smax_(props.numCells()*props.numPhases())
40  {
41  if (props.numPhases() != 2) {
42  OPM_THROW(std::runtime_error, "SimpleFluid2pWrapper requires 2 phases.");
43  }
44  const int num_cells = props.numCells();
45  std::vector<int> cells(num_cells);
46  for (int c = 0; c < num_cells; ++c) {
47  cells[c] = c;
48  }
49  props.satRange(num_cells, &cells[0], &smin_[0], &smax_[0]);
50  }
51 
52  inline double SimpleFluid2pWrappingProps::density(int phase) const
53  {
54  return props_.density()[phase];
55  }
56 
57  template <class Sat,
58  class Mob,
59  class DMob>
60  void SimpleFluid2pWrappingProps::mobility(int c, const Sat& s, Mob& mob, DMob& dmob) const
61  {
62  props_.relperm(1, &s[0], &c, &mob[0], &dmob[0]);
63  const double* mu = props_.viscosity();
64  mob[0] /= mu[0];
65  mob[1] /= mu[1];
66  // Recall that we use Fortran ordering for kr derivatives,
67  // therefore dmob[i*2 + j] is row j and column i of the
68  // matrix.
69  // Each row corresponds to a kr function, so which mu to
70  // divide by also depends on the row, j.
71  dmob[0*2 + 0] /= mu[0];
72  dmob[0*2 + 1] /= mu[1];
73  dmob[1*2 + 0] /= mu[0];
74  dmob[1*2 + 1] /= mu[1];
75  }
76 
77  template <class Sat,
78  class Pcap,
79  class DPcap>
80  void SimpleFluid2pWrappingProps::pc(int c, const Sat& s, Pcap& pcap, DPcap& dpcap) const
81  {
82  double pcow[2];
83  double dpcow[4];
84  props_.capPress(1, &s[0], &c, pcow, dpcow);
85  pcap = pcow[0];
86  assert(pcow[1] == 0.0);
87  dpcap = dpcow[0];
88  assert(dpcow[1] == 0.0);
89  assert(dpcow[2] == 0.0);
90  assert(dpcow[3] == 0.0);
91  }
92 
93  inline double SimpleFluid2pWrappingProps::s_min(int c) const
94  {
95  return smin_[2*c + 0];
96  }
97 
98  inline double SimpleFluid2pWrappingProps::s_max(int c) const
99  {
100  return smax_[2*c + 0];
101  }
102 
103 }
104 #endif // SIMPLEFLUID2PWRAPPINGPROPS_IMPL_HPP
Definition: IncompPropertiesInterface.hpp:35
virtual int numPhases() const =0
virtual void relperm(const int n, const double *s, const int *cells, double *kr, double *dkrds) const =0
void mobility(int c, const Sat &s, Mob &mob, DMob &dmob) const
Definition: SimpleFluid2pWrappingProps_impl.hpp:60
virtual void satRange(const int n, const int *cells, double *smin, double *smax) const =0
Definition: AnisotropicEikonal.hpp:43
virtual const double * viscosity() const =0
void pc(int c, const Sat &s, Pcap &pcap, DPcap &dpcap) const
Definition: SimpleFluid2pWrappingProps_impl.hpp:80
virtual const double * density() const =0
int numCells(const UnstructuredGrid &grid)
Get the number of cells of a grid.
double density(int phase) const
Definition: SimpleFluid2pWrappingProps_impl.hpp:52
double s_min(int c) const
Definition: SimpleFluid2pWrappingProps_impl.hpp:93
double s_max(int c) const
Definition: SimpleFluid2pWrappingProps_impl.hpp:98
SimpleFluid2pWrappingProps(const Opm::IncompPropertiesInterface &props)
Definition: SimpleFluid2pWrappingProps_impl.hpp:36
virtual int numCells() const =0
virtual void capPress(const int n, const double *s, const int *cells, double *pc, double *dpcds) const =0