MRFZoneList.C
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) 2012-2022 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 "MRFZoneList.H"
27 #include "volFields.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
33 (
34  const fvMesh& mesh,
35  const dictionary& dict
36 )
37 :
39  mesh_(mesh)
40 {
41  reset(dict);
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
46 
48 {}
49 
50 
51 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52 
54 {
55  label count = 0;
56  forAllConstIter(dictionary, dict, iter)
57  {
58  if (iter().isDict())
59  {
60  count++;
61  }
62  }
63 
64  this->setSize(count);
65  label i = 0;
66  forAllConstIter(dictionary, dict, iter)
67  {
68  if (iter().isDict())
69  {
70  const word& name = iter().keyword();
71  const dictionary& modelDict = iter().dict();
72 
73  Info<< " creating MRF zone: " << name << endl;
74 
75  this->set
76  (
77  i++,
78  new MRFZone(name, mesh_, modelDict)
79  );
80  }
81  }
82 }
83 
84 
86 {
87  bool allOk = true;
88  forAll(*this, i)
89  {
90  MRFZone& pm = this->operator[](i);
91  bool ok = pm.read(dict.subDict(pm.name()));
92  allOk = (allOk && ok);
93  }
94  return allOk;
95 }
96 
97 
99 {
100  forAll(*this, i)
101  {
102  os << nl;
103  this->operator[](i).writeData(os);
104  }
105 
106  return os.good();
107 }
108 
109 
111 (
112  const volVectorField& U,
113  volVectorField& ddtU
114 ) const
115 {
116  forAll(*this, i)
117  {
118  operator[](i).addCoriolis(U, ddtU);
119  }
120 }
121 
122 
124 {
125  forAll(*this, i)
126  {
127  operator[](i).addCoriolis(UEqn);
128  }
129 }
130 
131 
133 (
134  const volScalarField& rho,
135  fvVectorMatrix& UEqn
136 ) const
137 {
138  forAll(*this, i)
139  {
140  operator[](i).addCoriolis(rho, UEqn);
141  }
142 }
143 
144 
146 (
147  const volVectorField& U
148 ) const
149 {
151  (
153  (
154  "MRFZoneList:DDt",
155  U.mesh(),
157  )
158  );
159  volVectorField& DDt = tDDt.ref();
160 
161  forAll(*this, i)
162  {
163  operator[](i).addCoriolis(U, DDt);
164  }
165 
166  return tDDt;
167 }
168 
169 
171 (
172  const volScalarField& rho,
173  const volVectorField& U
174 ) const
175 {
176  return rho*DDt(U);
177 }
178 
179 
182 {
183  tmp<volVectorField> tcentrifugalAcceleration
184  (
186  (
187  "MRFZoneList:centrifugalAcceleration",
188  mesh_,
190  )
191  );
192  volVectorField& centrifugalAcceleration = tcentrifugalAcceleration.ref();
193 
194  forAll(*this, i)
195  {
196  operator[](i).addCentrifugalAcceleration(centrifugalAcceleration);
197  }
198 
199  return tcentrifugalAcceleration;
200 }
201 
202 
203 
205 {
206  forAll(*this, i)
207  {
208  operator[](i).makeRelative(U);
209  }
210 }
211 
212 
214 {
215  forAll(*this, i)
216  {
217  operator[](i).makeRelative(phi);
218  }
219 }
220 
221 
223 (
224  const tmp<surfaceScalarField>& tphi
225 ) const
226 {
227  if (size())
228  {
230  (
231  New
232  (
233  tphi,
234  "relative(" + tphi().name() + ')',
235  tphi().dimensions(),
236  true
237  )
238  );
239 
240  makeRelative(rphi.ref());
241 
242  tphi.clear();
243 
244  return rphi;
245  }
246  else
247  {
248  return tmp<surfaceScalarField>(tphi, true);
249  }
250 }
251 
252 
255 (
257 ) const
258 {
259  if (size())
260  {
261  tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
262 
263  forAll(*this, i)
264  {
265  operator[](i).makeRelative(rphi.ref());
266  }
267 
268  tphi.clear();
269 
270  return rphi;
271  }
272  else
273  {
274  return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
275  }
276 }
277 
278 
281 (
282  const tmp<Field<scalar>>& tphi,
283  const label patchi
284 ) const
285 {
286  if (size())
287  {
288  tmp<Field<scalar>> rphi(New(tphi, true));
289 
290  forAll(*this, i)
291  {
292  operator[](i).makeRelative(rphi.ref(), patchi);
293  }
294 
295  tphi.clear();
296 
297  return rphi;
298  }
299  else
300  {
301  return tmp<Field<scalar>>(tphi, true);
302  }
303 }
304 
305 
307 (
308  const surfaceScalarField& rho,
309  surfaceScalarField& phi
310 ) const
311 {
312  forAll(*this, i)
313  {
314  operator[](i).makeRelative(rho, phi);
315  }
316 }
317 
318 
320 {
321  forAll(*this, i)
322  {
323  operator[](i).makeAbsolute(U);
324  }
325 }
326 
327 
329 {
330  forAll(*this, i)
331  {
332  operator[](i).makeAbsolute(phi);
333  }
334 }
335 
336 
338 (
339  const tmp<surfaceScalarField>& tphi
340 ) const
341 {
342  if (size())
343  {
345  (
346  New
347  (
348  tphi,
349  "absolute(" + tphi().name() + ')',
350  tphi().dimensions(),
351  true
352  )
353  );
354 
355  makeAbsolute(rphi.ref());
356 
357  tphi.clear();
358 
359  return rphi;
360  }
361  else
362  {
363  return tmp<surfaceScalarField>(tphi, true);
364  }
365 }
366 
367 
369 (
370  const surfaceScalarField& rho,
371  surfaceScalarField& phi
372 ) const
373 {
374  forAll(*this, i)
375  {
376  operator[](i).makeAbsolute(rho, phi);
377  }
378 }
379 
380 
382 {
383  forAll(*this, i)
384  {
385  operator[](i).correctBoundaryVelocity(U);
386  }
387 }
388 
389 
391 {
392  if (mesh_.topoChanged())
393  {
394  forAll(*this, i)
395  {
396  operator[](i).update();
397  }
398  }
399 }
400 
401 
402 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
403 
404 Foam::Ostream& Foam::operator<<
405 (
406  Ostream& os,
407  const MRFZoneList& models
408 )
409 {
410  models.writeData(os);
411  return os;
412 }
413 
414 
415 // ************************************************************************* //
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:47
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:621
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:319
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:390
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:53
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:62
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, const Internal &, const PtrList< fvPatchField< vector >> &)
Return a temporary field constructed from name,.
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:96
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:330
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZoneList.H:64
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:1002
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis acceleration.
Definition: MRFZoneList.C:111
const dimensionSet dimTime
const dimensionSet & dimensions() const
Return dimensions.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:381
tmp< volVectorField > centrifugalAcceleration() const
Return the centrifugal acceleration.
Definition: MRFZoneList.C:181
const dimensionSet dimAcceleration
A class for handling words, derived from string.
Definition: word.H:59
const word & name() const
Return const access to the MRF region name.
Definition: MRFZoneI.H:26
bool topoChanged() const
Is mesh topology changing.
Definition: polyMesh.H:531
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:338
static const zero Zero
Definition: zero.H:97
tmp< volVectorField > DDt(const volVectorField &U) const
Return the Coriolis acceleration.
Definition: MRFZoneList.C:146
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
static const char nl
Definition: Ostream.H:260
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:98
List container for MRF zomes.
Definition: MRFZoneList.H:55
const Mesh & mesh() const
Return mesh.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:223
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:95
MRFZoneList(const fvMesh &mesh, const dictionary &dict)
Definition: MRFZoneList.C:33
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:53
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:204
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:85