stokesnccouplinglocalresidual.hh 15.9 KB
Newer Older
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
// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*****************************************************************************
 *   See the file COPYING for full copying permissions.                      *
 *                                                                           *
 *   This program 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.                                     *
 *                                                                           *
 *   This program 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 this program.  If not, see <http://www.gnu.org/licenses/>.   *
 *****************************************************************************/
/*!
 * \file
 *
 * \brief Element-wise calculation of the Jacobian matrix for problems
 *        using the coupled compositional Stokes box model.
 */

#ifndef DUMUX_STOKESNC_COUPLING_LOCAL_RESIDUAL_HH
#define DUMUX_STOKESNC_COUPLING_LOCAL_RESIDUAL_HH

29
30
#include <dune/common/deprecated.hh>

31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
#include <dumux/freeflow/stokesnc/stokesnclocalresidual.hh>
#include <dumux/freeflow/stokesnc/stokesncmodel.hh>

namespace Dumux
{
/*!
 * \ingroup ImplicitLocalResidual
 * \ingroup TwoPTwoCStokesTwoCModel
 * \ingroup TwoPTwoCZeroEqTwoCModel
 * \brief Element-wise calculation of the Jacobian matrix for problems
 *        using the coupled compositional Stokes box model.
 *        It is derived from the compositional Stokes box model.
 */
template<class TypeTag>
class StokesncCouplingLocalResidual : public StokesncLocalResidual<TypeTag>
{
47
    typedef StokesncLocalResidual<TypeTag> ParentType;
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

    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;

    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;

    enum {
        dim = GridView::dimension,
        dimWorld = GridView::dimensionworld,
        numEq = GET_PROP_VALUE(TypeTag, NumEq),
        numComponents = Indices::numComponents
    };
    enum {
        //indices of the equations
        massBalanceIdx = Indices::massBalanceIdx, //!< Index of the mass balance

        momentumXIdx = Indices::momentumXIdx, //!< Index of the x-component of the momentum balance
        momentumYIdx = Indices::momentumYIdx, //!< Index of the y-component of the momentum balance
        momentumZIdx = Indices::momentumZIdx, //!< Index of the z-component of the momentum balance
        lastMomentumIdx = Indices::lastMomentumIdx, //!< Index of the last component of the momentum balance
        transportEqIdx = Indices::transportEqIdx//!< Index of the transport equation
    };
    enum {
        //indices of phase and transported component
        phaseIdx = Indices::phaseIdx,
        transportCompIdx = Indices::transportCompIdx
    };
    enum {
        dimXIdx = Indices::dimXIdx, //!< Index for the first component of a vector
        dimYIdx = Indices::dimYIdx, //!< Index for the second component of a vector
        dimZIdx = Indices::dimZIdx //!< Index for the third component of a vector
    };

    typedef typename GridView::ctype CoordScalar;

    typedef Dune::FieldVector<CoordScalar, dimWorld> GlobalPosition;
    typedef Dune::FieldVector<Scalar, dim> DimVector;

    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;

    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) VolumeVariables;
    typedef typename GET_PROP_TYPE(TypeTag, FluxVariables) FluxVariables;
    typedef typename GET_PROP_TYPE(TypeTag, BoundaryTypes) BoundaryTypes;

