2pncvolumevariables.hh 19.4 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
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
// -*- 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 Contains the quantities which are constant within a
 *        finite volume in the two-phase, n-component model.
 */
#ifndef DUMUX_2PNC_VOLUME_VARIABLES_HH
#define DUMUX_2PNC_VOLUME_VARIABLES_HH

#include <iostream>
#include <vector>

#include <dumux/implicit/common/implicitmodel.hh>
#include <dumux/material/fluidstates/compositionalfluidstate.hh>
#include <dumux/common/math.hh>

#include "2pncproperties.hh"
#include "2pncindices.hh"
#include <dumux/material/constraintsolvers/computefromreferencephase2pnc.hh>
#include <dumux/material/constraintsolvers/miscible2pnccomposition.hh>

namespace Dumux
{

/*!
 * \ingroup TwoPNCModel
 * \ingroup ImplicitVolumeVariables
 * \brief Contains the quantities which are are constant within a
 *        finite volume in the two-phase, n-component model.
 */
template <class TypeTag>
class TwoPNCVolumeVariables : public ImplicitVolumeVariables<TypeTag>
{
    typedef ImplicitVolumeVariables<TypeTag> ParentType;
    typedef typename GET_PROP_TYPE(TypeTag, VolumeVariables) Implementation;

    typedef typename GET_PROP_TYPE(TypeTag, Scalar) Scalar;
    typedef typename GET_PROP_TYPE(TypeTag, Grid) Grid;
    typedef typename GET_PROP_TYPE(TypeTag, GridView) GridView;
    typedef typename GET_PROP_TYPE(TypeTag, Problem) Problem;
    typedef typename GET_PROP_TYPE(TypeTag, FVElementGeometry) FVElementGeometry;
    typedef typename GET_PROP_TYPE(TypeTag, PrimaryVariables) PrimaryVariables;
    typedef typename GET_PROP_TYPE(TypeTag, FluidSystem) FluidSystem;
    typedef typename GET_PROP_TYPE(TypeTag, MaterialLaw) MaterialLaw;
    typedef typename GET_PROP_TYPE(TypeTag, MaterialLawParams) MaterialLawParams;
    typedef typename GET_PROP_TYPE(TypeTag, Indices) Indices;
    enum 
    {
        dim = GridView::dimension,
        dimWorld=GridView::dimensionworld,

        numPhases = GET_PROP_VALUE(TypeTag, NumPhases),
        numComponents = GET_PROP_VALUE(TypeTag, NumComponents),
        numMajorComponents = GET_PROP_VALUE(TypeTag, NumMajorComponents),

        // formulations
        formulation = GET_PROP_VALUE(TypeTag, Formulation),
        plSg = TwoPNCFormulation::plSg,
        pgSl = TwoPNCFormulation::pgSl,

        // phase indices
        wPhaseIdx = FluidSystem::wPhaseIdx,
        nPhaseIdx = FluidSystem::nPhaseIdx,

        // component indices
        wCompIdx = FluidSystem::wCompIdx,
        nCompIdx = FluidSystem::nCompIdx,

        // phase presence enums
        nPhaseOnly = Indices::nPhaseOnly,
        wPhaseOnly = Indices::wPhaseOnly,
        bothPhases = Indices::bothPhases,

        // primary variable indices
        pressureIdx = Indices::pressureIdx,
        switchIdx = Indices::switchIdx,

    };

    typedef typename GridView::template Codim<0>::Entity Element;
    typedef typename Grid::ctype CoordScalar;
100
101
    typedef Dumux::miscible2pncComposition<Scalar, FluidSystem> miscible2pncComposition;
    typedef Dumux::computeFromReferencePhase2pnc<Scalar, FluidSystem> computeFromReferencePhase2pnc;
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
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
185
186
187
188
189
190
191
192
193
194
195
196
197

    enum { isBox = GET_PROP_VALUE(TypeTag, ImplicitIsBox) };
    enum { dofCodim = isBox ? dim : 0 };
public:
  
      //! The type of the object returned by the fluidState() method
      typedef CompositionalFluidState<Scalar, FluidSystem> FluidState;
    /*!
     * \copydoc ImplicitVolumeVariables::update
     */
    void update(const PrimaryVariables &primaryVariables,
                const Problem &problem,
                const Element &element,
                const FVElementGeometry &fvGeometry,
                int scvIdx,
                bool isOldSol)
    {
        ParentType::update(primaryVariables,
                           problem,
                           element,
                           fvGeometry,
                           scvIdx,
                           isOldSol);
        
        completeFluidState(primaryVariables, problem, element, fvGeometry, scvIdx, fluidState_, isOldSol);
            
        /////////////
        // calculate the remaining quantities
        /////////////
        const MaterialLawParams &materialParams = problem.spatialParams().materialLawParams(element, fvGeometry, scvIdx);
        
    // Second instance of a parameter cache.
        // Could be avoided if diffusion coefficients also
        // became part of the fluid state.
        typename FluidSystem::ParameterCache paramCache;
        paramCache.updateAll(fluidState_);
        for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
            {// relative permeabilities
                    Scalar kr;
                    if (phaseIdx == wPhaseIdx)
                        kr = MaterialLaw::krw(materialParams, saturation(wPhaseIdx));
                    else // ATTENTION: krn requires the liquid saturation
                        // as parameter!
                        kr = MaterialLaw::krn(materialParams, saturation(wPhaseIdx));
                        mobility_[phaseIdx] = kr / fluidState_.viscosity(phaseIdx);
                        Valgrind::CheckDefined(mobility_[phaseIdx]);
                    int compIIdx = phaseIdx;
                    for(int compIdx = 0; compIdx < numComponents; ++compIdx)
                    {
                    int compJIdx = compIdx;
                    // binary diffusion coefficents
                    diffCoeff_[phaseIdx][compIdx] = 0.0;
                    if(compIIdx!= compJIdx)
                    diffCoeff_[phaseIdx][compIdx] = FluidSystem::binaryDiffusionCoefficient(fluidState_,
                                                                                    paramCache,
                                                                                    phaseIdx,
                                                                                    compIIdx,
                                                                                    compJIdx);
                    Valgrind::CheckDefined(diffCoeff_[phaseIdx][compIdx]);
                }
            }

    // porosity 
    porosity_ = problem.spatialParams().porosity(element,
                                                        fvGeometry,
                                                        scvIdx);
    Valgrind::CheckDefined(porosity_);
    // energy related quantities not contained in the fluid state
            
    asImp_().updateEnergy_(primaryVariables, problem,element, fvGeometry, scvIdx, isOldSol);
    }

   /*!
    * \copydoc ImplicitModel::completeFluidState
    * \param isOldSol Specifies whether this is the previous solution or the current one
    */
    static void completeFluidState(const PrimaryVariables& primaryVariables,
                    const Problem& problem,
                    const Element& element,
                    const FVElementGeometry& fvGeometry,
                    int scvIdx,
                    FluidState& fluidState,
                    bool isOldSol = false)

