PiecewiseLinearTwoPhaseMaterialParams.hpp
Go to the documentation of this file.
1 // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2 // vi: set et ts=4 sw=4 sts=4:
3 /*
4  Copyright (C) 2009-2013 by Andreas Lauser
5 
6  This file is part of the Open Porous Media project (OPM).
7 
8  OPM is free software: you can redistribute it and/or modify
9  it under the terms of the GNU General Public License as published by
10  the Free Software Foundation, either version 2 of the License, or
11  (at your option) any later version.
12 
13  OPM is distributed in the hope that it will be useful,
14  but WITHOUT ANY WARRANTY; without even the implied warranty of
15  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
16  GNU General Public License for more details.
17 
18  You should have received a copy of the GNU General Public License
19  along with OPM. If not, see <http://www.gnu.org/licenses/>.
20 */
25 #ifndef OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_PARAMS_HPP
26 #define OPM_PIECEWISE_LINEAR_TWO_PHASE_MATERIAL_PARAMS_HPP
27 
28 #include <vector>
29 #include <cassert>
30 #include <cstddef>
31 
32 namespace Opm {
39 template<class TraitsT>
41 {
42  typedef typename TraitsT::Scalar Scalar;
43 
44 public:
45  typedef std::vector<Scalar> ValueVector;
46 
47  typedef TraitsT Traits;
48 
50  {
51 #ifndef NDEBUG
52  finalized_ = false;
53 #endif
54  }
55 
60  void finalize()
61  {
62 #ifndef NDEBUG
63  finalized_ = true;
64 #endif
65 
66  // revert the order of the sampling points if they were given
67  // in reverse direction.
68  if (SwPcwnSamples_.front() > SwPcwnSamples_.back())
69  swapOrder_(SwPcwnSamples_, pcwnSamples_);
70 
71  if (SwKrwSamples_.front() > SwKrwSamples_.back())
72  swapOrder_(SwKrwSamples_, krwSamples_);
73 
74 
75  if (SwKrnSamples_.front() > SwKrnSamples_.back())
76  swapOrder_(SwKrnSamples_, krnSamples_);
77 
78  }
79 
83  const ValueVector& SwKrwSamples() const
84  { assertFinalized_(); return SwKrwSamples_; }
85 
89  const ValueVector& SwKrnSamples() const
90  { assertFinalized_(); return SwKrnSamples_; }
91 
95  const ValueVector& SwPcwnSamples() const
96  { assertFinalized_(); return SwPcwnSamples_; }
97 
103  const ValueVector& pcnwSamples() const
104  { assertFinalized_(); return pcwnSamples_; }
105 
111  template <class Container>
112  void setPcnwSamples(const Container& SwValues, const Container& values)
113  {
114  assert(SwValues.size() == values.size());
115 
116  size_t n = SwValues.size();
117  SwPcwnSamples_.resize(n);
118  pcwnSamples_.resize(n);
119 
120  std::copy(SwValues.begin(), SwValues.end(), SwPcwnSamples_.begin());
121  std::copy(values.begin(), values.end(), pcwnSamples_.begin());
122  }
123 
130  const ValueVector& krwSamples() const
131  { assertFinalized_(); return krwSamples_; }
132 
139  template <class Container>
140  void setKrwSamples(const Container& SwValues, const Container& values)
141  {
142  assert(SwValues.size() == values.size());
143 
144  size_t n = SwValues.size();
145  SwKrwSamples_.resize(n);
146  krwSamples_.resize(n);
147 
148  std::copy(SwValues.begin(), SwValues.end(), SwKrwSamples_.begin());
149  std::copy(values.begin(), values.end(), krwSamples_.begin());
150  }
151 
158  const ValueVector& krnSamples() const
159  { assertFinalized_(); return krnSamples_; }
160 
167  template <class Container>
168  void setKrnSamples(const Container& SwValues, const Container& values)
169  {
170  assert(SwValues.size() == values.size());
171 
172  size_t n = SwValues.size();
173  SwKrnSamples_.resize(n);
174  krnSamples_.resize(n);
175 
176  std::copy(SwValues.begin(), SwValues.end(), SwKrnSamples_.begin());
177  std::copy(values.begin(), values.end(), krnSamples_.begin());
178  }
179 
180 private:
181 #ifndef NDEBUG
182  void assertFinalized_() const
183  { assert(finalized_); }
184 
185  bool finalized_;
186 #else
187  void assertFinalized_() const
188  { }
189 #endif
190 
191  void swapOrder_(ValueVector& swValues, ValueVector& values) const
192  {
193  if (swValues.front() > values.back()) {
194  for (unsigned origSampleIdx = 0;
195  origSampleIdx < swValues.size() / 2;
196  ++ origSampleIdx)
197  {
198  size_t newSampleIdx = swValues.size() - origSampleIdx - 1;
199 
200  std::swap(swValues[origSampleIdx], swValues[newSampleIdx]);
201  std::swap(values[origSampleIdx], values[newSampleIdx]);
202  }
203  }
204  }
205 
206  ValueVector SwPcwnSamples_;
207  ValueVector SwKrwSamples_;
208  ValueVector SwKrnSamples_;
209  ValueVector pcwnSamples_;
210  ValueVector krwSamples_;
211  ValueVector krnSamples_;
212 };
213 } // namespace Opm
214 
215 #endif
Definition: Air_Mesitylene.hpp:31
const ValueVector & krwSamples() const
Return the sampling points for the relative permeability curve of the wetting phase.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:130
void setKrwSamples(const Container &SwValues, const Container &values)
Set the sampling points for the relative permeability curve of the wetting phase. ...
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:140
void setPcnwSamples(const Container &SwValues, const Container &values)
Set the sampling points for the capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:112
const ValueVector & SwKrnSamples() const
Return the wetting-phase saturation values of all sampling points.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:89
Specification of the material parameters for a two-phase material law which uses a table and piecewis...
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:40
const ValueVector & SwPcwnSamples() const
Return the wetting-phase saturation values of all sampling points.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:95
std::vector< Scalar > ValueVector
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:45
void finalize()
Calculate all dependent quantities once the independent quantities of the parameter object have been ...
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:60
const ValueVector & SwKrwSamples() const
Return the wetting-phase saturation values of all sampling points.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:83
const ValueVector & pcnwSamples() const
Return the sampling points for the capillary pressure curve.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:103
const ValueVector & krnSamples() const
Return the sampling points for the relative permeability curve of the non-wetting phase...
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:158
void setKrnSamples(const Container &SwValues, const Container &values)
Set the sampling points for the relative permeability curve of the non-wetting phase.
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:168
TraitsT Traits
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:47
PiecewiseLinearTwoPhaseMaterialParams()
Definition: PiecewiseLinearTwoPhaseMaterialParams.hpp:49