28#ifndef EWOMS_DISCRETE_FRACTURE_PROBLEM_HH
29#define EWOMS_DISCRETE_FRACTURE_PROBLEM_HH
35#include <opm/material/common/Means.hpp>
37#include <dune/common/fvector.hh>
38#include <dune/common/fmatrix.hh>
47template<
class TypeTag>
58 enum { dimWorld = GridView::dimensionworld };
59 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
78 template <
class Context>
92 for (
unsigned i = 0; i < dimWorld; ++i)
93 for (
unsigned j = 0; j < dimWorld; ++j)
104 template <
class Context>
109 throw std::logic_error(
"Not implemented: Problem::fractureIntrinsicPermeability()");
120 template <
class Context>
125 throw std::logic_error(
"Not implemented: Problem::fracturePorosity()");
130 Implementation& asImp_()
131 {
return *
static_cast<Implementation *
>(
this); }
133 const Implementation& asImp_()
const
134 {
return *
static_cast<const Implementation *
>(
this); }
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:105
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:121
Simulator & simulator()
Returns Simulator object used by the simulation.
Definition fvbaseproblem.hh:682
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
Definition multiphasebaseproblem.hh:60
Defines the properties required for the immiscible multi-phase model which considers discrete fractur...
The base class for the problems of ECFV discretizations which deal with a multi-phase flow through a ...
This file contains a set of helper functions used by VFPProd / VFPInj.
Definition blackoilboundaryratevector.hh:37
constexpr auto getPropValue()
get the value data member of a property
Definition propertysystem.hh:242
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