All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PurePhaseModel.C
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 \*---------------------------------------------------------------------------*/
25 
26 #include "PurePhaseModel.H"
27 #include "phaseSystem.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class BasePhaseModel>
33 (
34  const phaseSystem& fluid,
35  const word& phaseName,
36  const bool referencePhase,
37  const label index
38 )
39 :
40  BasePhaseModel(fluid, phaseName, referencePhase, index)
41 {}
42 
43 
44 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
45 
46 template<class BasePhaseModel>
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
53 template<class BasePhaseModel>
55 {
56  return true;
57 }
58 
59 
60 template<class BasePhaseModel>
63 {
65  << "Cannot construct a species fraction equation for a pure phase"
66  << exit(FatalError);
67 
68  return tmp<fvScalarMatrix>();
69 }
70 
71 
72 template<class BasePhaseModel>
75 {
76  // Y_ has never been set, so we are returning an empty list
77 
78  return Y_;
79 }
80 
81 
82 template<class BasePhaseModel>
85 {
87  << "Cannot get a species fraction by name from a pure phase"
88  << exit(FatalError);
89 
90  return NullObjectRef<volScalarField>();
91 }
92 
93 
94 template<class BasePhaseModel>
97 {
99  << "Cannot access the species fractions of for a pure phase"
100  << exit(FatalError);
101 
102  return Y_;
103 }
104 
105 
106 template<class BasePhaseModel>
109 {
110  // Y_ has never been set, so we are returning an empty list
111 
112  return Y_;
113 }
114 
115 
116 template<class BasePhaseModel>
119 {
121  << "Cannot access the species fractions of for a pure phase"
122  << exit(FatalError);
123 
124  return Y_;
125 }
126 
127 
128 // ************************************************************************* //
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.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual UPtrList< volScalarField > & YActiveRef()
Access the active species mass fractions.
virtual PtrList< volScalarField > & YRef()
Access the species mass fractions.
virtual const PtrList< volScalarField > & Y() const
Return the species mass fractions.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
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
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