multiphasebaseextensivequantities.hh
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 */
26 #ifndef EWOMS_MULTI_PHASE_BASE_EXTENSIVE_QUANTITIES_HH
27 #define EWOMS_MULTI_PHASE_BASE_EXTENSIVE_QUANTITIES_HH
28 
30 
34 
35 #include <dune/common/fvector.hh>
36 
37 namespace Ewoms {
44 template <class TypeTag>
46  : public GET_PROP_TYPE(TypeTag, DiscExtensiveQuantities)
47  , public GET_PROP_TYPE(TypeTag, FluxModule)::FluxExtensiveQuantities
48 {
49  typedef typename GET_PROP_TYPE(TypeTag, DiscExtensiveQuantities) ParentType;
50 
51  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
52  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
53  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
54 
55  enum { numPhases = GET_PROP_VALUE(TypeTag, NumPhases) };
56 
57 
58  typedef typename GET_PROP_TYPE(TypeTag, FluxModule) FluxModule;
59  typedef typename FluxModule::FluxExtensiveQuantities FluxExtensiveQuantities;
60 
61 public:
65  static void registerParameters()
66  {
67  FluxModule::registerParameters();
68  }
69 
78  void update(const ElementContext &elemCtx, int scvfIdx, int timeIdx)
79  {
80  ParentType::update(elemCtx, scvfIdx, timeIdx);
81 
82  // compute the pressure potential gradients
83  FluxExtensiveQuantities::calculateGradients_(elemCtx, scvfIdx, timeIdx);
84 
85  // Check whether the pressure potential is in the same direction as the face
86  // normal or in the opposite one
87  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
88  if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
89  Valgrind::SetUndefined(upstreamScvIdx_[phaseIdx]);
90  Valgrind::SetUndefined(downstreamScvIdx_[phaseIdx]);
91  continue;
92  }
93 
94  upstreamScvIdx_[phaseIdx] = FluxExtensiveQuantities::upstreamIndex_(phaseIdx);
95  downstreamScvIdx_[phaseIdx] = FluxExtensiveQuantities::downstreamIndex_(phaseIdx);
96  }
97 
98  FluxExtensiveQuantities::calculateFluxes_(elemCtx, scvfIdx, timeIdx);
99  }
100 
101 
112  template <class Context, class FluidState>
113  void updateBoundary(const Context &context,
114  int bfIdx,
115  int timeIdx,
116  const FluidState &fluidState,
117  typename FluidSystem::ParameterCache &paramCache)
118  {
119  ParentType::updateBoundary(context, bfIdx, timeIdx, fluidState, paramCache);
120 
121  FluxExtensiveQuantities::calculateBoundaryGradients_(context.elementContext(),
122  bfIdx,
123  timeIdx,
124  fluidState,
125  paramCache);
126  FluxExtensiveQuantities::calculateBoundaryFluxes_(context.elementContext(),
127  bfIdx,
128  timeIdx);
129  }
130 
138  short upstreamIndex(int phaseIdx) const
139  { return upstreamScvIdx_[phaseIdx]; }
140 
148  short downstreamIndex(int phaseIdx) const
149  { return downstreamScvIdx_[phaseIdx]; }
150 
157  Scalar upstreamWeight(int phaseIdx) const
158  { return 1.0; }
159 
166  Scalar downstreamWeight(int phaseIdx) const
167  { return 1.0 - upstreamWeight(phaseIdx); }
168 
169 private:
170  short upstreamScvIdx_[numPhases];
171  short downstreamScvIdx_[numPhases];
172 };
173 
174 } // namespace Ewoms
175 
176 #endif
void update(const ElementContext &elemCtx, int scvfIdx, int timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: multiphasebaseextensivequantities.hh:78
Scalar downstreamWeight(int phaseIdx) const
Return the weight of the downstream control volume for a given phase as a function of the normal flux...
Definition: multiphasebaseextensivequantities.hh:166
#define GET_PROP_VALUE(TypeTag, PropTagName)
Access the value attribute of a property for a type tag.
Definition: propertysystem.hh:468
#define GET_PROP_TYPE(TypeTag, PropTagName)
Access the type attribute of a property for a type tag.
Definition: propertysystem.hh:485
This file provides the infrastructure to retrieve run-time parameters.
Definition: baseauxiliarymodule.hh:35
Provide the properties at a face which make sense indepentently of the conserved quantities.
This class calculates the pressure potential gradients and the filter velocities for multi-phase flow...
Definition: multiphasebaseextensivequantities.hh:45
This method contains all callback classes for quantities that are required by some extensive quantiti...
short downstreamIndex(int phaseIdx) const
Return the local index of the downstream control volume for a given phase as a function of the normal...
Definition: multiphasebaseextensivequantities.hh:148
static void registerParameters()
Register all run-time parameters for the extensive quantities.
Definition: multiphasebaseextensivequantities.hh:65
short upstreamIndex(int phaseIdx) const
Return the local index of the upstream control volume for a given phase as a function of the normal f...
Definition: multiphasebaseextensivequantities.hh:138
Defines the common properties required by the porous medium multi-phase models.
void updateBoundary(const Context &context, int bfIdx, int timeIdx, const FluidState &fluidState, typename FluidSystem::ParameterCache &paramCache)
Update the extensive quantities for a given boundary face.
Definition: multiphasebaseextensivequantities.hh:113
Scalar upstreamWeight(int phaseIdx) const
Return the weight of the upstream control volume for a given phase as a function of the normal flux...
Definition: multiphasebaseextensivequantities.hh:157