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-2015 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;
76  {
77  if (iter().isDict())
78  {
79  count++;
80  }
81  }
82 
83  this->setSize(count);
84  label i = 0;
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,
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(),
180  dimensionedVector("0", U.dimensions()/dimTime, vector::zero)
181  )
182  );
183  volVectorField& acceleration = tacceleration();
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 (
225 ) const
226 {
227  tmp<surfaceScalarField> rphi(phi.ptr());
228  makeRelative(rphi());
229  return rphi;
230 }
231 
232 
235 (
237 ) const
238 {
240 
241  forAll(*this, i)
242  {
243  operator[](i).makeRelative(rphi());
244  }
245 
246  return rphi;
247 }
248 
249 
251 (
252  const surfaceScalarField& rho,
254 ) const
255 {
256  forAll(*this, i)
257  {
258  operator[](i).makeRelative(rho, phi);
259  }
260 }
261 
262 
264 {
265  forAll(*this, i)
266  {
267  operator[](i).makeAbsolute(U);
268  }
269 }
270 
271 
273 {
274  forAll(*this, i)
275  {
276  operator[](i).makeAbsolute(phi);
277  }
278 }
279 
280 
282 (
284 ) const
285 {
286  tmp<surfaceScalarField> rphi(phi.ptr());
287  makeAbsolute(rphi());
288  return rphi;
289 }
290 
291 
293 (
294  const surfaceScalarField& rho,
296 ) const
297 {
298  forAll(*this, i)
299  {
300  operator[](i).makeAbsolute(rho, phi);
301  }
302 }
303 
304 
306 {
307  forAll(*this, i)
308  {
309  operator[](i).correctBoundaryVelocity(U);
310  }
311 }
312 
313 
315 (
316  const volVectorField& U,
318 ) const
319 {
321  (
322  relative(mesh_.Sf().boundaryField() & U.boundaryField())
323  );
324 
325  forAll(mesh_.boundary(), patchi)
326  {
327  if
328  (
329  isA<fixedValueFvsPatchScalarField>(phi.boundaryField()[patchi])
330  )
331  {
332  phi.boundaryField()[patchi] == phibf[patchi];
333  }
334  }
335 }
336 
337 
338 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
339 
340 Foam::Ostream& Foam::operator<<
341 (
342  Ostream& os,
343  const MRFZoneList& models
344 )
345 {
346  models.writeData(os);
347  return os;
348 }
349 
350 
351 // ************************************************************************* //
Foam::MRFZone
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition: MRFZone.H:65
volFields.H
Foam::MRFZoneList::addAcceleration
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
Definition: MRFZoneList.C:130
setSize
points setSize(newPointi)
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::MRFZoneList::MRFZoneList
MRFZoneList(const MRFZoneList &)
Disallow default bitwise copy construct.
Foam::MRFZoneList::writeData
bool writeData(Ostream &os) const
Write data to Ostream.
Definition: MRFZoneList.C:117
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::MRFZoneList
List container for MRF zomes.
Definition: MRFZoneList.H:55
Foam::MRFZoneList::reset
void reset(const dictionary &dict)
Reset the source list.
Definition: MRFZoneList.C:72
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::GeometricField::boundaryField
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
Definition: GeometricField.C:735
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::MRFZoneList::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZoneList.C:305
Foam::MRFZoneList::active
bool active(const bool warn=false) const
Return active status.
Definition: MRFZoneList.C:55
U
U
Definition: pEqn.H:46
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Foam::fvc::makeRelative
void makeRelative(surfaceScalarField &phi, const volVectorField &U)
Make the given flux relative.
Definition: fvcMeshPhi.C:74
Foam::MRFZone::read
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:564
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::label
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
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::MRFZoneList::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZoneList.C:204
Foam::MRFZoneList::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
Definition: MRFZoneList.C:282
Foam::MRFZoneList::DDt
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
Definition: MRFZoneList.C:165
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::PtrList< MRFZone >
dict
dictionary dict
Definition: searchingEngine.H:14
fixedValueFvsPatchFields.H
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::MRFZoneList::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
Definition: MRFZoneList.C:223
rho
rho
Definition: pEqn.H:3
Foam::MRFZoneList::makeAbsolute
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZoneList.C:263
MRFZoneList.H
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::MRFZoneList::correctBoundaryFlux
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
Definition: MRFZoneList.C:315
Foam::fvc::relative
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given absolute flux in relative form.
Definition: fvcMeshPhi.C:152
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::MRFZoneList::read
bool read(const dictionary &dict)
Read dictionary.
Definition: MRFZoneList.C:104
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::MRFZoneList::~MRFZoneList
~MRFZoneList()
Destructor.
Definition: MRFZoneList.C:49
Foam::fvc::makeAbsolute
void makeAbsolute(surfaceScalarField &phi, const volVectorField &U)
Make the given flux absolute.
Definition: fvcMeshPhi.C:113
Foam::IOstream::good
bool good() const
Return true if next operation might succeed.
Definition: IOstream.H:333
Foam::MRFZone::name
const word & name() const
Return const access to the MRF region name.
Definition: MRFZoneI.H:26
Foam::dictionary::subDict
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:631
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::fvc::DDt
tmp< GeometricField< Type, fvPatchField, volMesh > > DDt(const surfaceScalarField &phi, const GeometricField< Type, fvPatchField, volMesh > &psi)
Definition: fvcDDt.C:45