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-2022 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  forAll(patches, patchi)
91  {
92  const polyPatch& pp = patches[patchi];
93 
94  if (pp.coupled() || excludedPatchLabels_.contains(patchi))
95  {
96  forAll(pp, i)
97  {
98  label facei = pp.start()+i;
99 
100  if (zoneCell[own[facei]])
101  {
102  faceType[facei] = 2;
103  // ++nZoneFaces;
104  }
105  }
106  }
107  else if (!isA<emptyPolyPatch>(pp))
108  {
109  forAll(pp, i)
110  {
111  label facei = pp.start()+i;
112 
113  if (zoneCell[own[facei]])
114  {
115  faceType[facei] = 1;
116  // ++nZoneFaces;
117  }
118  }
119  }
120  }
121 
122  // Synchronize the faceType across processor patches
123  syncTools::syncFaceList(mesh_, faceType, maxEqOp<label>());
124 
125  // Now we have for faceType:
126  // 0 : face not in cellZone
127  // 1 : internal face or normal patch face
128  // 2 : coupled patch face or excluded patch face
129 
130  // Sort into lists per patch.
131 
132  internalFaces_.setSize(mesh_.nFaces());
133  label nInternal = 0;
134 
135  for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
136  {
137  if (faceType[facei] == 1)
138  {
139  internalFaces_[nInternal++] = facei;
140  }
141  }
142  internalFaces_.setSize(nInternal);
143 
144  labelList nIncludedFaces(patches.size(), Zero);
145  labelList nExcludedFaces(patches.size(), Zero);
146 
147  forAll(patches, patchi)
148  {
149  const polyPatch& pp = patches[patchi];
150 
151  forAll(pp, patchFacei)
152  {
153  label facei = pp.start() + patchFacei;
154 
155  if (faceType[facei] == 1)
156  {
157  nIncludedFaces[patchi]++;
158  }
159  else if (faceType[facei] == 2)
160  {
161  nExcludedFaces[patchi]++;
162  }
163  }
164  }
165 
166  includedFaces_.setSize(patches.size());
167  excludedFaces_.setSize(patches.size());
168  forAll(nIncludedFaces, patchi)
169  {
170  includedFaces_[patchi].setSize(nIncludedFaces[patchi]);
171  excludedFaces_[patchi].setSize(nExcludedFaces[patchi]);
172  }
173  nIncludedFaces = 0;
174  nExcludedFaces = 0;
175 
176  forAll(patches, patchi)
177  {
178  const polyPatch& pp = patches[patchi];
179 
180  forAll(pp, patchFacei)
181  {
182  label facei = pp.start() + patchFacei;
183 
184  if (faceType[facei] == 1)
185  {
186  includedFaces_[patchi][nIncludedFaces[patchi]++] = patchFacei;
187  }
188  else if (faceType[facei] == 2)
189  {
190  excludedFaces_[patchi][nExcludedFaces[patchi]++] = patchFacei;
191  }
192  }
193  }
194 
195 
196  if (debug)
197  {
198  faceSet internalFaces(mesh_, "internalFaces", internalFaces_);
199  Pout<< "Writing " << internalFaces.size()
200  << " internal faces in MRF zone to faceSet "
201  << internalFaces.name() << endl;
202  internalFaces.write();
203 
204  faceSet MRFFaces(mesh_, "includedFaces", 100);
205  forAll(includedFaces_, patchi)
206  {
207  forAll(includedFaces_[patchi], i)
208  {
209  label patchFacei = includedFaces_[patchi][i];
210  MRFFaces.insert(patches[patchi].start()+patchFacei);
211  }
212  }
213  Pout<< "Writing " << MRFFaces.size()
214  << " patch faces in MRF zone to faceSet "
215  << MRFFaces.name() << endl;
216  MRFFaces.write();
217 
218  faceSet excludedFaces(mesh_, "excludedFaces", 100);
219  forAll(excludedFaces_, patchi)
220  {
221  forAll(excludedFaces_[patchi], i)
222  {
223  label patchFacei = excludedFaces_[patchi][i];
224  excludedFaces.insert(patches[patchi].start()+patchFacei);
225  }
226  }
227  Pout<< "Writing " << excludedFaces.size()
228  << " faces in MRF zone with special handling to faceSet "
229  << excludedFaces.name() << endl;
230  excludedFaces.write();
231  }
232 }
233 
234 
235 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
236 
237 Foam::MRFZone::MRFZone
238 (
239  const word& name,
240  const fvMesh& mesh,
241  const dictionary& dict,
242  const word& cellZoneName
243 )
244 :
245  mesh_(mesh),
246  name_(name),
247  coeffs_(dict),
248  active_(true),
249  cellZoneName_(cellZoneName),
250  cellZoneID_(-1),
251  excludedPatchNames_(),
252  origin_(Zero),
253  axis_(Zero),
254  omega_(nullptr)
255 {
256  read(dict);
257 }
258 
259 
260 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
261 
263 {
264  return omega_->value(mesh_.time().timeOutputValue())*axis_;
265 }
266 
267 
269 (
270  const volVectorField& U,
271  volVectorField& ddtU
272 ) const
273 {
274  if (cellZoneID_ == -1)
275  {
276  return;
277  }
278 
279  const labelList& cells = mesh_.cellZones()[cellZoneID_];
280  vectorField& ddtUc = ddtU.primitiveFieldRef();
281  const vectorField& Uc = U;
282 
283  const vector Omega = this->Omega();
284 
285  forAll(cells, i)
286  {
287  label celli = cells[i];
288  ddtUc[celli] += (Omega ^ Uc[celli]);
289  }
290 }
291 
292 
293 void Foam::MRFZone::addCoriolis(fvVectorMatrix& UEqn, const bool rhs) const
294 {
295  if (cellZoneID_ == -1)
296  {
297  return;
298  }
299 
300  const labelList& cells = mesh_.cellZones()[cellZoneID_];
301  const scalarField& V = mesh_.V();
302  vectorField& Usource = UEqn.source();
303  const vectorField& U = UEqn.psi();
304 
305  const vector Omega = this->Omega();
306 
307  if (rhs)
308  {
309  forAll(cells, i)
310  {
311  label celli = cells[i];
312  Usource[celli] += V[celli]*(Omega ^ U[celli]);
313  }
314  }
315  else
316  {
317  forAll(cells, i)
318  {
319  label celli = cells[i];
320  Usource[celli] -= V[celli]*(Omega ^ U[celli]);
321  }
322  }
323 }
324 
325 
327 (
328  const volScalarField& rho,
330  const bool rhs
331 ) const
332 {
333  if (cellZoneID_ == -1)
334  {
335  return;
336  }
337 
338  const labelList& cells = mesh_.cellZones()[cellZoneID_];
339  const scalarField& V = mesh_.V();
340  vectorField& Usource = UEqn.source();
341  const vectorField& U = UEqn.psi();
342 
343  const vector Omega = this->Omega();
344 
345  if (rhs)
346  {
347  forAll(cells, i)
348  {
349  label celli = cells[i];
350  Usource[celli] += V[celli]*rho[celli]*(Omega ^ U[celli]);
351  }
352  }
353  else
354  {
355  forAll(cells, i)
356  {
357  label celli = cells[i];
358  Usource[celli] -= V[celli]*rho[celli]*(Omega ^ U[celli]);
359  }
360  }
361 }
362 
363 
365 {
366  if (cellZoneID_ == -1)
367  {
368  return;
369  }
370 
371  const volVectorField& C = mesh_.C();
372 
373  const vector Omega = this->Omega();
374 
375  const labelList& cells = mesh_.cellZones()[cellZoneID_];
376 
377  forAll(cells, i)
378  {
379  label celli = cells[i];
380  U[celli] -= (Omega ^ (C[celli] - origin_));
381  }
382 
383  // Included patches
384 
385  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
386 
387  forAll(includedFaces_, patchi)
388  {
389  forAll(includedFaces_[patchi], i)
390  {
391  label patchFacei = includedFaces_[patchi][i];
392  Ubf[patchi][patchFacei] = Zero;
393  }
394  }
395 
396  // Excluded patches
397  forAll(excludedFaces_, patchi)
398  {
399  forAll(excludedFaces_[patchi], i)
400  {
401  label patchFacei = excludedFaces_[patchi][i];
402  Ubf[patchi][patchFacei] -=
403  (Omega
404  ^ (C.boundaryField()[patchi][patchFacei] - origin_));
405  }
406  }
407 }
408 
411 {
412  makeRelativeRhoFlux(geometricOneField(), phi);
413 }
414 
417 {
418  makeRelativeRhoFlux(oneFieldField(), phi);
419 }
420 
421 
422 void Foam::MRFZone::makeRelative(Field<scalar>& phi, const label patchi) const
423 {
424  makeRelativeRhoFlux(oneField(), phi, patchi);
425 }
426 
427 
429 (
430  const surfaceScalarField& rho,
432 ) const
433 {
434  makeRelativeRhoFlux(rho, phi);
435 }
436 
437 
439 {
440  if (cellZoneID_ == -1)
441  {
442  return;
443  }
444 
445  const volVectorField& C = mesh_.C();
446 
447  const vector Omega = this->Omega();
448 
449  const labelList& cells = mesh_.cellZones()[cellZoneID_];
450 
451  forAll(cells, i)
452  {
453  label celli = cells[i];
454  U[celli] += (Omega ^ (C[celli] - origin_));
455  }
456 
457  // Included patches
458  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
459 
460  forAll(includedFaces_, patchi)
461  {
462  forAll(includedFaces_[patchi], i)
463  {
464  label patchFacei = includedFaces_[patchi][i];
465  Ubf[patchi][patchFacei] =
466  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
467  }
468  }
469 
470  // Excluded patches
471  forAll(excludedFaces_, patchi)
472  {
473  forAll(excludedFaces_[patchi], i)
474  {
475  label patchFacei = excludedFaces_[patchi][i];
476  Ubf[patchi][patchFacei] +=
477  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
478  }
479  }
480 }
481 
482 
484 {
485  makeAbsoluteRhoFlux(geometricOneField(), phi);
486 }
487 
488 
490 (
491  const surfaceScalarField& rho,
493 ) const
494 {
495  makeAbsoluteRhoFlux(rho, phi);
496 }
497 
498 
500 {
501  if (!active_)
502  {
503  return;
504  }
505 
506  const vector Omega = this->Omega();
507 
508  // Included patches
509  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
510 
511  forAll(includedFaces_, patchi)
512  {
513  const vectorField& patchC = mesh_.Cf().boundaryField()[patchi];
514 
515  vectorField pfld(Ubf[patchi]);
516 
517  forAll(includedFaces_[patchi], i)
518  {
519  label patchFacei = includedFaces_[patchi][i];
520 
521  pfld[patchFacei] = (Omega ^ (patchC[patchFacei] - origin_));
522  }
523 
524  Ubf[patchi] == pfld;
525  }
526 }
527 
528 
529 void Foam::MRFZone::writeData(Ostream& os) const
530 {
531  os << nl;
532  os.beginBlock(name_);
533 
534  os.writeEntry("active", active_);
535  os.writeEntry("cellZone", cellZoneName_);
536  os.writeEntry("origin", origin_);
537  os.writeEntry("axis", axis_);
538  omega_->writeData(os);
539 
540  if (excludedPatchNames_.size())
541  {
542  os.writeEntry("nonRotatingPatches", excludedPatchNames_);
543  }
544 
545  os.endBlock();
546 }
547 
548 
549 bool Foam::MRFZone::read(const dictionary& dict)
550 {
551  coeffs_ = dict;
552 
553  coeffs_.readIfPresent("active", active_);
554 
555  if (!active_)
556  {
557  cellZoneID_ = -1;
558  return true;
559  }
560 
561  coeffs_.readIfPresent("nonRotatingPatches", excludedPatchNames_);
562 
563  origin_ = coeffs_.get<vector>("origin");
564  axis_ = coeffs_.get<vector>("axis").normalise();
565  omega_.reset(Function1<scalar>::New("omega", coeffs_, &mesh_));
566 
567  const word oldCellZoneName = cellZoneName_;
568  if (cellZoneName_.empty())
569  {
570  coeffs_.readEntry("cellZone", cellZoneName_);
571  }
572  else
573  {
574  coeffs_.readIfPresent("cellZone", cellZoneName_);
575  }
576 
577  if (cellZoneID_ == -1 || oldCellZoneName != cellZoneName_)
578  {
579  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
580 
581  excludedPatchLabels_ =
582  mesh_.boundaryMesh().indices(excludedPatchNames_);
583 
584  if (!returnReduceOr(cellZoneID_ != -1))
585  {
587  << "cannot find MRF cellZone " << cellZoneName_
588  << exit(FatalError);
589  }
590 
591  setMRFFaces();
592  }
593 
594  return true;
595 }
596 
597 
599 {
600  if (mesh_.topoChanging())
601  {
602  setMRFFaces();
603  }
604 }
605 
606 
607 // ************************************************************************* //
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:542
Foam::surfaceFields.
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:262
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:357
dictionary dict
Graphite solid properties.
Definition: C.H:46
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:598
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1122
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:431
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void syncFaceList(const polyMesh &mesh, UList< T > &faceValues, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:432
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:321
GeometricBoundaryField< vector, fvPatchField, volMesh > Boundary
Type of boundary fields.
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZone.C:591
label nFaces() const noexcept
Number of mesh faces.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
virtual const fileName & name() const override
Get the name of the output serial stream. (eg, the name of the Fstream file name) ...
Definition: OSstream.H:128
bool contains(const T &val) const
True if the value is contained in the list.
Definition: UListI.H:300
A class representing the concept of a field of oneFields used to avoid unnecessary manipulations for ...
Definition: oneFieldField.H:48
A field of fields is a PtrList of fields with reference counting.
Definition: FieldField.H:51
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:485
void setSize(const label n)
Alias for resize()
Definition: List.H:316
dynamicFvMesh & mesh
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
const cellShapeList & cells
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:608
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual Ostream & endBlock()
Write end block group.
Definition: Ostream.C:108
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1116
virtual bool write(const token &tok)=0
Write token to stream or otherwise handle it.
label nInternalFaces() const noexcept
Number of internal faces.
Vector< scalar > vector
Definition: vector.H:57
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
int debug
Static debugging option.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
OBJstream os(runTime.globalPath()/outputName)
defineTypeNameAndDebug(combustionModel, 0)
volScalarField & C
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:492
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:255
Field< Type > & source() noexcept
Definition: fvMatrix.H:533
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition: oneField.H:47
Internal::FieldType & primitiveFieldRef(const bool updateAccessTime=true)
Return a reference to the internal field values.
U
Definition: pEqn.H:72
label nCells() const noexcept
Number of mesh cells.
fvVectorMatrix & UEqn
Definition: UEqn.H:13
fvMatrix< vector > fvVectorMatrix
Definition: fvMatricesFwd.H:40
const polyBoundaryMesh & patches
A special matrix type and solver, designed for finite volume solutions of scalar equations.
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:678
virtual Ostream & beginBlock(const keyType &kw)
Write begin block group with the given name.
Definition: Ostream.C:90
Field< vector > vectorField
Specialisation of Field<T> for vector.
List< label > labelList
A List of labels.
Definition: List.H:62
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:522
List< bool > boolList
A List of bools.
Definition: List.H:60
const Boundary & boundaryField() const noexcept
Return const-reference to the boundary field.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127