2cstokes2p2clocaloperator.hh 49.6 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 The local operator for the coupling of a two-component Stokes model
 *        and a two-phase two-component porous-medium model under isothermal conditions.
 */
#ifndef DUMUX_2CSTOKES_2P2C_LOCALOPERATOR_HH
#define DUMUX_2CSTOKES_2P2C_LOCALOPERATOR_HH

#include <iostream>

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

31
32
33
34
35
#include <dune/pdelab/multidomain/couplingutilities.hh>
#include <dune/pdelab/localoperator/pattern.hh>
#include <dune/pdelab/localoperator/idefault.hh>

#include <dumux/multidomain/common/multidomainproperties.hh>
Thomas Fetzer's avatar
Thomas Fetzer committed
36
#include <dumux/multidomain/2cstokes2p2c/2cstokes2p2cpropertydefaults.hh>
Thomas Fetzer's avatar
Thomas Fetzer committed
37
38
#include <dumux/freeflow/boundarylayermodel.hh>
#include <dumux/freeflow/masstransfermodel.hh>
39
40
41
#include <dumux/freeflow/stokesnc/stokesncmodel.hh>
#include <dumux/implicit/2p2c/2p2cmodel.hh>

Thomas Fetzer's avatar
Thomas Fetzer committed
42

43
44
45
namespace Dumux {

/*!
Thomas Fetzer's avatar
Thomas Fetzer committed
46
 * \ingroup TwoPTwoCStokesTwoCModel
Hao Wu's avatar
Hao Wu committed
47
 * \ingroup TwoPTwoCZeroEqTwoCModel
48
49
 * \brief The local operator for the coupling of a two-component Stokes model
 *        and a two-phase two-component porous-medium model under isothermal conditions.
Hao Wu's avatar
Hao Wu committed
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
 *
 * This model implements the coupling between a free-flow model
 * and a porous-medium flow model under isothermal conditions.
 * Here the coupling conditions for the individual balance are presented:
 *
 * The total mass balance equation:
 * \f[
 *  \left[
 *    \left( \varrho_\textrm{g} {\boldsymbol{v}}_\textrm{g} \right) \cdot \boldsymbol{n}
 *  \right]^\textrm{ff}
 *  = -\left[
 *      \left( \varrho_\textrm{g} \boldsymbol{v}_\textrm{g}
 *             + \varrho_\textrm{l} \boldsymbol{v}_\textrm{l} \right) \cdot \boldsymbol{n}
 *    \right]^\textrm{pm}
 * \f]
 * in which \f$n\f$ represents a vector normal to the interface pointing outside of
 * the specified subdomain.
 *
 * The momentum balance (tangential), which corresponds to the Beavers-Jospeh Saffman condition:
 * \f[
 *  \left[
 *   \left( {\boldsymbol{v}}_\textrm{g}
 *          + \frac{\sqrt{\left(\boldsymbol{K} \boldsymbol{t}_i \right) \cdot \boldsymbol{t}_i}}
 *                 {\alpha_\textrm{BJ} \mu_\textrm{g}} \boldsymbol{{\tau}}_\textrm{t} \boldsymbol{n}
 *          \right) \cdot \boldsymbol{t}_i
 *  \right]^\textrm{ff}
 *  = 0
 * \f]
 * with
 * \f$
 * \boldsymbol{{\tau}_\textrm{t}} = \left[ \mu_\textrm{g} + \mu_\textrm{g,t} \right]
 *                                  \nabla \left( \boldsymbol{v}_\textrm{g}
 *                                                + \boldsymbol{v}_\textrm{g}^\intercal \right)
 * \f$
 * in which the eddy viscosity \f$ \mu_\textrm{g,t} = 0 \f$ for the Stokes equation.
 *
 * The momentum balance (normal):
 * \f[
 *  \left[
 *    \left(
 *      \left\lbrace
 *        \varrho_\textrm{g} {\boldsymbol{v}}_\textrm{g} {\boldsymbol{v}}_\textrm{g}^\intercal
 *        - \boldsymbol{{\tau}}_\textrm{t}
 *        + {p}_\textrm{g} \boldsymbol{I}
 *      \right\rbrace \boldsymbol{n}
 *    \right) \cdot \boldsymbol{n}
 *  \right]^\textrm{ff}
 *  = p_\textrm{g}^\textrm{pm}
 * \f]
 *
 * The component mass balance equation (continuity of fluxes):
 * \f[
 *  \left[
 *    \left(
 *      \varrho_\textrm{g} {X}^\kappa_\textrm{g} {\boldsymbol{v}}_\textrm{g}
 *      - {\boldsymbol{j}}^\kappa_\textrm{g,ff,t,diff}
 *    \right) \cdot \boldsymbol{n}
 *  \right]^\textrm{ff}
 *  = -\left[
 *    \left(
 *      \varrho_\textrm{g} X^\kappa_\textrm{g} \boldsymbol{v}_\textrm{g}
 *      - \boldsymbol{j}^\kappa_\textrm{g,pm,diff}
 *      + \varrho_\textrm{l} \boldsymbol{v}_\textrm{l} X^\kappa_\textrm{l}
 *      - \boldsymbol{j}^\kappa_\textrm{l,pm,diff}
 *    \right) \cdot \boldsymbol{n}
 *  \right]^\textrm{pm}
 *  = 0
 * \f]
 * in which the diffusive fluxes \f$ j_\textrm{diff} \f$ are the diffusive fluxes as
 * they are implemented in the individual subdomain models.
 *
 * The component mass balance equation (continuity of mass/ mole fractions):
 * \f[
 *  \left[ {X}^{\kappa}_\textrm{g} \right]^\textrm{ff}
 *  = \left[ X^{\kappa}_\textrm{g} \right]^\textrm{pm}
 * \f]
 *
 * This is discretized by a fully-coupled vertex-centered finite volume
 * (box) scheme in space and by the implicit Euler method in time.
129
130
131
132
133
134
135
136
 */
template<class TypeTag>
class TwoCStokesTwoPTwoCLocalOperator :
        public Dune::PDELab::MultiDomain::CouplingOperatorDefaultFlags,
        public Dune::PDELab::MultiDomain::NumericalJacobianCoupling<TwoCStokesTwoPTwoCLocalOperator<TypeTag>>,
        public Dune::PDELab::MultiDomain::FullCouplingPattern,
        public Dune::PDELab::InstationaryLocalOperatorDefaultMethods<double>
{
137
public:
138
139
140
141
142
143
    typedef typename GET_PROP_TYPE(TypeTag, Problem) GlobalProblem;
    typedef typename GET_PROP_TYPE(TypeTag, MultiDomainCouplingLocalOperator) Implementation;

    typedef typename GET_PROP_TYPE(TypeTag, SubDomain1TypeTag) Stokes2cTypeTag;
    typedef typename GET_PROP_TYPE(TypeTag, SubDomain2TypeTag) TwoPTwoCTypeTag;

144
    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
145
    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, SpatialParams) SpatialParams;
146
147
148
149
150
151
152
153
154
155
156
157
158