    {
        Scalar t = Implementation::temperature_(primaryVariables, problem, element,
                                                fvGeometry, scvIdx);
        fluidState.setTemperature(t);
        
#if DUNE_VERSION_NEWER(DUNE_COMMON, 2, 4)
        int dofIdxGlobal = problem.model().dofMapper().subIndex(element, scvIdx, dofCodim);
#else
        int dofIdxGlobal = problem.model().dofMapper().map(element, scvIdx, dofCodim);
#endif        
        int phasePresence = problem.model().phasePresence(dofIdxGlobal, isOldSol);

198
        /////////////
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
        // set the saturations
        /////////////

    Scalar Sg;
        if (phasePresence == nPhaseOnly)
            Sg = 1.0;
        else if (phasePresence == wPhaseOnly) {
            Sg = 0.0;
        }
        else if (phasePresence == bothPhases) {
            if (formulation == plSg)
                Sg = primaryVariables[switchIdx];
            else if (formulation == pgSl)
                Sg = 1.0 - primaryVariables[switchIdx];
            else DUNE_THROW(Dune::InvalidStateException, "Formulation: " << formulation << " is invalid.");
        }   
    else DUNE_THROW(Dune::InvalidStateException, "phasePresence: " << phasePresence << " is invalid.");
        fluidState.setSaturation(nPhaseIdx, Sg);
        fluidState.setSaturation(wPhaseIdx, 1.0 - Sg);

219
        /////////////
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
        // set the pressures of the fluid phases
        /////////////

        // calculate capillary pressure
        const MaterialLawParams &materialParams
        = problem.spatialParams().materialLawParams(element, fvGeometry, scvIdx);
        Scalar pc = MaterialLaw::pc(materialParams, 1 - Sg);

        // extract the pressures
        if (formulation == plSg) {
            fluidState.setPressure(wPhaseIdx, primaryVariables[pressureIdx]);
            if (primaryVariables[pressureIdx] + pc < 0.0)
                 DUNE_THROW(Dumux::NumericalProblem,"Capillary pressure is too low");
            fluidState.setPressure(nPhaseIdx, primaryVariables[pressureIdx] + pc);
        }
        else if (formulation == pgSl) {
            fluidState.setPressure(nPhaseIdx, primaryVariables[pressureIdx]);
237
            // Here we check for (p_g - pc) in order to ensure that (p_l > 0)
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
            if (primaryVariables[pressureIdx] - pc < 0.0)
            {
                std::cout<< "p_g: "<< primaryVariables[pressureIdx]<<" Cap_press: "<< pc << std::endl;
                DUNE_THROW(Dumux::NumericalProblem,"Capillary pressure is too high");
            }
            fluidState.setPressure(wPhaseIdx, primaryVariables[pressureIdx] - pc);
        }
        else DUNE_THROW(Dune::InvalidStateException, "Formulation: " << formulation << " is invalid.");

        /////////////
        // calculate the phase compositions
        /////////////

    typename FluidSystem::ParameterCache paramCache;

        // now comes the tricky part: calculate phase composition
        if (phasePresence == bothPhases) {
            // both phases are present, phase composition results from
            // the gas <-> liquid equilibrium. This is
            // the job of the "MiscibleMultiPhaseComposition"
            // constraint solver

            // set the known mole fractions in the fluidState so that they
261
            // can be used by the miscible2pncComposition constraint solver
262
263
264
265
266
            for (int compIdx=numMajorComponents; compIdx<numComponents; ++compIdx)
            {
                fluidState.setMoleFraction(wPhaseIdx, compIdx, primaryVariables[compIdx]);
            }

267
            miscible2pncComposition::solve(fluidState,
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
                        paramCache,
                        wPhaseIdx,	//known phaseIdx
                        /*setViscosity=*/true,
                        /*setInternalEnergy=*/false);
        }
        else if (phasePresence == nPhaseOnly){

            Dune::FieldVector<Scalar, numComponents> moleFrac;


            moleFrac[wCompIdx] =  primaryVariables[switchIdx];

            for (int compIdx=numMajorComponents; compIdx<numComponents; ++compIdx)
                    moleFrac[compIdx] = primaryVariables[compIdx];


            Scalar sumMoleFracNotGas = 0;
            for (int compIdx=numMajorComponents; compIdx<numComponents; ++compIdx)
                    sumMoleFracNotGas+=moleFrac[compIdx];

            sumMoleFracNotGas += moleFrac[wCompIdx];
            moleFrac[nCompIdx] = 1 - sumMoleFracNotGas;


            // Set fluid state mole fractions
            for (int compIdx=0; compIdx<numComponents; ++compIdx)
                    fluidState.setMoleFraction(nPhaseIdx, compIdx, moleFrac[compIdx]);


            // calculate the composition of the remaining phases (as
            // well as the densities of all phases). this is the job
299
300
            // of the "computeFromReferencePhase2pnc" constraint solver
                computeFromReferencePhase2pnc::solve(fluidState,
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
                                                paramCache,
                                                nPhaseIdx,
                                                /*setViscosity=*/true,
                                                /*setInternalEnergy=*/false);

            }
        else if (phasePresence == wPhaseOnly){
        // only the liquid phase is present, i.e. liquid phase
        // composition is stored explicitly.
        // extract _mass_ fractions in the gas phase
            Dune::FieldVector<Scalar, numComponents> moleFrac;
            
            for (int compIdx=numMajorComponents; compIdx<numComponents; ++compIdx)
            {
                moleFrac[compIdx] = primaryVariables[compIdx];
            }
            moleFrac[nCompIdx] = primaryVariables[switchIdx];
            Scalar sumMoleFracNotWater = 0;
            for (int compIdx=numMajorComponents; compIdx<numComponents; ++compIdx)
            {
                    sumMoleFracNotWater+=moleFrac[compIdx];
            }
            sumMoleFracNotWater += moleFrac[nCompIdx];
            moleFrac[wCompIdx] = 1 -sumMoleFracNotWater;


            // convert mass to mole fractions and set the fluid state
            for (int compIdx=0; compIdx<numComponents; ++compIdx)
            {
                fluidState.setMoleFraction(wPhaseIdx, compIdx, moleFrac[compIdx]);
            }
            
            // calculate the composition of the remaining phases (as
            // well as the densities of all phases). this is the job
335
336
            // of the "computeFromReferencePhase2pnc" constraint solver
            computeFromReferencePhase2pnc::solve(fluidState,
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
372
373
374
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
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
                                                paramCache,
                                                wPhaseIdx,
                                                /*setViscosity=*/true,
                                                /*setInternalEnergy=*/false);
        }
        paramCache.updateAll(fluidState);
        for (int phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
        {
            Scalar rho = FluidSystem::density(fluidState, paramCache, phaseIdx);
            Scalar mu = FluidSystem::viscosity(fluidState, paramCache, phaseIdx);

            fluidState.setDensity(phaseIdx, rho);
            fluidState.setViscosity(phaseIdx, mu);
        }
    }
    
    /*!
     * \brief Returns the phase state for the control-volume.
     */
    const FluidState &fluidState() const
    { return fluidState_; }

    /*!
     * \brief Returns the saturation of a given phase within
     *        the control volume in \f$[-]\f$.
     *
     * \param phaseIdx The phase index
     */
    Scalar saturation(int phaseIdx) const
    { return fluidState_.saturation(phaseIdx); }

    /*!
     * \brief Returns the mass density of a given phase within the
     *        control volume.
     *
     * \param phaseIdx The phase index
     */
    Scalar density(int phaseIdx) const
    {
        if (phaseIdx < numPhases)
            return fluidState_.density(phaseIdx);
        
        else
            DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
    }

    /*!
     * \brief Returns the mass density of a given phase within the
     *        control volume.
     *
     * \param phaseIdx The phase index
     */
    Scalar molarDensity(int phaseIdx) const
    {
        if (phaseIdx < numPhases)
            return fluidState_.molarDensity(phaseIdx);

        else
            DUNE_THROW(Dune::InvalidStateException, "Invalid phase index " << phaseIdx);
    }

    /*!
     * \brief Returns the effective pressure of a given phase within
     *        the control volume.
     *
     * \param phaseIdx The phase index
     */
    Scalar pressure(int phaseIdx) const
    { 
        return fluidState_.pressure(phaseIdx); 
    }

    /*!
     * \brief Returns temperature inside the sub-control volume.
     *
     * Note that we assume thermodynamic equilibrium, i.e. the
     * temperature of the rock matrix and of all fluid phases are
     * identical.
     */
    Scalar temperature() const
    { return fluidState_.temperature(/*phaseIdx=*/0); }

    /*!
     * \brief Returns the effective mobility of a given phase within
     *        the control volume.
     *
     * \param phaseIdx The phase index
     */
    Scalar mobility(int phaseIdx) const
    {
        return mobility_[phaseIdx];
    }

    /*!
     * \brief Returns the effective capillary pressure within the control volume
     *        in \f$[kg/(m*s^2)=N/m^2=Pa]\f$.
     */
    Scalar capillaryPressure() const
    { return fluidState_.pressure(FluidSystem::nPhaseIdx) - fluidState_.pressure(FluidSystem::wPhaseIdx); }

    /*!
     * \brief Returns the average porosity within the control volume.
     */
    Scalar porosity() const
    { return porosity_; }


    /*!
     * \brief Returns the binary diffusion coefficients for a phase in \f$[m^2/s]\f$.
     */
    Scalar diffCoeff(int phaseIdx, int compIdx) const
    { return diffCoeff_[phaseIdx][compIdx]; }

protected:
  
    static Scalar temperature_(const PrimaryVariables &priVars,
                               const Problem& problem,
                               const Element &element,
                               const FVElementGeometry &fvGeometry,
                               int scvIdx)
    {
        return problem.temperatureAtPos(fvGeometry.subContVol[scvIdx].global);
    }
    
    template<class ParameterCache>
    static Scalar enthalpy_(const FluidState& fluidState,
                            const ParameterCache& paramCache,
                            int phaseIdx)
    {
        return 0;
    }

    /*!
        * \brief Called by update() to compute the energy related quantities
        */
    void updateEnergy_(const PrimaryVariables &priVars,
                        const Problem &problem,
                        const Element &element,
                        const FVElementGeometry &fvGeometry,
                        const int scvIdx,
                        bool isOldSol)
    { };

    Scalar porosity_;        //!< Effective porosity within the control volume
    Scalar mobility_[numPhases];  //!< Effective mobility within the control volume
    Scalar density_;
    FluidState fluidState_;
    Scalar theta_;
    Scalar InitialPorosity_;
    Scalar molWtPhase_[numPhases];
    Dune::FieldMatrix<Scalar, numPhases, numComponents> diffCoeff_;

private:
    Implementation &asImp_()
    { return *static_cast<Implementation*>(this); }

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

} // end namespace

#endif