2cstokes2p2clocaloperator.hh 30.5 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/version.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
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
 * \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.
 */
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>
{
 public:
    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;

    typedef typename GET_PROP_TYPE(TwoPTwoCTypeTag, FluidSystem) FluidSystem;

    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;

    // Multidomain Grid and Subgrid types
    typedef typename GET_PROP_TYPE(TypeTag, MultiDomainGrid) MDGrid;
    typedef typename MDGrid::Traits::template Codim<0>::Entity MDElement;
    typedef typename MDGrid::SubDomainGrid SDGrid;

    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;

    enum { dim = MDGrid::dimension };

    // FREE FLOW
    enum { numEq1 = GET_PROP_VALUE(Stokes2cTypeTag, NumEq) };
    enum { nPhaseIdx1 = Stokes2cIndices::phaseIdx };           	//!< Index of the free-flow phase of the fluidsystem
    enum { // equation indices in the Stokes domain
        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
    };
    enum { // indices of the components
        transportCompIdx1 = Stokes2cIndices::transportCompIdx, 	//!< Index of transported component
        phaseCompIdx1 = Stokes2cIndices::phaseCompIdx         	//!< Index of main component of the phase
    };

    // POROUS MEDIUM
    enum { numEq2 = GET_PROP_VALUE(TwoPTwoCTypeTag, NumEq) };
    enum { numPhases2 = GET_PROP_VALUE(TwoPTwoCTypeTag, NumPhases) };
    enum { // equation indices in the Darcy domain
        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)
    };
    enum { // component indices
    	wCompIdx2 = TwoPTwoCIndices::wCompIdx,          //!< Index of the liquids main component
        nCompIdx2 = TwoPTwoCIndices::nCompIdx           //!< Index of the main component of the gas
    };
    enum { // phase indices
        wPhaseIdx2 = TwoPTwoCIndices::wPhaseIdx,        //!< Index for the liquid phase
        nPhaseIdx2 = TwoPTwoCIndices::nPhaseIdx          //!< Index for the gas phase
    };

    typedef typename GET_PROP(TypeTag, ParameterTree) ParameterTree;
    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef Dune::FieldVector<Scalar, dim> DimVector;

    typedef typename Stokes2cGridView::template Codim<dim>::EntityPointer VertexPointer1;
    typedef typename TwoPTwoCGridView::template Codim<dim>::EntityPointer VertexPointer2;
    typedef typename MDGrid::Traits::template Codim<0>::EntityPointer MDElementPointer;

    TwoCStokesTwoPTwoCLocalOperator(GlobalProblem& globalProblem)
        : globalProblem_(globalProblem)
    {
Thomas Fetzer's avatar
Thomas Fetzer committed
136
137
        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
138

Thomas Fetzer's avatar
Thomas Fetzer committed
139
        if (blModel_ != 0)
Thomas Fetzer's avatar
Thomas Fetzer committed
140
            std::cout << "Using boundary layer model " << blModel_ << std::endl;
Thomas Fetzer's avatar
Thomas Fetzer committed
141
        if (massTransferModel_ != 0)
Thomas Fetzer's avatar
Thomas Fetzer committed
142
            std::cout << "Using mass transfer model " << massTransferModel_ << std::endl;
143
144
145
146
147
148
149
150
151
152
153
154
    }

    // multidomain flags
    static const bool doAlphaCoupling = true;
    static const bool doPatternCoupling = true;

    /*!
     * \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
155
     * \param lfsu_s local basis for the trial space of the Stokes domain
156
     * \param unknowns1 the unknowns vector of the Stokes element (formatted according to PDELab)
157
158
     * \param lfsv_s local basis for the test space of the Stokes domain
     * \param lfsu_n local basis for the trail space of the Darcy domain
159
     * \param unknowns2 the unknowns vector of the Darcy element (formatted according to PDELab)
160
     * \param lfsv_n local basis for the test space of the Darcy domain
161
162
163
164
165
166
     * \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>
167
168
169
170
    void alpha_coupling (const IntersectionGeom& intersectionGeometry,
                         const LFSU1& lfsu_s, const X& unknowns1, const LFSV1& lfsv_s,
                         const LFSU2& lfsu_n, const X& unknowns2, const LFSV2& lfsv_n,
                         RES& couplingRes1, RES& couplingRes2) const
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
    {
        const MDElementPointer mdElementPointer1 = intersectionGeometry.inside();
        const MDElementPointer mdElementPointer2 = intersectionGeometry.outside();

        const MDElement& mdElement1 = *mdElementPointer1;
        const MDElement& mdElement2 = *mdElementPointer2;

        // the subodmain elements
        const SDElement1& sdElement1 = *(globalProblem_.sdElementPointer1(mdElement1));
        const SDElement2& sdElement2 = *(globalProblem_.sdElementPointer2(mdElement2));

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

        // update fvElementGeometry and the element volume variables
        updateElemVolVars(lfsu_s, lfsu_n,
                          unknowns1, unknowns2,
                          sdElement1, sdElement2,
                          cParams);

        // first element
        const int faceIdx1 = intersectionGeometry.indexInInside();
193
194
        const Dune::ReferenceElement<typename MDGrid::ctype,dim>& referenceElement1 =
            Dune::ReferenceElements<typename MDGrid::ctype,dim>::general(mdElement1.type());
195
196
197
198
        const int numVerticesOfFace = referenceElement1.size(faceIdx1, 1, dim);

        // second element
        const int faceIdx2 = intersectionGeometry.indexInOutside();
199
200
        const Dune::ReferenceElement<typename MDGrid::ctype,dim>& referenceElement2 =
            Dune::ReferenceElements<typename MDGrid::ctype,dim>::general(mdElement2.type());
201
202
203
204
205
206
207
208
209
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

        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);

            globalProblem_.sdProblem1().boundaryTypes(cParams.boundaryTypes1, *vPtr1);
            globalProblem_.sdProblem2().boundaryTypes(cParams.boundaryTypes2, *vPtr2);

            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);

            asImp_()->evalCoupling12(lfsu_s, lfsu_n, // local function spaces
                                     vertInElem1, vertInElem2,
                                     sdElement1, sdElement2,
                                     boundaryVars1, boundaryVars2,
                                     cParams,
                                     couplingRes1, couplingRes2);
            asImp_()->evalCoupling21(lfsu_s, lfsu_n, // local function spaces
                                     vertInElem1, vertInElem2,
                                     sdElement1, sdElement2,
                                     boundaryVars1, boundaryVars2,
                                     cParams,
                                     couplingRes1, couplingRes2);
        }
    }

    /*!
     * \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.
     *
249
250
     * \param lfsu_s local basis for the trial space of the Stokes domain
     * \param lfsu_n local basis for the trial space of the Darcy domain
251
252
253
254
255
256
257
258
     * \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>
259
260
261
262
    void updateElemVolVars (const LFSU1& lfsu_s, const LFSU2& lfsu_n,
                            const X& unknowns1, const X& unknowns2,
                            const SDElement1& sdElement1, const SDElement2& sdElement2,
                            CParams &cParams) const
263
264
265
266
    {
        cParams.fvGeometry1.update(globalProblem_.sdGridView1(), sdElement1);
        cParams.fvGeometry2.update(globalProblem_.sdGridView2(), sdElement2);

267
268
269
270
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
        const int numVertsOfElem1 = sdElement1.subEntities(dim);
        const int numVertsOfElem2 = sdElement2.subEntities(dim);
#else
271
272
        const int numVertsOfElem1 = sdElement1.template count<dim>();
        const int numVertsOfElem2 = sdElement2.template count<dim>();
273
#endif
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

        //bring the local unknowns x_s into a form that can be passed to elemVolVarsCur.update()
        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)
                elementSol1[eqIdx1*numVertsOfElem1+idx] = unknowns1(lfsu_s.child(eqIdx1),idx);
            for (int eqIdx2=0; eqIdx2<numEq2; ++eqIdx2)
                elementSol2[eqIdx2*numVertsOfElem2+idx] = unknowns2(lfsu_n.child(eqIdx2),idx);
        }
#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);

    }

    /*!
     * \brief Evaluation of the coupling from Stokes (1 or s) to Darcy (2 or n).
     *
307
308
     * \param lfsu_s local basis for the trial space of the Stokes domain
     * \param lfsu_n local basis for the trial space of the Darcy domain
309
310
311
312
313
314
315
316
317
318
319
     * \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>
320
321
322
323
324
325
    void evalCoupling12(const LFSU1& lfsu_s, const LFSU2& lfsu_n,
                        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
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
    {
        const DimVector& globalPos1 = cParams.fvGeometry1.subContVol[vertInElem1].global;
        const DimVector& bfNormal1 = boundaryVars1.face().normal;
        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)
        if (cParams.boundaryTypes2.isCouplingInflow(massBalanceIdx2))
        {
            if (globalProblem_.sdProblem1().isCornerPoint(globalPos1))
            {
                couplingRes2.accumulate(lfsu_n.child(massBalanceIdx2), vertInElem2,
                                        -normalMassFlux);
            }
            else
            {
                couplingRes2.accumulate(lfsu_n.child(massBalanceIdx2), vertInElem2,
                                        globalProblem_.localResidual1().residual(vertInElem1)[massBalanceIdx1]);
            }
        }
        if (cParams.boundaryTypes2.isCouplingOutflow(massBalanceIdx2))
        {
            couplingRes2.accumulate(lfsu_n.child(massBalanceIdx2), vertInElem2,
                                    globalProblem_.localResidual1().residual(vertInElem1)[momentumYIdx1]
                                    -cParams.elemVolVarsCur1[vertInElem1].pressure());
        }

        if (cParams.boundaryTypes2.isCouplingInflow(contiWEqIdx2))
        {
            // only enter here, if a BOUNDARY LAYER MODEL is used for the computation of the diffusive fluxes
            if (blModel_)
            {
                const Scalar massFractionOut = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, FreeFlow, RefMassfrac);
                const Scalar M1 = FluidSystem::molarMass(transportCompIdx1);
                const Scalar M2 = FluidSystem::molarMass(phaseCompIdx1);
                const Scalar X2 = 1.0 - massFractionOut;
                const Scalar massToMoleDenominator = M2 + X2*(M1 - M2);
                const Scalar moleFractionOut = massFractionOut * M2 /massToMoleDenominator;

                const Scalar advectiveFlux =
                    normalMassFlux *
                    cParams.elemVolVarsCur1[vertInElem1].massFraction(transportCompIdx1);
                Scalar normalMoleFracGrad =
                    cParams.elemVolVarsCur1[vertInElem1].moleFraction(transportCompIdx1) -
                    moleFractionOut;

Thomas Fetzer's avatar
Thomas Fetzer committed
372
373
374
375
376
377
378
379
380
381
382
383
384
                const Scalar velocity = GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, FreeFlow, RefVelocity);
                // current position + additional virtual runup distance
                const Scalar distance = globalPos1[0] + GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, BoundaryLayer, Offset);
                const Scalar kinematicViscosity = cParams.elemVolVarsCur1[vertInElem2].kinematicViscosity();
                BoundaryLayerModel<TypeTag> boundaryLayerModel(velocity, distance, kinematicViscosity, blModel_);
                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));

                normalMoleFracGrad /= boundaryLayerModel.massBoundaryLayerThickness();
385
386
387
388

                const Scalar diffusiveFlux =
                    bfNormal1.two_norm() *
                    normalMoleFracGrad *
389
                    (boundaryVars1.diffusionCoeff(transportCompIdx1) + boundaryVars1.eddyDiffusivity()) *
390
391
392
393
394
395
396
397
398
399
400
401
                    boundaryVars1.molarDensity() *
                    FluidSystem::molarMass(transportCompIdx1);

                if (massTransferModel_ == 0)
                {
                    couplingRes2.accumulate(lfsu_n.child(contiWEqIdx2), vertInElem2,
                                            -(advectiveFlux - diffusiveFlux));
                }
                // transition from the mass transfer coefficient concept to the coupling via the local residual,
                // when saturations become small; only diffusive fluxes are scaled!
                else
                {
Thomas Fetzer's avatar
Thomas Fetzer committed
402
403
404
405
406
407
                    MassTransferModel<TypeTag> massTransferModel(cParams.elemVolVarsCur2[vertInElem2].saturation(wPhaseIdx2),
                                                    cParams.elemVolVarsCur2[vertInElem2].porosity(),
                                                    boundaryLayerModel.massBoundaryLayerThickness(),
                                                    massTransferModel_);
                    if (massTransferModel_ == 1)
                        massTransferModel.setMassTransferCoeff(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, MassTransfer, Coefficient));
Thomas Fetzer's avatar
Thomas Fetzer committed
408
409
                    if (massTransferModel_ == 2 || massTransferModel_ == 4)
                        massTransferModel.setCharPoreRadius(GET_RUNTIME_PARAM_FROM_GROUP(TypeTag, Scalar, MassTransfer, CharPoreRadius));
Thomas Fetzer's avatar
Thomas Fetzer committed
410
411
412
413
                    if (massTransferModel_ == 3)
                        massTransferModel.setCapillaryPressure(cParams.elemVolVarsCur2[vertInElem2].capillaryPressure());
                    const Scalar massTransferCoeff = massTransferModel.massTransferCoefficient();

414
                    if (massTransferCoeff > 1.0 || massTransferCoeff < 0.0)
415
                        std::cout << "MTC out of bounds, should be in between 0.0 and 1.0! >>> " << massTransferCoeff << std::endl;
416
417
418
419
420
421

                    if (globalProblem_.sdProblem1().isCornerPoint(globalPos1))
                    {
                        const Scalar diffusiveFluxAtCorner =
                            bfNormal1 *
                            boundaryVars1.moleFractionGrad(transportCompIdx1) *
422
                            (boundaryVars1.diffusionCoeff(transportCompIdx1) + boundaryVars1.eddyDiffusivity()) *
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
                            boundaryVars1.molarDensity() *
                            FluidSystem::molarMass(transportCompIdx1);

                        couplingRes2.accumulate(lfsu_n.child(contiWEqIdx2), vertInElem2,
                                                -massTransferCoeff*(advectiveFlux - diffusiveFlux) -
                                                (1.-massTransferCoeff)*(advectiveFlux - diffusiveFluxAtCorner));
                    }
                    else
                    {
                        couplingRes2.accumulate(lfsu_n.child(contiWEqIdx2), vertInElem2,
                                                -massTransferCoeff*(advectiveFlux - diffusiveFlux) +
                                                (1.-massTransferCoeff)*globalProblem_.localResidual1().residual(vertInElem1)[transportEqIdx1]);
                    }
                }
            }
            else
            {
                // compute fluxes explicitly at corner points - only quarter control volume
                if (globalProblem_.sdProblem1().isCornerPoint(globalPos1))
                {
                    const Scalar advectiveFlux =
                        normalMassFlux *
                        cParams.elemVolVarsCur1[vertInElem1].massFraction(transportCompIdx1);
                    const Scalar diffusiveFlux =
                        bfNormal1 *
                        boundaryVars1.moleFractionGrad(transportCompIdx1) *
449
                        (boundaryVars1.diffusionCoeff(transportCompIdx1) + boundaryVars1.eddyDiffusivity()) *
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
                        boundaryVars1.molarDensity() *
                        FluidSystem::molarMass(transportCompIdx1);

                    couplingRes2.accumulate(lfsu_n.child(contiWEqIdx2), vertInElem2,
                                            -(advectiveFlux - diffusiveFlux));
                }
                // coupling via the defect
                else
                {
                    // the component mass flux from the stokes domain
                    couplingRes2.accumulate(lfsu_n.child(contiWEqIdx2), vertInElem2,
                                            globalProblem_.localResidual1().residual(vertInElem1)[transportEqIdx1]);
                }
            }
        }
        if (cParams.boundaryTypes2.isCouplingOutflow(contiWEqIdx2))
466
        	std::cerr << "Upwind PM -> FF does not work for the transport equation for a 2-phase system!" << std::endl;
467
468
469
470
471
    }

    /*!
     * \brief Evaluation of the coupling from Darcy (2 or n) to Stokes (1 or s).
     *
472
473
     * \param lfsu_s local basis for the trial space of the Stokes domain
     * \param lfsu_n local basis for the trial space of the Darcy domain
474
475
476
477
478
479
480
481
482
483
484
     * \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>
485
486
487
488
489
490
    void evalCoupling21(const LFSU1& lfsu_s, const LFSU2& lfsu_n,
                        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
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
    {
        const DimVector& globalPos2 = cParams.fvGeometry2.subContVol[vertInElem2].global;
        DimVector normalFlux2(0.);

        // 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)
        if (cParams.boundaryTypes1.isCouplingOutflow(momentumYIdx1))
        {
            //p*A*n as NEUMANN condition for free flow (set, if B&J defined as NEUMANN condition)
            //pressure correction in stokeslocalresidual.hh
            couplingRes1.accumulate(lfsu_s.child(momentumYIdx1), vertInElem1,
                                    cParams.elemVolVarsCur2[vertInElem2].pressure(nPhaseIdx2) *
                                    boundaryVars2.face().area);
        }
        if (cParams.boundaryTypes1.isCouplingInflow(momentumYIdx1))
        {

            // v.n as Dirichlet condition for the Stokes domain
            // set residualStokes[momentumYIdx1] = vy in stokeslocalresidual.hh
            if (globalProblem_.sdProblem2().isCornerPoint(globalPos2))
            {
                couplingRes1.accumulate(lfsu_s.child(momentumYIdx1), vertInElem1,
                                        -((normalFlux2[nPhaseIdx2] + normalFlux2[wPhaseIdx2])
                                          / cParams.elemVolVarsCur1[vertInElem1].density()));
            }
            else
            {
                // v.n as DIRICHLET condition for the Stokes domain (negative sign!)
                couplingRes1.accumulate(lfsu_s.child(momentumYIdx1), vertInElem1,
                                        globalProblem_.localResidual2().residual(vertInElem2)[massBalanceIdx2]
                                        / cParams.elemVolVarsCur1[vertInElem1].density());
            }
        }

        //coupling residual is added to "real" residual
        //here each node is passed twice, hence only half of the dirichlet condition has to be set
        if (cParams.boundaryTypes1.isCouplingOutflow(transportEqIdx1))
        {
            // set residualStokes[transportEqIdx1] = x in stokes2clocalresidual.hh
            couplingRes1.accumulate(lfsu_s.child(transportEqIdx1), vertInElem1,
                                    -cParams.elemVolVarsCur2[vertInElem2].massFraction(nPhaseIdx2, wCompIdx2));
        }
        if (cParams.boundaryTypes1.isCouplingInflow(transportEqIdx1))
538
        	std::cerr << "Upwind PM -> FF does not work for the transport equation for a 2-phase system!" << std::endl;
539
540
541
542
543
544
545
546
547
548
549
550
551
    }

 protected:
    GlobalProblem& globalProblem() const
    { return globalProblem_; }

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

 private:
    /*!
552
553
     * \brief A struct that contains data of the FF and PM including boundary types,
     *        volume variables in both subdomains and geometric information
554
555
556
557
558
559
560
561
562
563
564
565
566
     */
    struct CParams
    {
        BoundaryTypes1 boundaryTypes1;
        BoundaryTypes2 boundaryTypes2;
        ElementVolumeVariables1 elemVolVarsPrev1;
        ElementVolumeVariables1 elemVolVarsCur1;
        ElementVolumeVariables2 elemVolVarsPrev2;
        ElementVolumeVariables2 elemVolVarsCur2;
        FVElementGeometry1 fvGeometry1;
        FVElementGeometry2 fvGeometry2;
    };

Thomas Fetzer's avatar
Thomas Fetzer committed
567
568
    unsigned int blModel_;
    unsigned int massTransferModel_;
569
570
571
572
573
574
575

    GlobalProblem& globalProblem_;
};

} // end namespace Dumux

#endif // DUMUX_2CSTOKES_2P2C_LOCALOPERATOR_HH