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  Copyright (C) 2014 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_PROBLEM_HH
27 #define EWOMS_DISCRETE_FRACTURE_PROBLEM_HH
28 
30 
32 
33 #include <opm/material/common/Means.hpp>
34 
35 #include <dune/common/fvector.hh>
36 #include <dune/common/fmatrix.hh>
37 
38 namespace Ewoms {
39 namespace Properties {
40 NEW_PROP_TAG(HeatConductionLawParams);
41 NEW_PROP_TAG(EnableGravity);
42 NEW_PROP_TAG(FluxModule);
43 }
44 
50 template<class TypeTag>
52  : public MultiPhaseBaseProblem<TypeTag>
53 {
55 
56  typedef typename GET_PROP_TYPE(TypeTag, Problem) Implementation;
57  typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
58  typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
59  typedef typename GET_PROP_TYPE(TypeTag, Simulator) Simulator;
60 
61  enum { dimWorld = GridView::dimensionworld };
62  typedef Dune::FieldMatrix<Scalar, dimWorld, dimWorld> DimMatrix;
63 
64 public:
69  : ParentType(simulator)
70  {}
71 
81  template <class Context>
82  void fractureFaceIntrinsicPermeability(DimMatrix &result,
83  const Context &context,
84  int localFaceIdx, int timeIdx) const
85  {
86  const auto &scvf = context.stencil(timeIdx).interiorFace(localFaceIdx);
87  int interiorElemIdx = scvf.interiorIndex();
88  int exteriorElemIdx = scvf.exteriorIndex();
89  const DimMatrix &K1 = asImp_().fractureIntrinsicPermeability(context, interiorElemIdx, timeIdx);
90  const DimMatrix &K2 = asImp_().fractureIntrinsicPermeability(context, exteriorElemIdx, timeIdx);
91 
92  // entry-wise harmonic mean. this is almost certainly wrong if
93  // you have off-main diagonal entries in your permeabilities!
94  for (int i = 0; i < dimWorld; ++i)
95  for (int j = 0; j < dimWorld; ++j)
96  result[i][j] = Opm::harmonicMean(K1[i][j], K2[i][j]);
97  }
106  template <class Context>
107  const DimMatrix &fractureIntrinsicPermeability(const Context &context,
108  int spaceIdx, int timeIdx) const
109  {
110  OPM_THROW(std::logic_error,
111  "Not implemented: Problem::fractureIntrinsicPermeability()");
112  }
113 
122  template <class Context>
123  Scalar fracturePorosity(const Context &context,
124  int spaceIdx, int timeIdx) const
125  {
126  OPM_THROW(std::logic_error,
127  "Not implemented: Problem::fracturePorosity()");
128  }
129 
130 private:
132  Implementation &asImp_()
133  { return *static_cast<Implementation *>(this); }
135  const Implementation &asImp_() const
136  { return *static_cast<const Implementation *>(this); }
137 };
138 
139 } // namespace Ewoms
140 
141 #endif
DiscreteFractureProblem(Simulator &simulator)
Definition: discretefractureproblem.hh:68
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition: fvbaseproblem.hh:504
NEW_PROP_TAG(Grid)
The type of the DUNE grid.
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: discretefractureproblem.hh:51
Manages the initializing and running of time dependent problems.
Definition: simulator.hh:73
Definition: baseauxiliarymodule.hh:35
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
void fractureFaceIntrinsicPermeability(DimMatrix &result, const Context &context, int localFaceIdx, int timeIdx) const
Returns the intrinsic permeability of a face due to a fracture.
Definition: discretefractureproblem.hh:82
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition: multiphasebaseproblem.hh:55
Defines the properties required for the immiscible multi-phase model which considers discrete fractur...
Scalar fracturePorosity(const Context &context, int spaceIdx, int timeIdx) const
Returns the porosity [] inside fractures for a given control volume.
Definition: discretefractureproblem.hh:123
const DimMatrix & fractureIntrinsicPermeability(const Context &context, int spaceIdx, int timeIdx) const
Returns the intrinsic permeability tensor at a given position due to a fracture. ...
Definition: discretefractureproblem.hh:107