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 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 2 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 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef EWOMS_MULTI_PHASE_BASE_EXTENSIVE_QUANTITIES_HH
29#define EWOMS_MULTI_PHASE_BASE_EXTENSIVE_QUANTITIES_HH
30
31#include <dune/common/fvector.hh>
32
33#include <opm/material/common/Valgrind.hpp>
34
39
40#include <array>
41
42namespace Opm {
49template <class TypeTag>
51 : public GetPropType<TypeTag, Properties::DiscExtensiveQuantities>
52 , public GetPropType<TypeTag, Properties::FluxModule>::FluxExtensiveQuantities
53{
57
58 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
59
61 using FluxExtensiveQuantities = typename FluxModule::FluxExtensiveQuantities;
62
63public:
67 static void registerParameters()
68 {
69 FluxModule::registerParameters();
70 }
71
80 void update(const ElementContext& elemCtx, unsigned scvfIdx, unsigned timeIdx)
81 {
82 ParentType::update(elemCtx, scvfIdx, timeIdx);
83
84 // compute the pressure potential gradients
85 FluxExtensiveQuantities::calculateGradients_(elemCtx, scvfIdx, timeIdx);
86
87 // Check whether the pressure potential is in the same direction as the face
88 // normal or in the opposite one
89 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
90 if (!elemCtx.model().phaseIsConsidered(phaseIdx)) {
91 Valgrind::SetUndefined(upstreamScvIdx_[phaseIdx]);
92 Valgrind::SetUndefined(downstreamScvIdx_[phaseIdx]);
93 continue;
94 }
95
96 upstreamScvIdx_[phaseIdx] = FluxExtensiveQuantities::upstreamIndex_(phaseIdx);
97 downstreamScvIdx_[phaseIdx] = FluxExtensiveQuantities::downstreamIndex_(phaseIdx);
98 }
99
100 FluxExtensiveQuantities::calculateFluxes_(elemCtx, scvfIdx, timeIdx);
101 }
102
103
113 template <class Context, class FluidState>
114 void updateBoundary(const Context& context,
115 unsigned bfIdx,
116 unsigned timeIdx,
117 const FluidState& fluidState)
118 {
119 ParentType::updateBoundary(context, bfIdx, timeIdx, fluidState);
120
121 FluxExtensiveQuantities::calculateBoundaryGradients_(context.elementContext(),
122 bfIdx,
123 timeIdx,
124 fluidState);
125 FluxExtensiveQuantities::calculateBoundaryFluxes_(context.elementContext(),
126 bfIdx,
127 timeIdx);
128 }
129
137 short upstreamIndex(unsigned phaseIdx) const
138 { return upstreamScvIdx_[phaseIdx]; }
139
147 short downstreamIndex(unsigned phaseIdx) const
148 { return downstreamScvIdx_[phaseIdx]; }
149
156 Scalar upstreamWeight(unsigned) const
157 { return 1.0; }
158
165 Scalar downstreamWeight(unsigned phaseIdx) const
166 { return 1.0 - upstreamWeight(phaseIdx); }
167
168private:
169 std::array<short, numPhases> upstreamScvIdx_{};
170 std::array<short, numPhases> downstreamScvIdx_{};
171};
172
173} // namespace Opm
174
175#endif
This class calculates the pressure potential gradients and the filter velocities for multi-phase flow...
Definition: multiphasebaseextensivequantities.hh:53
void updateBoundary(const Context &context, unsigned bfIdx, unsigned timeIdx, const FluidState &fluidState)
Update the extensive quantities for a given boundary face.
Definition: multiphasebaseextensivequantities.hh:114
Scalar downstreamWeight(unsigned phaseIdx) const
Return the weight of the downstream control volume for a given phase as a function of the normal flux...
Definition: multiphasebaseextensivequantities.hh:165
Scalar upstreamWeight(unsigned) const
Return the weight of the upstream control volume for a given phase as a function of the normal flux.
Definition: multiphasebaseextensivequantities.hh:156
short upstreamIndex(unsigned 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:137
void update(const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: multiphasebaseextensivequantities.hh:80
short downstreamIndex(unsigned phaseIdx) const
Return the local index of the downstream control volume for a given phase as a function of the normal...
Definition: multiphasebaseextensivequantities.hh:147
static void registerParameters()
Register all run-time parameters for the extensive quantities.
Definition: multiphasebaseextensivequantities.hh:67
Defines the common properties required by the porous medium multi-phase models.
Definition: blackoilbioeffectsmodules.hh:43
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:233
This file provides the infrastructure to retrieve run-time parameters.
This method contains all callback classes for quantities that are required by some extensive quantiti...