ExplicitArraysSatDerivativesFluidState.hpp
Go to the documentation of this file.
1 /*
2  Copyright 2015 Andreas Lauser
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_EXPLICIT_ARRAYS_SAT_DERIVATIVES_FLUID_STATE_HEADER_INCLUDED
21 #define OPM_EXPLICIT_ARRAYS_SAT_DERIVATIVES_FLUID_STATE_HEADER_INCLUDED
22 
23 #include <opm/material/localad/Evaluation.hpp>
24 #include <opm/material/localad/Math.hpp>
25 
27 
28 #include <array>
29 
30 namespace Opm
31 {
32 
33 // this class does not need to be implemented. its only purpose is to make compiler
34 // messages for local-AD framework more explicit (because the class name of the
35 // Evaluation will include a tag name so that it is clear which derivatives are handled).
37 {};
38 
48 {
49 public:
51  enum { numComponents = 3 };
52 
53  typedef Opm::LocalAd::Evaluation<double, SaturationDerivativesTag, numPhases> Evaluation;
54  typedef Evaluation Scalar;
55 
57  : phaseUsage_(phaseUsage)
58  {
59  globalSaturationArray_ = 0;
60 
61  // initialize the evaluation objects for the saturations
62  for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
63  saturation_[phaseIdx] = Evaluation::createVariable(0.0, phaseIdx);
64  }
65  }
66 
73  void setIndex(unsigned arrayIdx)
74  {
75  int np = phaseUsage_.num_phases;
76 
77  // copy the saturations values from the global value. the derivatives do not need
78  // to be modified for these...
79  for (int phaseIdx = 0; phaseIdx < numPhases; ++ phaseIdx) {
80  if (!phaseUsage_.phase_used[phaseIdx]) {
81  saturation_[phaseIdx].value = 0.0;
82  }
83  else {
84  saturation_[phaseIdx].value = globalSaturationArray_[np*arrayIdx + phaseUsage_.phase_pos[phaseIdx]];
85  }
86  }
87  }
88 
97  void setSaturationArray(const double* saturations)
98  { globalSaturationArray_ = saturations; }
99 
103  const Evaluation& saturation(int phaseIdx) const
104  { return saturation_[phaseIdx]; }
105 
106  // TODO (?) temperature, pressure, composition, etc
107 
108 private:
109  const PhaseUsage phaseUsage_;
110  const double* globalSaturationArray_;
111 
112  std::array<Evaluation, numPhases> saturation_;
113 };
114 
115 } // namespace Opm
116 
117 #endif
Definition: AnisotropicEikonal.hpp:43
void setIndex(unsigned arrayIdx)
Sets the currently used array index.
Definition: ExplicitArraysSatDerivativesFluidState.hpp:73
Evaluation Scalar
Definition: ExplicitArraysSatDerivativesFluidState.hpp:54
int phase_used[MaxNumPhases]
Definition: BlackoilPhases.hpp:39
int phase_pos[MaxNumPhases]
Definition: BlackoilPhases.hpp:40
Definition: ExplicitArraysSatDerivativesFluidState.hpp:51
int num_phases
Definition: BlackoilPhases.hpp:38
void setSaturationArray(const double *saturations)
Set the array containing the phase saturations.
Definition: ExplicitArraysSatDerivativesFluidState.hpp:97
ExplicitArraysSatDerivativesFluidState(const PhaseUsage &phaseUsage)
Definition: ExplicitArraysSatDerivativesFluidState.hpp:56
This is a fluid state which translates global arrays and translates them to a subset of the fluid sta...
Definition: ExplicitArraysSatDerivativesFluidState.hpp:47
Opm::LocalAd::Evaluation< double, SaturationDerivativesTag, numPhases > Evaluation
Definition: ExplicitArraysSatDerivativesFluidState.hpp:53
static const int MaxNumPhases
Definition: BlackoilPhases.hpp:30
Definition: ExplicitArraysSatDerivativesFluidState.hpp:36
Definition: ExplicitArraysSatDerivativesFluidState.hpp:50
Definition: BlackoilPhases.hpp:36
const Evaluation & saturation(int phaseIdx) const
Returns the saturation of a phase for the current cell index.
Definition: ExplicitArraysSatDerivativesFluidState.hpp:103