File: discretefractureproblem.hh

package info (click to toggle)
opm-simulators 2024.10%2Bds-6
  • links: PTS, VCS
  • area: main
  • in suites: trixie
  • size: 19,416 kB
  • sloc: cpp: 165,337; sh: 1,285; lisp: 1,108; python: 355; makefile: 24; awk: 10
file content (139 lines) | stat: -rw-r--r-- 5,385 bytes parent folder | download | duplicates (2)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
  This file is part of the Open Porous Media project (OPM).

  OPM is free software: you can redistribute it and/or modify
  it under the terms of the GNU General Public License as published by
  the Free Software Foundation, either version 2 of the License, or
  (at your option) any later version.

  OPM is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.

  You should have received a copy of the GNU General Public License
  along with OPM.  If not, see <http://www.gnu.org/licenses/>.

  Consult the COPYING file in the top-level source directory of this
  module for the precise wording of the license and the list of
  copyright holders.
*/
/*!
 * \file
 *
 * \copydoc Opm::DiscreteFractureProblem
 */
#ifndef EWOMS_DISCRETE_FRACTURE_PROBLEM_HH
#define EWOMS_DISCRETE_FRACTURE_PROBLEM_HH

#include "discretefractureproperties.hh"

#include <opm/models/common/multiphasebaseproblem.hh>

#include <opm/material/common/Means.hpp>

#include <dune/common/fvector.hh>
#include <dune/common/fmatrix.hh>

namespace Opm {

/*!
 * \ingroup DiscreteFractureModel
 * \brief The base class for the problems of ECFV discretizations which deal
 *        with a multi-phase flow through a porous medium.
 */
template<class TypeTag>
class DiscreteFractureProblem
    : public MultiPhaseBaseProblem<TypeTag>
{
    using ParentType = Opm::MultiPhaseBaseProblem<TypeTag>;

    using Implementation = GetPropType<TypeTag, Properties::Problem>;
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using GridView = GetPropType<TypeTag, Properties::GridView>;
    using Simulator = GetPropType<TypeTag, Properties::Simulator>;

    enum { dimWorld = GridView::dimensionworld };
    using DimMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;

public:
    /*!
     * \copydoc Problem::FvBaseProblem(Simulator& )
     */
    DiscreteFractureProblem(Simulator& simulator)
        : ParentType(simulator)
    {}

    /*!
     * \brief Returns the intrinsic permeability of a face due to a fracture.
     *
     * This method is specific to the finite volume discretizations. If left unspecified,
     * it calls the intrinsicPermeability() methods for the face's interior and exterior
     * finite volume and averages them harmonically. Note that if this function is
     * defined, the intrinsicPermeability() method does not need to be defined by the
     * problem (if a finite-volume discretization is used).
     */
    template <class Context>
    void fractureFaceIntrinsicPermeability(DimMatrix& result,
                                           const Context& context,
                                           unsigned localFaceIdx,
                                           unsigned timeIdx) const
    {
        const auto& scvf = context.stencil(timeIdx).interiorFace(localFaceIdx);
        unsigned interiorElemIdx = scvf.interiorIndex();
        unsigned exteriorElemIdx = scvf.exteriorIndex();
        const DimMatrix& K1 = asImp_().fractureIntrinsicPermeability(context, interiorElemIdx, timeIdx);
        const DimMatrix& K2 = asImp_().fractureIntrinsicPermeability(context, exteriorElemIdx, timeIdx);

        // entry-wise harmonic mean. this is almost certainly wrong if
        // you have off-main diagonal entries in your permeabilities!
        for (unsigned i = 0; i < dimWorld; ++i)
            for (unsigned j = 0; j < dimWorld; ++j)
                result[i][j] = Opm::harmonicMean(K1[i][j], K2[i][j]);
    }
    /*!
     * \brief Returns the intrinsic permeability tensor \f$[m^2]\f$ at a given position due to a fracture
     *
     * \param context Reference to the object which represents the
     *                current execution context.
     * \param spaceIdx The local index of spatial entity defined by the context
     * \param timeIdx The index used by the time discretization.
     */
    template <class Context>
    const DimMatrix& fractureIntrinsicPermeability(const Context&,
                                                   unsigned,
                                                   unsigned) const
    {
        throw std::logic_error("Not implemented: Problem::fractureIntrinsicPermeability()");
    }

    /*!
     * \brief Returns the porosity [] inside fractures for a given control volume.
     *
     * \param context Reference to the object which represents the
     *                current execution context.
     * \param spaceIdx The local index of spatial entity defined by the context
     * \param timeIdx The index used by the time discretization.
     */
    template <class Context>
    Scalar fracturePorosity(const Context&,
                            unsigned,
                            unsigned) const
    {
        throw std::logic_error("Not implemented: Problem::fracturePorosity()");
    }

private:
    //! Returns the implementation of the problem (i.e. static polymorphism)
    Implementation& asImp_()
    { return *static_cast<Implementation *>(this); }
    //! \copydoc asImp_()
    const Implementation& asImp_() const
    { return *static_cast<const Implementation *>(this); }
};

} // namespace Opm

#endif