    typedef typename GET_PROP_TYPE(Stokes2cTypeTag, ElementVolumeVariables) ElementVolumeVariables1;
    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, ElementVolumeVariables) ElementVolumeVariables2;

    typedef typename GET_PROP_TYPE(Stokes2cTypeTag, FluxVariables) BoundaryVariables1;
    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, FluxVariables) BoundaryVariables2;

    typedef typename GET_PROP_TYPE(Stokes2cTypeTag, BoundaryTypes) BoundaryTypes1;
    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, BoundaryTypes) BoundaryTypes2;

    typedef typename GET_PROP_TYPE(Stokes2cTypeTag, FVElementGeometry) FVElementGeometry1;
    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, FVElementGeometry) FVElementGeometry2;

Thomas Fetzer's avatar
Thomas Fetzer committed
159
    // Multidomain Grid types
160
161
162
163
164
165
166
167
168
169
170
    typedef typename GET_PROP_TYPE(TypeTag, MultiDomainGrid) MDGrid;
    typedef typename MDGrid::Traits::template Codim<0>::Entity MDElement;

    typedef typename GET_PROP_TYPE(Stokes2cTypeTag, GridView) Stokes2cGridView;
    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, GridView) TwoPTwoCGridView;
    typedef typename Stokes2cGridView::template Codim<0>::Entity SDElement1;
    typedef typename TwoPTwoCGridView::template Codim<0>::Entity SDElement2;

    typedef typename GET_PROP_TYPE(Stokes2cTypeTag, Indices) Stokes2cIndices;
    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, Indices) TwoPTwoCIndices;

171
172
173
174
    enum {
        dim = MDGrid::dimension,
        dimWorld = MDGrid::dimensionworld
    };
175

176
    // Stokes
177
    enum { numEq1 = GET_PROP_VALUE(Stokes2cTypeTag, NumEq) };
178
    enum { numComponents1 = Stokes2cIndices::numComponents };
179
180
181
182
183
184
185
    enum { // equation indices
        momentumXIdx1 = Stokes2cIndices::momentumXIdx,         //!< Index of the x-component of the momentum balance
        momentumYIdx1 = Stokes2cIndices::momentumYIdx,         //!< Index of the y-component of the momentum balance
        momentumZIdx1 = Stokes2cIndices::momentumZIdx,         //!< Index of the z-component of the momentum balance
        lastMomentumIdx1 = Stokes2cIndices::lastMomentumIdx,   //!< Index of the last component of the momentum balance
        massBalanceIdx1 = Stokes2cIndices::massBalanceIdx,     //!< Index of the mass balance
        transportEqIdx1 = Stokes2cIndices::transportEqIdx      //!< Index of the transport equation
186
    };
187
188
189
    enum { // component indices
        transportCompIdx1 = Stokes2cIndices::transportCompIdx, //!< Index of transported component
        phaseCompIdx1 = Stokes2cIndices::phaseCompIdx    //!< Index of main component of the phase
190
191
    };

192
    // Darcy
193
194
    enum { numEq2 = GET_PROP_VALUE(TwoPTwoCTypeTag, NumEq) };
    enum { numPhases2 = GET_PROP_VALUE(TwoPTwoCTypeTag, NumPhases) };
195
196
197
    enum { // equation indices
        contiWEqIdx2 = TwoPTwoCIndices::contiWEqIdx,     //!< Index of the continuity equation for water component
        massBalanceIdx2 = TwoPTwoCIndices::contiNEqIdx   //!< Index of the total mass balance (if one comopnent balance is replaced)
198
199
    };
    enum { // component indices
200
201
        wCompIdx2 = TwoPTwoCIndices::wCompIdx,           //!< Index of the liquids main component
        nCompIdx2 = TwoPTwoCIndices::nCompIdx            //!< Index of the main component of the gas
202
203
    };
    enum { // phase indices
204
        wPhaseIdx2 = TwoPTwoCIndices::wPhaseIdx,         //!< Index for the liquid phase
205
206
207
        nPhaseIdx2 = TwoPTwoCIndices::nPhaseIdx          //!< Index for the gas phase
    };

208
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
209
210
211
    typedef typename MDGrid::ctype CoordScalar;
    typedef Dune::FieldVector<CoordScalar, dimWorld> GlobalPosition;

212
213
214
    typedef typename Stokes2cGridView::template Codim<dim>::EntityPointer VertexPointer1;
    typedef typename TwoPTwoCGridView::template Codim<dim>::EntityPointer VertexPointer2;

215
216
217
218
219
    // multidomain flags
    static const bool doAlphaCoupling = true;
    static const bool doPatternCoupling = true;

public:
Thomas Fetzer's avatar
Thomas Fetzer committed
220
    //! \brief The constructor
221
222
223
    TwoCStokesTwoPTwoCLocalOperator(GlobalProblem& globalProblem)
        : globalProblem_(globalProblem)
    {
Thomas Fetzer's avatar
Thomas Fetzer committed
224
225
226
        static_assert(GET_PROP_VALUE(Stokes2cTypeTag, UseMoles) == GET_PROP_VALUE(TwoPTwoCTypeTag, UseMoles),
                      "The coupling conditions is only implemented for same formulations (mass or mole) in both subdomains.");

Thomas Fetzer's avatar
Thomas Fetzer committed
227
228
        blModel_ = GET_PARAM_FROM_GROUP(TypeTag, int, BoundaryLayer, Model);
        massTransferModel_ = GET_PARAM_FROM_GROUP(TypeTag, int, MassTransfer, Model);
Thomas Fetzer's avatar
Thomas Fetzer committed
229

Thomas Fetzer's avatar
Thomas Fetzer committed
230
        if (blModel_ != 0)
Thomas Fetzer's avatar
Thomas Fetzer committed
231
            std::cout << "Using boundary layer model " << blModel_ << std::endl;
Thomas Fetzer's avatar
Thomas Fetzer committed
232
        if (massTransferModel_ != 0)
Thomas Fetzer's avatar
Thomas Fetzer committed
233
            std::cout << "Using mass transfer model " << massTransferModel_ << std::endl;
234
235
236
237
238
239
240
241
    }

    /*!
     * \brief Do the coupling. The unknowns are transferred from dune-multidomain.
     *        Based on them, a coupling residual is calculated and added at the
     *        respective positions in the matrix.
     *
     * \param intersectionGeometry the geometry of the intersection
242
     * \param lfsu1 local basis for the trial space of the Stokes domain
243
     * \param unknowns1 the unknowns vector of the Stokes element (formatted according to PDELab)
244
245
     * \param lfsv1 local basis for the test space of the Stokes domain
     * \param lfsu2 local basis for the trail space of the Darcy domain
246
     * \param unknowns2 the unknowns vector of the Darcy element (formatted according to PDELab)
247
     * \param lfsv2 local basis for the test space of the Darcy domain
248
249
250
251
252
     * \param couplingRes1 the coupling residual from the Stokes domain
     * \param couplingRes2 the coupling residual from the Darcy domain
     */
    template<typename IntersectionGeom, typename LFSU1, typename LFSU2,
             typename X, typename LFSV1, typename LFSV2,typename RES>
