File: discretefracturelocalresidual.hh

package info (click to toggle)
opm-models 2022.10%2Bds-4
  • links: PTS, VCS
  • area: main
  • in suites: bookworm
  • size: 6,216 kB
  • sloc: cpp: 37,910; ansic: 1,897; sh: 277; xml: 182; makefile: 10
file content (159 lines) | stat: -rw-r--r-- 6,091 bytes parent folder | download | duplicates (3)
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
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
// -*- 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::DiscreteFractureLocalResidual
 */
#ifndef EWOMS_DISCRETE_FRACTURE_LOCAL_RESIDUAL_BASE_HH
#define EWOMS_DISCRETE_FRACTURE_LOCAL_RESIDUAL_BASE_HH

#include <opm/models/immiscible/immisciblelocalresidual.hh>

namespace Opm {

/*!
 * \ingroup DiscreteFractureModel
 *
 * \brief Calculates the local residual of the discrete fracture
 *        immiscible multi-phase model.
 */
template <class TypeTag>
class DiscreteFractureLocalResidual : public ImmiscibleLocalResidual<TypeTag>
{
    using ParentType = ImmiscibleLocalResidual<TypeTag>;

    using ElementContext = GetPropType<TypeTag, Properties::ElementContext>;
    using Indices = GetPropType<TypeTag, Properties::Indices>;
    using EqVector = GetPropType<TypeTag, Properties::EqVector>;
    using RateVector = GetPropType<TypeTag, Properties::RateVector>;
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;

    enum { conti0EqIdx = Indices::conti0EqIdx };
    enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
    enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };

    using EnergyModule = Opm::EnergyModule<TypeTag, enableEnergy>;

public:
    /*!
     * \brief Adds the amount all conservation quantities (e.g. phase
     *        mass) within a single fluid phase
     *
     * \copydetails Doxygen::storageParam
     * \copydetails Doxygen::dofCtxParams
     * \copydetails Doxygen::phaseIdxParam
     */
    void addPhaseStorage(EqVector& storage,
                         const ElementContext& elemCtx,
                         unsigned dofIdx,
                         unsigned timeIdx,
                         unsigned phaseIdx) const
    {
        EqVector phaseStorage(0.0);
        ParentType::addPhaseStorage(phaseStorage, elemCtx, dofIdx, timeIdx, phaseIdx);

        const auto& problem = elemCtx.problem();
        const auto& fractureMapper = problem.fractureMapper();
        unsigned globalIdx = elemCtx.globalSpaceIndex(dofIdx, timeIdx);

        if (!fractureMapper.isFractureVertex(globalIdx)) {
            // don't do anything in addition to the immiscible model for degrees of
            // freedom that do not feature fractures
            storage += phaseStorage;
            return;
        }

        const auto& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
        const auto& scv = elemCtx.stencil(timeIdx).subControlVolume(dofIdx);

        // reduce the matrix storage by the fracture volume
        phaseStorage *= 1 - intQuants.fractureVolume()/scv.volume();

        // add the storage term inside the fractures
        const auto& fsFracture = intQuants.fractureFluidState();

        phaseStorage[conti0EqIdx + phaseIdx] +=
            intQuants.fracturePorosity()*
            fsFracture.saturation(phaseIdx) *
            fsFracture.density(phaseIdx) *
            intQuants.fractureVolume()/scv.volume();

        EnergyModule::addFracturePhaseStorage(phaseStorage, intQuants, scv,
                                              phaseIdx);

        // add the result to the overall storage term
        storage += phaseStorage;
    }

    /*!
     * \copydoc FvBaseLocalResidual::computeFlux
     */
    void computeFlux(RateVector& flux,
                     const ElementContext& elemCtx,
                     unsigned scvfIdx,
                     unsigned timeIdx) const
    {
        ParentType::computeFlux(flux, elemCtx, scvfIdx, timeIdx);

        const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);

        unsigned i = extQuants.interiorIndex();
        unsigned j = extQuants.exteriorIndex();
        unsigned I = elemCtx.globalSpaceIndex(i, timeIdx);
        unsigned J = elemCtx.globalSpaceIndex(j, timeIdx);
        const auto& fractureMapper = elemCtx.problem().fractureMapper();
        if (!fractureMapper.isFractureEdge(I, J))
            // do nothing if the edge from i to j is not part of a
            // fracture
            return;

        const auto& scvf = elemCtx.stencil(timeIdx).interiorFace(scvfIdx);
        Scalar scvfArea = scvf.area();

        // advective mass fluxes of all phases
        for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
            if (!elemCtx.model().phaseIsConsidered(phaseIdx))
                continue;

            // reduce the matrix mass flux by the width of the scv
            // face that is occupied by the fracture. As usual, the
            // fracture is shared between two SCVs, so the its width
            // needs to be divided by two.
            flux[conti0EqIdx + phaseIdx] *=
                1 - extQuants.fractureWidth() / (2 * scvfArea);

            // intensive quantities of the upstream and the downstream DOFs
            unsigned upIdx = static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
            const auto& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
            flux[conti0EqIdx + phaseIdx] +=
                extQuants.fractureVolumeFlux(phaseIdx) * up.fractureFluidState().density(phaseIdx);
        }

        EnergyModule::handleFractureFlux(flux, elemCtx, scvfIdx, timeIdx);
    }
};

} // namespace Opm

#endif