Calculates the chemical equilibrium from the component fugacities in a phase.
More...
#include <CompositionFromFugacities.hpp>
|
template<class FluidState > |
static void | guessInitial (FluidState &fluidState, ParameterCache &, unsigned phaseIdx, const ComponentVector &) |
| Guess an initial value for the composition of the phase. More...
|
|
template<class FluidState > |
static void | solve (FluidState &fluidState, ParameterCache ¶mCache, unsigned phaseIdx, const ComponentVector &targetFug) |
| Calculates the chemical equilibrium from the component fugacities in a phase. More...
|
|
|
template<class FluidState > |
static void | solveIdealMix_ (FluidState &fluidState, ParameterCache ¶mCache, unsigned phaseIdx, const ComponentVector &fugacities) |
|
template<class FluidState > |
static Scalar | linearize_ (Dune::FieldMatrix< Evaluation, numComponents, numComponents > &J, Dune::FieldVector< Evaluation, numComponents > &defect, FluidState &fluidState, ParameterCache ¶mCache, unsigned phaseIdx, const ComponentVector &targetFug) |
|
template<class FluidState > |
static Scalar | update_ (FluidState &fluidState, ParameterCache ¶mCache, Dune::FieldVector< Evaluation, numComponents > &x, Dune::FieldVector< Evaluation, numComponents > &, unsigned phaseIdx, const Dune::FieldVector< Evaluation, numComponents > &targetFug) |
|
template<class FluidState > |
static Scalar | calculateDefect_ (const FluidState ¶ms, unsigned phaseIdx, const ComponentVector &targetFug) |
|
template<class Scalar, class FluidSystem, class Evaluation = Scalar>
class Opm::CompositionFromFugacities< Scalar, FluidSystem, Evaluation >
Calculates the chemical equilibrium from the component fugacities in a phase.
template<class Scalar , class FluidSystem , class Evaluation = Scalar>
template<class Scalar , class FluidSystem , class Evaluation = Scalar>
template<class FluidState >
template<class Scalar , class FluidSystem , class Evaluation = Scalar>
template<class FluidState >
Guess an initial value for the composition of the phase.
template<class Scalar , class FluidSystem , class Evaluation = Scalar>
template<class FluidState >
static Scalar Opm::CompositionFromFugacities< Scalar, FluidSystem, Evaluation >::linearize_ |
( |
Dune::FieldMatrix< Evaluation, numComponents, numComponents > & |
J, |
|
|
Dune::FieldVector< Evaluation, numComponents > & |
defect, |
|
|
FluidState & |
fluidState, |
|
|
ParameterCache & |
paramCache, |
|
|
unsigned |
phaseIdx, |
|
|
const ComponentVector & |
targetFug |
|
) |
| |
|
inlinestaticprotected |
template<class Scalar , class FluidSystem , class Evaluation = Scalar>
template<class FluidState >
Calculates the chemical equilibrium from the component fugacities in a phase.
The phase's fugacities must already be set.
References Valgrind::CheckDefined(), Opm::CompositionFromFugacities< Scalar, FluidSystem, Evaluation >::linearize_(), Opm::CompositionFromFugacities< Scalar, FluidSystem, Evaluation >::solveIdealMix_(), and Opm::CompositionFromFugacities< Scalar, FluidSystem, Evaluation >::update_().
Referenced by Opm::ComputeFromReferencePhase< Scalar, FluidSystem, Evaluation >::solve().
template<class Scalar , class FluidSystem , class Evaluation = Scalar>
template<class FluidState >
template<class Scalar , class FluidSystem , class Evaluation = Scalar>
template<class FluidState >
static Scalar Opm::CompositionFromFugacities< Scalar, FluidSystem, Evaluation >::update_ |
( |
FluidState & |
fluidState, |
|
|
ParameterCache & |
paramCache, |
|
|
Dune::FieldVector< Evaluation, numComponents > & |
x, |
|
|
Dune::FieldVector< Evaluation, numComponents > & |
, |
|
|
unsigned |
phaseIdx, |
|
|
const Dune::FieldVector< Evaluation, numComponents > & |
targetFug |
|
) |
| |
|
inlinestaticprotected |
The documentation for this class was generated from the following file:
|