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-2018 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  active(true);
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
48 
50 {}
51 
52 
53 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
54 
55 bool Foam::MRFZoneList::active(const bool warn) const
56 {
57  bool a = false;
58  forAll(*this, i)
59  {
60  a = a || this->operator[](i).active();
61  }
62 
63  if (warn && this->size() && !a)
64  {
65  Info<< " No MRF zones active" << endl;
66  }
67 
68  return a;
69 }
70 
71 
73 {
74  label count = 0;
75  forAllConstIter(dictionary, dict, iter)
76  {
77  if (iter().isDict())
78  {
79  count++;
80  }
81  }
82 
83  this->setSize(count);
84  label i = 0;
85  forAllConstIter(dictionary, dict, iter)
86  {
87  if (iter().isDict())
88  {
89  const word& name = iter().keyword();
90  const dictionary& modelDict = iter().dict();
91 
92  Info<< " creating MRF zone: " << name << endl;
93 
94  this->set
95  (
96  i++,
97  new MRFZone(name, mesh_, modelDict)
98  );
99  }
100  }
101 }
102 
103 
105 {
106  bool allOk = true;
107  forAll(*this, i)
108  {
109  MRFZone& pm = this->operator[](i);
110  bool ok = pm.read(dict.subDict(pm.name()));
111  allOk = (allOk && ok);
112  }
113  return allOk;
114 }
115 
116 
118 {
119  forAll(*this, i)
120  {
121  os << nl;
122  this->operator[](i).writeData(os);
123  }
124 
125  return os.good();
126 }
127 
128 
130 (
131  const volVectorField& U,
132  volVectorField& ddtU
133 ) const
134 {
135  forAll(*this, i)
136  {
137  operator[](i).addCoriolis(U, ddtU);
138  }
139 }
140 
141 
143 {
144  forAll(*this, i)
145  {
146  operator[](i).addCoriolis(UEqn);
147  }
148 }
149 
150 
152 (
153  const volScalarField& rho,
154  fvVectorMatrix& UEqn
155 ) const
156 {
157  forAll(*this, i)
158  {
159  operator[](i).addCoriolis(rho, UEqn);
160  }
161 }
162 
163 
165 (
166  const volVectorField& U
167 ) const
168 {
169  tmp<volVectorField> tacceleration
170  (
172  (
173  "MRFZoneList:acceleration",
174  U.mesh(),
176  )
177  );
178  volVectorField& acceleration = tacceleration.ref();
179 
180  forAll(*this, i)
181  {
182  operator[](i).addCoriolis(U, acceleration);
183  }
184 
185  return tacceleration;
186 }
187 
188 
190 (
191  const volScalarField& rho,
192  const volVectorField& U
193 ) const
194 {
195  return rho*DDt(U);
196 }
197 
198 
200 {
201  forAll(*this, i)
202  {
203  operator[](i).makeRelative(U);
204  }
205 }
206 
207 
209 {
210  forAll(*this, i)
211  {
212  operator[](i).makeRelative(phi);
213  }
214 }
215 
216 
218 (
219  const tmp<surfaceScalarField>& tphi
220 ) const
221 {
222  if (size())
223  {
225  (
226  New
227  (
228  tphi,
229  "relative(" + tphi().name() + ')',
230  tphi().dimensions(),
231  true
232  )
233  );
234 
235  makeRelative(rphi.ref());
236 
237  tphi.clear();
238 
239  return rphi;
240  }
241  else
242  {
243  return tmp<surfaceScalarField>(tphi, true);
244  }
245 }
246 
247 
250 (
252 ) const
253 {
254  if (size())
255  {
256  tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
257 
258  forAll(*this, i)
259  {
260  operator[](i).makeRelative(rphi.ref());
261  }
262 
263  tphi.clear();
264 
265  return rphi;
266  }
267  else
268  {
269  return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
270  }
271 }
272 
273 
276 (
277  const tmp<Field<scalar>>& tphi,
278  const label patchi
279 ) const
280 {
281  if (size())
282  {
283  tmp<Field<scalar>> rphi(New(tphi, true));
284 
285  forAll(*this, i)
286  {
287  operator[](i).makeRelative(rphi.ref(), patchi);
288  }
289 
290  tphi.clear();
291 
292  return rphi;
293  }
294  else
295  {
296  return tmp<Field<scalar>>(tphi, true);
297  }
298 }
299 
300 
302 (
303  const surfaceScalarField& rho,
304  surfaceScalarField& phi
305 ) const
306 {
307  forAll(*this, i)
308  {
309  operator[](i).makeRelative(rho, phi);
310  }
311 }
312 
313 
315 {
316  forAll(*this, i)
317  {
318  operator[](i).makeAbsolute(U);
319  }
320 }
321 
322 
324 {
325  forAll(*this, i)
326  {
327  operator[](i).makeAbsolute(phi);
328  }
329 }
330 
331 
333 (
334  const tmp<surfaceScalarField>& tphi
335 ) const
336 {
337  if (size())
338  {
340  (
341  New
342  (
343  tphi,
344  "absolute(" + tphi().name() + ')',
345  tphi().dimensions(),
346  true
347  )
348  );
349 
350  makeAbsolute(rphi.ref());
351 
352  tphi.clear();
353 
354  return rphi;
355  }
356  else
357  {
358  return tmp<surfaceScalarField>(tphi, true);
359  }
360 }
361 
362 
364 (
365  const surfaceScalarField& rho,
366  surfaceScalarField& phi
367 ) const
368 {
369  forAll(*this, i)
370  {
371  operator[](i).makeAbsolute(rho, phi);
372  }
373 }
374 
375 
377 {
378  forAll(*this, i)
379  {
380  operator[](i).correctBoundaryVelocity(U);
381  }
382 }
383 
384 
386 (
387  const volVectorField& U,
388  surfaceScalarField& phi
389 ) const
390 {
392  (
394  );
395 
396 
397  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
398 
400  {
401  if
402  (
403  isA<fixedValueFvsPatchScalarField>(phibf[patchi])
404  )
405  {
406  phibf[patchi] == Uf[patchi];
407  }
408  }
409 }
410 
411 
413 {
414  if (mesh_.topoChanging())
415  {
416  forAll(*this, i)
417  {
418  operator[](i).update();
419  }
420  }
421 }
422 
423 
424 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
425 
426 Foam::Ostream& Foam::operator<<
427 (
428  Ostream& os,
429  const MRFZoneList& models
430 )
431 {
432  models.writeData(os);
433  return os;
434 }
435 
436 
437 // ************************************************************************* //
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:49
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:591
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:314
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
const surfaceVectorField & Sf() const
Return cell face area vectors.
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZoneList.C:412
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:72
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:230
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:130
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
const Boundary & boundaryField() const
Return const-reference to the boundary field.
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:174
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
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
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:934
const dimensionSet & dimensions() const
Return dimensions.
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:55
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:376
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
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:386
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:333
static const zero Zero
Definition: zero.H:97
autoPtr< surfaceVectorField > Uf
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:165
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:117
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:218
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
bool topoChanging() const
Is mesh topology changing.
Definition: polyMesh.H:526
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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:199
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:104
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:540