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-2023 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 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
27 
28 namespace Foam
29 {
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 template<class Type>
35 {
36  public:
37 
39  Type t;
40 
42  {}
43 
45  :
46  wl(is),
47  t(is)
48  {}
49 };
50 
51 template<class Type>
53 {
54  return is >> wlat.wl >> wlat.t;
55 }
56 
57 template<class Type>
59 {
60  return os << wlat.wl << wlat.t;
61 }
62 
63 template<class Type>
64 inline bool operator==
65 (
66  const wordListAndType<Type>& a,
68 )
69 {
70  return a.wl == b.wl && a.t == b.t;
71 }
72 
73 template<class Type>
74 inline bool operator!=
75 (
76  const wordListAndType<Type>& a,
78 )
79 {
80  return !(a == b);
81 }
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 } // End namespace Foam
86 
87 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
88 
90 {
91  return mesh_;
92 }
93 
94 
96 {
97  return pimple_;
98 }
99 
100 
103 {
104  return phaseModels_;
105 }
106 
107 
110 {
111  return phaseModels_;
112 }
113 
114 
117 {
118  return movingPhaseModels_;
119 }
120 
121 
124 {
125  return movingPhaseModels_;
126 }
127 
128 
131 {
132  return stationaryPhaseModels_;
133 }
134 
135 
138 {
139  return stationaryPhaseModels_;
140 }
141 
142 
145 {
146  return anisothermalPhaseModels_;
147 }
148 
149 
152 {
153  return anisothermalPhaseModels_;
154 }
155 
156 
159 {
160  return multicomponentPhaseModels_;
161 }
162 
163 
166 {
167  return multicomponentPhaseModels_;
168 }
169 
170 
172 (
173  const phaseModel& phase
174 ) const
175 {
176  if (phaseModels_.size() != 2)
177  {
179  << "Call from a two-phase model in a multi-phase system."
180  << exit(FatalError);
181  }
182 
183  if (&phase == &phaseModels_[0])
184  {
185  return phaseModels_[1];
186  }
187  else
188  {
189  return phaseModels_[0];
190  }
191 }
192 
193 
195 {
196  return phi_;
197 }
198 
199 
201 {
202  return phi_;
203 }
204 
205 
207 {
208  return dpdt_;
209 }
210 
211 
213 {
214  return dpdt_;
215 }
216 
217 
219 {
220  return MRF_;
221 }
222 
223 
225 {
226  return Foam::fvModels::New(mesh);
227 }
228 
229 
231 {
232  return Foam::fvModels::New(mesh_);
233 }
234 
235 
237 {
238  return Foam::fvConstraints::New(mesh);
239 }
240 
241 
243 {
244  return Foam::fvConstraints::New(mesh_);
245 }
246 
247 
249 {
250  return deltaN_;
251 }
252 
253 
254 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
255 
256 namespace Foam
257 {
258 
259 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
260 
261 template<class GeoField>
262 inline void addField
263 (
264  const label phasei,
265  const word& name,
266  tmp<GeoField> field,
267  PtrList<GeoField>& fieldList
268 )
269 {
270  if (fieldList.set(phasei))
271  {
272  fieldList[phasei] += field;
273  }
274  else
275  {
276  fieldList.set
277  (
278  phasei,
279  new GeoField
280  (
281  name, // IOobject::groupName(name, group.name()),
282  field
283  )
284  );
285  }
286 }
287 
288 
289 template<class GeoField, class Group>
290 inline void addField
291 (
292  const Group& group,
293  const word& name,
294  tmp<GeoField> field,
295  PtrList<GeoField>& fieldList
296 )
297 {
298  if (fieldList.set(group.index()))
299  {
300  fieldList[group.index()] += field;
301  }
302  else
303  {
304  fieldList.set
305  (
306  group.index(),
307  new GeoField
308  (
309  IOobject::groupName(name, group.name()),
310  field
311  )
312  );
313  }
314 }
315 
316 
317 template<class GeoField, class Group>
318 inline void addField
319 (
320  const Group& group,
321  const word& name,
322  const GeoField& field,
323  PtrList<GeoField>& fieldList
324 )
325 {
326  addField(group, name, tmp<GeoField>(field), fieldList);
327 }
328 
329 
330 template<class GeoField, class Group>
331 inline void addField
332 (
333  const Group& group,
334  const word& name,
335  tmp<GeoField> field,
336  HashPtrTable<GeoField>& fieldTable
337 )
338 {
339  if (fieldTable.found(group.name()))
340  {
341  *fieldTable[group.name()] += field;
342  }
343  else
344  {
345  fieldTable.set
346  (
347  group.name(),
348  new GeoField
349  (
350  IOobject::groupName(name, group.name()),
351  field
352  )
353  );
354  }
355 }
356 
357 
358 template<class GeoField, class Group>
359 inline void addField
360 (
361  const Group& group,
362  const word& name,
363  const GeoField& field,
364  HashPtrTable<GeoField>& fieldTable
365 )
366 {
367  addField(group, name, tmp<GeoField>(field), fieldTable);
368 }
369 
370 
371 template<class Type, template<class> class PatchField, class GeoMesh>
373 (
376 )
377 {
379 
380  forAll(Afs, i)
381  {
382  Afs.set(i, As[i][i]*fs[i]);
383 
384  forAll(Afs, j)
385  {
386  if (j != i)
387  {
388  Afs[i] += As[i][j]*fs[j];
389  }
390  }
391  }
392 
393  return Afs;
394 }
395 
396 
397 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
398 
399 } // End namespace Foam
400 
401 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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:100
Finite volume constraints.
Definition: fvConstraints.H:67
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:99
Finite volume models.
Definition: fvModels.H:65
const IOMRFZoneList & MRF() const
Return MRF zones.
Definition: phaseSystemI.H:218
const phaseModel & otherPhase(const phaseModel &phase) const
Return the phase not given as an argument in a two-phase system.
Definition: phaseSystemI.H:172
const fvMesh & mesh_
Reference to the mesh.
Definition: phaseSystem.H:128
const dimensionedScalar & deltaN() const
Stabilisation for normalisation of the interface normal.
Definition: phaseSystemI.H:248
const Foam::fvModels & fvModels() const
Access the fvModels.
Definition: phaseSystemI.H:230
const surfaceScalarField & phi() const
Return the mixture flux.
Definition: phaseSystemI.H:194
const volScalarField & dpdt() const
Return the rate of change of the pressure.
Definition: phaseSystemI.H:206
const phaseModelPartialList & multicomponentPhases() const
Return the models for phases that have multiple species.
Definition: phaseSystemI.H:158
const phaseModelPartialList & movingPhases() const
Return the models for phases that are moving.
Definition: phaseSystemI.H:116
const phaseModelList & phases() const
Return the phase models.
Definition: phaseSystemI.H:102
const phaseModelPartialList & anisothermalPhases() const
Return the models for phases that have variable temperature.
Definition: phaseSystemI.H:144
const pimpleNoLoopControl & pimple() const
Return pimpleNoLoopControl.
Definition: phaseSystemI.H:95
const phaseModelPartialList & stationaryPhases() const
Return the models for phases that are stationary.
Definition: phaseSystemI.H:130
const fvMesh & mesh() const
Return the mesh.
Definition: phaseSystemI.H:89
const Foam::fvConstraints & fvConstraints() const
Access the fvConstraints.
Definition: phaseSystemI.H:242
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:44
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
volScalarField & b
Definition: createFields.H:25
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:263
Istream & operator>>(Istream &, pistonPointEdgeData &)
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError