phaseSystemI.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) 2014-2025 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 "phaseSystem.H"
27 
28 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
29 
30 namespace Foam
31 {
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 template<class Type>
37 {
38  public:
39 
41  Type t;
42 
44  {}
45 
47  :
48  wl(is),
49  t(is)
50  {}
51 };
52 
53 template<class Type>
55 {
56  return is >> wlat.wl >> wlat.t;
57 }
58 
59 template<class Type>
61 {
62  return os << wlat.wl << wlat.t;
63 }
64 
65 template<class Type>
66 inline bool operator==
67 (
68  const wordListAndType<Type>& a,
70 )
71 {
72  return a.wl == b.wl && a.t == b.t;
73 }
74 
75 template<class Type>
76 inline bool operator!=
77 (
78  const wordListAndType<Type>& a,
80 )
81 {
82  return !(a == b);
83 }
84 
85 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
86 
87 } // End namespace Foam
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 {
93  return mesh_;
94 }
95 
96 
98 {
99  return pimple_;
100 }
101 
102 
105 {
106  return phaseModels_;
107 }
108 
109 
112 {
113  return phaseModels_;
114 }
115 
116 
119 {
120  return movingPhaseModels_;
121 }
122 
123 
126 {
127  return movingPhaseModels_;
128 }
129 
130 
133 {
134  return stationaryPhaseModels_;
135 }
136 
137 
140 {
141  return stationaryPhaseModels_;
142 }
143 
144 
147 {
148  return thermalPhaseModels_;
149 }
150 
151 
154 {
155  return thermalPhaseModels_;
156 }
157 
158 
161 {
162  return multicomponentPhaseModels_;
163 }
164 
165 
168 {
169  return multicomponentPhaseModels_;
170 }
171 
172 
174 (
175  const phaseModel& phase
176 ) const
177 {
178  if (phaseModels_.size() != 2)
179  {
181  << "Call from a two-phase model in a multi-phase system."
182  << exit(FatalError);
183  }
184 
185  if (&phase == &phaseModels_[0])
186  {
187  return phaseModels_[1];
188  }
189  else
190  {
191  return phaseModels_[0];
192  }
193 }
194 
195 
197 {
198  return phi_;
199 }
200 
201 
203 {
204  return phi_;
205 }
206 
207 
209 {
210  return dpdt_;
211 }
212 
213 
215 {
216  return dpdt_;
217 }
218 
219 
221 {
222  return MRF_;
223 }
224 
225 
227 {
228  return Foam::fvModels::New(mesh);
229 }
230 
231 
233 {
234  return Foam::fvModels::New(mesh_);
235 }
236 
237 
239 {
241 }
242 
243 
245 {
246  return Foam::fvConstraints::New(mesh_);
247 }
248 
249 
251 {
252  return deltaN_;
253 }
254 
255 
256 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
257 
258 namespace Foam
259 {
260 
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263 template<class GeoField>
264 inline void addField
265 (
266  const label phasei,
267  const word& name,
268  tmp<GeoField> field,
269  PtrList<GeoField>& fieldList
270 )
271 {
272  if (fieldList.set(phasei))
273  {
274  fieldList[phasei] += field;
275  }
276  else
277  {
278  fieldList.set
279  (
280  phasei,
281  new GeoField
282  (
283  name, // IOobject::groupName(name, group.name()),
284  field
285  )
286  );
287  }
288 }
289 
290 
291 template<class GeoField, class Group>
292 inline void addField
293 (
294  const Group& group,
295  const word& name,
296  tmp<GeoField> field,
297  PtrList<GeoField>& fieldList
298 )
299 {
300  if (fieldList.set(group.index()))
301  {
302  fieldList[group.index()] += field;
303  }
304  else
305  {
306  fieldList.set
307  (
308  group.index(),
309  new GeoField
310  (
311  IOobject::groupName(name, group.name()),
312  field
313  )
314  );
315  }
316 }
317 
318 
319 template<class GeoField, class Group>
320 inline void addField
321 (
322  const Group& group,
323  const word& name,
324  const GeoField& field,
325  PtrList<GeoField>& fieldList
326 )
327 {
328  addField(group, name, tmp<GeoField>(field), fieldList);
329 }
330 
331 
332 template<class GeoField, class Group>
333 inline void addField
334 (
335  const Group& group,
336  const word& name,
337  tmp<GeoField> field,
338  HashPtrTable<GeoField>& fieldTable
339 )
340 {
341  if (fieldTable.found(group.name()))
342  {
343  *fieldTable[group.name()] += field;
344  }
345  else
346  {
347  fieldTable.set
348  (
349  group.name(),
350  new GeoField
351  (
352  IOobject::groupName(name, group.name()),
353  field
354  )
355  );
356  }
357 }
358 
359 
360 template<class GeoField, class Group>
361 inline void addField
362 (
363  const Group& group,
364  const word& name,
365  const GeoField& field,
366  HashPtrTable<GeoField>& fieldTable
367 )
368 {
369  addField(group, name, tmp<GeoField>(field), fieldTable);
370 }
371 
372 
373 template<class Type, class GeoMesh>
375 (
378 )
379 {
381 
382  forAll(Afs, i)
383  {
384  Afs.set(i, As[i][i]*fs[i]);
385 
386  forAll(Afs, j)
387  {
388  if (j != i)
389  {
390  Afs[i] += As[i][j]*fs[j];
391  }
392  }
393  }
394 
395  return Afs;
396 }
397 
398 
399 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
400 
401 } // End namespace Foam
402 
403 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
Generic GeometricField class.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
bool set(const Key &, const T &newElmt)
Set a new hashedEntry, overwriting existing entries.
Definition: HashTableI.H:91
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
List of MRF zones with IO functionality. MRF zones are specified by a list of dictionary entries,...
Definition: IOMRFZoneList.H:68
static word groupName(Name name, const word &group)
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Finite volume constraints.
Definition: fvConstraints.H:67
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:96
Finite volume models.
Definition: fvModels.H:65
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:220
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:174
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:156
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
Definition: phaseSystemI.H:250
const Foam::fvModels & fvModels() const
Access the fvModels.
Definition: phaseSystemI.H:232
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:196
const volScalarField & dpdt() const
Return the rate of change of the pressure.
Definition: phaseSystemI.H:208
const phaseModelPartialList & multicomponentPhases() const
Return the models for phases that have multiple species.
Definition: phaseSystemI.H:160
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:118
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:104
const pimpleNoLoopControl & pimple() const
Return pimpleNoLoopControl.
Definition: phaseSystemI.H:97
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:132
const phaseModelPartialList & thermalPhases() const
Return the models for phases that have variable temperature.
Definition: phaseSystemI.H:146
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:91
const Foam::fvConstraints & fvConstraints() const
Access the fvConstraints.
Definition: phaseSystemI.H:244
Pimple no-loop control class. Implements various option flags, but leaves loop controls to the deriva...
A class for managing temporary objects.
Definition: tmp.H:55
wordListAndType(Istream &is)
Definition: phaseSystemI.H:46
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField & b
Definition: createFields.H:27
const char *const group
Group name for atomic constants.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:265
Istream & operator>>(Istream &, pistonPointEdgeData &)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError