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-2021 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 {
150  tmp<volVectorField> tacceleration
151  (
153  (
154  "MRFZoneList:acceleration",
155  U.mesh(),
157  )
158  );
159  volVectorField& acceleration = tacceleration.ref();
160 
161  forAll(*this, i)
162  {
163  operator[](i).addCoriolis(U, acceleration);
164  }
165 
166  return tacceleration;
167 }
168 
169 
171 (
172  const volScalarField& rho,
173  const volVectorField& U
174 ) const
175 {
176  return rho*DDt(U);
177 }
178 
179 
181 {
182  forAll(*this, i)
183  {
184  operator[](i).makeRelative(U);
185  }
186 }
187 
188 
190 {
191  forAll(*this, i)
192  {
193  operator[](i).makeRelative(phi);
194  }
195 }
196 
197 
199 (
200  const tmp<surfaceScalarField>& tphi
201 ) const
202 {
203  if (size())
204  {
206  (
207  New
208  (
209  tphi,
210  "relative(" + tphi().name() + ')',
211  tphi().dimensions(),
212  true
213  )
214  );
215 
216  makeRelative(rphi.ref());
217 
218  tphi.clear();
219 
220  return rphi;
221  }
222  else
223  {
224  return tmp<surfaceScalarField>(tphi, true);
225  }
226 }
227 
228 
231 (
233 ) const
234 {
235  if (size())
236  {
237  tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
238 
239  forAll(*this, i)
240  {
241  operator[](i).makeRelative(rphi.ref());
242  }
243 
244  tphi.clear();
245 
246  return rphi;
247  }
248  else
249  {
250  return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
251  }
252 }
253 
254 
257 (
258  const tmp<Field<scalar>>& tphi,
259  const label patchi
260 ) const
261 {
262  if (size())
263  {
264  tmp<Field<scalar>> rphi(New(tphi, true));
265 
266  forAll(*this, i)
267  {
268  operator[](i).makeRelative(rphi.ref(), patchi);
269  }
270 
271  tphi.clear();
272 
273  return rphi;
274  }
275  else
276  {
277  return tmp<Field<scalar>>(tphi, true);
278  }
279 }
280 
281 
283 (
284  const surfaceScalarField& rho,
285  surfaceScalarField& phi
286 ) const
287 {
288  forAll(*this, i)
289  {
290  operator[](i).makeRelative(rho, phi);
291  }
292 }
293 
294 
296 {
297  forAll(*this, i)
298  {
299  operator[](i).makeAbsolute(U);
300  }
301 }
302 
303 
305 {
306  forAll(*this, i)
307  {
308  operator[](i).makeAbsolute(phi);
309  }
310 }
311 
312 
314 (
315  const tmp<surfaceScalarField>& tphi
316 ) const
317 {
318  if (size())
319  {
321  (
322  New
323  (
324  tphi,
325  "absolute(" + tphi().name() + ')',
326  tphi().dimensions(),
327  true
328  )
329  );
330 
331  makeAbsolute(rphi.ref());
332 
333  tphi.clear();
334 
335  return rphi;
336  }
337  else
338  {
339  return tmp<surfaceScalarField>(tphi, true);
340  }
341 }
342 
343 
345 (
346  const surfaceScalarField& rho,
347  surfaceScalarField& phi
348 ) const
349 {
350  forAll(*this, i)
351  {
352  operator[](i).makeAbsolute(rho, phi);
353  }
354 }
355 
356 
358 {
359  forAll(*this, i)
360  {
361  operator[](i).correctBoundaryVelocity(U);
362  }
363 }
364 
365 
367 {
368  if (mesh_.topoChanging())
369  {
370  forAll(*this, i)
371  {
372  operator[](i).update();
373  }
374  }
375 }
376 
377 
378 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
379 
380 Foam::Ostream& Foam::operator<<
381 (
382  Ostream& os,
383  const MRFZoneList& models
384 )
385 {
386  models.writeData(os);
387  return os;
388 }
389 
390 
391 // ************************************************************************* //
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:47
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:577
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:295
#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:366
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:53
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
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:111
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:65
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
Generic field type.
Definition: FieldField.H:51
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:982
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:357
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
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:314
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame 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
Internal & ref()
Return a reference to the dimensioned internal field.
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:199
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:78
bool topoChanging() const
Is mesh topology changing.
Definition: polyMesh.H:526
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:180
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:85