253
254
255
256
    void alpha_coupling(const IntersectionGeom& intersectionGeometry,
                        const LFSU1& lfsu1, const X& unknowns1, const LFSV1& lfsv1,
                        const LFSU2& lfsu2, const X& unknowns2, const LFSV2& lfsv2,
                        RES& couplingRes1, RES& couplingRes2) const
257
    {
258
259
        const MDElement& mdElement1 = intersectionGeometry.inside();
        const MDElement& mdElement2 = intersectionGeometry.outside();
260
261

        // the subodmain elements
262
263
        const SDElement1& sdElement1 = globalProblem_.sdElementPointer1(mdElement1);
        const SDElement2& sdElement2 = globalProblem_.sdElementPointer2(mdElement2);
264
265
266
267
268

        // a container for the parameters on each side of the coupling interface (see below)
        CParams cParams;

        // update fvElementGeometry and the element volume variables
269
        updateElemVolVars(lfsu1, lfsu2,
270
271
272
273
274
275
                          unknowns1, unknowns2,
                          sdElement1, sdElement2,
                          cParams);

        // first element
        const int faceIdx1 = intersectionGeometry.indexInInside();
276
277
        const Dune::ReferenceElement<typename MDGrid::ctype,dim>& referenceElement1 =
            Dune::ReferenceElements<typename MDGrid::ctype,dim>::general(mdElement1.type());
278
279
280
281
        const int numVerticesOfFace = referenceElement1.size(faceIdx1, 1, dim);

        // second element
        const int faceIdx2 = intersectionGeometry.indexInOutside();
282
283
        const Dune::ReferenceElement<typename MDGrid::ctype,dim>& referenceElement2 =
            Dune::ReferenceElements<typename MDGrid::ctype,dim>::general(mdElement2.type());
284
285
286
287
288
289
290
291
292
293
294
295
296

        for (int vertexInFace = 0; vertexInFace < numVerticesOfFace; ++vertexInFace)
        {
            const int vertInElem1 = referenceElement1.subEntity(faceIdx1, 1, vertexInFace, dim);
            const int vertInElem2 = referenceElement2.subEntity(faceIdx2, 1, vertexInFace, dim);

            const int boundaryFaceIdx1 = cParams.fvGeometry1.boundaryFaceIndex(faceIdx1, vertexInFace);
            const int boundaryFaceIdx2 = cParams.fvGeometry2.boundaryFaceIndex(faceIdx2, vertexInFace);

            // obtain the boundary types
            const VertexPointer1 vPtr1 = sdElement1.template subEntity<dim>(vertInElem1);
            const VertexPointer2 vPtr2 = sdElement2.template subEntity<dim>(vertInElem2);

297
298
            globalProblem_.sdProblem1().boundaryTypes(cParams.boundaryTypes1, vPtr1);
            globalProblem_.sdProblem2().boundaryTypes(cParams.boundaryTypes2, vPtr2);
299
300
301
302
303
304
305
306
307
308
309
310
311
312

            const BoundaryVariables1 boundaryVars1(globalProblem_.sdProblem1(),
                                                   sdElement1,
                                                   cParams.fvGeometry1,
                                                   boundaryFaceIdx1,
                                                   cParams.elemVolVarsCur1,
                                                   /*onBoundary=*/true);
            const BoundaryVariables2 boundaryVars2(globalProblem_.sdProblem2(),
                                                   sdElement2,
                                                   cParams.fvGeometry2,
                                                   boundaryFaceIdx2,
                                                   cParams.elemVolVarsCur2,
                                                   /*onBoundary=*/true);

313
314
315
316
317
318
            asImp_()->evalCoupling(lfsu1, lfsu2,
                                   vertInElem1, vertInElem2,
                                   sdElement1, sdElement2,
                                   boundaryVars1, boundaryVars2,
                                   cParams,
                                   couplingRes1, couplingRes2);
319
320
321
322
323
324
325
        }
    }

    /*!
     * \brief Update the volume variables of the element and extract the unknowns from dune-pdelab vectors
     *        and bring them into a form which fits to dumux.
     *
326
327
     * \param lfsu1 local basis for the trial space of the Stokes domain
     * \param lfsu2 local basis for the trial space of the Darcy domain
328
329
330
331
332
333
334
     * \param unknowns1 the unknowns vector of the Stokes element (formatted according to PDELab)
     * \param unknowns2 the unknowns vector of the Darcy element (formatted according to PDELab)
     * \param sdElement1 the element in the Stokes domain
     * \param sdElement2 the element in the Darcy domain
     * \param cParams a parameter container
     */
    template<typename LFSU1, typename LFSU2, typename X, typename CParams>
335
336
337
338
    void updateElemVolVars(const LFSU1& lfsu1, const LFSU2& lfsu2,
                           const X& unknowns1, const X& unknowns2,
                           const SDElement1& sdElement1, const SDElement2& sdElement2,
                           CParams &cParams) const
339
340
341
342
    {
        cParams.fvGeometry1.update(globalProblem_.sdGridView1(), sdElement1);
        cParams.fvGeometry2.update(globalProblem_.sdGridView2(), sdElement2);

343
344
        const int numVertsOfElem1 = sdElement1.subEntities(dim);
        const int numVertsOfElem2 = sdElement2.subEntities(dim);
345

346
        // bring the local unknowns x1 into a form that can be passed to elemVolVarsCur.update()
347
348
349
350
351
352
353
354
        Dune::BlockVector<Dune::FieldVector<Scalar,1>> elementSol1(0.);
        Dune::BlockVector<Dune::FieldVector<Scalar,1>> elementSol2(0.);
        elementSol1.resize(unknowns1.size());
        elementSol2.resize(unknowns2.size());

        for (int idx=0; idx<numVertsOfElem1; ++idx)
        {
            for (int eqIdx1=0; eqIdx1<numEq1; ++eqIdx1)
355
                elementSol1[eqIdx1*numVertsOfElem1+idx] = unknowns1(lfsu1.child(eqIdx1),idx);
356
            for (int eqIdx2=0; eqIdx2<numEq2; ++eqIdx2)
357
                elementSol2[eqIdx2*numVertsOfElem2+idx] = unknowns2(lfsu2.child(eqIdx2),idx);
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
        }
#if HAVE_VALGRIND
        for (unsigned int i = 0; i < elementSol1.size(); i++)
            Valgrind::CheckDefined(elementSol1[i]);
        for (unsigned int i = 0; i < elementSol2.size(); i++)
            Valgrind::CheckDefined(elementSol2[i]);
#endif // HAVE_VALGRIND


        // evaluate the local residual with the PDELab solution
        globalProblem_.localResidual1().evalPDELab(sdElement1, cParams.fvGeometry1, elementSol1,
                                                   cParams.elemVolVarsPrev1, cParams.elemVolVarsCur1);
        globalProblem_.localResidual2().evalPDELab(sdElement2, cParams.fvGeometry2, elementSol2,
                                                   cParams.elemVolVarsPrev2, cParams.elemVolVarsCur2);

    }

