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 
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  // Synchronize 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 
236 Foam::MRFZone::MRFZone
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  (
294  "MRFZone"
295  "("
296  "const word&, "
297  "const fvMesh&, "
298  "const dictionary&, "
299  "const word&"
300  ")"
301  )
302  << "cannot find MRF cellZone " << cellZoneName_
303  << exit(FatalError);
304  }
305 
306  setMRFFaces();
307  }
308 }
309 
310 
311 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
312 
314 {
315  return omega_->value(mesh_.time().timeOutputValue())*axis_;
316 }
317 
318 
320 (
321  const volVectorField& U,
322  volVectorField& ddtU
323 ) const
324 {
325  if (cellZoneID_ == -1)
326  {
327  return;
328  }
329 
330  const labelList& cells = mesh_.cellZones()[cellZoneID_];
331  vectorField& ddtUc = ddtU.internalField();
332  const vectorField& Uc = U.internalField();
333 
334  const vector Omega = this->Omega();
335 
336  forAll(cells, i)
337  {
338  label celli = cells[i];
339  ddtUc[celli] += (Omega ^ Uc[celli]);
340  }
341 }
342 
343 
344 void Foam::MRFZone::addCoriolis(fvVectorMatrix& UEqn, const bool rhs) const
345 {
346  if (cellZoneID_ == -1)
347  {
348  return;
349  }
350 
351  const labelList& cells = mesh_.cellZones()[cellZoneID_];
352  const scalarField& V = mesh_.V();
353  vectorField& Usource = UEqn.source();
354  const vectorField& U = UEqn.psi();
355 
356  const vector Omega = this->Omega();
357 
358  if (rhs)
359  {
360  forAll(cells, i)
361  {
362  label celli = cells[i];
363  Usource[celli] += V[celli]*(Omega ^ U[celli]);
364  }
365  }
366  else
367  {
368  forAll(cells, i)
369  {
370  label celli = cells[i];
371  Usource[celli] -= V[celli]*(Omega ^ U[celli]);
372  }
373  }
374 }
375 
376 
378 (
379  const volScalarField& rho,
380  fvVectorMatrix& UEqn,
381  const bool rhs
382 ) const
383 {
384  if (cellZoneID_ == -1)
385  {
386  return;
387  }
388 
389  const labelList& cells = mesh_.cellZones()[cellZoneID_];
390  const scalarField& V = mesh_.V();
391  vectorField& Usource = UEqn.source();
392  const vectorField& U = UEqn.psi();
393 
394  const vector Omega = this->Omega();
395 
396  if (rhs)
397  {
398  forAll(cells, i)
399  {
400  label celli = cells[i];
401  Usource[celli] += V[celli]*rho[celli]*(Omega ^ U[celli]);
402  }
403  }
404  else
405  {
406  forAll(cells, i)
407  {
408  label celli = cells[i];
409  Usource[celli] -= V[celli]*rho[celli]*(Omega ^ U[celli]);
410  }
411  }
412 }
413 
414 
416 {
417  const volVectorField& C = mesh_.C();
418 
419  const vector Omega = this->Omega();
420 
421  const labelList& cells = mesh_.cellZones()[cellZoneID_];
422 
423  forAll(cells, i)
424  {
425  label celli = cells[i];
426  U[celli] -= (Omega ^ (C[celli] - origin_));
427  }
428 
429  // Included patches
430  forAll(includedFaces_, patchi)
431  {
432  forAll(includedFaces_[patchi], i)
433  {
434  label patchFacei = includedFaces_[patchi][i];
435  U.boundaryField()[patchi][patchFacei] = vector::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  U.boundaryField()[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 
466 (
467  const surfaceScalarField& rho,
468  surfaceScalarField& phi
469 ) const
470 {
471  makeRelativeRhoFlux(rho, phi);
472 }
473 
474 
476 {
477  const volVectorField& C = mesh_.C();
478 
479  const vector Omega = this->Omega();
480 
481  const labelList& cells = mesh_.cellZones()[cellZoneID_];
482 
483  forAll(cells, i)
484  {
485  label celli = cells[i];
486  U[celli] += (Omega ^ (C[celli] - origin_));
487  }
488 
489  // Included patches
490  forAll(includedFaces_, patchi)
491  {
492  forAll(includedFaces_[patchi], i)
493  {
494  label patchFacei = includedFaces_[patchi][i];
495  U.boundaryField()[patchi][patchFacei] =
496  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
497  }
498  }
499 
500  // Excluded patches
501  forAll(excludedFaces_, patchi)
502  {
503  forAll(excludedFaces_[patchi], i)
504  {
505  label patchFacei = excludedFaces_[patchi][i];
506  U.boundaryField()[patchi][patchFacei] +=
507  (Omega ^ (C.boundaryField()[patchi][patchFacei] - origin_));
508  }
509  }
510 }
511 
512 
514 {
515  makeAbsoluteRhoFlux(geometricOneField(), phi);
516 }
517 
518 
520 (
521  const surfaceScalarField& rho,
522  surfaceScalarField& phi
523 ) const
524 {
525  makeAbsoluteRhoFlux(rho, phi);
526 }
527 
528 
530 {
531  const vector Omega = this->Omega();
532 
533  // Included patches
534  forAll(includedFaces_, patchi)
535  {
536  const vectorField& patchC = mesh_.Cf().boundaryField()[patchi];
537 
538  vectorField pfld(U.boundaryField()[patchi]);
539 
540  forAll(includedFaces_[patchi], i)
541  {
542  label patchFacei = includedFaces_[patchi][i];
543 
544  pfld[patchFacei] = (Omega ^ (patchC[patchFacei] - origin_));
545  }
546 
547  U.boundaryField()[patchi] == pfld;
548  }
549 }
550 
551 
553 {
554  os << nl;
555  os.write(name_) << nl;
556  os << token::BEGIN_BLOCK << incrIndent << nl;
557  os.writeKeyword("active") << active_ << token::END_STATEMENT << nl;
558  os.writeKeyword("cellZone") << cellZoneName_ << token::END_STATEMENT << nl;
559  os.writeKeyword("origin") << origin_ << token::END_STATEMENT << nl;
560  os.writeKeyword("axis") << axis_ << token::END_STATEMENT << nl;
561  omega_->writeData(os);
562 
563  if (excludedPatchNames_.size())
564  {
565  os.writeKeyword("nonRotatingPatches") << excludedPatchNames_
566  << token::END_STATEMENT << nl;
567  }
568 
569  os << decrIndent << token::END_BLOCK << nl;
570 }
571 
572 
574 {
575  coeffs_ = dict;
576 
577  active_ = coeffs_.lookupOrDefault("active", true);
578  coeffs_.lookup("cellZone") >> cellZoneName_;
579  cellZoneID_ = mesh_.cellZones().findZoneID(cellZoneName_);
580 
581  return true;
582 }
583 
584 
585 // ************************************************************************* //
void reduce(const List< UPstream::commsStruct > &comms, T &Value, const BinaryOp &bop, const int tag, const label comm)
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::surfaceFields.
Generic field type.
Definition: FieldField.H:51
vector Omega() const
Return the current Omega vector.
Definition: MRFZone.C:313
dimensioned< scalar > mag(const dimensioned< Type > &)
GeometricBoundaryField & boundaryField()
Return reference to GeometricBoundaryField.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
C()
Construct null.
Definition: C.C:41
void correctBoundaryVelocity(volVectorField &U) const
Correct the boundary velocity for the rotation of the MRF region.
Definition: MRFZone.C:529
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
bool read(const dictionary &dict)
Read MRF dictionary.
Definition: MRFZone.C:573
A class for handling words, derived from string.
Definition: word.H:59
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
InternalField & internalField()
Return internal field.
void makeRelative(volVectorField &U) const
Make the given absolute velocity relative within the MRF region.
Definition: MRFZone.C:415
Graphite solid properties.
Definition: C.H:57
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:68
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
static void syncFaceList(const polyMesh &mesh, UList< T > &l, const CombineOp &cop)
Synchronize values on all mesh faces.
Definition: syncTools.H:380
Namespace for OpenFOAM.
void addCoriolis(const volVectorField &U, volVectorField &ddtU) const
Add the Coriolis force contribution to the acceleration field.
Definition: MRFZone.C:320
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:39
Ostream & incrIndent(Ostream &os)
Increment the indent level.
Definition: Ostream.H:228
dictionary dict
const GeometricField< Type, fvPatchField, volMesh > & psi() const
Definition: fvMatrix.H:281
static const char nl
Definition: Ostream.H:260
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Ostream & decrIndent(Ostream &os)
Decrement the indent level.
Definition: Ostream.H:235
void makeAbsolute(volVectorField &U) const
Make the given relative velocity absolute within the MRF region.
Definition: MRFZone.C:475
A class representing the concept of a field of oneFields used to avoid unnecessary manipulations for ...
Definition: oneFieldField.H:50
void writeData(Ostream &os) const
Write.
Definition: MRFZone.C:552
#define forAll(list, i)
Definition: UList.H:421
virtual Ostream & write(const token &)=0
Write next token to stream.
label patchi
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
Ostream & writeKeyword(const keyType &)
Write the keyword followed by an appropriate indentation.
Definition: Ostream.C:59
#define FatalErrorIn(functionName)
Report an error message using Foam::FatalError.
Definition: error.H:314
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:210
const cellShapeList & cells
error FatalError
T lookupOrDefault(const word &, const T &, bool recursive=false, bool patternMatch=true) const
Find and return a T,.
List< label > labelList
A List of labels.
Definition: labelList.H:56
static const Vector zero
Definition: Vector.H:80
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:53
Field< Type > & source()
Definition: fvMatrix.H:291
Top level data entry class for use in dictionaries. Provides a mechanism to specify a variable as a c...
Definition: DataEntry.H:52
static const word null
An empty word.
Definition: word.H:77
defineTypeNameAndDebug(combustionModel, 0)
A special matrix type and solver, designed for finite volume solutions of scalar equations.
prefixOSstream Pout(cout,"Pout")
Definition: IOstreams.H:53
List< bool > boolList
Bool container classes.
Definition: boolList.H:50