28#ifndef EWOMS_FLASH_LOCAL_RESIDUAL_HH 
   29#define EWOMS_FLASH_LOCAL_RESIDUAL_HH 
   35#include <opm/material/common/Valgrind.hpp> 
   44template <
class TypeTag>
 
   54    enum { numEq = getPropValue<TypeTag, Properties::NumEq>() };
 
   55    enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
 
   56    enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
 
   57    enum { conti0EqIdx = Indices::conti0EqIdx };
 
   59    enum { enableDiffusion = getPropValue<TypeTag, Properties::EnableDiffusion>() };
 
   62    enum { enableEnergy = getPropValue<TypeTag, Properties::EnableEnergy>() };
 
   65    using Toolbox = Opm::MathToolbox<Evaluation>;
 
   71    template <
class LhsEval>
 
   73                         const ElementContext& elemCtx,
 
   76                         unsigned phaseIdx)
 const 
   78        const IntensiveQuantities& intQuants = elemCtx.intensiveQuantities(dofIdx, timeIdx);
 
   79        const auto& fs = intQuants.fluidState();
 
   82        for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
 
   83            unsigned eqIdx = conti0EqIdx + compIdx;
 
   85                Toolbox::template decay<LhsEval>(fs.molarity(phaseIdx, compIdx))
 
   86                * Toolbox::template decay<LhsEval>(fs.saturation(phaseIdx))
 
   87                * Toolbox::template decay<LhsEval>(intQuants.porosity());
 
   90        EnergyModule::addPhaseStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx), phaseIdx);
 
   96    template <
class LhsEval>
 
   98                        const ElementContext& elemCtx,
 
  100                        unsigned timeIdx)
 const 
  103        for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx)
 
  106        EnergyModule::addSolidEnergyStorage(storage, elemCtx.intensiveQuantities(dofIdx, timeIdx));
 
  113                     const ElementContext& elemCtx,
 
  115                     unsigned timeIdx)
 const 
  119        Opm::Valgrind::CheckDefined(flux);
 
  122        Opm::Valgrind::CheckDefined(flux);
 
  129                          const ElementContext& elemCtx,
 
  131                          unsigned timeIdx)
 const 
  133        const auto& extQuants = elemCtx.extensiveQuantities(scvfIdx, timeIdx);
 
  135        unsigned focusDofIdx = elemCtx.focusDofIndex();
 
  136        for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  138            unsigned upIdx = 
static_cast<unsigned>(extQuants.upstreamIndex(phaseIdx));
 
  139            const IntensiveQuantities& up = elemCtx.intensiveQuantities(upIdx, timeIdx);
 
  144            if (upIdx == focusDofIdx) {
 
  146                    up.fluidState().molarDensity(phaseIdx)
 
  147                    * extQuants.volumeFlux(phaseIdx);
 
  149                for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
 
  150                    flux[conti0EqIdx + compIdx] +=
 
  151                        tmp*up.fluidState().moleFraction(phaseIdx, compIdx);
 
  156                    Toolbox::value(up.fluidState().molarDensity(phaseIdx))
 
  157                    * extQuants.volumeFlux(phaseIdx);
 
  159                for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
 
  160                    flux[conti0EqIdx + compIdx] +=
 
  161                        tmp*Toolbox::value(up.fluidState().moleFraction(phaseIdx, compIdx));
 
  166        EnergyModule::addAdvectiveFlux(flux, elemCtx, scvfIdx, timeIdx);
 
  173                          const ElementContext& elemCtx,
 
  175                          unsigned timeIdx)
 const 
  177        DiffusionModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
 
  178        EnergyModule::addDiffusiveFlux(flux, elemCtx, scvfIdx, timeIdx);
 
  185                       const ElementContext& elemCtx,
 
  187                       unsigned timeIdx)
 const 
  189        Opm::Valgrind::SetUndefined(source);
 
  190        elemCtx.problem().source(source, elemCtx, dofIdx, timeIdx);
 
  191        Opm::Valgrind::CheckDefined(source);
 
Provides the auxiliary methods required for consideration of the diffusion equation.
Definition: diffusionmodule.hh:48
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Calculates the local residual of the compositional multi-phase model based on flash calculations.
Definition: flash/flashlocalresidual.hh:46
void computeStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Evaluate the amount all conservation quantities (e.g. phase mass) within a finite sub-control volume.
Definition: flash/flashlocalresidual.hh:97
void addDiffusiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Adds the diffusive flux at a given flux integration point.
Definition: flash/flashlocalresidual.hh:172
void addPhaseStorage(Dune::FieldVector< LhsEval, numEq > &storage, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx, unsigned phaseIdx) const
Adds the amount all conservation quantities (e.g. phase mass) within a single fluid phase.
Definition: flash/flashlocalresidual.hh:72
void computeSource(RateVector &source, const ElementContext &elemCtx, unsigned dofIdx, unsigned timeIdx) const
Calculate the source term of the equation.
Definition: flash/flashlocalresidual.hh:184
void addAdvectiveFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Add the advective mass flux at a given flux integration point.
Definition: flash/flashlocalresidual.hh:128
void computeFlux(RateVector &flux, const ElementContext &elemCtx, unsigned scvfIdx, unsigned timeIdx) const
Evaluates the total mass flux of all conservation quantities over a face of a sub-control volume.
Definition: flash/flashlocalresidual.hh:112
Classes required for molecular diffusion.
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:37
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:235