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 | Copyright (C) 2012-2016 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 
32 Foam::MRFZoneList::MRFZoneList
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  (
171  new volVectorField
172  (
173  IOobject
174  (
175  "MRFZoneList:acceleration",
176  U.mesh().time().timeName(),
177  U.mesh()
178  ),
179  U.mesh(),
181  )
182  );
183  volVectorField& acceleration = tacceleration.ref();
184 
185  forAll(*this, i)
186  {
187  operator[](i).addCoriolis(U, acceleration);
188  }
189 
190  return tacceleration;
191 }
192 
193 
195 (
196  const volScalarField& rho,
197  const volVectorField& U
198 ) const
199 {
200  return rho*DDt(U);
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  const volVectorField& U,
393  surfaceScalarField& phi
394 ) const
395 {
397  (
399  );
400 
401 
402  surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
403 
405  {
406  if
407  (
408  isA<fixedValueFvsPatchScalarField>(phibf[patchi])
409  )
410  {
411  phibf[patchi] == Uf[patchi];
412  }
413  }
414 }
415 
416 
417 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
418 
419 Foam::Ostream& Foam::operator<<
420 (
421  Ostream& os,
422  const MRFZoneList& models
423 )
424 {
425  models.writeData(os);
426  return os;
427 }
428 
429 
430 // ************************************************************************* //
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:49
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:319
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
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
const surfaceVectorField & Sf() const
Return cell face area vectors.
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:65
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:633
const word & name() const
Return const access to the MRF region name.
Definition: MRFZoneI.H:26
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
Uf
Definition: pEqn.H:67
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:223
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:381
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZoneList.H:75
Generic field type.
Definition: FieldField.H:51
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:391
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
const Boundary & boundaryField() const
Return const-reference to the boundary field.
A class for handling words, derived from string.
Definition: word.H:59
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:71
static const zero Zero
Definition: zero.H:91
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:55
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
const dimensionSet & dimensions() const
Return dimensions.
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
static const char nl
Definition: Ostream.H:262
List container for MRF zomes.
Definition: MRFZoneList.H:55
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< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:165
const T & operator[](const label) const
Return element const reference.
Definition: UPtrListI.H:103
label patchi
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
const Mesh & mesh() const
Return mesh.
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:338
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
A class for managing temporary objects.
Definition: PtrList.H:54
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
const fvBoundaryMesh & boundary() const
Return reference to boundary mesh.
Definition: fvMesh.C:545
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:104
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:130
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:117
label size() const
Return the number of elements in the UPtrList.
Definition: UPtrListI.H:29