ReactingParcelI.H
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class ParcelType>
32 inline
34 :
35  ParcelType::constantProperties(),
36  pMin_(this->dict_, 0.0),
37  constantVolume_(this->dict_, false),
38  volUpdateType_(this->dict_, mUndefined)
39 {}
40 
41 
42 template<class ParcelType>
44 (
45  const constantProperties& cp
46 )
47 :
48  ParcelType::constantProperties(cp),
49  pMin_(cp.pMin_),
50  constantVolume_(cp.constantVolume_),
51  volUpdateType_(cp.volUpdateType_)
52 {}
53 
54 
55 template<class ParcelType>
57 (
58  const dictionary& parentDict
59 )
60 :
61  ParcelType::constantProperties(parentDict),
62  pMin_(this->dict_, "pMin", 1000.0),
63  constantVolume_(this->dict_, "constantVolume", false),
64  volUpdateType_(this->dict_, "volumeUpdateMethod")
65 {
66  word updateName;
67 
68  // If constantVolume found use it
69  if (this->dict_.found("constantVolume"))
70  {
71  volUpdateType_.setValue(mUndefined);
72  }
73  else if (this->dict_.readIfPresent("volumeUpdateMethod", updateName))
74  {
75  // Manual version of Enum<volumeUpdateType>
76  if (updateName == "constantRho")
77  {
78  volUpdateType_.setValue(mConstRho);
79  }
80  else if (updateName == "constantVolume")
81  {
82  volUpdateType_.setValue(mConstVol);
83  }
84  else if (updateName == "updateRhoAndVol")
85  {
86  volUpdateType_.setValue(mUpdateRhoAndVol);
87  }
88  else
89  {
90  // FatalIOErrorInLookup
91 
92  FatalIOErrorInFunction(this->dict_)
93  << "Unknown volumeUpdateMethod type " << updateName
94  << "\n\nValid volumeUpdateMethod types :\n"
95  << "(constantRho constantVolume updateRhoAndVol)" << nl
96  << exit(FatalIOError);
97  }
98  }
99  else
100  {
101  constantVolume_.setValue(false);
102  }
103 }
104 
105 
106 template<class ParcelType>
108 (
109  const polyMesh& mesh,
110  const barycentric& coordinates,
111  const label celli,
112  const label tetFacei,
113  const label tetPti
114 )
115 :
116  ParcelType(mesh, coordinates, celli, tetFacei, tetPti),
117  mass0_(0.0),
118  Y_(0)
119 {}
120 
121 
122 template<class ParcelType>
124 (
125  const polyMesh& mesh,
126  const vector& position,
127  const label celli
128 )
129 :
130  ParcelType(mesh, position, celli),
131  mass0_(0.0),
132  Y_(0)
133 {}
134 
135 
136 template<class ParcelType>
138 (
139  const polyMesh& mesh,
140  const barycentric& coordinates,
141  const label celli,
142  const label tetFacei,
143  const label tetPti,
144  const label typeId,
145  const scalar nParticle0,
146  const scalar d0,
147  const scalar dTarget0,
148  const vector& U0,
149  const vector& f0,
150  const vector& angularMomentum0,
151  const vector& torque0,
152  const scalarField& Y0,
153  const constantProperties& constProps
154 )
155 :
156  ParcelType
157  (
158  mesh,
159  coordinates,
160  celli,
161  tetFacei,
162  tetPti,
163  typeId,
164  nParticle0,
165  d0,
166  dTarget0,
167  U0,
168  f0,
169  angularMomentum0,
170  torque0,
171  constProps
172  ),
173  mass0_(0.0),
174  Y_(Y0)
175 {
176  // Set initial parcel mass
177  mass0_ = this->mass();
178 }
179 
180 
181 // * * * * * * * * * constantProperties Member Functions * * * * * * * * * * //
182 
183 template<class ParcelType>
184 inline Foam::scalar
186 {
187  return pMin_.value();
188 }
189 
190 
191 template<class ParcelType>
192 inline bool
194 {
195  return constantVolume_.value();
196 }
197 
198 
199 template<class ParcelType>
200 inline Foam::label
202 {
203  return volUpdateType_.value();
204 }
205 
206 
207 // * * * * * * * * * * ThermoParcel Member Functions * * * * * * * * * * * * //
208 
209 template<class ParcelType>
210 inline Foam::scalar Foam::ReactingParcel<ParcelType>::mass0() const
211 {
212  return mass0_;
213 }
214 
215 
216 template<class ParcelType>
218 {
219  return Y_;
220 }
221 
222 
223 template<class ParcelType>
225 {
226  return Y_;
227 }
228 
229 
230 template<class ParcelType>
231 inline const Foam::scalarField&
233 {
234  return Y_;
235 }
236 
237 
238 template<class ParcelType>
239 inline const Foam::scalarField&
241 {
242  return Y_;
243 }
244 
245 
246 template<class ParcelType>
248 {
249  return mass0_;
250 }
251 
252 
253 template<class ParcelType>
255 {
256  return Y_;
257 }
258 
259 
260 // ************************************************************************* //
Foam::ReactingParcel::mass0
scalar mass0() const
Return const access to initial mass [kg].
Definition: ReactingParcelI.H:210
Foam::ReactingParcel::YSolid
const scalarField & YSolid() const
Return const access to mass fractions of solids.
Definition: ReactingParcelI.H:240
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::FatalIOError
IOerror FatalIOError
Foam::ReactingParcel::Y
const scalarField & Y() const
Return const access to mass fractions of mixture [].
Definition: ReactingParcelI.H:217
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
Foam::ReactingParcel::YLiquid
const scalarField & YLiquid() const
Return const access to mass fractions of liquids.
Definition: ReactingParcelI.H:232
Foam::Field< scalar >
coordinates
PtrList< coordinateSystem > coordinates(solidRegions.size())
Foam::Barycentric< scalar >
Foam::ReactingParcel::constantProperties::volUpdateType
label volUpdateType() const
Return const access to the constant volume flag.
Definition: ReactingParcelI.H:201
Foam::ReactingParcel::YGas
const scalarField & YGas() const
Return const access to mass fractions of gases.
Definition: ReactingParcelI.H:224
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::ReactingParcel::ReactingParcel
ReactingParcel(const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
Construct from mesh, coordinates and topology.
Definition: ReactingParcelI.H:108
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::ReactingParcel::mass0_
scalar mass0_
Initial mass [kg].
Definition: ReactingParcel.H:198
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::ReactingParcel::constantProperties::constantVolume
bool constantVolume() const
Return const access to the constant volume flag.
Definition: ReactingParcelI.H:193
Foam::cp
bool cp(const fileName &src, const fileName &dst, const bool followLink=true)
Copy the source to the destination (recursively if necessary).
Definition: MSwindows.C:802
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::ReactingParcel::constantProperties
Class to hold reacting parcel constant properties.
Definition: ReactingParcel.H:83
Foam::ReactingParcel::constantProperties::constantProperties
constantProperties()
Null constructor.
Definition: ReactingParcelI.H:33
Foam::ReactingParcel::constantProperties::pMin
scalar pMin() const
Return const access to the minimum pressure.
Definition: ReactingParcelI.H:185
Foam::ReactingParcel::Y_
scalarField Y_
Mass fractions of mixture [].
Definition: ReactingParcel.H:201
FatalIOErrorInFunction
#define FatalIOErrorInFunction(ios)
Report an error message using Foam::FatalIOError.
Definition: error.H:473
Y0
scalarList Y0(nSpecie, Zero)