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