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-2024 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"
31 #include "surfaceFields.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 Foam::MRFZoneList::MRFZoneList
37 (
38  const fvMesh& mesh,
39  const dictionary& dict
40 )
41 :
42  PtrList<MRFZone>(),
43  mesh_(mesh)
44 {
45  reset(dict);
46 
47  active(true);
48 }
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
53 bool Foam::MRFZoneList::active(const bool warn) const
54 {
55  bool a = false;
56  forAll(*this, i)
57  {
58  a = a || this->operator[](i).active();
59  }
60 
61  if (warn && this->size() && !a)
62  {
63  Info<< " No MRF zones active" << endl;
64  }
65 
66  return a;
67 }
68 
69 
70 void Foam::MRFZoneList::reset(const dictionary& dict)
71 {
72  label count = 0;
73  for (const entry& dEntry : dict)
74  {
75  if (dEntry.isDict())
76  {
77  ++count;
78  }
79  }
80 
81  this->resize(count);
82 
83  count = 0;
84  for (const entry& dEntry : dict)
85  {
86  if (dEntry.isDict())
87  {
88  const word& name = dEntry.keyword();
89  const dictionary& modelDict = dEntry.dict();
90 
91  Info<< " creating MRF zone: " << name << endl;
92 
93  this->set
94  (
95  count++,
96  new MRFZone(name, mesh_, modelDict)
97  );
98  }
99  }
100 }
101 
102 
104 (
105  const word& name
106 ) const
107 {
108  DynamicList<word> names;
109  for (const auto& mrf: *this)
110  {
111  if (mrf.name() == name)
112  {
113  return mrf;
114  }
115 
116  names.append(mrf.name());
117  }
118 
120  << "Unable to find MRFZone " << name
121  << ". Available zones are: " << names
122  << exit(FatalError);
123 
124  return first();
125 }
126 
127 
129 {
130  bool allOk = true;
131  for (auto& mrf: *this)
132  {
133  bool ok = mrf.read(dict.subDict(mrf.name()));
134  allOk = (allOk && ok);
135  }
136  return allOk;
137 }
138 
139 
140 bool Foam::MRFZoneList::writeData(Ostream& os) const
141 {
142  for (const auto& mrf: *this)
143  {
144  os << nl;
145  mrf.writeData(os);
146  }
147 
148  return os.good();
149 }
150 
151 
153 (
154  const volVectorField& U,
155  volVectorField& ddtU
156 ) const
157 {
158  for (const auto& mrf: *this)
159  {
160  mrf.addCoriolis(U, ddtU);
161  }
162 }
163 
164 
166 {
167  for (const auto& mrf: *this)
168  {
169  mrf.addCoriolis(UEqn);
170  }
171 }
172 
173 
175 (
176  const volScalarField& rho,
178 ) const
179 {
180  for (const auto& mrf: *this)
181  {
182  mrf.addCoriolis(rho, UEqn);
183  }
184 }
185 
186 
188 (
189  const volVectorField& U
190 ) const
191 {
192  auto tacceleration = volVectorField::New
193  (
194  IOobject::scopedName("MRFZoneList", "acceleration"),
196  U.mesh(),
197  dimensionedVector(U.dimensions()/dimTime, Zero)
198  );
199  auto& acceleration = tacceleration.ref();
200 
201  for (const auto& mrf: *this)
202  {
203  mrf.addCoriolis(U, acceleration);
204  }
205 
206  return tacceleration;
207 }
208 
209 
211 (
212  const volScalarField& rho,
214 ) const
215 {
216  return rho*DDt(U);
217 }
218 
219 
221 {
222  auto tphi = surfaceScalarField::New
223  (
224  "phiMRF",
226  mesh_,
228  );
229  auto& phi = tphi.ref();
230 
231  for (const auto& mrf : *this)
232  {
233  mrf.makeAbsolute(phi);
234  }
235 
236  return tphi;
237 }
238 
239 
241 {
242  for (const auto& mrf: *this)
243  {
244  mrf.makeRelative(U);
245  }
246 }
247 
248 
250 {
251  for (const auto& mrf: *this)
252  {
253  mrf.makeRelative(phi);
254  }
255 }
256 
257 
259 (
260  const tmp<surfaceScalarField>& tphi
261 ) const
262 {
263  if (size())
264  {
265  tmp<surfaceScalarField> rphi
266  (
267  New
268  (
269  tphi,
270  "relative(" + tphi().name() + ')',
271  tphi().dimensions(),
272  true
273  )
274  );
275 
276  makeRelative(rphi.ref());
277 
278  tphi.clear();
279 
280  return rphi;
281  }
283  return tmp<surfaceScalarField>(tphi, true);
284 }
285 
286 
289 (
291 ) const
292 {
293  if (size())
294  {
295  tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
296 
297  for (const auto& mrf: *this)
298  {
299  mrf.makeRelative(rphi.ref());
300  }
301 
302  tphi.clear();
303 
304  return rphi;
305  }
307  return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
308 }
309 
310 
313 (
314  const tmp<Field<scalar>>& tphi,
315  const label patchi
316 ) const
317 {
318  if (size())
319  {
320  tmp<Field<scalar>> rphi(New(tphi, true));
321 
322  for (const auto& mrf: *this)
323  {
324  mrf.makeRelative(rphi.ref(), patchi);
325  }
326 
327  tphi.clear();
328 
329  return rphi;
330  }
331 
332  return tmp<Field<scalar>>(tphi, true);
333 }
334 
335 
337 (
338  const surfaceScalarField& rho,
340 ) const
341 {
342  for (const auto& mrf: *this)
343  {
344  mrf.makeRelative(rho, phi);
345  }
346 }
347 
348 
350 {
351  for (const auto& mrf: *this)
352  {
353  mrf.makeAbsolute(U);
354  }
355 }
356 
357 
359 {
360  for (const auto& mrf: *this)
361  {
362  mrf.makeAbsolute(phi);
363  }
364 }
365 
366 
368 (
369  const tmp<surfaceScalarField>& tphi
370 ) const
371 {
372  if (size())
373  {
374  tmp<surfaceScalarField> rphi
375  (
376  New
377  (
378  tphi,
379  "absolute(" + tphi().name() + ')',
380  tphi().dimensions(),
381  true
382  )
383  );
384 
385  makeAbsolute(rphi.ref());
386 
387  tphi.clear();
388 
389  return rphi;
390  }
391 
392  return tmp<surfaceScalarField>(tphi, true);
393 }
394 
395 
397 (
398  const surfaceScalarField& rho,
400 ) const
401 {
402  for (const auto& mrf: *this)
403  {
404  mrf.makeAbsolute(rho, phi);
405  }
406 }
407 
408 
410 {
411  for (const auto& mrf: *this)
412  {
413  mrf.correctBoundaryVelocity(U);
414  }
415 }
416 
417 
419 (
420  const volVectorField& U,
422 ) const
423 {
424  FieldField<fvsPatchField, scalar> Uf
425  (
426  relative(mesh_.Sf().boundaryField() & U.boundaryField())
427  );
428 
429  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
430 
431  forAll(mesh_.boundary(), patchi)
432  {
433  if (isA<fixedValueFvsPatchScalarField>(phibf[patchi]))
434  {
435  phibf[patchi] == Uf[patchi];
436  }
437  }
438 }
439 
440 
442 {
443  if (mesh_.topoChanging())
444  {
445  for (auto& mrf: *this)
446  {
447  mrf.update();
448  }
449  }
450 }
451 
452 
453 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
454 
455 Foam::Ostream& Foam::operator<<
456 (
457  Ostream& os,
458  const MRFZoneList& models
459 )
460 {
461  models.writeData(os);
462  return os;
463 }
464 
465 
466 // ************************************************************************* //
Foam::surfaceFields.
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:342
dictionary dict
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:434
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:63
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:40
List< word > names(const UPtrList< T > &list, const UnaryMatchPredicate &matcher)
List of names generated by calling name() for each list item and filtered for matches.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:146
patchWriters resize(patchIds.size())
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:521
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:64
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
Type of boundary fields.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
const MRFZone & getFromName(const word &name) const
Return the MRF with a given name.
Definition: MRFZoneList.C:97
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:46
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:402
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:412
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:361
autoPtr< surfaceVectorField > Uf
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:181
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
tmp< surfaceScalarField > phi() const
Return the MRF absolute flux.
Definition: MRFZoneList.C:213
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:133
OBJstream os(runTime.globalPath()/outputName)
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:252
U
Definition: pEqn.H:72
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
bool good() const noexcept
True if next operation might succeed.
Definition: IOstream.H:281
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:109
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:70
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:233
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:121
Do not request registration (bool: false)
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:57
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition: fvcMeshPhi.C:148
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const dimensionSet dimVelocity