28#ifndef EWOMS_DISCRETE_FRACTURE_PROBLEM_HH
29#define EWOMS_DISCRETE_FRACTURE_PROBLEM_HH
31#include <dune/common/fmatrix.hh>
32#include <dune/common/fvector.hh>
34#include <opm/material/common/Means.hpp>
47template<
class TypeTag>
58 enum { dimWorld = GridView::dimensionworld };
59 using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
78 template <
class Context>
80 const Context& context,
81 unsigned localFaceIdx,
82 unsigned timeIdx)
const
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);
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]);
106 template <
class Context>
111 throw std::logic_error(
"Not implemented: Problem::fractureIntrinsicPermeability()");
122 template <
class Context>
127 throw std::logic_error(
"Not implemented: Problem::fracturePorosity()");
132 Implementation& asImp_()
133 {
return *
static_cast<Implementation *
>(
this); }
136 const Implementation& asImp_()
const
137 {
return *
static_cast<const Implementation *
>(
this); }
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:694
MultiPhaseBaseProblem(Simulator &simulator)
Definition multiphasebaseproblem.hh:91
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 blackoilbioeffectsmodules.hh:45
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