MRFZoneList.C
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) 2012-2017 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "MRFZoneList.H"
29 #include "volFields.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 Foam::MRFZoneList::MRFZoneList
35 (
36  const fvMesh& mesh,
37  const dictionary& dict
38 )
39 :
41  mesh_(mesh)
42 {
43  reset(dict);
44 
45  active(true);
46 }
47 
48 
49 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
50 
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
57 bool Foam::MRFZoneList::active(const bool warn) const
58 {
59  bool a = false;
60  forAll(*this, i)
61  {
62  a = a || this->operator[](i).active();
63  }
64 
65  if (warn && this->size() && !a)
66  {
67  Info<< " No MRF zones active" << endl;
68  }
69 
70  return a;
71 }
72 
73 
75 {
76  label count = 0;
77  for (const entry& dEntry : dict)
78  {
79  if (dEntry.isDict())
80  {
81  ++count;
82  }
83  }
84 
85  this->resize(count);
86 
87  count = 0;
88  for (const entry& dEntry : dict)
89  {
90  if (dEntry.isDict())
91  {
92  const word& name = dEntry.keyword();
93  const dictionary& modelDict = dEntry.dict();
94 
95  Info<< " creating MRF zone: " << name << endl;
96 
97  this->set
98  (
99  count++,
100  new MRFZone(name, mesh_, modelDict)
101  );
102  }
103  }
104 }
105 
106 
108 {
109  bool allOk = true;
110  forAll(*this, i)
111  {
112  MRFZone& pm = this->operator[](i);
113  bool ok = pm.read(dict.subDict(pm.name()));
114  allOk = (allOk && ok);
115  }
116  return allOk;
117 }
118 
119 
121 {
122  forAll(*this, i)
123  {
124  os << nl;
125  this->operator[](i).writeData(os);
126  }
127 
128  return os.good();
129 }
130 
131 
133 (
134  const volVectorField& U,
135  volVectorField& ddtU
136 ) const
137 {
138  forAll(*this, i)
139  {
140  operator[](i).addCoriolis(U, ddtU);
141  }
142 }
143 
144 
146 {
147  forAll(*this, i)
148  {
149  operator[](i).addCoriolis(UEqn);
150  }
151 }
152 
153 
155 (
156  const volScalarField& rho,
158 ) const
159 {
160  forAll(*this, i)
161  {
162  operator[](i).addCoriolis(rho, UEqn);
163  }
164 }
165 
166 
168 (
169  const volVectorField& U
170 ) const
171 {
172  tmp<volVectorField> tacceleration
173  (
174  new volVectorField
175  (
176  IOobject
177  (
178  "MRFZoneList:acceleration",
179  U.mesh().time().timeName(),
180  U.mesh()
181  ),
182  U.mesh(),
183  dimensionedVector(U.dimensions()/dimTime, Zero)
184  )
185  );
186  volVectorField& acceleration = tacceleration.ref();
187 
188  forAll(*this, i)
189  {
190  operator[](i).addCoriolis(U, acceleration);
191  }
192 
193  return tacceleration;
194 }
195 
196 
198 (
199  const volScalarField& rho,
200  const volVectorField& U
201 ) const
202 {
203  return rho*DDt(U);
204 }
205 
206 
208 {
209  forAll(*this, i)
210  {
211  operator[](i).makeRelative(U);
212  }
213 }
214 
215 
217 {
218  forAll(*this, i)
219  {
220  operator[](i).makeRelative(phi);
221  }
222 }
223 
224 
226 (
227  const tmp<surfaceScalarField>& tphi
228 ) const
229 {
230  if (size())
231  {
233  (
234  New
235  (
236  tphi,
237  "relative(" + tphi().name() + ')',
238  tphi().dimensions(),
239  true
240  )
241  );
242 
243  makeRelative(rphi.ref());
244 
245  tphi.clear();
246 
247  return rphi;
248  }
249  else
250  {
251  return tmp<surfaceScalarField>(tphi, true);
252  }
253 }
254 
255 
258 (
260 ) const
261 {
262  if (size())
263  {
264  tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
265 
266  forAll(*this, i)
267  {
268  operator[](i).makeRelative(rphi.ref());
269  }
270 
271  tphi.clear();
272 
273  return rphi;
274  }
275  else
276  {
277  return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
278  }
279 }
280 
281 
284 (
285  const tmp<Field<scalar>>& tphi,
286  const label patchi
287 ) const
288 {
289  if (size())
290  {
291  tmp<Field<scalar>> rphi(New(tphi, true));
292 
293  forAll(*this, i)
294  {
295  operator[](i).makeRelative(rphi.ref(), patchi);
296  }
297 
298  tphi.clear();
299 
300  return rphi;
301  }
302  else
303  {
304  return tmp<Field<scalar>>(tphi, true);
305  }
306 }
307 
308 
310 (
311  const surfaceScalarField& rho,
313 ) const
314 {
315  forAll(*this, i)
316  {
317  operator[](i).makeRelative(rho, phi);
318  }
319 }
320 
321 
323 {
324  forAll(*this, i)
325  {
326  operator[](i).makeAbsolute(U);
327  }
328 }
329 
330 
332 {
333  forAll(*this, i)
334  {
335  operator[](i).makeAbsolute(phi);
336  }
337 }
338 
339 
341 (
342  const tmp<surfaceScalarField>& tphi
343 ) const
344 {
345  if (size())
346  {
348  (
349  New
350  (
351  tphi,
352  "absolute(" + tphi().name() + ')',
353  tphi().dimensions(),
354  true
355  )
356  );
357 
358  makeAbsolute(rphi.ref());
359 
360  tphi.clear();
361 
362  return rphi;
363  }
364  else
365  {
366  return tmp<surfaceScalarField>(tphi, true);
367  }
368 }
369 
370 
372 (
373  const surfaceScalarField& rho,
375 ) const
376 {
377  forAll(*this, i)
378  {
379  operator[](i).makeAbsolute(rho, phi);
380  }
381 }
382 
383 
385 {
386  forAll(*this, i)
387  {
388  operator[](i).correctBoundaryVelocity(U);
389  }
390 }
391 
392 
394 (
395  const volVectorField& U,
397 ) const
398 {
400  (
401  relative(mesh_.Sf().boundaryField() & U.boundaryField())
402  );
403 
404 
405  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
406 
407  forAll(mesh_.boundary(), patchi)
408  {
409  if
410  (
411  isA<fixedValueFvsPatchScalarField>(phibf[patchi])
412  )
413  {
414  phibf[patchi] == Uf[patchi];
415  }
416  }
417 }
418 
419 
421 {
422  if (mesh_.topoChanging())
423  {
424  forAll(*this, i)
425  {
426  operator[](i).update();
427  }
428  }
429 }
430 
431 
432 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
433 
434 Foam::Ostream& Foam::operator<<
435 (
436  Ostream& os,
437  const MRFZoneList& models
438 )
439 {
440  models.writeData(os);
441  return os;
442 }
443 
444 
445 // ************************************************************************* //
Foam::MRFZone
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:67
Foam::entry
A keyword and a list of tokens is an 'entry'.
Definition: entry.H:67
volFields.H
Foam::MRFZoneList::addAcceleration
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:133
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::BitOps::set
void set(List< bool > &bools, const labelRange &range)
Set the specified range 'on' in a boolList.
Definition: BitOps.C:37
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:291
Foam::MRFZoneList::writeData
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:120
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::MRFZoneList
List container for MRF zomes.
Definition: MRFZoneList.H:57
Foam::MRFZoneList::reset
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:74
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
rho
rho
Definition: readInitialConditions.H:88
Foam::MRFZoneList::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:384
Foam::MRFZoneList::active
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:57
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:77
Foam::dimensionedVector
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Definition: dimensionedVector.H:50
Foam::MRFZone::read
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:592
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:54
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:228
resize
patchWriters resize(patchIds.size())
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (uses stdout - output is on the master only)
Foam::MRFZoneList::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:207
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::MRFZoneList::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:341
Foam::MRFZoneList::DDt
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:168
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:528
Foam::PtrList< MRFZone >
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dict
dictionary dict
Definition: searchingEngine.H:14
fixedValueFvsPatchFields.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
Foam::MRFZoneList::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:226
Foam::New
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
Global function forwards to reuseTmpDimensionedField::New.
Definition: DimensionedFieldReuseFunctions.H:105
U
U
Definition: pEqn.H:72
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::BitOps::count
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition: BitOps.H:77
Foam::GeometricField::ref
Internal & ref(const bool updateAccessTime=true)
Return a reference to the dimensioned internal field.
Definition: GeometricField.C:749
Foam::GeometricField::boundaryFieldRef
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
Definition: GeometricField.C:783
Foam::MRFZoneList::makeAbsolute
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:322
MRFZoneList.H
Foam::MRFZoneList::correctBoundaryFlux
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:394
Foam::MRFZoneList::update
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:420
Foam::fvc::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition: fvcMeshPhi.C:155
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::MRFZoneList::read
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:107
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::MRFZoneList::~MRFZoneList
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:51
Foam::fvc::makeAbsolute
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:116
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:224
Foam::MRFZone::name
const word & name() const
Return const access to the MRF region name.
Definition: MRFZoneI.H:28
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fvc::DDt
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:47