375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
    /*!
     * \brief Evaluation of the coupling between the Stokes (1) and Darcy (2).
     *
     * Dirichlet-like and Neumann-like conditions for the respective domain are evaluated.
     *
     * \param lfsu1 local basis for the trial space of the Stokes domain
     * \param lfsu2 local basis for the trial space of the Darcy domain
     * \param vertInElem1 local vertex index in element1
     * \param vertInElem2 local vertex index in element2
     * \param sdElement1 the element in the Stokes domain
     * \param sdElement2 the element in the Darcy domain
     * \param boundaryVars1 the boundary variables at the interface of the Stokes domain
     * \param boundaryVars2 the boundary variables at the interface of the Darcy domain
     * \param cParams a parameter container
     * \param couplingRes1 the coupling residual from the Stokes domain
     * \param couplingRes2 the coupling residual from the Darcy domain
     */
    template<typename LFSU1, typename LFSU2, typename RES1, typename RES2, typename CParams>
    void evalCoupling(const LFSU1& lfsu1, const LFSU2& lfsu2,
                      const int vertInElem1, const int vertInElem2,
                      const SDElement1& sdElement1, const SDElement2& sdElement2,
                      const BoundaryVariables1& boundaryVars1, const BoundaryVariables2& boundaryVars2,
                      const CParams &cParams,
                      RES1& couplingRes1, RES2& couplingRes2) const
    {
        const GlobalPosition& globalPos1 = cParams.fvGeometry1.subContVol[vertInElem1].global;
        const GlobalPosition& globalPos2 = cParams.fvGeometry2.subContVol[vertInElem2].global;

        const GlobalPosition& bfNormal1 = boundaryVars1.face().normal;
        const Scalar normalMassFlux1 = boundaryVars1.normalVelocity()
                                       * cParams.elemVolVarsCur1[vertInElem1].density();

        // MASS Balance
        // Neumann-like conditions
        if (cParams.boundaryTypes2.isCouplingNeumann(massBalanceIdx2))
        {
            static_assert(!GET_PROP_VALUE(TwoPTwoCTypeTag, UseMoles),
                          "This coupling condition is only implemented for mass fraction formulation.");

            if (globalProblem_.sdProblem1().isCornerPoint(globalPos1))
            {
                couplingRes2.accumulate(lfsu2.child(massBalanceIdx2), vertInElem2,
                                        -normalMassFlux1);
            }
            else
            {
                couplingRes2.accumulate(lfsu2.child(massBalanceIdx2), vertInElem2,
                                        globalProblem_.localResidual1().residual(vertInElem1)[massBalanceIdx1]);
            }
        }

        // Dirichlet-like
        if (cParams.boundaryTypes2.isCouplingDirichlet(massBalanceIdx2))
        {
            couplingRes2.accumulate(lfsu2.child(massBalanceIdx2), vertInElem2,
                                    globalProblem_.localResidual1().residual(vertInElem1)[momentumYIdx1]
                                    -cParams.elemVolVarsCur1[vertInElem1].pressure());
        }


        // MOMENTUM_X Balance
436
        SpatialParams spatialParams = globalProblem_.sdProblem2().spatialParams();
437
438
        Scalar beaversJosephCoeff = spatialParams.beaversJosephCoeffAtPos(globalPos1);
        assert(beaversJosephCoeff > 0);
439
        beaversJosephCoeff /= std::sqrt(spatialParams.intrinsicPermeability(sdElement2, cParams.fvGeometry2, vertInElem2));
440

441
        // Neumann-like conditions
442
443
444
445
446
447
448
449
450
        if (cParams.boundaryTypes1.isCouplingNeumann(momentumXIdx1))
        {
            const Scalar normalComp = boundaryVars1.velocity()*bfNormal1;
            GlobalPosition normalV = bfNormal1;
            normalV *= normalComp; // v*n*n

            // v_tau = v - (v.n)n
            const GlobalPosition tangentialV = boundaryVars1.velocity() - normalV;

451
            // Implementation as Neumann-like condition: (v.n)n
452
453
454
455
456
457
458
459
460
461
462
            for (int dimIdx=0; dimIdx < dim; ++dimIdx)
            {
                couplingRes1.accumulate(lfsu1.child(momentumXIdx1), vertInElem1,
                                        beaversJosephCoeff
                                        * boundaryVars1.face().area
                                        * tangentialV[dimIdx]
                                        * (boundaryVars1.dynamicViscosity()
                                          + boundaryVars1.dynamicEddyViscosity()));
            }
        }

463
464
465
466
467
468
469
470
471
472
473
474
475
476
        // Dirichlet-like conditions
        if (cParams.boundaryTypes1.isCouplingDirichlet(momentumXIdx1))
        {
            // NOTE: This boundary condition is not implemented anymore because curPrimaryVars_ is protected
            DUNE_THROW(Dune::NotImplemented, "The boundary condition isCouplingNeumann(momentumXIdx1) on the Stokes side is not implemented anymore.");

            // tangential component: vx = sqrt K /alpha * (grad v n(unity))t
            // GlobalPosition tangentialVelGrad(0);
            // boundaryVars1.velocityGrad().umv(elementUnitNormal, tangentialVelGrad);
            // tangentialVelGrad /= -beaversJosephCoeff; // was - before
            // this->residual_[vertInElem1][momentumXIdx1] =
            //        tangentialVelGrad[momentumXIdx1] - globalProblem_.localResidual1().curPriVars_(vertInElem1)[momentumXIdx1]);
        }

477
478
479
480
481

        // MOMENTUM_Y Balance
        // Neumann-like conditions
        if (cParams.boundaryTypes1.isCouplingNeumann(momentumYIdx1))
        {
482
483
484
485
486
487
488
489
490
491
492
            // p*A*n as condition for free flow
            // pressure correction is done in stokeslocalresidual.hh
            couplingRes1.accumulate(lfsu1.child(momentumYIdx1), vertInElem1,
                                    cParams.elemVolVarsCur2[vertInElem2].pressure(nPhaseIdx2) *
                                    boundaryVars2.face().area);
        }

        // Dirichlet-like conditions
        if (cParams.boundaryTypes1.isCouplingDirichlet(momentumYIdx1))
        {
            // v.n as Dirichlet-like condition for the Stokes domain
493
494
495
496
497
498
499
500
501
502
503
504
505
506
            if (globalProblem_.sdProblem2().isCornerPoint(globalPos2))
            {
                Scalar sumNormalPhaseFluxes = 0.0;
                for (int phaseIdx=0; phaseIdx<numPhases2; ++phaseIdx)
                {
                    sumNormalPhaseFluxes -= boundaryVars2.volumeFlux(phaseIdx)
                                            * cParams.elemVolVarsCur2[vertInElem2].density(phaseIdx);
                }
                couplingRes1.accumulate(lfsu1.child(momentumYIdx1), vertInElem1,
                                        -sumNormalPhaseFluxes
                                        / cParams.elemVolVarsCur1[vertInElem1].density());
            }
            else
            {
507
                // set residualStokes[momentumYIdx1] = v_y in stokesnccouplinglocalresidual.hh
508
509
510
511
512
513
514
515
516
517
518
                couplingRes1.accumulate(lfsu1.child(momentumYIdx1), vertInElem1,
                                        globalProblem_.localResidual2().residual(vertInElem2)[massBalanceIdx2]
                                        / cParams.elemVolVarsCur1[vertInElem1].density());
            }
        }


        // COMPONENT Balance
        // Neumann-like conditions
        if (cParams.boundaryTypes1.isCouplingNeumann(transportEqIdx1))
        {
519
520
            DUNE_THROW(Dune::NotImplemented, "The boundary condition isCouplingNeumann(transportEqIdx1) is not implemented \
                                              for the Stokes side for multicomponent systems.");
521
522
523
524
525
526
527
        }
        if (cParams.boundaryTypes2.isCouplingNeumann(contiWEqIdx2))
        {
            // only enter here, if a boundary layer model is used for the computation of the diffusive fluxes
            if (blModel_)
            {
                Scalar diffusiveFlux = bfNormal1.two_norm()
528
                                       * evalBoundaryLayerConcentrationGradient(cParams, vertInElem1)
529
530
531
532
533
534
535
536
                                       * (boundaryVars1.diffusionCoeff(transportCompIdx1)
                                         + boundaryVars1.eddyDiffusivity())
                                       * boundaryVars1.molarDensity()
                                       * FluidSystem::molarMass(transportCompIdx1);

                Scalar advectiveFlux = normalMassFlux1
                                       * cParams.elemVolVarsCur1[vertInElem1].massFraction(transportCompIdx1);

537
                const Scalar massTransferCoeff = evalMassTransferCoefficient(cParams, vertInElem1, vertInElem2);
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
                // TODO: unify this behavior with the one in the non-isothermal LOP
                if (globalProblem_.sdProblem1().isCornerPoint(globalPos1) && massTransferModel_)
                {
                    Scalar diffusiveFluxAtCorner = bfNormal1
                                                   * boundaryVars1.moleFractionGrad(transportCompIdx1)
                                                   * (boundaryVars1.diffusionCoeff(transportCompIdx1)
                                                      + boundaryVars1.eddyDiffusivity())
                                                   * boundaryVars1.molarDensity()
                                                   * FluidSystem::molarMass(transportCompIdx1);

                    couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
                                            -massTransferCoeff*(advectiveFlux - diffusiveFlux) -
                                            (1.-massTransferCoeff)*(advectiveFlux - diffusiveFluxAtCorner));
                }
                else
                {
                    couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
                                            -massTransferCoeff*(advectiveFlux - diffusiveFlux) +
                                            (1.-massTransferCoeff)*globalProblem_.localResidual1().residual(vertInElem1)[transportEqIdx1]);
                }
            }
            else if (globalProblem_.sdProblem1().isCornerPoint(globalPos1))
            {
                static_assert(!GET_PROP_VALUE(TwoPTwoCTypeTag, UseMoles),
                              "This coupling condition is only implemented for mass fraction formulation.");

                Scalar advectiveFlux = normalMassFlux1
                                       * cParams.elemVolVarsCur1[vertInElem1].massFraction(transportCompIdx1);
                Scalar diffusiveFlux = bfNormal1
                                       * boundaryVars1.moleFractionGrad(transportCompIdx1)
                                       * (boundaryVars1.diffusionCoeff(transportCompIdx1)
                                          + boundaryVars1.eddyDiffusivity())
                                       * boundaryVars1.molarDensity()
                                       * FluidSystem::molarMass(transportCompIdx1);

                couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
                                        -(advectiveFlux - diffusiveFlux));
            }
            else
            {
                static_assert(GET_PROP_VALUE(Stokes2cTypeTag, UseMoles) == GET_PROP_VALUE(TwoPTwoCTypeTag, UseMoles),
                              "This coupling condition is not implemented for different formulations (mass/mole) in the subdomains.");

                // the component mass flux from the stokes domain
                couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
                                        globalProblem_.localResidual1().residual(vertInElem1)[transportEqIdx1]);
            }
        }

