discretefracturelocalresidual.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_DISCRETE_FRACTURE_LOCAL_RESIDUAL_BASE_HH
29#define EWOMS_DISCRETE_FRACTURE_LOCAL_RESIDUAL_BASE_HH
30
32
33namespace Opm {
34
41template <class TypeTag>
43{
45
51
52 enum { conti0EqIdx = Indices::conti0EqIdx };
53 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
54 enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
55
57
58public:
67 void addPhaseStorage(EqVector& storage,
68 const ElementContext& elemCtx,
69 unsigned dofIdx,
70 unsigned timeIdx,
71 unsigned phaseIdx) const
72 {
73 EqVector phaseStorage(0.0);
74 ParentType::addPhaseStorage(phaseStorage, elemCtx, dofIdx, timeIdx, phaseIdx);
75
76 const auto& problem = elemCtx.problem();
77 const auto& fractureMapper = problem.fractureMapper();
78 unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);
79
80 if (!fractureMapper.isFractureVertex(globalIdx)) {
81 // don't do anything in addition to the immiscible model for degrees of
82 // freedom that do not feature fractures
83 storage += phaseStorage;
84 return;
85 }
86
87 const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
88 const auto& scv = elemCtx.stencil(timeIdx).subControlVolume(dofIdx);
89
90 // reduce the matrix storage by the fracture volume
91 phaseStorage *= 1 - intQuants.fractureVolume()/scv.volume();
92
93 // add the storage term inside the fractures
94 const auto& fsFracture = intQuants.fractureFluidState();
95
96 phaseStorage[conti0EqIdx + phaseIdx] +=
97 intQuants.fracturePorosity()*
98 fsFracture.saturation(phaseIdx) *
99 fsFracture.density(phaseIdx) *
100 intQuants.fractureVolume()/scv.volume();
101
102 EnergyModule::addFracturePhaseStorage(phaseStorage, intQuants, scv,
103 phaseIdx);
104
105 // add the result to the overall storage term
106 storage += phaseStorage;
107 }
108
112 void computeFlux(RateVector& flux,
113 const ElementContext& elemCtx,
114 unsigned scvfIdx,
115 unsigned timeIdx) const
116 {
117 ParentType::computeFlux(flux, elemCtx, scvfIdx, timeIdx);
118
119 const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
120
121 unsigned i = extQuants.interiorIndex();
122 unsigned j = extQuants.exteriorIndex();
123 unsigned I = elemCtx.globalSpaceIndex(i, timeIdx);
124 unsigned J = elemCtx.globalSpaceIndex(j, timeIdx);
125 const auto& fractureMapper = elemCtx.problem().fractureMapper();
126 if (!fractureMapper.isFractureEdge(I, J))
127 // do nothing if the edge from i to j is not part of a
128 // fracture
129 return;
130
131 const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
132 Scalar scvfArea = scvf.area();
133
134 // advective mass fluxes of all phases
135 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
136 if (!elemCtx.model().phaseIsConsidered(phaseIdx))
137 continue;
138
139 // reduce the matrix mass flux by the width of the scv
140 // face that is occupied by the fracture. As usual, the
141 // fracture is shared between two SCVs, so the its width
142 // needs to be divided by two.
143 flux[conti0EqIdx + phaseIdx] *=
144 1 - extQuants.fractureWidth() / (2 * scvfArea);
145
146 // intensive quantities of the upstream and the downstream DOFs
147 unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
148 const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
149 flux[conti0EqIdx + phaseIdx] +=
150 extQuants.fractureVolumeFlux(phaseIdx) * up.fractureFluidState().density(phaseIdx);
151 }
152
153 EnergyModule::handleFractureFlux(flux, elemCtx, scvfIdx, timeIdx);
154 }
155};
156
157} // namespace Opm
158
159#endif
Calculates the local residual of the discrete fracture immiscible multi-phase model.
Definition: discretefracturelocalresidual.hh:43
void addPhaseStorage(EqVector &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase.
Definition: discretefracturelocalresidual.hh:67
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: discretefracturelocalresidual.hh:112
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Calculates the local residual of the immiscible multi-phase model.
Definition: immisciblelocalresidual.hh:46
void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase.
Definition: immisciblelocalresidual.hh:76
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: immisciblelocalresidual.hh:114
Definition: blackoilboundaryratevector.hh:37
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:235