ptflash/flashprimaryvariables.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*
4 This file is part of the Open Porous Media project (OPM).
5
6 OPM is free software: you can redistribute it and/or modify
7 it under the terms of the GNU General Public License as published by
8 the Free Software Foundation, either version 2 of the License, or
9 (at your option) any later version.
10
11 OPM is distributed in the hope that it will be useful,
12 but WITHOUT ANY WARRANTY; without even the implied warranty of
13 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14 GNU General Public License for more details.
15
16 You should have received a copy of the GNU General Public License
17 along with OPM. If not, see <http://www.gnu.org/licenses/>.
18
19 Consult the COPYING file in the top-level source directory of this
20 module for the precise wording of the license and the list of
21 copyright holders.
22*/
28#ifndef OPM_PTFLASH_PRIMARY_VARIABLES_HH
29#define OPM_PTFLASH_PRIMARY_VARIABLES_HH
30
31#include "flashindices.hh"
32
35
36#include <opm/material/constraintsolvers/NcpFlash.hpp>
37#include <opm/material/fluidstates/CompositionalFluidState.hpp>
38#include <opm/material/common/Valgrind.hpp>
39
40#include <dune/common/fvector.hh>
41
42#include <iostream>
43
44namespace Opm {
45
55template <class TypeTag>
56class FlashPrimaryVariables : public FvBasePrimaryVariables<TypeTag>
57{
58 using ParentType = FvBasePrimaryVariables<TypeTag>;
59
60 using Scalar = GetPropType<TypeTag, Properties::Scalar>;
61 using Evaluation = GetPropType<TypeTag, Properties::Evaluation>;
62 using MaterialLawParams = GetPropType<TypeTag, Properties::MaterialLawParams>;
63 using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>;
64 using Indices = GetPropType<TypeTag, Properties::Indices>;
65
66 // primary variable indices
67 enum { z0Idx = Indices::z0Idx };
68 enum { pressure0Idx = Indices::pressure0Idx };
69
70 enum { numPhases = getPropValue<TypeTag, Properties::NumPhases>() };
71 enum { numComponents = getPropValue<TypeTag, Properties::NumComponents>() };
72
73 using Toolbox = typename Opm::MathToolbox<Evaluation>;
75
76public:
77 FlashPrimaryVariables() : ParentType()
78 { Opm::Valgrind::SetDefined(*this); }
79
83 FlashPrimaryVariables(Scalar value) : ParentType(value)
84 {
85 Opm::Valgrind::CheckDefined(value);
86 Opm::Valgrind::SetDefined(*this);
87 }
88
95
99 template <class FluidState>
100 void assignMassConservative(const FluidState& fluidState,
101 const MaterialLawParams&,
102 bool = false)
103 {
104 // there is no difference between naive and mass conservative
105 // assignment in the flash model. (we only need the total
106 // concentrations.)
107 assignNaive(fluidState);
108 }
109
113 template <class FluidState>
114 void assignNaive(const FluidState& fluidState)
115 {
116 // reset everything
117 (*this) = 0.0;
118
119 // assign the phase temperatures. this is out-sourced to
120 // the energy module
121 EnergyModule::setPriVarTemperatures(*this, fluidState);
122
123 // determine the component fractions
124 Dune::FieldVector<Scalar, numComponents> z(0.0);
125 Scalar sumMoles = 0.0;
126 for (unsigned phaseIdx = 0; phaseIdx < numPhases; ++phaseIdx) {
127 for (unsigned compIdx = 0; compIdx < numComponents; ++compIdx) {
128 Scalar tmp = Opm::getValue(fluidState.molarity(phaseIdx, compIdx) * fluidState.saturation(phaseIdx));
129 z[compIdx] += Opm::max(tmp, 1e-8);
130 sumMoles += tmp;
131 }
132 }
133 z /= sumMoles;
134
135 for (int i = 0; i < numComponents - 1; ++i)
136 (*this)[z0Idx + i] = z[i];
137
138 (*this)[pressure0Idx] = Opm::getValue(fluidState.pressure(0));
139 }
140
146 void print(std::ostream& os = std::cout) const
147 {
148 os << "(p_" << FluidSystem::phaseName(0) << " = "
149 << this->operator[](pressure0Idx);
150 for (unsigned compIdx = 0; compIdx < numComponents - 2; ++compIdx) {
151 os << ", z_" << FluidSystem::componentName(compIdx) << " = "
152 << this->operator[](z0Idx + compIdx);
153 }
154 os << ")" << std::flush;
155 }
156};
157
158} // namespace Opm
159
160#endif
Provides the auxiliary methods required for consideration of the energy equation.
Definition: energymodule.hh:50
Represents the primary variables used by the compositional flow model based on flash calculations.
Definition: flash/flashprimaryvariables.hh:57
void print(std::ostream &os=std::cout) const
Prints the names of the primary variables and their values.
Definition: ptflash/flashprimaryvariables.hh:146
void assignNaive(const FluidState &fluidState)
Directly retrieve the primary variables from an arbitrary fluid state.
Definition: flash/flashprimaryvariables.hh:113
FlashPrimaryVariables(Scalar value)
Constructor with assignment from scalar.
Definition: ptflash/flashprimaryvariables.hh:83
FlashPrimaryVariables()
Definition: ptflash/flashprimaryvariables.hh:77
FlashPrimaryVariables(const FlashPrimaryVariables &value)=default
void assignMassConservative(const FluidState &fluidState, const MaterialLawParams &, bool=false)
Set the primary variables from an arbitrary fluid state in a mass conservative way.
Definition: ptflash/flashprimaryvariables.hh:100
FlashPrimaryVariables & operator=(const FlashPrimaryVariables &value)=default
Contains the classes required to consider energy as a conservation quantity in a multi-phase module.
Definition: blackoilboundaryratevector.hh:37