587
        // Dirichlet-like conditions
588
589
590
591
592
        if (cParams.boundaryTypes1.isCouplingDirichlet(transportEqIdx1))
        {
            static_assert(!GET_PROP_VALUE(Stokes2cTypeTag, UseMoles),
                          "This coupling condition is only implemented for mass fraction formulation.");

593
594
            // set residualStokes[transportEqIdx1] = x in stokesnccouplinglocalresidual.hh
            // coupling residual is added to "real" residual
595
596
597
598
599
            couplingRes1.accumulate(lfsu1.child(transportEqIdx1), vertInElem1,
                                    -cParams.elemVolVarsCur2[vertInElem2].massFraction(nPhaseIdx2, wCompIdx2));
        }
        if (cParams.boundaryTypes2.isCouplingDirichlet(contiWEqIdx2))
        {
600
601
            DUNE_THROW(Dune::NotImplemented, "The boundary condition isCouplingDirichlet(contiWEqIdx2) is not implemented \
                                              for the Darcy side for multicomponent systems.");
602
603
604
        }
    }

605
    /*!
606
     * \brief Evaluation of the coupling from Stokes (1) to Darcy (2).
607
     *
608
609
     * \param lfsu1 local basis for the trial space of the Stokes domain
     * \param lfsu2 local basis for the trial space of the Darcy domain
610
611
612
613
614
615
616
617
618
619
620
     * \param vertInElem1 local vertex index in element1
     * \param vertInElem2 local vertex index in element2
     * \param sdElement1 the element in the Stokes domain
     * \param sdElement2 the element in the Darcy domain
     * \param boundaryVars1 the boundary variables at the interface of the Stokes domain
     * \param boundaryVars2 the boundary variables at the interface of the Darcy domain
     * \param cParams a parameter container
     * \param couplingRes1 the coupling residual from the Stokes domain
     * \param couplingRes2 the coupling residual from the Darcy domain
     */
    template<typename LFSU1, typename LFSU2, typename RES1, typename RES2, typename CParams>
