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-2021 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_(true),
251  cellZoneName_(cellZoneName),
252  cellZoneID_(-1),
253  excludedPatchNames_(wordRes()),
254  origin_(Zero),
255  axis_(Zero),
256  omega_(nullptr)
257 {
258  read(dict);
259 }
260 
261 
262 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
263 
265 {
266  return omega_->value(mesh_.time().timeOutputValue())*axis_;
267 }
268 
269 
271 (
272  const volVectorField& U,
273  volVectorField& ddtU
274 ) const
275 {
276  if (cellZoneID_ == -1)
277  {
278  return;
279  }
280 
281  const labelList& cells = mesh_.cellZones()[cellZoneID_];
282  vectorField& ddtUc = ddtU.primitiveFieldRef();
283  const vectorField& Uc = U;
284 
285  const vector Omega = this->Omega();
286 
287  forAll(cells, i)
288  {
289  label celli = cells[i];
290  ddtUc[celli] += (Omega ^ Uc[celli]);
291  }
292 }
293 
294 
295 void Foam::MRFZone::addCoriolis(fvVectorMatrix& UEqn, const bool rhs) const
296 {
297  if (cellZoneID_ == -1)
298  {
299  return;
300  }
301 
302  const labelList& cells = mesh_.cellZones()[cellZoneID_];
303  const scalarField& V = mesh_.V();
304  vectorField& Usource = UEqn.source();
305  const vectorField& U = UEqn.psi();
306 
307  const vector Omega = this->Omega();
308 
309  if (rhs)
310  {
311  forAll(cells, i)
312  {
313  label celli = cells[i];
314  Usource[celli] += V[celli]*(Omega ^ U[celli]);
315  }
316  }
317  else
318  {
319  forAll(cells, i)
320  {
321  label celli = cells[i];
322  Usource[celli] -= V[celli]*(Omega ^ U[celli]);
323  }
324  }
325 }
326 
327 
329 (
330  const volScalarField& rho,
332  const bool rhs
333 ) const
334 {
335  if (cellZoneID_ == -1)
336  {
337  return;
338  }
339 
340  const labelList& cells = mesh_.cellZones()[cellZoneID_];
341  const scalarField& V = mesh_.V();
342  vectorField& Usource = UEqn.source();
343  const vectorField& U = UEqn.psi();
344 
345  const vector Omega = this->Omega();
346 
347  if (rhs)
348  {
349  forAll(cells, i)
350  {
351  label celli = cells[i];
352  Usource[celli] += V[celli]*rho[celli]*(Omega ^ U[celli]);
353  }
354  }
355  else
356  {
357  forAll(cells, i)
358  {
359  label celli = cells[i];
360  Usource[celli] -= V[celli]*rho[celli]*(Omega ^ U[celli]);
361  }
362  }
363 }
364 
365 
367 {
368  if (cellZoneID_ == -1)
369  {
370  return;
371  }
372 
373  const volVectorField& C = mesh_.C();
374 
375  const vector Omega = this->Omega();
376 
377  const labelList& cells = mesh_.cellZones()[cellZoneID_];
378 
379  forAll(cells, i)
380  {
381  label celli = cells[i];
382  U[celli] -= (Omega ^ (C[celli] - origin_));
383  }
384 
385  // Included patches
386 
387  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
388 
389  forAll(includedFaces_, patchi)
390  {
391  forAll(includedFaces_[patchi], i)
392  {
393  label patchFacei = includedFaces_[patchi][i];
394  Ubf[patchi][patchFacei] = Zero;
395  }
396  }
397 
398  // Excluded patches
399  forAll(excludedFaces_, patchi)
400  {
401  forAll(excludedFaces_[patchi], i)
402  {
403  label patchFacei = excludedFaces_[patchi][i];
404  Ubf[patchi][patchFacei] -=
405  (Omega
406  ^ (C.boundaryField()[patchi][patchFacei] - origin_));
407  }
408  }
409 }
410 
411 
413 {
414  makeRelativeRhoFlux(geometricOneField(), phi);
415 }
416 
417 
419 {
420  makeRelativeRhoFlux(oneFieldField(), phi);
421 }
422 
423 
424 void Foam::MRFZone::makeRelative(Field<scalar>& phi, const label patchi) const
425 {
426  makeRelativeRhoFlux(oneField(), phi, patchi);
427 }
428 
429 
431 (
432  const surfaceScalarField& rho,
434 ) const
435 {
436  makeRelativeRhoFlux(rho, phi);
437 }
438 
439 
441 {
442  if (cellZoneID_ == -1)
443  {
444  return;
445  }
446 
447  const volVectorField& C = mesh_.C();
448 
449  const vector Omega = this->Omega();
450 
451  const labelList& cells = mesh_.cellZones()[cellZoneID_];
452 
453  forAll(cells, i)
454  {
455  label celli = cells[i];
456  U[celli] += (Omega ^ (C[celli] - origin_));
457  }
458 
459  // Included patches
460  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
461 
462  forAll(includedFaces_, patchi)
463  {
464  forAll(includedFaces_[patchi], i)
465  {
466  label patchFacei = includedFaces_[patchi][i];
467  Ubf[patchi][patchFacei] =
468  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
469  }
470  }
471 
472  // Excluded patches
473  forAll(excludedFaces_, patchi)
474  {
475  forAll(excludedFaces_[patchi], i)
476  {
477  label patchFacei = excludedFaces_[patchi][i];
478  Ubf[patchi][patchFacei] +=
479  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
480  }
481  }
482 }
483 
484 
486 {
487  makeAbsoluteRhoFlux(geometricOneField(), phi);
488 }
489 
490 
492 (
493  const surfaceScalarField& rho,
495 ) const
496 {
497  makeAbsoluteRhoFlux(rho, phi);
498 }
499 
500 
502 {
503  if (!active_)
504  {
505  return;
506  }
507 
508  const vector Omega = this->Omega();
509 
510  // Included patches
511  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
512 
513  forAll(includedFaces_, patchi)
514  {
515  const vectorField& patchC = mesh_.Cf().boundaryField()[patchi];
516 
517  vectorField pfld(Ubf[patchi]);
518 
519  forAll(includedFaces_[patchi], i)
520  {
521  label patchFacei = includedFaces_[patchi][i];
522 
523  pfld[patchFacei] = (Omega ^ (patchC[patchFacei] - origin_));
524  }
525 
526  Ubf[patchi] == pfld;
527  }
528 }
529 
530 
532 {
533  os << nl;
534  os.beginBlock(name_);
535 
536  os.writeEntry("active", active_);
537  os.writeEntry("cellZone", cellZoneName_);
538  os.writeEntry("origin", origin_);
539  os.writeEntry("axis", axis_);
540  omega_->writeData(os);
541 
542  if (excludedPatchNames_.size())
543  {
544  os.writeEntry("nonRotatingPatches", excludedPatchNames_);
545  }
546 
547  os.endBlock();
548 }
549 
550 
552 {
553  coeffs_ = dict;
554 
555  coeffs_.readIfPresent("active", active_);
556 
557  if (!active_)
558  {
559  cellZoneID_ = -1;
560  return true;
561  }
562 
563  coeffs_.readIfPresent("nonRotatingPatches", excludedPatchNames_);
564 
565  origin_ = coeffs_.get<vector>("origin");
566  axis_ = coeffs_.get<vector>("axis").normalise();
567  omega_.reset(Function1<scalar>::New("omega", coeffs_, &mesh_));
568 
569  const word oldCellZoneName = cellZoneName_;
570  if (cellZoneName_ == word::null)
571  {
572  coeffs_.readEntry("cellZone", cellZoneName_);
573  }
574  else
575  {
576  coeffs_.readIfPresent("cellZone", cellZoneName_);
577  }
578 
579  if (cellZoneID_ == -1 || oldCellZoneName != cellZoneName_)
580  {
581  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
582 
583  const labelHashSet excludedPatchSet
584  (
585  mesh_.boundaryMesh().patchSet(excludedPatchNames_)
586  );
587 
588  excludedPatchLabels_.setSize(excludedPatchSet.size());
589 
590  label i = 0;
591  for (const label patchi : excludedPatchSet)
592  {
593  excludedPatchLabels_[i++] = patchi;
594  }
595 
596  bool cellZoneFound = (cellZoneID_ != -1);
597 
598  reduce(cellZoneFound, orOp<bool>());
599 
600  if (!cellZoneFound)
601  {
603  << "cannot find MRF cellZone " << cellZoneName_
604  << exit(FatalError);
605  }
606 
607  setMRFFaces();
608  }
609 
610  return true;
611 }
612 
613 
615 {
616  if (mesh_.topoChanging())
617  {
618  setMRFFaces();
619  }
620 }
621 
622 
623 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:67
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:65
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
MRFZone.H
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
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:440
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:65
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:369
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::HashSet< label, Hash< label > >
rho
rho
Definition: readInitialConditions.H:88
syncTools.H
Foam::Function1
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: propellerInfo.H:291
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:551
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:396
Foam::primitiveMesh::nCells
label nCells() const noexcept
Number of mesh cells.
Definition: primitiveMeshI.H:96
Foam::Field< vector >
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const noexcept
Return cell zone mesh.
Definition: polyMesh.H:492
Foam::OSstream::name
virtual const fileName & name() const
Get the name of the stream.
Definition: OSstream.H:107
Foam::List::setSize
void setSize(const label n)
Alias for resize()
Definition: List.H:222
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::fvMatrix::psi
const GeometricField< Type, fvPatchField, volMesh > & psi(const label i=0) const
Return psi.
Definition: fvMatrix.H:399
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:123
os
OBJstream os(runTime.globalPath()/outputName)
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:85
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::MRFZone::writeData
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:531
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:264
Foam::MRFZone::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:501
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:366
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
geometricOneField.H
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Vector< scalar >
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const noexcept
Number of internal faces.
Definition: primitiveMeshI.H:78
Foam::List< label >
Foam::fvMatrix::source
Field< Type > & source()
Definition: fvMatrix.H:445
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:80
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:68
Foam::Ostream::writeEntry
Ostream & writeEntry(const keyType &key, const T &value)
Write a keyword/value entry.
Definition: Ostream.H:236
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
UEqn
fvVectorMatrix & UEqn
Definition: UEqn.H:13
Foam::primitiveMesh::nFaces
label nFaces() const noexcept
Number of mesh faces.
Definition: primitiveMeshI.H:90
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 of labels, uses label hasher.
Definition: HashSet.H:85
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:614
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1113
Foam::dictionary::readIfPresent
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryTemplates.C:405
Foam::MRFZone::addCoriolis
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:271