25 #ifndef OPM_COMPUTE_FROM_REFERENCE_PHASE_HPP 
   26 #define OPM_COMPUTE_FROM_REFERENCE_PHASE_HPP 
   30 #include <opm/common/Exceptions.hpp> 
   31 #include <opm/common/ErrorMacros.hpp> 
   34 #include <dune/common/fvector.hh> 
   63 template <
class Scalar, 
class Flu
idSystem, 
class Evaluation = Scalar>
 
   66     enum { numPhases = FluidSystem::numPhases };
 
   67     enum { numComponents = FluidSystem::numComponents };
 
   69     typedef Dune::FieldVector<Evaluation, numComponents> ComponentVector;
 
  106     template <
class Flu
idState, 
class ParameterCache>
 
  107     static void solve(FluidState &fluidState,
 
  108                       ParameterCache ¶mCache,
 
  109                       unsigned refPhaseIdx,
 
  117         paramCache.updatePhase(fluidState, refPhaseIdx);
 
  118         fluidState.setDensity(refPhaseIdx,
 
  119                               FluidSystem::density(fluidState,
 
  124             fluidState.setEnthalpy(refPhaseIdx,
 
  125                                    FluidSystem::enthalpy(fluidState,
 
  130             fluidState.setViscosity(refPhaseIdx,
 
  131                                     FluidSystem::viscosity(fluidState,
 
  136         for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
 
  137             fluidState.setFugacityCoefficient(refPhaseIdx,
 
  139                                               FluidSystem::fugacityCoefficient(fluidState,
 
  146         for (
unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
 
  147             if (phaseIdx == refPhaseIdx)
 
  150             ComponentVector fugVec;
 
  151             for (
unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
 
  152                 const auto& fug = fluidState.fugacity(refPhaseIdx, compIdx);
 
  153                 fugVec[compIdx] = FsToolbox::template toLhs<Evaluation>(fug);
 
  159                 fluidState.setViscosity(phaseIdx,
 
  160                                         FluidSystem::viscosity(fluidState,
 
  165                 fluidState.setEnthalpy(phaseIdx,
 
  166                                        FluidSystem::enthalpy(fluidState,
 
Calculates the chemical equilibrium from the component fugacities in a phase. 
Definition: CompositionFromFugacities.hpp:47
Definition: Air_Mesitylene.hpp:31
Some templates to wrap the valgrind client request macros. 
static void solve(FluidState &fluidState, ParameterCache ¶mCache, unsigned phaseIdx, const ComponentVector &targetFug)
Calculates the chemical equilibrium from the component fugacities in a phase. 
Definition: CompositionFromFugacities.hpp:84
Calculates the chemical equilibrium from the component fugacities in a phase. 
Computes all quantities of a generic fluid state if a reference phase has been specified. 
Definition: ComputeFromReferencePhase.hpp:64
static void solve(FluidState &fluidState, ParameterCache ¶mCache, unsigned refPhaseIdx, bool setViscosity, bool setEnthalpy)
Computes all quantities of a generic fluid state if a reference phase has been specified. 
Definition: ComputeFromReferencePhase.hpp:107