Loading...
Searching...
No Matches
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 | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2012-2017 OpenFOAM Foundation
9 Copyright (C) 2021-2024 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "MRFZoneList.H"
30#include "volFields.H"
31#include "surfaceFields.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36Foam::MRFZoneList::MRFZoneList
37(
38 const fvMesh& mesh,
39 const dictionary& dict
40)
41:
43 mesh_(mesh)
44{
45 reset(dict);
47 active(true);
48}
49
50
51// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
52
53bool Foam::MRFZoneList::active(const bool warn) const
54{
55 bool a = false;
56 forAll(*this, i)
57 {
58 a = a || this->operator[](i).active();
59 }
60
61 if (warn && this->size() && !a)
62 {
63 Info<< " No MRF zones active" << endl;
64 }
65
66 return a;
67}
68
69
71{
72 label count = 0;
73 for (const entry& dEntry : dict)
74 {
75 if (dEntry.isDict())
76 {
77 ++count;
78 }
79 }
80
81 this->resize(count);
82
83 count = 0;
84 for (const entry& dEntry : dict)
85 {
86 if (dEntry.isDict())
87 {
88 const word& name = dEntry.keyword();
89 const dictionary& modelDict = dEntry.dict();
90
91 Info<< " creating MRF zone: " << name << endl;
92
93 this->set
94 (
95 count++,
96 new MRFZone(name, mesh_, modelDict)
97 );
98 }
99 }
100}
101
102
103const Foam::MRFZone& Foam::MRFZoneList::getFromName
104(
105 const word& name
106) const
107{
109 for (const auto& mrf: *this)
110 {
111 if (mrf.name() == name)
112 {
113 return mrf;
114 }
115
116 names.append(mrf.name());
117 }
118
120 << "Unable to find MRFZone " << name
121 << ". Available zones are: " << names
122 << exit(FatalError);
123
124 return first();
125}
126
127
129{
130 bool allOk = true;
131 for (auto& mrf: *this)
132 {
133 bool ok = mrf.read(dict.subDict(mrf.name()));
134 allOk = (allOk && ok);
135 }
136 return allOk;
137}
138
139
140bool Foam::MRFZoneList::writeData(Ostream& os) const
141{
142 for (const auto& mrf: *this)
143 {
144 os << nl;
145 mrf.writeData(os);
146 }
147
148 return os.good();
149}
150
151
153(
154 const volVectorField& U,
155 volVectorField& ddtU
156) const
157{
158 for (const auto& mrf: *this)
159 {
160 mrf.addCoriolis(U, ddtU);
161 }
162}
163
164
166{
167 for (const auto& mrf: *this)
168 {
169 mrf.addCoriolis(UEqn);
170 }
171}
172
173
175(
176 const volScalarField& rho,
178) const
179{
180 for (const auto& mrf: *this)
181 {
182 mrf.addCoriolis(rho, UEqn);
183 }
184}
185
186
188(
189 const volVectorField& U
190) const
191{
192 auto tacceleration = volVectorField::New
193 (
194 IOobject::scopedName("MRFZoneList", "acceleration"),
196 U.mesh(),
197 dimensionedVector(U.dimensions()/dimTime, Zero)
198 );
199 auto& acceleration = tacceleration.ref();
200
201 for (const auto& mrf: *this)
202 {
203 mrf.addCoriolis(U, acceleration);
204 }
205
206 return tacceleration;
207}
208
209
211(
212 const volScalarField& rho,
214) const
215{
216 return rho*DDt(U);
217}
218
219
221{
222 auto tphi = surfaceScalarField::New
223 (
224 "phiMRF",
226 mesh_,
228 );
229 auto& phi = tphi.ref();
230
231 for (const auto& mrf : *this)
232 {
233 mrf.makeAbsolute(phi);
234 }
235
236 return tphi;
237}
238
239
241{
242 for (const auto& mrf: *this)
243 {
244 mrf.makeRelative(U);
245 }
246}
247
248
250{
251 for (const auto& mrf: *this)
252 {
253 mrf.makeRelative(phi);
254 }
255}
256
257
259(
260 const tmp<surfaceScalarField>& tphi
261) const
262{
263 if (size())
264 {
265 tmp<surfaceScalarField> rphi
266 (
267 New
268 (
269 tphi,
270 "relative(" + tphi().name() + ')',
271 tphi().dimensions(),
272 true
273 )
274 );
275
276 makeRelative(rphi.ref());
277
278 tphi.clear();
279
280 return rphi;
282
283 return tmp<surfaceScalarField>(tphi, true);
284}
285
286
289(
290 const tmp<FieldField<fvsPatchField, scalar>>& tphi
291) const
292{
293 if (size())
294 {
295 tmp<FieldField<fvsPatchField, scalar>> rphi(New(tphi, true));
296
297 for (const auto& mrf: *this)
298 {
299 mrf.makeRelative(rphi.ref());
300 }
301
302 tphi.clear();
303
304 return rphi;
306
307 return tmp<FieldField<fvsPatchField, scalar>>(tphi, true);
308}
309
310
313(
314 const tmp<Field<scalar>>& tphi,
315 const label patchi
316) const
317{
318 if (size())
319 {
320 tmp<Field<scalar>> rphi(New(tphi, true));
321
322 for (const auto& mrf: *this)
323 {
324 mrf.makeRelative(rphi.ref(), patchi);
325 }
326
327 tphi.clear();
328
329 return rphi;
330 }
331
332 return tmp<Field<scalar>>(tphi, true);
333}
334
335
337(
338 const surfaceScalarField& rho,
340) const
341{
342 for (const auto& mrf: *this)
343 {
344 mrf.makeRelative(rho, phi);
345 }
346}
347
348
350{
351 for (const auto& mrf: *this)
352 {
353 mrf.makeAbsolute(U);
354 }
355}
356
357
359{
360 for (const auto& mrf: *this)
361 {
362 mrf.makeAbsolute(phi);
363 }
364}
365
366
368(
369 const tmp<surfaceScalarField>& tphi
370) const
371{
372 if (size())
373 {
374 tmp<surfaceScalarField> rphi
375 (
376 New
377 (
378 tphi,
379 "absolute(" + tphi().name() + ')',
380 tphi().dimensions(),
381 true
382 )
383 );
384
385 makeAbsolute(rphi.ref());
386
387 tphi.clear();
388
389 return rphi;
390 }
391
392 return tmp<surfaceScalarField>(tphi, true);
393}
394
395
397(
398 const surfaceScalarField& rho,
400) const
401{
402 for (const auto& mrf: *this)
403 {
404 mrf.makeAbsolute(rho, phi);
405 }
406}
407
408
410{
411 for (const auto& mrf: *this)
412 {
413 mrf.correctBoundaryVelocity(U);
414 }
415}
416
417
419(
420 const volVectorField& U,
422) const
423{
424 FieldField<fvsPatchField, scalar> Uf
425 (
426 relative(mesh_.Sf().boundaryField() & U.boundaryField())
427 );
428
429 surfaceScalarField::Boundary& phibf = phi.boundaryFieldRef();
430
431 forAll(mesh_.boundary(), patchi)
432 {
433 if (isA<fixedValueFvsPatchScalarField>(phibf[patchi]))
435 phibf[patchi] == Uf[patchi];
436 }
437 }
438}
439
440
442{
443 if (mesh_.topoChanging())
444 {
445 for (auto& mrf: *this)
446 {
447 mrf.update();
449 }
450}
451
452
453// * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
454
455Foam::Ostream& Foam::operator<<
456(
457 Ostream& os,
458 const MRFZoneList& models
459)
460{
461 models.writeData(os);
462 return os;
463}
464
465
466// ************************************************************************* //
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
A field of fields is a PtrList of fields with reference counting.
Definition FieldField.H:77
Generic templated field type that is much like a Foam::List except that it is expected to hold numeri...
Definition Field.H:172
static tmp< GeometricField< vector, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< vector >::calculatedType())
GeometricBoundaryField< scalar, fvsPatchField, surfaceMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
List container for MRF zomes.
Definition MRFZoneList.H:57
tmp< surfaceScalarField > relative(const tmp< surfaceScalarField > &phi) const
Return the given absolute flux relative within the MRF region.
void reset(const dictionary &dict)
Reset the source list.
Definition MRFZoneList.C:63
const fvMesh & mesh_
Reference to the mesh database.
Definition MRFZoneList.H:78
bool writeData(Ostream &os) const
Write data to Ostream.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &phi) const
Return the given relative flux absolute within the MRF region.
bool read(const dictionary &dict)
Read dictionary.
void addAcceleration(const volVectorField &U, volVectorField &ddtU) const
Add the frame acceleration.
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
const MRFZone & getFromName(const word &name) const
Return the MRF with a given name.
Definition MRFZoneList.C:97
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
void update()
Update MRFZone faces if the mesh topology changes.
tmp< surfaceScalarField > phi() const
Return the MRF absolute flux.
tmp< volVectorField > DDt(const volVectorField &U) const
Return the frame acceleration.
bool active(const bool warn=false) const
Return active status.
Definition MRFZoneList.C:46
void correctBoundaryFlux(const volVectorField &U, surfaceScalarField &phi) const
Correct the boundary flux for the rotation of the MRF region.
MRF zone definition based on cell zone and parameters obtained from a control dictionary constructed ...
Definition MRFZone.H:65
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition Ostream.H:59
const MRFZone * set(const label i) const
Definition PtrList.H:171
constexpr PtrList() noexcept
Definition PtrListI.H:29
const MRFZone & operator[](const label i) const
Definition UPtrListI.H:289
label size() const noexcept
Definition UPtrListI.H:106
label count() const noexcept
Definition UPtrList.H:890
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
A class for managing temporary objects.
Definition tmp.H:75
void clear() const noexcept
If object pointer points to valid object: delete object and set pointer to nullptr.
Definition tmpI.H:289
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition tmpI.H:235
A class for handling words, derived from Foam::string.
Definition word.H:66
U
Definition pEqn.H:72
MRF makeRelative(fvc::interpolate(rho), phiHbyA)
fvVectorMatrix & UEqn
Definition UEqn.H:13
patchWriters resize(patchIds.size())
dynamicFvMesh & mesh
autoPtr< surfaceVectorField > Uf
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
OBJstream os(runTime.globalPath()/outputName)
auto & name
auto & names
void set(List< bool > &bools, const labelUList &locations)
Set the listed locations (assign 'true').
Definition BitOps.C:35
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of 'true' entries.
Definition BitOps.H:73
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh > > &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere).
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
const Type * isA(const U &obj)
Attempt dynamic_cast to Type.
Definition typeInfo.H:87
fvMatrix< vector > fvVectorMatrix
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< vector > dimensionedVector
Dimensioned vector obtained from generic dimensioned type.
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
Foam::surfaceFields.