    static const bool useMoles = GET_PROP_VALUE(TypeTag, UseMoles);

public:
    /*!
97
98
99
100
     * \brief Implementation of the boundary evaluation for the Stokes model
     *
     * Evaluate one part of the Dirichlet-like coupling conditions for a single
     * sub-control volume face; rest is done in the local coupling operator
101
102
103
     */
    void evalBoundary_()
    {
104
        ParentType::evalBoundary_();
105
106
107

        typedef Dune::ReferenceElements<Scalar, dim> ReferenceElements;
        typedef Dune::ReferenceElement<Scalar, dim> ReferenceElement;
108
109
110
        const ReferenceElement &refElement = ReferenceElements::general(this->element_().geometry().type());

        // loop over vertices of the element
111
        for (int scvIdx = 0; scvIdx < this->fvGeometry_().numScv; scvIdx++)
112
113
        {
            // consider only SCVs on the boundary
114
            if (this->fvGeometry_().subContVol[scvIdx].inner)
115
116
                continue;

117
118
119
            const BoundaryTypes &bcTypes = this->bcTypes_(scvIdx);

            // evaluate boundary conditions for the intersections of the current element
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
            for (const auto& intersection : Dune::intersections(this->gridView_(), this->element_()))
            {
                // handle only intersections on the boundary
                if (!intersection.boundary())
                    continue;

                // assemble the boundary for all vertices of the current face
                const int fIdx = intersection.indexInInside();
                const int numFaceVertices = refElement.size(fIdx, 1, dim);

                // loop over the single vertices on the current face
                for (int faceVertexIdx = 0; faceVertexIdx < numFaceVertices; ++faceVertexIdx)
                {
                    // only evaluate, if we consider the same face vertex as in the outer
                    // loop over the element vertices
135
                    if (refElement.subEntity(fIdx, 1, faceVertexIdx, dim) != scvIdx)
136
137
138
139
140
141
142
143
144
                        continue;

                    const int boundaryFaceIdx = this->fvGeometry_().boundaryFaceIndex(fIdx, faceVertexIdx);
                    const FluxVariables boundaryVars(this->problem_(),
                                                     this->element_(),
                                                     this->fvGeometry_(),
                                                     boundaryFaceIdx,
                                                     this->curVolVars_(),
                                                     true);
145
146
147
                    const VolumeVariables &volVars = this->curVolVars_()[scvIdx];

                    // set velocity normal to the interface
148
                    if (bcTypes.isCouplingDirichlet(momentumYIdx))
149
150
151
152
153
154
155
156
157
                    {
                        this->residual_[scvIdx][momentumYIdx] = volVars.velocity()
                                                                * boundaryVars.face().normal
                                                                / boundaryVars.face().normal.two_norm();
                        Valgrind::CheckDefined(this->residual_[scvIdx][momentumYIdx]);
                    }

                    // add pressure correction - required for pressure coupling,
                    // if p.n comes from the pm
158
                    if (bcTypes.isCouplingNeumann(momentumYIdx) || bcTypes.isCouplingMortar(momentumYIdx))
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
                    {
                        DimVector pressureCorrection(intersection.centerUnitOuterNormal());
                        pressureCorrection *= volVars.pressure();
                        this->residual_[scvIdx][momentumYIdx] += pressureCorrection[momentumYIdx]
                                                                 * boundaryVars.face().area;
                        Valgrind::CheckDefined(this->residual_[scvIdx][momentumYIdx]);
                    }

                    // set mole or mass fraction for the transported components
                    for (int compIdx = 0; compIdx < numComponents; compIdx++)
                    {
                        int eqIdx = dim + compIdx; // TODO: ist das so richtig
                        if (eqIdx != massBalanceIdx)
                        {
                            if (bcTypes.isCouplingDirichlet(eqIdx))
                            {
                                if(useMoles)
                                    this->residual_[scvIdx][eqIdx] = volVars.moleFraction(compIdx);
                                else
                                    this->residual_[scvIdx][eqIdx] = volVars.massFraction(compIdx);
                            }
                        }
                    }
                }
            }
        }
185
186
    }

187
188
189
    /*!
     * \brief Removes the stabilization for the Stokes model.
     */
190
191
192
    void evalBoundaryPDELab_()
    {
        // loop over vertices of the element
193
        for (int scvIdx = 0; scvIdx < this->fvGeometry_().numScv; scvIdx++)
194
195
        {
            // consider only SCVs on the boundary
196
            if (this->fvGeometry_().subContVol[scvIdx].inner)
197
198
                continue;

199
200
            this->removeStabilizationAtBoundary_(scvIdx);
        }
201
202
203
204
205
206
207
208
    }

protected:
    /*!
     * \brief Evaluate one part of the Dirichlet-like coupling conditions for a single
     *        sub-control volume face; rest is done in the local coupling operator
     */
    template <class IntersectionIterator>
209
    DUNE_DEPRECATED_MSG("evalCouplingVertex_ is deprecated. Its functionality is now included in evalBoundary_.")
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
    void evalCouplingVertex_(const IntersectionIterator &isIt,
                             const int scvIdx,
                             const int boundaryFaceIdx,
                             const FluxVariables &boundaryVars)
    {
        const VolumeVariables &volVars = this->curVolVars_()[scvIdx];
        const BoundaryTypes &bcTypes = this->bcTypes_(scvIdx);

        // set velocity normal to the interface
        if (bcTypes.isCouplingNeumann(momentumYIdx))
            this->residual_[scvIdx][momentumYIdx] =
                    volVars.velocity() *
                    boundaryVars.face().normal /
                    boundaryVars.face().normal.two_norm();
        Valgrind::CheckDefined(this->residual_[scvIdx][momentumYIdx]);

        // add pressure correction - required for pressure coupling,
        // if p.n comes from the pm
        if (bcTypes.isCouplingDirichlet(momentumYIdx) || bcTypes.isCouplingMortar(momentumYIdx))
        {
            DimVector pressureCorrection(isIt->centerUnitOuterNormal());
            pressureCorrection *= volVars.pressure();
            this->residual_[scvIdx][momentumYIdx] += pressureCorrection[momentumYIdx]*
                    boundaryVars.face().area;
        }

        // set mole or mass fraction for the transported components
        for (int compIdx = 0; compIdx < numComponents; compIdx++)
        {
            int eqIdx = dim + compIdx;
            if (eqIdx != massBalanceIdx) {
                if (bcTypes.isCouplingDirichlet(eqIdx))
                {
                    if(useMoles)
                        this->residual_[scvIdx][eqIdx] = volVars.moleFraction(compIdx);
                    else
                        this->residual_[scvIdx][eqIdx] = volVars.massFraction(compIdx);
                    Valgrind::CheckDefined(this->residual_[scvIdx][eqIdx]);
                }
            }
        }
    }

