MRFZone.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2021 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 
45 void Foam::MRFZone::setMRFFaces()
46 {
47  const polyBoundaryMesh& patches = mesh_.boundaryMesh();
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  // Synchronise the faceType across processor patches
122  syncTools::syncFaceList(mesh_, faceType, maxEqOp<label>());
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 
131  internalFaces_.setSize(mesh_.nFaces());
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 
146  forAll(patches, patchi)
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 
165  includedFaces_.setSize(patches.size());
166  excludedFaces_.setSize(patches.size());
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 
175  forAll(patches, patchi)
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);
204  forAll(includedFaces_, patchi)
205  {
206  forAll(includedFaces_[patchi], i)
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);
218  forAll(excludedFaces_, patchi)
219  {
220  forAll(excludedFaces_[patchi], i)
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  cellZoneName_(cellZoneName),
248  cellZoneID_(),
249  excludedPatchNames_
250  (
251  wordReList(coeffs_.lookupOrDefault("nonRotatingPatches", wordReList()))
252  ),
253  origin_(coeffs_.lookup("origin")),
254  axis_(coeffs_.lookup("axis")),
255  omega_(Function1<scalar>::New("omega", coeffs_))
256 {
257  if (cellZoneName_ == word::null)
258  {
259  coeffs_.lookup("cellZone") >> cellZoneName_;
260  }
261 
262  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
263 
264  axis_ = axis_/mag(axis_);
265 
266  const labelHashSet excludedPatchSet
267  (
268  mesh_.boundaryMesh().patchSet(excludedPatchNames_)
269  );
270 
271  excludedPatchLabels_.setSize(excludedPatchSet.size());
272 
273  label i = 0;
274  forAllConstIter(labelHashSet, excludedPatchSet, iter)
275  {
276  excludedPatchLabels_[i++] = iter.key();
277  }
278 
279  bool cellZoneFound = (cellZoneID_ != -1);
280 
281  reduce(cellZoneFound, orOp<bool>());
282 
283  if (!cellZoneFound)
284  {
286  << "cannot find MRF cellZone " << cellZoneName_
287  << exit(FatalError);
288  }
289 
290  setMRFFaces();
291 }
292 
293 
294 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
295 
297 {
298  return omega_->value(mesh_.time().timeOutputValue())*axis_;
299 }
300 
301 
303 (
304  const volVectorField& U,
305  volVectorField& ddtU
306 ) const
307 {
308  if (cellZoneID_ == -1)
309  {
310  return;
311  }
312 
313  const labelList& cells = mesh_.cellZones()[cellZoneID_];
314  vectorField& ddtUc = ddtU.primitiveFieldRef();
315  const vectorField& Uc = U;
316 
317  const vector Omega = this->Omega();
318 
319  forAll(cells, i)
320  {
321  label celli = cells[i];
322  ddtUc[celli] += (Omega ^ Uc[celli]);
323  }
324 }
325 
326 
327 void Foam::MRFZone::addCoriolis(fvVectorMatrix& UEqn, const bool rhs) const
328 {
329  if (cellZoneID_ == -1)
330  {
331  return;
332  }
333 
334  const labelList& cells = mesh_.cellZones()[cellZoneID_];
335  const scalarField& V = mesh_.V();
336  vectorField& Usource = UEqn.source();
337  const vectorField& U = UEqn.psi();
338 
339  const vector Omega = this->Omega();
340 
341  if (rhs)
342  {
343  forAll(cells, i)
344  {
345  label celli = cells[i];
346  Usource[celli] += V[celli]*(Omega ^ U[celli]);
347  }
348  }
349  else
350  {
351  forAll(cells, i)
352  {
353  label celli = cells[i];
354  Usource[celli] -= V[celli]*(Omega ^ U[celli]);
355  }
356  }
357 }
358 
359 
361 (
362  const volScalarField& rho,
363  fvVectorMatrix& UEqn,
364  const bool rhs
365 ) const
366 {
367  if (cellZoneID_ == -1)
368  {
369  return;
370  }
371 
372  const labelList& cells = mesh_.cellZones()[cellZoneID_];
373  const scalarField& V = mesh_.V();
374  vectorField& Usource = UEqn.source();
375  const vectorField& U = UEqn.psi();
376 
377  const vector Omega = this->Omega();
378 
379  if (rhs)
380  {
381  forAll(cells, i)
382  {
383  label celli = cells[i];
384  Usource[celli] += V[celli]*rho[celli]*(Omega ^ U[celli]);
385  }
386  }
387  else
388  {
389  forAll(cells, i)
390  {
391  label celli = cells[i];
392  Usource[celli] -= V[celli]*rho[celli]*(Omega ^ U[celli]);
393  }
394  }
395 }
396 
397 
399 {
400  if (cellZoneID_ == -1)
401  {
402  return;
403  }
404 
405  const volVectorField& C = mesh_.C();
406 
407  const vector Omega = this->Omega();
408 
409  const labelList& cells = mesh_.cellZones()[cellZoneID_];
410 
411  forAll(cells, i)
412  {
413  label celli = cells[i];
414  U[celli] -= (Omega ^ (C[celli] - origin_));
415  }
416 
417  // Included patches
418 
419  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
420 
421  forAll(includedFaces_, patchi)
422  {
423  forAll(includedFaces_[patchi], i)
424  {
425  label patchFacei = includedFaces_[patchi][i];
426  Ubf[patchi][patchFacei] = 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  Ubf[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 
456 void Foam::MRFZone::makeRelative(Field<scalar>& phi, const label patchi) const
457 {
458  makeRelativeRhoFlux(oneField(), phi, patchi);
459 }
460 
461 
463 (
464  const surfaceScalarField& rho,
465  surfaceScalarField& phi
466 ) const
467 {
468  makeRelativeRhoFlux(rho, phi);
469 }
470 
471 
473 {
474  if (cellZoneID_ == -1)
475  {
476  return;
477  }
478 
479  const volVectorField& C = mesh_.C();
480 
481  const vector Omega = this->Omega();
482 
483  const labelList& cells = mesh_.cellZones()[cellZoneID_];
484 
485  forAll(cells, i)
486  {
487  label celli = cells[i];
488  U[celli] += (Omega ^ (C[celli] - origin_));
489  }
490 
491  // Included patches
492  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
493 
494  forAll(includedFaces_, patchi)
495  {
496  forAll(includedFaces_[patchi], i)
497  {
498  label patchFacei = includedFaces_[patchi][i];
499  Ubf[patchi][patchFacei] =
500  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
501  }
502  }
503 
504  // Excluded patches
505  forAll(excludedFaces_, patchi)
506  {
507  forAll(excludedFaces_[patchi], i)
508  {
509  label patchFacei = excludedFaces_[patchi][i];
510  Ubf[patchi][patchFacei] +=
511  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
512  }
513  }
514 }
515 
516 
518 {
519  makeAbsoluteRhoFlux(geometricOneField(), phi);
520 }
521 
522 
524 (
525  const surfaceScalarField& rho,
526  surfaceScalarField& phi
527 ) const
528 {
529  makeAbsoluteRhoFlux(rho, phi);
530 }
531 
532 
534 {
535  const vector Omega = this->Omega();
536 
537  // Included patches
538  volVectorField::Boundary& Ubf = U.boundaryFieldRef();
539 
540  forAll(includedFaces_, patchi)
541  {
542  const vectorField& patchC = mesh_.Cf().boundaryField()[patchi];
543 
544  vectorField pfld(Ubf[patchi]);
545 
546  forAll(includedFaces_[patchi], i)
547  {
548  label patchFacei = includedFaces_[patchi][i];
549 
550  pfld[patchFacei] = (Omega ^ (patchC[patchFacei] - origin_));
551  }
552 
553  Ubf[patchi] == pfld;
554  }
555 }
556 
557 
559 {
560  os << nl;
561  os.write(name_) << nl;
562  os << token::BEGIN_BLOCK << incrIndent << nl;
563  writeEntry(os, "cellZone", cellZoneName_);
564  writeEntry(os, "origin", origin_);
565  writeEntry(os, "axis", axis_);
566  omega_->write(os);
567 
568  if (excludedPatchNames_.size())
569  {
570  writeEntry(os, "nonRotatingPatches", excludedPatchNames_);
571  }
572 
573  os << decrIndent << token::END_BLOCK << nl;
574 }
575 
576 
578 {
579  coeffs_ = dict;
580 
581  coeffs_.lookup("cellZone") >> cellZoneName_;
582  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
583 
584  return true;
585 }
586 
587 
589 {
590  if (mesh_.topoChanging())
591  {
592  setMRFFaces();
593  }
594 }
595 
596 
597 // ************************************************************************* //
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:577
Foam::surfaceFields.
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:303
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:82
Run-time selectable general function of one variable.
Definition: Function1.H:52
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:398
dictionary dict
Graphite solid properties.
Definition: C.H:48
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
virtual Ostream & write(const char)=0
Write character.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
MRFZone(const word &name, const fvMesh &mesh, const dictionary &dict, const word &cellZoneName=word::null)
Construct from fvMesh.
Definition: MRFZone.C:237
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
const Boundary & boundaryField() const
Return const-reference to the boundary field.
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:282
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:472
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
void update()
Update MRFZone faces if the mesh topology changes.
Definition: MRFZone.C:588
A class representing the concept of a field of oneFields used to avoid unnecessary manipulations for ...
Definition: oneFieldField.H:50
Generic field type.
Definition: FieldField.H:51
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
C()
Construct null.
Definition: C.C:40
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
const cellShapeList & cells
A class for handling words, derived from string.
Definition: word.H:59
static const word null
An empty word.
Definition: word.H:77
A special matrix type and solver, designed for finite volume solutions of scalar equations. Face addressing is used to make all matrix assembly and solution loops vectorise.
Definition: fvPatchField.H:72
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const zero Zero
Definition: zero.H:97
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:54
Internal::FieldType & primitiveFieldRef()
Return a reference to the internal field.
static const char nl
Definition: Ostream.H:260
Field< Type > & source()
Definition: fvMatrix.H:292
defineTypeNameAndDebug(combustionModel, 0)
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
void writeEntry(Ostream &os, const HashTable< T, Key, Hash > &ht)
Definition: HashTableIO.C:96
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:533
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:296
A class representing the concept of a field of 1 used to avoid unnecessary manipulations for objects ...
Definition: oneField.H:50
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronise values on all mesh faces.
Definition: syncTools.H:381
label patchi
U
Definition: pEqn.H:72
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
A special matrix type and solver, designed for finite volume solutions of scalar equations.
dimensioned< scalar > mag(const dimensioned< Type > &)
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:558
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:844