621
    DUNE_DEPRECATED_MSG("evalCoupling12() is deprecated. Use evalCoupling() instead.")
622
    void evalCoupling12(const LFSU1& lfsu1, const LFSU2& lfsu2,
623
624
625
626
627
                        const int vertInElem1, const int vertInElem2,
                        const SDElement1& sdElement1, const SDElement2& sdElement2,
                        const BoundaryVariables1& boundaryVars1, const BoundaryVariables2& boundaryVars2,
                        const CParams &cParams,
                        RES1& couplingRes1, RES2& couplingRes2) const
628
    {
629
630
        const GlobalPosition& globalPos1 = cParams.fvGeometry1.subContVol[vertInElem1].global;
        const GlobalPosition& bfNormal1 = boundaryVars1.face().normal;
Thomas Fetzer's avatar
Thomas Fetzer committed
631

632
633
634
635
        const Scalar normalMassFlux = boundaryVars1.normalVelocity() *
            cParams.elemVolVarsCur1[vertInElem1].density();

        //rho*v*n as NEUMANN condition for porous medium (set, if B&J defined as NEUMANN condition)
636
        if (cParams.boundaryTypes2.isCouplingNeumann(massBalanceIdx2))
637
        {
Thomas Fetzer's avatar
Thomas Fetzer committed
638
639
640
            static_assert(!GET_PROP_VALUE(Stokes2cTypeTag, UseMoles),
                          "This coupling condition is only implemented for mass fraction formulation.");

641
642
            if (globalProblem_.sdProblem1().isCornerPoint(globalPos1))
            {
643
                couplingRes2.accumulate(lfsu2.child(massBalanceIdx2), vertInElem2,
644
645
646
647
                                        -normalMassFlux);
            }
            else
            {
648
                couplingRes2.accumulate(lfsu2.child(massBalanceIdx2), vertInElem2,
649
650
651
                                        globalProblem_.localResidual1().residual(vertInElem1)[massBalanceIdx1]);
            }
        }
652
        if (cParams.boundaryTypes2.isCouplingDirichlet(massBalanceIdx2))
653
        {
654
            couplingRes2.accumulate(lfsu2.child(massBalanceIdx2), vertInElem2,
655
656
657
658
                                    globalProblem_.localResidual1().residual(vertInElem1)[momentumYIdx1]
                                    -cParams.elemVolVarsCur1[vertInElem1].pressure());
        }

659
        if (cParams.boundaryTypes2.isCouplingNeumann(contiWEqIdx2))
660
        {
661
            // only enter here, if a boundary layer model is used for the computation of the diffusive fluxes
662
663
664
            if (blModel_)
            {
                const Scalar diffusiveFlux =
665
                    bfNormal1.two_norm()
666
                    * evalBoundaryLayerConcentrationGradient(cParams, vertInElem1)
667
668
669
670
671
672
                    * (boundaryVars1.diffusionCoeff(transportCompIdx1)
                       + boundaryVars1.eddyDiffusivity())
                    * boundaryVars1.molarDensity()
                    * FluidSystem::molarMass(transportCompIdx1);

                Scalar advectiveFlux = normalMassFlux * cParams.elemVolVarsCur1[vertInElem1].massFraction(transportCompIdx1);
673

674
675
                const Scalar massTransferCoeff = evalMassTransferCoefficient(cParams, vertInElem1, vertInElem2);

676
                if (globalProblem_.sdProblem1().isCornerPoint(globalPos1) && massTransferModel_)
677
                {
678
                    const Scalar diffusiveFluxAtCorner =
679
680
                        bfNormal1 *
                        boundaryVars1.moleFractionGrad(transportCompIdx1) *
681
                        (boundaryVars1.diffusionCoeff(transportCompIdx1) + boundaryVars1.eddyDiffusivity()) *
682
683
684
                        boundaryVars1.molarDensity() *
                        FluidSystem::molarMass(transportCompIdx1);

685
                    couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
686
687
                                            -massTransferCoeff*(advectiveFlux - diffusiveFlux) -
                                            (1.-massTransferCoeff)*(advectiveFlux - diffusiveFluxAtCorner));
688
689
690
                }
                else
                {
691
                    couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
692
693
                                            -massTransferCoeff*(advectiveFlux - diffusiveFlux) +
                                            (1.-massTransferCoeff)*globalProblem_.localResidual1().residual(vertInElem1)[transportEqIdx1]);
694
695
                }
            }
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
            else if (globalProblem_.sdProblem1().isCornerPoint(globalPos1))
            {
                static_assert(!GET_PROP_VALUE(TwoPTwoCTypeTag, UseMoles),
                              "This coupling condition is only implemented for mass fraction formulation.");

                const Scalar advectiveFlux =
                    normalMassFlux *
                    cParams.elemVolVarsCur1[vertInElem1].massFraction(transportCompIdx1);
                const Scalar diffusiveFlux =
                    bfNormal1 *
                    boundaryVars1.moleFractionGrad(transportCompIdx1) *
                    (boundaryVars1.diffusionCoeff(transportCompIdx1) + boundaryVars1.eddyDiffusivity()) *
                    boundaryVars1.molarDensity() *
                    FluidSystem::molarMass(transportCompIdx1);

711
                couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
712
713
714
715
716
                                        -(advectiveFlux - diffusiveFlux));
            }
            else
            {
                static_assert(GET_PROP_VALUE(Stokes2cTypeTag, UseMoles) == GET_PROP_VALUE(TwoPTwoCTypeTag, UseMoles),
717
                              "This coupling condition is not implemented dor different formulations (mass/mole) in the subdomains.");
718
719

                // the component mass flux from the stokes domain
720
                couplingRes2.accumulate(lfsu2.child(contiWEqIdx2), vertInElem2,
721
722
                                        globalProblem_.localResidual1().residual(vertInElem1)[transportEqIdx1]);
            }
723
        }
724
        if (cParams.boundaryTypes2.isCouplingDirichlet(contiWEqIdx2))
