1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 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 198 199 200 201 202 203 204
  
     | 
    
      // -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
// vi: set et ts=4 sw=4 sts=4:
/*
  This file is part of the Open Porous Media project (OPM).
  OPM is free software: you can redistribute it and/or modify
  it under the terms of the GNU General Public License as published by
  the Free Software Foundation, either version 2 of the License, or
  (at your option) any later version.
  OPM is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.
  You should have received a copy of the GNU General Public License
  along with OPM.  If not, see <http://www.gnu.org/licenses/>.
  Consult the COPYING file in the top-level source directory of this
  module for the precise wording of the license and the list of
  copyright holders.
*/
/*!
 * \file
 *
 * \copydoc Opm::PvsBoundaryRateVector
 */
#ifndef EWOMS_PVS_BOUNDARY_RATE_VECTOR_HH
#define EWOMS_PVS_BOUNDARY_RATE_VECTOR_HH
#include "pvsproperties.hh"
#include <opm/models/common/energymodule.hh>
#include <opm/material/common/Valgrind.hpp>
namespace Opm {
/*!
 * \ingroup PvsModel
 *
 * \brief Implements a rate vector on the boundary for the fully
 *        implicit compositional multi-phase primary variable
 *        switching compositional model.
 */
template <class TypeTag>
class PvsBoundaryRateVector : public GetPropType<TypeTag, Properties::RateVector>
{
    using ParentType = GetPropType<TypeTag, Properties::RateVector>;
    using ExtensiveQuantities = GetPropType<TypeTag, Properties::ExtensiveQuantities>;
    using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
    using Scalar = GetPropType<TypeTag, Properties::Scalar>;
    using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
    using Indices = GetPropType<TypeTag, Properties::Indices>;
    enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
    enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
    enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
    enum { conti0EqIdx = Indices::conti0EqIdx };
    enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
    using EnergyModule = Opm::EnergyModule<TypeTag, enableEnergy>;
    using Toolbox = Opm::MathToolbox<Evaluation>;
public:
    PvsBoundaryRateVector()
        : ParentType()
    {}
    /*!
     * \copydoc
     * ImmiscibleBoundaryRateVector::ImmiscibleBoundaryRateVector(Scalar)
     */
    PvsBoundaryRateVector(const Evaluation& value)
        : ParentType(value)
    {}
    /*!
     * \copydoc ImmiscibleBoundaryRateVector::ImmiscibleBoundaryRateVector(const
     * ImmiscibleBoundaryRateVector& )
     */
    PvsBoundaryRateVector(const PvsBoundaryRateVector& value) = default;
    PvsBoundaryRateVector& operator=(const PvsBoundaryRateVector& value) = default;
    /*!
     * \copydoc ImmiscibleBoundaryRateVector::setFreeFlow
     */
    template <class Context, class FluidState>
    void setFreeFlow(const Context& context, unsigned bfIdx, unsigned timeIdx, const FluidState& fluidState)
    {
        ExtensiveQuantities extQuants;
        extQuants.updateBoundary(context, bfIdx, timeIdx, fluidState);
        const auto& insideIntQuants = context.intensiveQuantities(bfIdx, timeIdx);
        unsigned focusDofIdx = context.focusDofIndex();
        unsigned interiorDofIdx = context.interiorScvIndex(bfIdx, timeIdx);
        ////////
        // advective fluxes of all components in all phases
        ////////
        (*this) = Evaluation(0.0);
        for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
            Evaluation density;
            if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
                if (focusDofIdx == interiorDofIdx)
                    density = fluidState.density(phaseIdx);
                else
                    density = Opm::getValue(fluidState.density(phaseIdx));
            }
            else if (focusDofIdx == interiorDofIdx)
                density = insideIntQuants.fluidState().density(phaseIdx);
            else
                density = Opm::getValue(insideIntQuants.fluidState().density(phaseIdx));
            for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
                Evaluation molarity;
                if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
                    if (focusDofIdx == interiorDofIdx)
                        molarity = fluidState.molarity(phaseIdx, compIdx);
                    else
                        molarity = Opm::getValue(fluidState.molarity(phaseIdx, compIdx));
                }
                else if (focusDofIdx == interiorDofIdx)
                    molarity = insideIntQuants.fluidState().molarity(phaseIdx, compIdx);
                else
                    molarity = Opm::getValue(insideIntQuants.fluidState().molarity(phaseIdx, compIdx));
                // add advective flux of current component in current
                // phase
                (*this)[conti0EqIdx + compIdx] += extQuants.volumeFlux(phaseIdx)*molarity;
            }
            if (enableEnergy) {
                Evaluation specificEnthalpy;
                if (fluidState.pressure(phaseIdx) > insideIntQuants.fluidState().pressure(phaseIdx)) {
                    if (focusDofIdx == interiorDofIdx)
                        specificEnthalpy = fluidState.enthalpy(phaseIdx);
                    else
                        specificEnthalpy = Opm::getValue(fluidState.enthalpy(phaseIdx));
                }
                else if (focusDofIdx == interiorDofIdx)
                    specificEnthalpy = insideIntQuants.fluidState().enthalpy(phaseIdx);
                else
                    specificEnthalpy = Opm::getValue(insideIntQuants.fluidState().enthalpy(phaseIdx));
                Evaluation enthalpyRate = density*extQuants.volumeFlux(phaseIdx)*specificEnthalpy;
                EnergyModule::addToEnthalpyRate(*this, enthalpyRate);
            }
        }
        if (enableEnergy)
            // heat conduction
            EnergyModule::addToEnthalpyRate(*this, EnergyModule::thermalConductionRate(extQuants));
#ifndef NDEBUG
        for (unsigned i = 0; i < numEq; ++i)
            Opm::Valgrind::CheckDefined((*this)[i]);
#endif
    }
    /*!
     * \copydoc ImmiscibleBoundaryRateVector::setInFlow
     */
    template <class Context, class FluidState>
    void setInFlow(const Context& context,
                   unsigned bfIdx,
                   unsigned timeIdx,
                   const FluidState& fluidState)
    {
        this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
        // we only allow fluxes in the direction opposite to the outer unit normal
        for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            Evaluation& val = this->operator[](eqIdx);
            val = Toolbox::min(0.0, val);
        }
    }
    /*!
     * \copydoc ImmiscibleBoundaryRateVector::setOutFlow
     */
    template <class Context, class FluidState>
    void setOutFlow(const Context& context,
                    unsigned bfIdx,
                    unsigned timeIdx,
                    const FluidState& fluidState)
    {
        this->setFreeFlow(context, bfIdx, timeIdx, fluidState);
        // we only allow fluxes in the same direction as the outer unit normal
        for (unsigned eqIdx = 0; eqIdx < numEq; ++eqIdx) {
            Evaluation& val = this->operator[](eqIdx);
            val = Toolbox::max(0.0, val);
        }
    }
    /*!
     * \copydoc ImmiscibleBoundaryRateVector::setNoFlow
     */
    void setNoFlow()
    { (*this) = Evaluation(0.0); }
};
} // namespace Opm
#endif
 
     |