discretefractureproblem.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_PROBLEM_HH
29#define EWOMS_DISCRETE_FRACTURE_PROBLEM_HH
30
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
33
34#include <opm/material/common/Means.hpp>
35
37
38#include <stdexcept>
39
40namespace Opm {
41
47template<class TypeTag>
49 : public MultiPhaseBaseProblem<TypeTag>
50{
52
53 using Implementation = GetPropType<TypeTag, Properties::Problem>;
57
58 enum { dimWorld = GridView::dimensionworld };
59 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
60
61public:
67 {}
68
78 template <class Context>
79 void fractureFaceIntrinsicPermeability(DimMatrix& result,
80 const Context& context,
81 unsigned localFaceIdx,
82 unsigned timeIdx) const
83 {
84 const auto& scvf = context.stencil(timeIdx).interiorFace(localFaceIdx);
85 const unsigned interiorElemIdx = scvf.interiorIndex();
86 const unsigned exteriorElemIdx = scvf.exteriorIndex();
87 const DimMatrix& K1 = asImp_().fractureIntrinsicPermeability(context, interiorElemIdx, timeIdx);
88 const DimMatrix& K2 = asImp_().fractureIntrinsicPermeability(context, exteriorElemIdx, timeIdx);
89
90 // entry-wise harmonic mean. this is almost certainly wrong if
91 // you have off-main diagonal entries in your permeabilities!
92 for (unsigned i = 0; i < dimWorld; ++i) {
93 for (unsigned j = 0; j < dimWorld; ++j) {
94 result[i][j] = harmonicMean(K1[i][j], K2[i][j]);
95 }
96 }
97 }
106 template <class Context>
107 const DimMatrix& fractureIntrinsicPermeability(const Context&,
108 unsigned,
109 unsigned) const
110 {
111 throw std::logic_error("Not implemented: Problem::fractureIntrinsicPermeability()");
112 }
113
122 template <class Context>
123 Scalar fracturePorosity(const Context&,
124 unsigned,
125 unsigned) const
126 {
127 throw std::logic_error("Not implemented: Problem::fracturePorosity()");
128 }
129
130private:
132 Implementation& asImp_()
133 { return *static_cast<Implementation *>(this); }
134
136 const Implementation& asImp_() const
137 { return *static_cast<const Implementation *>(this); }
138};
139
140} // namespace Opm
141
142#endif
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: discretefractureproblem.hh:50
const DimMatrix & fractureIntrinsicPermeability(const Context &, unsigned, unsigned) const
Returns the intrinsic permeability tensor at a given position due to a fracture.
Definition: discretefractureproblem.hh:107
void fractureFaceIntrinsicPermeability(DimMatrix &result, const Context &context, unsigned localFaceIdx, unsigned timeIdx) const
Returns the intrinsic permeability of a face due to a fracture.
Definition: discretefractureproblem.hh:79
DiscreteFractureProblem(Simulator &simulator)
Definition: discretefractureproblem.hh:65
Scalar fracturePorosity(const Context &, unsigned, unsigned) const
Returns the porosity [] inside fractures for a given control volume.
Definition: discretefractureproblem.hh:123
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition: fvbaseproblem.hh:677
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:65
Definition: blackoilboundaryratevector.hh:39
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