discretefractureextensivequantities.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_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
27 #define EWOMS_DISCRETE_FRACTURE_EXTENSIVE_QUANTITIES_HH
28 
30 
31 #include <dune/common/fvector.hh>
32 #include <dune/common/fmatrix.hh>
33 
34 namespace Ewoms {
35 
42 template <class TypeTag>
44 {
46 
47  typedef typename GET_PROP_TYPE(TypeTag, ElementContext) ElementContext;
48  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
49  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
50  typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
51 
52  enum { dimWorld = GridView::dimensionworld };
53  enum { numPhases = FluidSystem::numPhases };
54 
55  typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
56  typedef Dune::FieldVector<Scalar, dimWorld> DimVector;
57 
58 public:
62  void update(const ElementContext &elemCtx, int scvfIdx, int timeIdx)
63  {
64  ParentType::update(elemCtx, scvfIdx, timeIdx);
65 
66  const auto &extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
67  const auto &stencil = elemCtx.stencil(timeIdx);
68  const auto &scvf = stencil.interiorFace(scvfIdx);
69  int insideScvIdx = scvf.interiorIndex();
70  int outsideScvIdx = scvf.exteriorIndex();
71 
72  int globalI = elemCtx.globalSpaceIndex(insideScvIdx, timeIdx);
73  int globalJ = elemCtx.globalSpaceIndex(outsideScvIdx, timeIdx);
74  const auto &fractureMapper = elemCtx.problem().fractureMapper();
75  if (!fractureMapper.isFractureEdge(globalI, globalJ))
76  // do nothing if no fracture goes though the current edge
77  return;
78 
79  // average the intrinsic permeability of the fracture
80  elemCtx.problem().fractureFaceIntrinsicPermeability(fractureIntrinsicPermeability_,
81  elemCtx, scvfIdx, timeIdx);
82 
83  auto distDirection = elemCtx.pos(outsideScvIdx, timeIdx);
84  distDirection -= elemCtx.pos(insideScvIdx, timeIdx);
85  distDirection /= distDirection.two_norm();
86 
87  const auto &problem = elemCtx.problem();
88  fractureWidth_ = problem.fractureWidth(elemCtx, insideScvIdx,
89  outsideScvIdx, timeIdx);
90  assert(fractureWidth_ < scvf.area());
91 
92  for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
93  const auto &pGrad = extQuants.potentialGrad(phaseIdx);
94 
95  int upstreamIdx = extQuants.upstreamIndex(phaseIdx);
96  const auto &up = elemCtx.intensiveQuantities(upstreamIdx, timeIdx);
97 
98  // multiply with the fracture mobility of the upstream vertex
99  fractureIntrinsicPermeability_.mv(pGrad,
100  fractureFilterVelocity_[phaseIdx]);
101  fractureFilterVelocity_[phaseIdx] *= -up.fractureMobility(phaseIdx);
102 
103  // divide the volume flux by two. This is required because
104  // a fracture is always shared by two sub-control-volume
105  // faces.
106  fractureVolumeFlux_[phaseIdx] = 0;
107  for (int dimIdx = 0; dimIdx < dimWorld; ++dimIdx)
108  fractureVolumeFlux_[phaseIdx] +=
109  (fractureFilterVelocity_[phaseIdx][dimIdx] * distDirection[dimIdx])
110  * (fractureWidth_ / 2.0) / scvf.area();
111  }
112  }
113 
114 public:
115  const DimMatrix &fractureIntrinsicPermeability() const
116  { return fractureIntrinsicPermeability_; }
117 
118  Scalar fractureVolumeFlux(int phaseIdx) const
119  { return fractureVolumeFlux_[phaseIdx]; }
120 
121  Scalar fractureWidth() const
122  { return fractureWidth_; }
123 
124  const DimVector &fractureFilterVelocity(int phaseIdx) const
125  { return fractureFilterVelocity_[phaseIdx]; }
126 
127 private:
128  DimMatrix fractureIntrinsicPermeability_;
129  DimVector fractureFilterVelocity_[numPhases];
130  Scalar fractureVolumeFlux_[numPhases];
131  Scalar fractureWidth_;
132 };
133 
134 } // namespace Ewoms
135 
136 #endif
This class provides the data all quantities that are required to calculate the fluxes of the fluid ph...
void update(const ElementContext &elemCtx, int scvfIdx, int timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: immiscibleextensivequantities.hh:71
const DimMatrix & fractureIntrinsicPermeability() const
Definition: discretefractureextensivequantities.hh:115
const DimVector & fractureFilterVelocity(int phaseIdx) const
Definition: discretefractureextensivequantities.hh:124
This class expresses all intensive quantities of the discrete fracture model.
Definition: discretefractureextensivequantities.hh:43
This class provides the data all quantities that are required to calculate the fluxes of the fluid ph...
Definition: immiscibleextensivequantities.hh:47
Definition: baseauxiliarymodule.hh:35
Scalar fractureVolumeFlux(int phaseIdx) const
Definition: discretefractureextensivequantities.hh:118
void update(const ElementContext &elemCtx, int scvfIdx, int timeIdx)
Update the extensive quantities for a given sub-control-volume-face.
Definition: discretefractureextensivequantities.hh:62
Scalar fractureWidth() const
Definition: discretefractureextensivequantities.hh:121