PurePhaseModel.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2015-2020 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 Class
25  Foam::PurePhaseModel
26 
27 Description
28  Class which represents pure phases, i.e. without any species. Returns an
29  empty list of mass fractions.
30 
31 SourceFiles
32  PurePhaseModel.C
33 
34 \*---------------------------------------------------------------------------*/
35 
36 #ifndef PurePhaseModel_H
37 #define PurePhaseModel_H
38 
39 #include "phaseModel.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
45 
46 /*---------------------------------------------------------------------------*\
47  Class PurePhaseModel Declaration
48 \*---------------------------------------------------------------------------*/
49 
50 template<class BasePhaseModel>
51 class PurePhaseModel
52 :
53  public BasePhaseModel
54 {
55 protected:
56 
57  // Protected data
58 
59  //- Empty mass fraction field list
61 
62 
63 public:
64 
65  // Constructors
66 
68  (
69  const phaseSystem& fluid,
70  const word& phaseName,
71  const bool referencePhase,
72  const label index
73  );
74 
75 
76  //- Destructor
77  virtual ~PurePhaseModel();
78 
79 
80  // Member Functions
81 
82  // Species
83 
84  //- Return whether the phase is pure (i.e., not multi-component)
85  virtual bool pure() const;
86 
87  //- Return the species fraction equation
89 
90  //- Return the species mass fractions
91  virtual const PtrList<volScalarField>& Y() const;
92 
93  //- Return a species mass fraction by name
94  virtual const volScalarField& Y(const word& name) const;
95 
96  //- Access the species mass fractions
97  virtual PtrList<volScalarField>& YRef();
98 
99  //- Return the active species mass fractions
100  virtual const UPtrList<volScalarField>& YActive() const;
101 
102  //- Access the active species mass fractions
104 };
105 
106 
107 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
108 
109 } // End namespace Foam
110 
111 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
112 
113 #ifdef NoRepository
114  #include "PurePhaseModel.C"
115 #endif
116 
117 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
118 
119 #endif
120 
121 // ************************************************************************* //
virtual tmp< fvScalarMatrix > YiEqn(volScalarField &Yi)
Return the species fraction equation.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
virtual const UPtrList< volScalarField > & YActive() const
Return the active species mass fractions.
PtrList< volScalarField > Y_
Empty mass fraction field list.
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
Generic GeometricField class.
virtual const PtrList< volScalarField > & Y() const
Return the species mass fractions.
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:67
A class for handling words, derived from string.
Definition: word.H:59
phaseSystem & fluid
Definition: createFields.H:11
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:54
virtual bool pure() const
Return whether the phase is pure (i.e., not multi-component)
virtual ~PurePhaseModel()
Destructor.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Class which represents pure phases, i.e. without any species. Returns an empty list of mass fractions...
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
PurePhaseModel(const phaseSystem &fluid, const word &phaseName, const bool referencePhase, const label index)
A class for managing temporary objects.
Definition: PtrList.H:53
Namespace for OpenFOAM.