MRFZone.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2018-2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 "MRFZone.H"
30 #include "fvMesh.H"
31 #include "volFields.H"
32 #include "surfaceFields.H"
33 #include "fvMatrices.H"
34 #include "faceSet.H"
35 #include "geometricOneField.H"
36 #include "syncTools.H"
37 
38 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42  defineTypeNameAndDebug(MRFZone, 0);
43 }
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
48 void Foam::MRFZone::setMRFFaces()
49 {
50  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
51 
52  // Type per face:
53  // 0:not in zone
54  // 1:moving with frame
55  // 2:other
56  labelList faceType(mesh_.nFaces(), Zero);
57 
58  // Determine faces in cell zone
59  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
60  // (without constructing cells)
61 
62  const labelList& own = mesh_.faceOwner();
63  const labelList& nei = mesh_.faceNeighbour();
64 
65  // Cells in zone
66  boolList zoneCell(mesh_.nCells(), false);
67 
68  if (cellZoneID_ != -1)
69  {
70  const labelList& cellLabels = mesh_.cellZones()[cellZoneID_];
71  forAll(cellLabels, i)
72  {
73  zoneCell[cellLabels[i]] = true;
74  }
75  }
76 
77 
78  label nZoneFaces = 0;
79 
80  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
81  {
82  if (zoneCell[own[facei]] || zoneCell[nei[facei]])
83  {
84  faceType[facei] = 1;
85  nZoneFaces++;
86  }
87  }
88 
89 
90  labelHashSet excludedPatches(excludedPatchLabels_);
91 
92  forAll(patches, patchi)
93  {
94  const polyPatch& pp = patches[patchi];
95 
96  if (pp.coupled() || excludedPatches.found(patchi))
97  {
98  forAll(pp, i)
99  {
100  label facei = pp.start()+i;
101 
102  if (zoneCell[own[facei]])
103  {
104  faceType[facei] = 2;
105  nZoneFaces++;
106  }
107  }
108  }
109  else if (!isA<emptyPolyPatch>(pp))
110  {
111  forAll(pp, i)
112  {
113  label facei = pp.start()+i;
114 
115  if (zoneCell[own[facei]])
116  {
117  faceType[facei] = 1;
118  nZoneFaces++;
119  }
120  }
121  }
122  }
123 
124  // Synchronize the faceType across processor patches
125  syncTools::syncFaceList(mesh_, faceType, maxEqOp<label>());
126 
127  // Now we have for faceType:
128  // 0 : face not in cellZone
129  // 1 : internal face or normal patch face
130  // 2 : coupled patch face or excluded patch face
131 
132  // Sort into lists per patch.
133 
134  internalFaces_.setSize(mesh_.nFaces());
135  label nInternal = 0;
136 
137  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
138  {
139  if (faceType[facei] == 1)
140  {
141  internalFaces_[nInternal++] = facei;
142  }
143  }
144  internalFaces_.setSize(nInternal);
145 
146  labelList nIncludedFaces(patches.size(), Zero);
147  labelList nExcludedFaces(patches.size(), Zero);
148 
149  forAll(patches, patchi)
150  {
151  const polyPatch& pp = patches[patchi];
152 
153  forAll(pp, patchFacei)
154  {
155  label facei = pp.start() + patchFacei;
156 
157  if (faceType[facei] == 1)
158  {
159  nIncludedFaces[patchi]++;
160  }
161  else if (faceType[facei] == 2)
162  {
163  nExcludedFaces[patchi]++;
164  }
165  }
166  }
167 
168  includedFaces_.setSize(patches.size());
169  excludedFaces_.setSize(patches.size());
170  forAll(nIncludedFaces, patchi)
171  {
172  includedFaces_[patchi].setSize(nIncludedFaces[patchi]);
173  excludedFaces_[patchi].setSize(nExcludedFaces[patchi]);
174  }
175  nIncludedFaces = 0;
176  nExcludedFaces = 0;
177 
178  forAll(patches, patchi)
179  {
180  const polyPatch& pp = patches[patchi];
181 
182  forAll(pp, patchFacei)
183  {
184  label facei = pp.start() + patchFacei;
185 
186  if (faceType[facei] == 1)
187  {
188  includedFaces_[patchi][nIncludedFaces[patchi]++] = patchFacei;
189  }
190  else if (faceType[facei] == 2)
191  {
192  excludedFaces_[patchi][nExcludedFaces[patchi]++] = patchFacei;
193  }
194  }
195  }
196 
197 
198  if (debug)
199  {
200  faceSet internalFaces(mesh_, "internalFaces", internalFaces_);
201  Pout<< "Writing " << internalFaces.size()
202  << " internal faces in MRF zone to faceSet "
203  << internalFaces.name() << endl;
204  internalFaces.write();
205 
206  faceSet MRFFaces(mesh_, "includedFaces", 100);
207  forAll(includedFaces_, patchi)
208  {
209  forAll(includedFaces_[patchi], i)
210  {
211  label patchFacei = includedFaces_[patchi][i];
212  MRFFaces.insert(patches[patchi].start()+patchFacei);
213  }
214  }
215  Pout<< "Writing " << MRFFaces.size()
216  << " patch faces in MRF zone to faceSet "
217  << MRFFaces.name() << endl;
218  MRFFaces.write();
219 
220  faceSet excludedFaces(mesh_, "excludedFaces", 100);
221  forAll(excludedFaces_, patchi)
222  {
223  forAll(excludedFaces_[patchi], i)
224  {
225  label patchFacei = excludedFaces_[patchi][i];
226  excludedFaces.insert(patches[patchi].start()+patchFacei);
227  }
228  }
229  Pout<< "Writing " << excludedFaces.size()
230  << " faces in MRF zone with special handling to faceSet "
231  << excludedFaces.name() << endl;
232  excludedFaces.write();
233  }
234 }
235 
236 
237 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
238 
239 Foam::MRFZone::MRFZone
240 (
241  const word& name,
242  const fvMesh& mesh,
243  const dictionary& dict,
244  const word& cellZoneName
245 )
246 :
247  mesh_(mesh),
248  name_(name),
249  coeffs_(dict),
250  active_(coeffs_.getOrDefault("active", true)),
251  cellZoneName_(cellZoneName),
252  cellZoneID_(),
253  excludedPatchNames_
254  (
255  coeffs_.getOrDefault<wordRes>("nonRotatingPatches", wordRes())
256  ),
257  origin_(coeffs_.get<vector>("origin")),
258  axis_(coeffs_.get<vector>("axis").normalise()),
259  omega_(Function1<scalar>::New("omega", coeffs_))
260 {
261  if (cellZoneName_ == word::null)
262  {
263  coeffs_.readEntry("cellZone", cellZoneName_);
264  }
265 
266  if (!active_)
267  {
268  cellZoneID_ = -1;
269  }
270  else
271  {
272  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
273 
274  const labelHashSet excludedPatchSet
275  (
276  mesh_.boundaryMesh().patchSet(excludedPatchNames_)
277  );
278 
279  excludedPatchLabels_.setSize(excludedPatchSet.size());
280 
281  label i = 0;
282  for (const label patchi : excludedPatchSet)
283  {
284  excludedPatchLabels_[i++] = patchi;
285  }
286 
287  bool cellZoneFound = (cellZoneID_ != -1);
288 
289  reduce(cellZoneFound, orOp<bool>());
290 
291  if (!cellZoneFound)
292  {
294  << "cannot find MRF cellZone " << cellZoneName_
295  << exit(FatalError);
296  }
297 
298  setMRFFaces();
299  }
300 }
301 
302 
303 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
304 
306 {
307  return omega_->value(mesh_.time().timeOutputValue())*axis_;
308 }
309 
310 
312 (
313  const volVectorField& U,
314  volVectorField& ddtU
315 ) const
316 {
317  if (cellZoneID_ == -1)
318  {
319  return;
320  }
321 
322  const labelList& cells = mesh_.cellZones()[cellZoneID_];
323  vectorField& ddtUc = ddtU.primitiveFieldRef();
324  const vectorField& Uc = U;
325 
326  const vector Omega = this->Omega();
327 
328  forAll(cells, i)
329  {
330  label celli = cells[i];
331  ddtUc[celli] += (Omega ^ Uc[celli]);
332  }
333 }
334 
335 
336 void Foam::MRFZone::addCoriolis(fvVectorMatrix& UEqn, const bool rhs) const
337 {
338  if (cellZoneID_ == -1)
339  {
340  return;
341  }
342 
343  const labelList& cells = mesh_.cellZones()[cellZoneID_];
344  const scalarField& V = mesh_.V();
345  vectorField& Usource = UEqn.source();
346  const vectorField& U = UEqn.psi();
347 
348  const vector Omega = this->Omega();
349 
350  if (rhs)
351  {
352  forAll(cells, i)
353  {
354  label celli = cells[i];
355  Usource[celli] += V[celli]*(Omega ^ U[celli]);
356  }
357  }
358  else
359  {
360  forAll(cells, i)
361  {
362  label celli = cells[i];
363  Usource[celli] -= V[celli]*(Omega ^ U[celli]);
364  }
365  }
366 }
367 
368 
370 (
371  const volScalarField& rho,
373  const bool rhs
374 ) const
375 {
376  if (cellZoneID_ == -1)
377  {
378  return;
379  }
380 
381  const labelList& cells = mesh_.cellZones()[cellZoneID_];
382  const scalarField& V = mesh_.V();
383  vectorField& Usource = UEqn.source();
384  const vectorField& U = UEqn.psi();
385 
386  const vector Omega = this->Omega();
387 
388  if (rhs)
389  {
390  forAll(cells, i)
391  {
392  label celli = cells[i];
393  Usource[celli] += V[celli]*rho[celli]*(Omega ^ U[celli]);
394  }
395  }
396  else
397  {
398  forAll(cells, i)
399  {
400  label celli = cells[i];
401  Usource[celli] -= V[celli]*rho[celli]*(Omega ^ U[celli]);
402  }
403  }
404 }
405 
406 
408 {
409  if (cellZoneID_ == -1)
410  {
411  return;
412  }
413 
414  const volVectorField& C = mesh_.C();
415 
416  const vector Omega = this->Omega();
417 
418  const labelList& cells = mesh_.cellZones()[cellZoneID_];
419 
420  forAll(cells, i)
421  {
422  label celli = cells[i];
423  U[celli] -= (Omega ^ (C[celli] - origin_));
424  }
425 
426  // Included patches
427 
428  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
429 
430  forAll(includedFaces_, patchi)
431  {
432  forAll(includedFaces_[patchi], i)
433  {
434  label patchFacei = includedFaces_[patchi][i];
435  Ubf[patchi][patchFacei] = Zero;
436  }
437  }
438 
439  // Excluded patches
440  forAll(excludedFaces_, patchi)
441  {
442  forAll(excludedFaces_[patchi], i)
443  {
444  label patchFacei = excludedFaces_[patchi][i];
445  Ubf[patchi][patchFacei] -=
446  (Omega
447  ^ (C.boundaryField()[patchi][patchFacei] - origin_));
448  }
449  }
450 }
451 
452 
454 {
455  makeRelativeRhoFlux(geometricOneField(), phi);
456 }
457 
458 
460 {
461  makeRelativeRhoFlux(oneFieldField(), phi);
462 }
463 
464 
465 void Foam::MRFZone::makeRelative(Field<scalar>& phi, const label patchi) const
466 {
467  makeRelativeRhoFlux(oneField(), phi, patchi);
468 }
469 
470 
472 (
473  const surfaceScalarField& rho,
475 ) const
476 {
477  makeRelativeRhoFlux(rho, phi);
478 }
479 
480 
482 {
483  if (cellZoneID_ == -1)
484  {
485  return;
486  }
487 
488  const volVectorField& C = mesh_.C();
489 
490  const vector Omega = this->Omega();
491 
492  const labelList& cells = mesh_.cellZones()[cellZoneID_];
493 
494  forAll(cells, i)
495  {
496  label celli = cells[i];
497  U[celli] += (Omega ^ (C[celli] - origin_));
498  }
499 
500  // Included patches
501  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
502 
503  forAll(includedFaces_, patchi)
504  {
505  forAll(includedFaces_[patchi], i)
506  {
507  label patchFacei = includedFaces_[patchi][i];
508  Ubf[patchi][patchFacei] =
509  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
510  }
511  }
512 
513  // Excluded patches
514  forAll(excludedFaces_, patchi)
515  {
516  forAll(excludedFaces_[patchi], i)
517  {
518  label patchFacei = excludedFaces_[patchi][i];
519  Ubf[patchi][patchFacei] +=
520  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
521  }
522  }
523 }
524 
525 
527 {
528  makeAbsoluteRhoFlux(geometricOneField(), phi);
529 }
530 
531 
533 (
534  const surfaceScalarField& rho,
536 ) const
537 {
538  makeAbsoluteRhoFlux(rho, phi);
539 }
540 
541 
543 {
544  if (!active_)
545  {
546  return;
547  }
548 
549  const vector Omega = this->Omega();
550 
551  // Included patches
552  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
553 
554  forAll(includedFaces_, patchi)
555  {
556  const vectorField& patchC = mesh_.Cf().boundaryField()[patchi];
557 
558  vectorField pfld(Ubf[patchi]);
559 
560  forAll(includedFaces_[patchi], i)
561  {
562  label patchFacei = includedFaces_[patchi][i];
563 
564  pfld[patchFacei] = (Omega ^ (patchC[patchFacei] - origin_));
565  }
566 
567  Ubf[patchi] == pfld;
568  }
569 }
570 
571 
573 {
574  os << nl;
575  os.beginBlock(name_);
576 
577  os.writeEntry("active", active_);
578  os.writeEntry("cellZone", cellZoneName_);
579  os.writeEntry("origin", origin_);
580  os.writeEntry("axis", axis_);
581  omega_->writeData(os);
582 
583  if (excludedPatchNames_.size())
584  {
585  os.writeEntry("nonRotatingPatches", excludedPatchNames_);
586  }
587 
588  os.endBlock();
589 }
590 
591 
593 {
594  coeffs_ = dict;
595 
596  active_ = coeffs_.getOrDefault("active", true);
597  coeffs_.readEntry("cellZone", cellZoneName_);
598  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
599 
600  return true;
601 }
602 
603 
605 {
606  if (mesh_.topoChanging())
607  {
608  setMRFFaces();
609  }
610 }
611 
612 
613 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:71
volFields.H
Foam::oneFieldField
A class representing the concept of a field of oneFields used to avoid unnecessary manipulations for ...
Definition: oneFieldField.H:53
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::FieldField
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:53
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::primitiveMesh::nFaces
label nFaces() const
Number of mesh faces.
Definition: primitiveMeshI.H:90
MRFZone.H
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:55
Foam::MRFZone::makeAbsolute
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:481
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:69
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:350
surfaceFields.H
Foam::surfaceFields.
Foam::Ostream::beginBlock
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:91
Foam::Pout
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Foam::Vector::normalise
Vector< Cmpt > & normalise()
Normalise the vector by its magnitude.
Definition: VectorI.H:123
Foam::HashSet< label, Hash< label > >
rho
rho
Definition: readInitialConditions.H:88
syncTools.H
Foam::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:287
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: Function1.H:86
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::MRFZone::read
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:592
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:51
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:390
Foam::primitiveMesh::nCells
label nCells() const
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field< vector >
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:107
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1107
faceSet.H
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::Ostream::endBlock
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:109
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::C::C
C()
Construct null.
Definition: C.C:43
Foam::Ostream::write
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:83
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::MRFZone::writeData
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:572
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::GeometricField::primitiveFieldRef
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field.
Definition: GeometricField.C:766
Foam::MRFZone::Omega
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:305
Foam::MRFZone::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:542
U
U
Definition: pEqn.H:72
Foam::MRFZone::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:407
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:381
geometricOneField.H
Foam::nl
constexpr char nl
Definition: Ostream.H:385
Foam::Vector< scalar >
Foam::List< label >
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:297
Foam::oneField
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition: oneField.H:53
Foam::wordRes
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:51
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:76
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:232
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::C
Graphite solid properties.
Definition: C.H:50
Foam::labelHashSet
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys and label hasher.
Definition: HashSet.H:409
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::dictionary::getOrDefault
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:122
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::orOp
Definition: ops.H:234
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::MRFZone::update
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZone.C:604
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::MRFZone::addCoriolis
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:312