    template <class IntersectionIterator>
254
    DUNE_DEPRECATED_MSG("evalBeaversJoseph_ is deprecated. Its functionality is now included in the LOP function evalCoupling21().")
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
    void evalBeaversJoseph_(const IntersectionIterator &isIt,
                            const int scvIdx,
                            const int boundaryFaceIdx,
                            const FluxVariables &boundaryVars)
    {
        const BoundaryTypes &bcTypes = this->bcTypes_(scvIdx);

        const GlobalPosition &globalPos = this->fvGeometry_().boundaryFace[boundaryFaceIdx].ipGlobal;
        Scalar beaversJosephCoeff = this->problem_().beaversJosephCoeffAtPos(globalPos);

        const Scalar Kxx = this->problem_().permeability(this->element_(),
                                                         this->fvGeometry_(),
                                                         scvIdx);
        beaversJosephCoeff /= std::sqrt(Kxx);
        const DimVector& elementUnitNormal = isIt->centerUnitOuterNormal();

        // implementation as NEUMANN condition /////////////////////////////////////////////
        // (v.n)n
        if (bcTypes.isCouplingDirichlet(momentumXIdx))
        {
            const Scalar normalComp = boundaryVars.velocity()*elementUnitNormal;
            DimVector normalV = elementUnitNormal;
            normalV *= normalComp; // v*n*n

            // v_tau = v - (v.n)n
            const DimVector tangentialV = boundaryVars.velocity() - normalV;
            const Scalar boundaryFaceArea = boundaryVars.face().area;

            for (int dimIdx=0; dimIdx < dim; ++dimIdx)
                this->residual_[scvIdx][dimIdx] += beaversJosephCoeff *
                                        boundaryFaceArea *
                                        tangentialV[dimIdx] *
                                        (boundaryVars.dynamicViscosity()
                                          + boundaryVars.dynamicEddyViscosity());

            ////////////////////////////////////////////////////////////////////////////////////
            //normal component has only to be set if no coupling conditions are defined
            //set NEUMANN flux (set equal to pressure in problem)
//             PrimaryVariables priVars(0.0);
//             this->problem_().neumann(priVars, this->element_(), this->fvGeometry_(),
//                             intersection, scvIdx, boundaryFaceIdx);
//             for (int dimIdx=0; dimIdx < dim; ++dimIdx)
//                 this->residual_[scvIdx][dimIdx] += priVars[dimIdx]*
//                                                    boundaryFaceArea;
        }
        if (bcTypes.isCouplingNeumann(momentumXIdx))
        {
            assert(beaversJosephCoeff > 0);
            ///////////////////////////////////////////////////////////////////////////////////////////
            //IMPLEMENTATION AS DIRICHLET CONDITION
            //tangential componment: vx = sqrt K /alpha * (grad v n(unity))t
            DimVector tangentialVelGrad(0);
            boundaryVars.velocityGrad().umv(elementUnitNormal, tangentialVelGrad);
            tangentialVelGrad /= -beaversJosephCoeff; // was - before
            this->residual_[scvIdx][momentumXIdx] =
                    tangentialVelGrad[momentumXIdx] - this->curPriVars_(scvIdx)[momentumXIdx];

            ////////////////////////////////////////////////////////////////////////////////////
            //for testing Beavers and Joseph without adjacent porous medium set vy to 0
//                Scalar normalVel(0);
//                this->residual_[scvIdx][momentumYIdx] = this->curPrimaryVars_(scvIdx)[momentumYIdx] - normalVel;
            ////////////////////////////////////////////////////////////////////////////////////
        }
    }

    //! \brief Return true, if at least one equation on the boundary has a  coupling condition
    DUNE_DEPRECATED_MSG("boundaryHasCoupling_() is deprecated")
    bool boundaryHasCoupling_(const BoundaryTypes& bcTypes) const
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
            if (bcTypes.isCouplingInflow(eqIdx) || bcTypes.isCouplingOutflow(eqIdx))
                return true;
        return false;
    }

    //! \brief Return true, if at least one equation on the boundary has a mortar coupling condition
    DUNE_DEPRECATED_MSG("boundaryHasMortarCoupling_() is deprecated")
    bool boundaryHasMortarCoupling_(const BoundaryTypes& bcTypes) const
    {
        for (int eqIdx = 0; eqIdx < numEq; ++eqIdx)
            if (bcTypes.isMortarCoupling(eqIdx))
                return true;
        return false;
    }
};

} // namespace Dumux

#endif // DUMUX_STOKESNC_COUPLING_LOCAL_RESIDUAL_HH