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 | Copyright (C) 2011-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 "MRFZone.H"
27 #include "fvMesh.H"
28 #include "volFields.H"
29 #include "surfaceFields.H"
30 #include "fvMatrices.H"
31 #include "faceSet.H"
32 #include "geometricOneField.H"
33 #include "syncTools.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(MRFZone, 0);
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
46 {
48 
49  // Type per face:
50  // 0:not in zone
51  // 1:moving with frame
52  // 2:other
53  labelList faceType(mesh_.nFaces(), 0);
54 
55  // Determine faces in cell zone
56  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
57  // (without constructing cells)
58 
59  const labelList& own = mesh_.faceOwner();
60  const labelList& nei = mesh_.faceNeighbour();
61 
62  // Cells in zone
63  boolList zoneCell(mesh_.nCells(), false);
64 
65  if (cellZoneID_ != -1)
66  {
67  const labelList& cellLabels = mesh_.cellZones()[cellZoneID_];
68  forAll(cellLabels, i)
69  {
70  zoneCell[cellLabels[i]] = true;
71  }
72  }
73 
74 
75  label nZoneFaces = 0;
76 
77  for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
78  {
79  if (zoneCell[own[faceI]] || zoneCell[nei[faceI]])
80  {
81  faceType[faceI] = 1;
82  nZoneFaces++;
83  }
84  }
85 
86 
87  labelHashSet excludedPatches(excludedPatchLabels_);
88 
89  forAll(patches, patchI)
90  {
91  const polyPatch& pp = patches[patchI];
92 
93  if (pp.coupled() || excludedPatches.found(patchI))
94  {
95  forAll(pp, i)
96  {
97  label faceI = pp.start()+i;
98 
99  if (zoneCell[own[faceI]])
100  {
101  faceType[faceI] = 2;
102  nZoneFaces++;
103  }
104  }
105  }
106  else if (!isA<emptyPolyPatch>(pp))
107  {
108  forAll(pp, i)
109  {
110  label faceI = pp.start()+i;
111 
112  if (zoneCell[own[faceI]])
113  {
114  faceType[faceI] = 1;
115  nZoneFaces++;
116  }
117  }
118  }
119  }
120 
121  // Synchronize the faceType across processor patches
123 
124  // Now we have for faceType:
125  // 0 : face not in cellZone
126  // 1 : internal face or normal patch face
127  // 2 : coupled patch face or excluded patch face
128 
129  // Sort into lists per patch.
130 
132  label nInternal = 0;
133 
134  for (label faceI = 0; faceI < mesh_.nInternalFaces(); faceI++)
135  {
136  if (faceType[faceI] == 1)
137  {
138  internalFaces_[nInternal++] = faceI;
139  }
140  }
141  internalFaces_.setSize(nInternal);
142 
143  labelList nIncludedFaces(patches.size(), 0);
144  labelList nExcludedFaces(patches.size(), 0);
145 
147  {
148  const polyPatch& pp = patches[patchi];
149 
150  forAll(pp, patchFacei)
151  {
152  label faceI = pp.start() + patchFacei;
153 
154  if (faceType[faceI] == 1)
155  {
156  nIncludedFaces[patchi]++;
157  }
158  else if (faceType[faceI] == 2)
159  {
160  nExcludedFaces[patchi]++;
161  }
162  }
163  }
164 
167  forAll(nIncludedFaces, patchi)
168  {
169  includedFaces_[patchi].setSize(nIncludedFaces[patchi]);
170  excludedFaces_[patchi].setSize(nExcludedFaces[patchi]);
171  }
172  nIncludedFaces = 0;
173  nExcludedFaces = 0;
174 
176  {
177  const polyPatch& pp = patches[patchi];
178 
179  forAll(pp, patchFacei)
180  {
181  label faceI = pp.start() + patchFacei;
182 
183  if (faceType[faceI] == 1)
184  {
185  includedFaces_[patchi][nIncludedFaces[patchi]++] = patchFacei;
186  }
187  else if (faceType[faceI] == 2)
188  {
189  excludedFaces_[patchi][nExcludedFaces[patchi]++] = patchFacei;
190  }
191  }
192  }
193 
194 
195  if (debug)
196  {
197  faceSet internalFaces(mesh_, "internalFaces", internalFaces_);
198  Pout<< "Writing " << internalFaces.size()
199  << " internal faces in MRF zone to faceSet "
200  << internalFaces.name() << endl;
201  internalFaces.write();
202 
203  faceSet MRFFaces(mesh_, "includedFaces", 100);
205  {
207  {
208  label patchFacei = includedFaces_[patchi][i];
209  MRFFaces.insert(patches[patchi].start()+patchFacei);
210  }
211  }
212  Pout<< "Writing " << MRFFaces.size()
213  << " patch faces in MRF zone to faceSet "
214  << MRFFaces.name() << endl;
215  MRFFaces.write();
216 
217  faceSet excludedFaces(mesh_, "excludedFaces", 100);
219  {
221  {
222  label patchFacei = excludedFaces_[patchi][i];
223  excludedFaces.insert(patches[patchi].start()+patchFacei);
224  }
225  }
226  Pout<< "Writing " << excludedFaces.size()
227  << " faces in MRF zone with special handling to faceSet "
228  << excludedFaces.name() << endl;
229  excludedFaces.write();
230  }
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
235 
237 (
238  const word& name,
239  const fvMesh& mesh,
240  const dictionary& dict,
241  const word& cellZoneName
242 )
243 :
244  mesh_(mesh),
245  name_(name),
246  coeffs_(dict),
247  active_(coeffs_.lookupOrDefault("active", true)),
248  cellZoneName_(cellZoneName),
249  cellZoneID_(),
250  excludedPatchNames_
251  (
252  wordReList(coeffs_.lookupOrDefault("nonRotatingPatches", wordReList()))
253  ),
254  origin_(coeffs_.lookup("origin")),
255  axis_(coeffs_.lookup("axis")),
256  omega_(DataEntry<scalar>::New("omega", coeffs_))
257 {
258  if (cellZoneName_ == word::null)
259  {
260  coeffs_.lookup("cellZone") >> cellZoneName_;
261  }
262 
263  if (!active_)
264  {
265  cellZoneID_ = -1;
266  }
267  else
268  {
269  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
270 
271  axis_ = axis_/mag(axis_);
272 
273  const labelHashSet excludedPatchSet
274  (
275  mesh_.boundaryMesh().patchSet(excludedPatchNames_)
276  );
277 
278  excludedPatchLabels_.setSize(excludedPatchSet.size());
279 
280  label i = 0;
281  forAllConstIter(labelHashSet, excludedPatchSet, iter)
282  {
283  excludedPatchLabels_[i++] = iter.key();
284  }
285 
286  bool cellZoneFound = (cellZoneID_ != -1);
287 
288  reduce(cellZoneFound, orOp<bool>());
289 
290  if (!cellZoneFound)
291  {
293  << "cannot find MRF cellZone " << cellZoneName_
294  << exit(FatalError);
295  }
296 
297  setMRFFaces();
298  }
299 }
300 
301 
302 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
303 
305 {
306  return omega_->value(mesh_.time().timeOutputValue())*axis_;
307 }
308 
309 
311 (
312  const volVectorField& U,
313  volVectorField& ddtU
314 ) const
315 {
316  if (cellZoneID_ == -1)
317  {
318  return;
319  }
320 
321  const labelList& cells = mesh_.cellZones()[cellZoneID_];
322  vectorField& ddtUc = ddtU.internalField();
323  const vectorField& Uc = U.internalField();
324 
325  const vector Omega = this->Omega();
326 
327  forAll(cells, i)
328  {
329  label celli = cells[i];
330  ddtUc[celli] += (Omega ^ Uc[celli]);
331  }
332 }
333 
334 
335 void Foam::MRFZone::addCoriolis(fvVectorMatrix& UEqn, const bool rhs) const
336 {
337  if (cellZoneID_ == -1)
338  {
339  return;
340  }
341 
342  const labelList& cells = mesh_.cellZones()[cellZoneID_];
343  const scalarField& V = mesh_.V();
344  vectorField& Usource = UEqn.source();
345  const vectorField& U = UEqn.psi();
346 
347  const vector Omega = this->Omega();
348 
349  if (rhs)
350  {
351  forAll(cells, i)
352  {
353  label celli = cells[i];
354  Usource[celli] += V[celli]*(Omega ^ U[celli]);
355  }
356  }
357  else
358  {
359  forAll(cells, i)
360  {
361  label celli = cells[i];
362  Usource[celli] -= V[celli]*(Omega ^ U[celli]);
363  }
364  }
365 }
366 
367 
369 (
370  const volScalarField& rho,
372  const bool rhs
373 ) const
374 {
375  if (cellZoneID_ == -1)
376  {
377  return;
378  }
379 
380  const labelList& cells = mesh_.cellZones()[cellZoneID_];
381  const scalarField& V = mesh_.V();
382  vectorField& Usource = UEqn.source();
383  const vectorField& U = UEqn.psi();
384 
385  const vector Omega = this->Omega();
386 
387  if (rhs)
388  {
389  forAll(cells, i)
390  {
391  label celli = cells[i];
392  Usource[celli] += V[celli]*rho[celli]*(Omega ^ U[celli]);
393  }
394  }
395  else
396  {
397  forAll(cells, i)
398  {
399  label celli = cells[i];
400  Usource[celli] -= V[celli]*rho[celli]*(Omega ^ U[celli]);
401  }
402  }
403 }
404 
405 
407 {
408  const volVectorField& C = mesh_.C();
409 
410  const vector Omega = this->Omega();
411 
412  const labelList& cells = mesh_.cellZones()[cellZoneID_];
413 
414  forAll(cells, i)
415  {
416  label celli = cells[i];
417  U[celli] -= (Omega ^ (C[celli] - origin_));
418  }
419 
420  // Included patches
421  forAll(includedFaces_, patchi)
422  {
423  forAll(includedFaces_[patchi], i)
424  {
425  label patchFacei = includedFaces_[patchi][i];
426  U.boundaryField()[patchi][patchFacei] = vector::zero;
427  }
428  }
429 
430  // Excluded patches
431  forAll(excludedFaces_, patchi)
432  {
433  forAll(excludedFaces_[patchi], i)
434  {
435  label patchFacei = excludedFaces_[patchi][i];
436  U.boundaryField()[patchi][patchFacei] -=
437  (Omega
438  ^ (C.boundaryField()[patchi][patchFacei] - origin_));
439  }
440  }
441 }
442 
443 
445 {
446  makeRelativeRhoFlux(geometricOneField(), phi);
447 }
448 
449 
451 {
452  makeRelativeRhoFlux(oneFieldField(), phi);
453 }
454 
455 
457 (
458  const surfaceScalarField& rho,
460 ) const
461 {
462  makeRelativeRhoFlux(rho, phi);
463 }
464 
465 
467 {
468  const volVectorField& C = mesh_.C();
469 
470  const vector Omega = this->Omega();
471 
472  const labelList& cells = mesh_.cellZones()[cellZoneID_];
473 
474  forAll(cells, i)
475  {
476  label celli = cells[i];
477  U[celli] += (Omega ^ (C[celli] - origin_));
478  }
479 
480  // Included patches
481  forAll(includedFaces_, patchi)
482  {
483  forAll(includedFaces_[patchi], i)
484  {
485  label patchFacei = includedFaces_[patchi][i];
486  U.boundaryField()[patchi][patchFacei] =
487  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
488  }
489  }
490 
491  // Excluded patches
492  forAll(excludedFaces_, patchi)
493  {
494  forAll(excludedFaces_[patchi], i)
495  {
496  label patchFacei = excludedFaces_[patchi][i];
497  U.boundaryField()[patchi][patchFacei] +=
498  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
499  }
500  }
501 }
502 
503 
505 {
506  makeAbsoluteRhoFlux(geometricOneField(), phi);
507 }
508 
509 
511 (
512  const surfaceScalarField& rho,
514 ) const
515 {
516  makeAbsoluteRhoFlux(rho, phi);
517 }
518 
519 
521 {
522  const vector Omega = this->Omega();
523 
524  // Included patches
525  forAll(includedFaces_, patchi)
526  {
527  const vectorField& patchC = mesh_.Cf().boundaryField()[patchi];
528 
529  vectorField pfld(U.boundaryField()[patchi]);
530 
531  forAll(includedFaces_[patchi], i)
532  {
533  label patchFacei = includedFaces_[patchi][i];
534 
535  pfld[patchFacei] = (Omega ^ (patchC[patchFacei] - origin_));
536  }
537 
538  U.boundaryField()[patchi] == pfld;
539  }
540 }
541 
542 
544 {
545  os << nl;
546  os.write(name_) << nl;
547  os << token::BEGIN_BLOCK << incrIndent << nl;
548  os.writeKeyword("active") << active_ << token::END_STATEMENT << nl;
549  os.writeKeyword("cellZone") << cellZoneName_ << token::END_STATEMENT << nl;
550  os.writeKeyword("origin") << origin_ << token::END_STATEMENT << nl;
551  os.writeKeyword("axis") << axis_ << token::END_STATEMENT << nl;
552  omega_->writeData(os);
553 
554  if (excludedPatchNames_.size())
555  {
556  os.writeKeyword("nonRotatingPatches") << excludedPatchNames_
557  << token::END_STATEMENT << nl;
558  }
559 
560  os << decrIndent << token::END_BLOCK << nl;
561 }
562 
563 
565 {
566  coeffs_ = dict;
567 
568  active_ = coeffs_.lookupOrDefault("active", true);
569  coeffs_.lookup("cellZone") >> cellZoneName_;
570  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
571 
572  return true;
573 }
574 
575 
576 // ************************************************************************* //
volFields.H
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::Vector< scalar >::zero
static const Vector zero
Definition: Vector.H:80
Foam::oneFieldField
A class representing the concept of a field of oneFields used to avoid unnecessary manipulations for ...
Definition: oneFieldField.H:50
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::FieldField
Generic field type.
Definition: FieldField.H:51
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::MRFZone::MRFZone
MRFZone(const MRFZone &)
Disallow default bitwise copy construct.
MRFZone.H
Foam::polyPatch::coupled
virtual bool coupled() const
Return true if this patch is geometrically coupled (i.e. faces and.
Definition: polyPatch.H:322
Foam::MRFZone::setMRFFaces
void setMRFFaces()
Divide faces in frame according to patch.
Definition: MRFZone.C:45
Foam::wordReList
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
Foam::maxEqOp
Definition: ops.H:77
Foam::geometricOneField
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Definition: geometricOneField.H:52
Foam::MRFZone::internalFaces_
labelList internalFaces_
Internal faces that are part of MRF.
Definition: MRFZone.H:92
Foam::dictionary::lookupOrDefault
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
Definition: dictionaryTemplates.C:33
Foam::MRFZone::makeAbsolute
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:466
Foam::polyMesh::cellZones
const cellZoneMesh & cellZones() const
Return cell zone mesh.
Definition: polyMesh.H:469
Foam::regIOobject::write
virtual bool write() const
Write using setting from DB.
Definition: regIOobjectWrite.C:126
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:421
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
surfaceFields.H
Foam::surfaceFields.
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
Foam::HashSet< label, Hash< label > >
Foam::incrIndent
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
syncTools.H
U
U
Definition: pEqn.H:46
forAllConstIter
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
fvMatrices.H
A special matrix type and solver, designed for finite volume solutions of scalar equations.
Foam::primitiveMesh::nCells
label nCells() const
Definition: primitiveMeshI.H:64
Foam::MRFZone::read
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:564
Foam::reduce
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Definition: PstreamReduceOps.H:43
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::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::orOp
Definition: ops.H:178
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:66
Foam::polyMesh::faceOwner
virtual const labelList & faceOwner() const
Return face owner.
Definition: polyMesh.C:1017
faceSet.H
UEqn
tmp< fvVectorMatrix > UEqn(fvm::div(phi, U)+MRF.DDt(U)+turbulence->divDevReff(U)==fvOptions(U))
Foam::GeometricField::internalField
InternalField & internalField()
Return internal field.
Definition: GeometricField.C:724
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:137
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
Foam::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::C::C
C()
Construct null.
Definition: C.C:41
Foam::HashTable::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
fvMesh.H
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::decrIndent
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
Foam::polyPatch::start
label start() const
Return start label of this patch in the polyMesh face list.
Definition: polyPatch.H:312
Foam::MRFZone::writeData
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:543
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::MRFZone::cellZoneID_
label cellZoneID_
Cell zone ID.
Definition: MRFZone.H:85
Foam::MRFZone::Omega
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:304
rho
rho
Definition: pEqn.H:3
Foam::MRFZone::correctBoundaryVelocity
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:520
Foam::MRFZone::excludedPatchLabels_
labelList excludedPatchLabels_
Definition: MRFZone.H:89
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::MRFZone::makeRelative
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:406
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
geometricOneField.H
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Foam::token::END_BLOCK
@ END_BLOCK
Definition: token.H:105
Foam::Vector< scalar >
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::Ostream::writeKeyword
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
patchi
label patchi
Definition: getPatchFieldScalar.H:1
Foam::MRFZone::excludedFaces_
labelListList excludedFaces_
Excluded faces (per patch) that do not move with the MRF.
Definition: MRFZone.H:98
Foam::MRFZone::mesh_
const fvMesh & mesh_
Reference to the mesh database.
Definition: MRFZone.H:70
Foam::MRFZone::includedFaces_
labelListList includedFaces_
Outside faces (per patch) that move with the MRF.
Definition: MRFZone.H:95
Foam::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::word::null
static const word null
An empty word.
Definition: word.H:77
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::token::BEGIN_BLOCK
@ BEGIN_BLOCK
Definition: token.H:104
patches
patches[0]
Definition: createSingleCellMesh.H:36
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
Foam::syncTools::syncFaceList
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:381
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::C
Graphite solid properties.
Definition: C.H:57
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::DataEntry
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::polyMesh::faceNeighbour
virtual const labelList & faceNeighbour() const
Return face neighbour.
Definition: polyMesh.C:1023
Foam::MRFZone::addCoriolis
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:311