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  Copyright (C) 2021 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 #include "MRFZoneList.H"
30 #include "volFields.H"
32 
33 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
34 
35 Foam::MRFZoneList::MRFZoneList
36 (
37  const fvMesh& mesh,
38  const dictionary& dict
39 )
40 :
42  mesh_(mesh)
43 {
44  reset(dict);
45 
46  active(true);
47 }
48 
49 
50 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
51 
52 bool Foam::MRFZoneList::active(const bool warn) const
53 {
54  bool a = false;
55  forAll(*this, i)
56  {
57  a = a || this->operator[](i).active();
58  }
59 
60  if (warn && this->size() && !a)
61  {
62  Info<< " No MRF zones active" << endl;
63  }
64 
65  return a;
66 }
67 
68 
70 {
71  label count = 0;
72  for (const entry& dEntry : dict)
73  {
74  if (dEntry.isDict())
75  {
76  ++count;
77  }
78  }
79 
80  this->resize(count);
81 
82  count = 0;
83  for (const entry& dEntry : dict)
84  {
85  if (dEntry.isDict())
86  {
87  const word& name = dEntry.keyword();
88  const dictionary& modelDict = dEntry.dict();
89 
90  Info<< " creating MRF zone: " << name << endl;
91 
92  this->set
93  (
94  count++,
95  new MRFZone(name, mesh_, modelDict)
96  );
97  }
98  }
99 }
100 
101 
103 (
104  const word& name
105 ) const
106 {
108  for (const auto& mrf: *this)
109  {
110  if (mrf.name() == name)
111  {
112  return mrf;
113  }
114 
115  names.append(mrf.name());
116  }
117 
119  << "Unable to find MRFZone " << name
120  << ". Available zones are: " << names
121  << exit(FatalError);
122 
123  return first();
124 }
125 
126 
128 {
129  bool allOk = true;
130  for (auto& mrf: *this)
131  {
132  bool ok = mrf.read(dict.subDict(mrf.name()));
133  allOk = (allOk && ok);
134  }
135  return allOk;
136 }
137 
138 
140 {
141  for (const auto& mrf: *this)
142  {
143  os << nl;
144  mrf.writeData(os);
145  }
146 
147  return os.good();
148 }
149 
150 
152 (
153  const volVectorField& U,
154  volVectorField& ddtU
155 ) const
156 {
157  for (const auto& mrf: *this)
158  {
159  mrf.addCoriolis(U, ddtU);
160  }
161 }
162 
163 
165 {
166  for (const auto& mrf: *this)
167  {
168  mrf.addCoriolis(UEqn);
169  }
170 }
171 
172 
174 (
175  const volScalarField& rho,
177 ) const
178 {
179  for (const auto& mrf: *this)
180  {
181  mrf.addCoriolis(rho, UEqn);
182  }
183 }
184 
185 
187 (
188  const volVectorField& U
189 ) const
190 {
191  auto tacceleration =
193  (
194  IOobject
195  (
196  "MRFZoneList:acceleration",
197  U.mesh().time().timeName(),
198  U.mesh()
199  ),
200  U.mesh(),
201  dimensionedVector(U.dimensions()/dimTime, Zero)
202  );
203  auto& acceleration = tacceleration.ref();
204 
205  for (const auto& mrf: *this)
206  {
207  mrf.addCoriolis(U, acceleration);
208  }
209 
210  return tacceleration;
211 }
212 
213 
215 (
216  const volScalarField& rho,
217  const volVectorField& U
218 ) const
219 {
220  return rho*DDt(U);
221 }
222 
223 
225 {
226  for (const auto& mrf: *this)
227  {
228  mrf.makeRelative(U);
229  }
230 }
231 
232 
234 {
235  for (const auto& mrf: *this)
236  {
237  mrf.makeRelative(phi);
238  }
239 }
240 
241 
243 (
244  const tmp<surfaceScalarField>& tphi
245 ) const
246 {
247  if (size())
248  {
250  (
251  New
252  (
253  tphi,
254  "relative(" + tphi().name() + ')',
255  tphi().dimensions(),
256  true
257  )
258  );
259 
260  makeRelative(rphi.ref());
261 
262  tphi.clear();
263 
264  return rphi;
265  }
266  else
267  {
268  return tmp<surfaceScalarField>(tphi, true);
269  }
270 }
271 
272 
275 (
277 ) const
278 {
279  if (size())
280  {
281  tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
282 
283  for (const auto& mrf: *this)
284  {
285  mrf.makeRelative(rphi.ref());
286  }
287 
288  tphi.clear();
289 
290  return rphi;
291  }
292  else
293  {
294  return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
295  }
296 }
297 
298 
301 (
302  const tmp<Field<scalar>>& tphi,
303  const label patchi
304 ) const
305 {
306  if (size())
307  {
308  tmp<Field<scalar>> rphi(New(tphi, true));
309 
310  for (const auto& mrf: *this)
311  {
312  mrf.makeRelative(rphi.ref(), patchi);
313  }
314 
315  tphi.clear();
316 
317  return rphi;
318  }
319  else
320  {
321  return tmp<Field<scalar>>(tphi, true);
322  }
323 }
324 
325 
327 (
328  const surfaceScalarField& rho,
330 ) const
331 {
332  for (const auto& mrf: *this)
333  {
334  mrf.makeRelative(rho, phi);
335  }
336 }
337 
338 
340 {
341  for (const auto& mrf: *this)
342  {
343  mrf.makeAbsolute(U);
344  }
345 }
346 
347 
349 {
350  for (const auto& mrf: *this)
351  {
352  mrf.makeAbsolute(phi);
353  }
354 }
355 
356 
358 (
359  const tmp<surfaceScalarField>& tphi
360 ) const
361 {
362  if (size())
363  {
365  (
366  New
367  (
368  tphi,
369  "absolute(" + tphi().name() + ')',
370  tphi().dimensions(),
371  true
372  )
373  );
374 
375  makeAbsolute(rphi.ref());
376 
377  tphi.clear();
378 
379  return rphi;
380  }
381  else
382  {
383  return tmp<surfaceScalarField>(tphi, true);
384  }
385 }
386 
387 
389 (
390  const surfaceScalarField& rho,
392 ) const
393 {
394  for (const auto& mrf: *this)
395  {
396  mrf.makeAbsolute(rho, phi);
397  }
398 }
399 
400 
402 {
403  for (const auto& mrf: *this)
404  {
405  mrf.correctBoundaryVelocity(U);
406  }
407 }
408 
409 
411 (
412  const volVectorField& U,
414 ) const
415 {
417  (
418  relative(mesh_.Sf().boundaryField() & U.boundaryField())
419  );
420 
421  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
422 
423  forAll(mesh_.boundary(), patchi)
424  {
425  if (isA<fixedValueFvsPatchScalarField>(phibf[patchi]))
426  {
427  phibf[patchi] == Uf[patchi];
428  }
429  }
430 }
431 
432 
434 {
435  if (mesh_.topoChanging())
436  {
437  for (auto& mrf: *this)
438  {
439  mrf.update();
440  }
441  }
442 }
443 
444 
445 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
446 
447 Foam::Ostream& Foam::operator<<
448 (
449  Ostream& os,
450  const MRFZoneList& models
451 )
452 {
453  models.writeData(os);
454  return os;
455 }
456 
457 
458 // ************************************************************************* //
Foam::MRFZone
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:68
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:152
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
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:65
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:287
Foam::MRFZoneList::writeData
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:139
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:58
Foam::DynamicList< word >
Foam::MRFZoneList::reset
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:69
Uf
autoPtr< surfaceVectorField > Uf
Definition: createUfIfPresent.H:33
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:175
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
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:401
Foam::MRFZoneList::active
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:52
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::IOstream::good
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:233
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::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
resize
patchWriters resize(patchIds.size())
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::MRFZoneList::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:224
Foam::MRFZoneList::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:358
Foam::MRFZoneList::DDt
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:187
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:460
Foam::PtrList< MRFZone >
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dict
dictionary dict
Definition: searchingEngine.H:14
fixedValueFvsPatchFields.H
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:85
Foam::MRFZoneList::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:243
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
reset
meshPtr reset(new Foam::fvMesh(Foam::IOobject(regionName, runTime.timeName(), runTime, Foam::IOobject::MUST_READ), false))
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
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::nl
constexpr char nl
Definition: Ostream.H:404
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::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:339
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:411
Foam::MRFZoneList::update
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:433
Foam::tmp::New
static tmp< T > New(Args &&... args)
Construct tmp of T with forwarding arguments.
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:68
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::MRFZoneList::read
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:127
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::MRFZoneList::getFromName
const MRFZone & getFromName(const word &name) const
Return the MRF with a given name.
Definition: MRFZoneList.C:103
Foam::fvc::makeAbsolute
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:116
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::PtrListOps::names
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
Foam::fvc::DDt
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:47