725
        {
726
            DUNE_THROW(Dune::NotImplemented, "The boundary condition isCouplingDirichlet(contiWEqIdx2) is not implemented for the Darcy side.");
727
        }
728
729
730
    }

    /*!
731
     * \brief Evaluation of the coupling from Darcy (2) to Stokes (1).
732
     *
733
734
     * \param lfsu1 local basis for the trial space of the Stokes domain
     * \param lfsu2 local basis for the trial space of the Darcy domain
735
736
737
738
739
740
741
742
743
744
745
     * \param vertInElem1 local vertex index in element1
     * \param vertInElem2 local vertex index in element2
     * \param sdElement1 the element in the Stokes domain
     * \param sdElement2 the element in the Darcy domain
     * \param boundaryVars1 the boundary variables at the interface of the Stokes domain
     * \param boundaryVars2 the boundary variables at the interface of the Darcy domain
     * \param cParams a parameter container
     * \param couplingRes1 the coupling residual from the Stokes domain
     * \param couplingRes2 the coupling residual from the Darcy domain
     */
    template<typename LFSU1, typename LFSU2, typename RES1, typename RES2, typename CParams>
746
    DUNE_DEPRECATED_MSG("evalCoupling21() is deprecated. Use evalCoupling() instead.")
747
    void evalCoupling21(const LFSU1& lfsu1, const LFSU2& lfsu2,
748
749
750
751
752
                        const int vertInElem1, const int vertInElem2,
                        const SDElement1& sdElement1, const SDElement2& sdElement2,
                        const BoundaryVariables1& boundaryVars1, const BoundaryVariables2& boundaryVars2,
                        const CParams &cParams,
                        RES1& couplingRes1, RES2& couplingRes2) const
753
    {
754
755
        const GlobalPosition& globalPos2 = cParams.fvGeometry2.subContVol[vertInElem2].global;
        GlobalPosition normalFlux2(0.);
756
757
758
759
760
761
762

        // velocity*normal*area*rho
        for (int phaseIdx=0; phaseIdx<numPhases2; ++phaseIdx)
            normalFlux2[phaseIdx] = -boundaryVars2.volumeFlux(phaseIdx)*
                cParams.elemVolVarsCur2[vertInElem2].density(phaseIdx);

        //p*n as NEUMANN condition for free flow (set, if B&J defined as NEUMANN condition)
763
        if (cParams.boundaryTypes1.isCouplingDirichlet(momentumYIdx1))
764
765
766
        {
            //p*A*n as NEUMANN condition for free flow (set, if B&J defined as NEUMANN condition)
            //pressure correction in stokeslocalresidual.hh
767
            couplingRes1.accumulate(lfsu1.child(momentumYIdx1), vertInElem1,
768
769
770
                                    cParams.elemVolVarsCur2[vertInElem2].pressure(nPhaseIdx2) *
                                    boundaryVars2.face().area);
        }
771
        if (cParams.boundaryTypes1.isCouplingNeumann(momentumYIdx1))
772
        {
773
            // v.n as Dirichlet-like condition for the Stokes domain
774
775
776
            // set residualStokes[momentumYIdx1] = vy in stokeslocalresidual.hh
            if (globalProblem_.sdProblem2().isCornerPoint(globalPos2))
            {
777
                couplingRes1.accumulate(lfsu1.child(momentumYIdx1), vertInElem1,
778
779
780
781
782
783
                                        -((normalFlux2[nPhaseIdx2] + normalFlux2[wPhaseIdx2])
                                          / cParams.elemVolVarsCur1[vertInElem1].density()));
            }
            else
            {
                // v.n as DIRICHLET condition for the Stokes domain (negative sign!)
784
                couplingRes1.accumulate(lfsu1.child(momentumYIdx1), vertInElem1,
785
786
787
788
789
                                        globalProblem_.localResidual2().residual(vertInElem2)[massBalanceIdx2]
                                        / cParams.elemVolVarsCur1[vertInElem1].density());
            }
        }

790
791
        typedef typename GET_PROP_TYPE(Stokes2cTypeTag, SpatialParams) SpatialParams1;
        SpatialParams1 spatialParams = globalProblem_.sdProblem1().spatialParams();
792
793
794
795
796
797
798
799
        const GlobalPosition& globalPos = cParams.fvGeometry1.subContVol[vertInElem1].global;
        Scalar beaversJosephCoeff = spatialParams.beaversJosephCoeffAtPos(globalPos);
        assert(beaversJosephCoeff > 0);

        const Scalar Kxx = spatialParams.intrinsicPermeability(sdElement1, cParams.fvGeometry1,
                                                               vertInElem1);

        beaversJosephCoeff /= std::sqrt(Kxx);
800
        const GlobalPosition& elementUnitNormal = boundaryVars1.face().normal;
801

802
        // Implementation as Neumann-like condition: (v.n)n
803
804
805
        if (cParams.boundaryTypes1.isCouplingDirichlet(momentumXIdx1))
        {
            const Scalar normalComp = boundaryVars1.velocity()*elementUnitNormal;
806
            GlobalPosition normalV = elementUnitNormal;
807
808
809
            normalV *= normalComp; // v*n*n

            // v_tau = v - (v.n)n
810
            const GlobalPosition tangentialV = boundaryVars1.velocity() - normalV;
811
812
813
814
            const Scalar boundaryFaceArea = boundaryVars1.face().area;

            for (int dimIdx=0; dimIdx < dim; ++dimIdx)
            {
815
                couplingRes1.accumulate(lfsu1.child(momentumXIdx1), vertInElem1,
816
817
818
819
820
821
822
                                        beaversJosephCoeff
                                        * boundaryFaceArea
                                        * tangentialV[dimIdx]
                                        * (boundaryVars1.dynamicViscosity()
                                          + boundaryVars1.dynamicEddyViscosity()));
            }
        }
823
        // Implementation as Dirichlet-like condition
824
825
826
        // tangential component: vx = sqrt K /alpha * (grad v n(unity))t
        if (cParams.boundaryTypes1.isCouplingNeumann(momentumXIdx1))
        {
827
            DUNE_THROW(Dune::NotImplemented, "The boundary conditionisCouplingNeumann(momentumXIdx1) on the Stokes side is not implemented anymore.");
828
829
830
            // GlobalPosition tangentialVelGrad(0);
            // boundaryVars1.velocityGrad().umv(elementUnitNormal, tangentialVelGrad);
            // tangentialVelGrad /= -beaversJosephCoeff; // was - before
831
            // couplingRes1.accumulate(lfsu1.child(momentumXIdx1), vertInElem1,
832
833
834
835
            // this->residual_[vertInElem1][momentumXIdx1] =
            //        tangentialVelGrad[momentumXIdx1] - globalProblem_.localResidual1().curPriVars_(vertInElem1)[momentumXIdx1]);
        }

836
837
        //coupling residual is added to "real" residual
        //here each node is passed twice, hence only half of the dirichlet condition has to be set
838
        if (cParams.boundaryTypes1.isCouplingDirichlet(transportEqIdx1))
839
840
        {
            // set residualStokes[transportEqIdx1] = x in stokes2clocalresidual.hh
Thomas Fetzer's avatar
Thomas Fetzer committed
841
842
843
844
845


            static_assert(!GET_PROP_VALUE(Stokes2cTypeTag, UseMoles),
                          "This coupling condition is only implemented for mass fraction formulation.");

846
            couplingRes1.accumulate(lfsu1.child(transportEqIdx1), vertInElem1,
847
848
                                    -cParams.elemVolVarsCur2[vertInElem2].massFraction(nPhaseIdx2, wCompIdx2));
        }
849
        if (cParams.boundaryTypes1.isCouplingNeumann(transportEqIdx1))
850
        {
851
            DUNE_THROW(Dune::NotImplemented, "The boundary condition isCouplingNeumann(transportEqIdx1) is not implemented for the Stokes side.");
852
        }
853
854
    }

855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
    /*!
     * \brief Returns a BoundaryLayerModel object
     *
     * This function is reused in Child LocalOperators and used for extracting
     * the respective boundary layer thickness.<br>
     * \todo This function could be moved to a more model specific place, because
     *       of its runtime parameters.
     *
     * \param cParams a parameter container
     * \param scvIdx1 The local index of the sub-control volume of the Stokes domain
     */
    template<typename CParams>
    BoundaryLayerModel<TypeTag> evalBoundaryLayerModel(CParams cParams, const int scvIdx1) const
    {
        // current position + additional virtual runup distance
        const Scalar distance = cParams.fvGeometry1.subContVol[scvIdx1].global[0]
                                + GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, BoundaryLayer, Offset);
872
873
874
875
        BoundaryLayerModel<TypeTag> boundaryLayerModel(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, FreeFlow, RefVelocity),
                                                       distance,
                                                       cParams.elemVolVarsCur1[scvIdx1].kinematicViscosity(),
                                                       blModel_);
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
        if (blModel_ == 1)
            boundaryLayerModel.setConstThickness(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, BoundaryLayer, ConstThickness));
        if (blModel_ >= 4)
            boundaryLayerModel.setYPlus(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, BoundaryLayer, YPlus));
        if (blModel_ >= 5)
            boundaryLayerModel.setRoughnessLength(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, BoundaryLayer, RoughnessLength));
        if (blModel_ == 7)
            boundaryLayerModel.setHydraulicDiameter(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, BoundaryLayer, HydraulicDiameter));

        return boundaryLayerModel;
    }

    /*!
     * \brief Returns the concentration gradient through the boundary layer
     *
     * \todo This function could be moved to a more model specific place, because
     *       of its runtime parameters.
     *
     * \param cParams a parameter container
     * \param scvIdx1 The local index of the sub-control volume of the Stokes domain
     */
    template<typename CParams>
    Scalar evalBoundaryLayerConcentrationGradient(CParams cParams, const int scvIdx1) const
    {
900
901
        static_assert(numComponents1 == 2,
                      "This coupling condition is only implemented for two components.");
902
903
904
905
906
907
908
909
910
        Scalar massFractionOut = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, FreeFlow, RefMassfrac);
        Scalar M1 = FluidSystem::molarMass(transportCompIdx1);
        Scalar M2 = FluidSystem::molarMass(phaseCompIdx1);
        Scalar X2 = 1.0 - massFractionOut;
        Scalar massToMoleDenominator = M2 + X2*(M1 - M2);
        Scalar moleFractionOut = massFractionOut * M2 /massToMoleDenominator;

        Scalar normalMoleFracGrad = cParams.elemVolVarsCur1[scvIdx1].moleFraction(transportCompIdx1)
                                    - moleFractionOut;
911
        return normalMoleFracGrad / evalBoundaryLayerModel(cParams, scvIdx1).massBoundaryLayerThickness();
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
    }

    /*!
     * \brief Returns the mass transfer coefficient
     *
     * This function is reused in Child LocalOperators.
     * \todo This function could be moved to a more model specific place, because
     *       of its runtime parameters.
     *
     * \param cParams a parameter container
     * \param scvIdx1 The local index of the sub-control volume of the Stokes domain
     * \param scvIdx1 The local index of the sub-control volume of the Darcy domain
     */
    template<typename CParams>
    Scalar evalMassTransferCoefficient(CParams cParams, const int scvIdx1, const int scvIdx2) const
    {
        MassTransferModel<TypeTag> massTransferModel(cParams.elemVolVarsCur2[scvIdx2].saturation(wPhaseIdx2),
                                                     cParams.elemVolVarsCur2[scvIdx2].porosity(),
930
                                                     evalBoundaryLayerModel(cParams, scvIdx1).massBoundaryLayerThickness(),
931
932
933
934
935
936
937
938
                                                     massTransferModel_);
        if (massTransferModel_ == 1)
            massTransferModel.setMassTransferCoeff(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, MassTransfer, Coefficient));
        if (massTransferModel_ == 2 || massTransferModel_ == 4)
            massTransferModel.setCharPoreRadius(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, MassTransfer, CharPoreRadius));
        if (massTransferModel_ == 3)
            massTransferModel.setCapillaryPressure(cParams.elemVolVarsCur2[scvIdx2].capillaryPressure());

939
        return massTransferModel.massTransferCoefficient();
940
941
    }

942
943
944
945
946
947
948
949
950
 protected:
    GlobalProblem& globalProblem() const
    { return globalProblem_; }

    Implementation *asImp_()
    { return static_cast<Implementation *> (this); }
    const Implementation *asImp_() const
    { return static_cast<const Implementation *> (this); }

951
952
953
    unsigned int blModel_;
    unsigned int massTransferModel_;

954
955
 private:
    /*!
956
957
     * \brief A struct that contains data of the FF and PM including boundary types,
     *        volume variables in both subdomains and geometric information
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
     */
    struct CParams
    {
        BoundaryTypes1 boundaryTypes1;
        BoundaryTypes2 boundaryTypes2;
        ElementVolumeVariables1 elemVolVarsPrev1;
        ElementVolumeVariables1 elemVolVarsCur1;
        ElementVolumeVariables2 elemVolVarsPrev2;
        ElementVolumeVariables2 elemVolVarsCur2;
        FVElementGeometry1 fvGeometry1;
        FVElementGeometry2 fvGeometry2;
    };

    GlobalProblem& globalProblem_;
};

} // end namespace Dumux

#endif // DUMUX_2CSTOKES_2P2C_LOCALOPERATOR_HH