faceZone.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 "faceZone.H"
28 #include "faceZoneMesh.H"
29 #include "polyMesh.H"
30 #include "primitiveMesh.H"
31 #include "demandDrivenData.H"
32 #include "mapPolyMesh.H"
33 #include "syncTools.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39  defineTypeNameAndDebug(faceZone, 0);
40  defineRunTimeSelectionTable(faceZone, dictionary);
41  addToRunTimeSelectionTable(faceZone, faceZone, dictionary);
42 }
43 
44 const char* const Foam::faceZone::labelsName = "faceLabels";
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 {
50  if (debug)
51  {
52  Info<< "void faceZone::calcFaceZonePatch() const : "
53  << "Calculating primitive patch"
54  << endl;
55  }
56 
57  if (patchPtr_)
58  {
60  << "primitive face zone patch already calculated"
61  << abort(FatalError);
62  }
63 
64  patchPtr_ =
66  (
67  faceList(size()),
68  zoneMesh().mesh().points()
69  );
70 
71  primitiveFacePatch& patch = *patchPtr_;
72 
73  const faceList& f = zoneMesh().mesh().faces();
74 
75  const labelList& addr = *this;
76  const boolList& flip = flipMap();
77 
78  forAll(addr, faceI)
79  {
80  if (flip[faceI])
81  {
82  patch[faceI] = f[addr[faceI]].reverseFace();
83  }
84  else
85  {
86  patch[faceI] = f[addr[faceI]];
87  }
88  }
89 
90  if (debug)
91  {
92  Info<< "void faceZone::calcFaceZonePatch() const : "
93  << "Finished calculating primitive patch"
94  << endl;
95  }
96 }
97 
98 
100 {
101  if (debug)
102  {
103  Info<< "void Foam::faceZone::calcCellLayers() const : "
104  << "calculating master cells"
105  << endl;
106  }
107 
108  // It is an error to attempt to recalculate edgeCells
109  // if the pointer is already set
110  if (masterCellsPtr_ || slaveCellsPtr_)
111  {
113  << "cell layers already calculated"
114  << abort(FatalError);
115  }
116  else
117  {
118  // Go through all the faces in the master zone. Choose the
119  // master or slave cell based on the face flip
120 
121  const labelList& own = zoneMesh().mesh().faceOwner();
122  const labelList& nei = zoneMesh().mesh().faceNeighbour();
123 
124  const labelList& mf = *this;
125 
126  const boolList& faceFlip = flipMap();
127 
128  masterCellsPtr_ = new labelList(mf.size());
129  labelList& mc = *masterCellsPtr_;
130 
131  slaveCellsPtr_ = new labelList(mf.size());
132  labelList& sc = *slaveCellsPtr_;
133 
134  forAll(mf, faceI)
135  {
136  label ownCellI = own[mf[faceI]];
137  label neiCellI =
138  (
139  zoneMesh().mesh().isInternalFace(mf[faceI])
140  ? nei[mf[faceI]]
141  : -1
142  );
143 
144  if (!faceFlip[faceI])
145  {
146  // Face is oriented correctly, no flip needed
147  mc[faceI] = neiCellI;
148  sc[faceI] = ownCellI;
149  }
150  else
151  {
152  mc[faceI] = ownCellI;
153  sc[faceI] = neiCellI;
154  }
155  }
156  //Info<< "masterCells: " << mc << endl;
157  //Info<< "slaveCells: " << sc << endl;
158  }
159 }
160 
161 
163 {
164  if (size() != flipMap_.size())
165  {
167  << "Size of addressing: " << size()
168  << " size of flip map: " << flipMap_.size()
169  << abort(FatalError);
170  }
171 
172  const labelList& mf = *this;
173 
174  // Note: nFaces, nCells might not be set yet on mesh so use owner size
175  const label nFaces = zoneMesh().mesh().faceOwner().size();
176 
177  bool hasWarned = false;
178  forAll(mf, i)
179  {
180  if (!hasWarned && (mf[i] < 0 || mf[i] >= nFaces))
181  {
183  << "Illegal face index " << mf[i] << " outside range 0.."
184  << nFaces-1 << endl;
185  hasWarned = true;
186  }
187  }
188 }
189 
190 
191 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
192 
193 // Construct from components
195 (
196  const word& name,
197  const labelUList& addr,
198  const boolList& fm,
199  const label index,
200  const faceZoneMesh& zm
201 )
202 :
203  zone(name, addr, index),
204  flipMap_(fm),
205  zoneMesh_(zm),
206  patchPtr_(NULL),
207  masterCellsPtr_(NULL),
208  slaveCellsPtr_(NULL),
209  mePtr_(NULL)
210 {
211  checkAddressing();
212 }
213 
214 
216 (
217  const word& name,
218  const Xfer<labelList>& addr,
219  const Xfer<boolList>& fm,
220  const label index,
221  const faceZoneMesh& zm
222 )
223 :
224  zone(name, addr, index),
225  flipMap_(fm),
226  zoneMesh_(zm),
227  patchPtr_(NULL),
228  masterCellsPtr_(NULL),
229  slaveCellsPtr_(NULL),
230  mePtr_(NULL)
231 {
232  checkAddressing();
233 }
234 
235 
237 (
238  const word& name,
239  const dictionary& dict,
240  const label index,
241  const faceZoneMesh& zm
242 )
243 :
244  zone(name, dict, this->labelsName, index),
245  flipMap_(dict.lookup("flipMap")),
246  zoneMesh_(zm),
247  patchPtr_(NULL),
248  masterCellsPtr_(NULL),
249  slaveCellsPtr_(NULL),
250  mePtr_(NULL)
251 {
252  checkAddressing();
253 }
254 
255 
257 (
258  const faceZone& fz,
259  const labelUList& addr,
260  const boolList& fm,
261  const label index,
262  const faceZoneMesh& zm
263 )
264 :
265  zone(fz, addr, index),
266  flipMap_(fm),
267  zoneMesh_(zm),
268  patchPtr_(NULL),
269  masterCellsPtr_(NULL),
270  slaveCellsPtr_(NULL),
271  mePtr_(NULL)
272 {
273  checkAddressing();
274 }
275 
276 
278 (
279  const faceZone& fz,
280  const Xfer<labelList>& addr,
281  const Xfer<boolList>& fm,
282  const label index,
283  const faceZoneMesh& zm
284 )
285 :
286  zone(fz, addr, index),
287  flipMap_(fm),
288  zoneMesh_(zm),
289  patchPtr_(NULL),
290  masterCellsPtr_(NULL),
291  slaveCellsPtr_(NULL),
292  mePtr_(NULL)
293 {
294  checkAddressing();
295 }
296 
297 
298 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
299 
301 {
302  clearAddressing();
303 }
304 
305 
306 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
307 
309 {
310  return zoneMesh_;
311 }
312 
313 
314 Foam::label Foam::faceZone::whichFace(const label globalFaceID) const
315 {
316  return zone::localID(globalFaceID);
317 }
318 
319 
321 {
322  if (!patchPtr_)
323  {
324  calcFaceZonePatch();
325  }
326 
327  return *patchPtr_;
328 }
329 
330 
332 {
333  if (!masterCellsPtr_)
334  {
335  calcCellLayers();
336  }
337 
338  return *masterCellsPtr_;
339 }
340 
341 
343 {
344  if (!slaveCellsPtr_)
345  {
346  calcCellLayers();
347  }
348 
349  return *slaveCellsPtr_;
350 }
351 
352 
354 {
355  if (!mePtr_)
356  {
357  //labelList faceCells(size());
358  //
359  //const labelList& own = zoneMesh().mesh().faceOwner();
360  //
361  //const labelList& faceLabels = *this;
362  //
363  //forAll(faceCells, faceI)
364  //{
365  // faceCells[faceI] = own[faceLabels[faceI]];
366  //}
367  //
368  //mePtr_ =
369  // new labelList
370  // (
371  // operator()().meshEdges
372  // (
373  // zoneMesh().mesh().edges(),
374  // zoneMesh().mesh().cellEdges(),
375  // faceCells
376  // )
377  // );
378 
379  mePtr_ =
380  new labelList
381  (
382  operator()().meshEdges
383  (
384  zoneMesh().mesh().edges(),
385  zoneMesh().mesh().pointEdges()
386  )
387  );
388  }
389 
390  return *mePtr_;
391 }
392 
393 
395 {
397 
398  deleteDemandDrivenData(patchPtr_);
399 
400  deleteDemandDrivenData(masterCellsPtr_);
401  deleteDemandDrivenData(slaveCellsPtr_);
402 
403  deleteDemandDrivenData(mePtr_);
404 }
405 
406 
408 (
409  const labelUList& addr,
410  const boolList& flipMap
411 )
412 {
413  clearAddressing();
414  labelList::operator=(addr);
415  flipMap_ = flipMap;
416 }
417 
418 
420 {
421  clearAddressing();
422 
423  labelList newAddressing(size());
424  boolList newFlipMap(flipMap_.size());
425  label nFaces = 0;
426 
427  const labelList& faceMap = mpm.reverseFaceMap();
428 
429  forAll(*this, i)
430  {
431  const label faceI = operator[](i);
432 
433  if (faceMap[faceI] >= 0)
434  {
435  newAddressing[nFaces] = faceMap[faceI];
436  newFlipMap[nFaces] = flipMap_[i]; // Keep flip map.
437  nFaces++;
438  }
439  }
440 
441  newAddressing.setSize(nFaces);
442  newFlipMap.setSize(nFaces);
443 
444  transfer(newAddressing);
445  flipMap_.transfer(newFlipMap);
446 }
447 
448 
449 bool Foam::faceZone::checkDefinition(const bool report) const
450 {
451  return zone::checkDefinition(zoneMesh().mesh().faces().size(), report);
452 }
453 
454 
455 bool Foam::faceZone::checkParallelSync(const bool report) const
456 {
457  const polyMesh& mesh = zoneMesh().mesh();
458  const polyBoundaryMesh& bm = mesh.boundaryMesh();
459 
460  bool hasError = false;
461 
462 
463  // Check that zone faces are synced
464  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
465 
466  {
467  boolList neiZoneFace(mesh.nFaces()-mesh.nInternalFaces(), false);
468  boolList neiZoneFlip(mesh.nFaces()-mesh.nInternalFaces(), false);
469  forAll(*this, i)
470  {
471  const label faceI = operator[](i);
472 
473  if (!mesh.isInternalFace(faceI))
474  {
475  neiZoneFace[faceI-mesh.nInternalFaces()] = true;
476  neiZoneFlip[faceI-mesh.nInternalFaces()] = flipMap()[i];
477  }
478  }
479  boolList myZoneFace(neiZoneFace);
481  boolList myZoneFlip(neiZoneFlip);
483 
484  forAll(*this, i)
485  {
486  const label faceI = operator[](i);
487  const label patchI = bm.whichPatch(faceI);
488 
489  if (patchI != -1 && bm[patchI].coupled())
490  {
491  const label bFaceI = faceI-mesh.nInternalFaces();
492 
493  // Check face in zone on both sides
494  if (myZoneFace[bFaceI] != neiZoneFace[bFaceI])
495  {
496  hasError = true;
497 
498  if (report)
499  {
500  Pout<< " ***Problem with faceZone " << index()
501  << " named " << name()
502  << ". Face " << faceI
503  << " on coupled patch "
504  << bm[patchI].name()
505  << " is not consistent with its coupled neighbour."
506  << endl;
507  }
508  else
509  {
510  // w/o report - can stop checking now
511  break;
512  }
513  }
514  else if (myZoneFlip[bFaceI] == neiZoneFlip[bFaceI])
515  {
516  // Flip state should be opposite.
517  hasError = true;
518 
519  if (report)
520  {
521  Pout<< " ***Problem with faceZone " << index()
522  << " named " << name()
523  << ". Face " << faceI
524  << " on coupled patch "
525  << bm[patchI].name()
526  << " does not have consistent flipMap"
527  << " across coupled faces."
528  << endl;
529  }
530  else
531  {
532  // w/o report - can stop checking now
533  break;
534  }
535  }
536  }
537  }
538  }
539 
540  return returnReduce(hasError, orOp<bool>());
541 }
542 
543 
545 {
546  if (patchPtr_)
547  {
548  patchPtr_->movePoints(p);
549  }
550 }
551 
553 {
554  os << nl << name()
555  << nl << static_cast<const labelList&>(*this)
556  << nl << flipMap();
557 }
558 
559 
561 {
562  os << nl << name() << nl << token::BEGIN_BLOCK << nl
563  << " type " << type() << token::END_STATEMENT << nl;
564 
565  writeEntry(this->labelsName, os);
566  flipMap().writeEntry("flipMap", os);
567 
568  os << token::END_BLOCK << endl;
569 }
570 
571 
572 // * * * * * * * * * * * * * * * Ostream Operator * * * * * * * * * * * * * //
573 
575 {
576  zn.write(os);
577  os.check("Ostream& operator<<(Ostream&, const faceZone&");
578  return os;
579 }
580 
581 
582 // ************************************************************************* //
Foam::faceZone::clearAddressing
virtual void clearAddressing()
Clear addressing.
Definition: faceZone.C:394
Foam::token::END_STATEMENT
@ END_STATEMENT
Definition: token.H:99
Foam::faceZone::labelsName
static const char *const labelsName
The name associated with the zone-labels dictionary entry.
Definition: faceZone.H:126
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::faceZone::write
virtual void write(Ostream &) const
Write.
Definition: faceZone.C:552
p
p
Definition: pEqn.H:62
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::faceZone::checkAddressing
void checkAddressing() const
Check addressing.
Definition: faceZone.C:162
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::returnReduce
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
Definition: PstreamReduceOps.H:86
Foam::polyBoundaryMesh
Foam::polyBoundaryMesh.
Definition: polyBoundaryMesh.H:60
Foam::faceZone::meshEdges
const labelList & meshEdges() const
Return global edge index for local edges.
Definition: faceZone.C:353
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::ZoneMesh::mesh
const MeshType & mesh() const
Return the mesh reference.
Definition: ZoneMesh.H:127
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::defineRunTimeSelectionTable
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
Foam::faceZone::movePoints
virtual void movePoints(const pointField &)
Correct patch after moving points.
Definition: faceZone.C:544
mapPolyMesh.H
Foam::dictionary::lookup
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:449
Foam::zone
Base class for zones.
Definition: zone.H:57
Foam::faceZone::faceZone
faceZone(const faceZone &)
Disallow default bitwise copy construct.
Foam::faceZone::resetAddressing
virtual void resetAddressing(const labelUList &, const boolList &)
Reset addressing and flip map (clearing demand-driven data)
Definition: faceZone.C:408
Foam::zone::clearAddressing
virtual void clearAddressing()
Clear addressing.
Definition: zone.C:187
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
polyMesh.H
Foam::List::operator=
void operator=(const UList< T > &)
Assignment from UList operator. Takes linear time.
syncTools.H
Foam::faceZone::~faceZone
virtual ~faceZone()
Destructor.
Definition: faceZone.C:300
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::primitiveFacePatch
PrimitivePatch< face, List, const pointField & > primitiveFacePatch
Foam::primitiveFacePatch.
Definition: primitiveFacePatch.H:45
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
Foam::faceZone::patchPtr_
primitiveFacePatch * patchPtr_
Primitive patch made out of correctly flipped faces.
Definition: faceZone.H:88
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::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::orOp
Definition: ops.H:178
Foam::nl
static const char nl
Definition: Ostream.H:260
Foam::Info
messageStream Info
Foam::faceZone::checkDefinition
virtual bool checkDefinition(const bool report=false) const
Check zone definition. Return true if in error.
Definition: faceZone.C:449
Foam::faceZone::slaveCells
const labelList & slaveCells() const
Return labels of slave cells.
Definition: faceZone.C:342
Foam::ZoneMesh< faceZone, polyMesh >
Foam::faceZone::zoneMesh
const faceZoneMesh & zoneMesh() const
Return zoneMesh reference.
Definition: faceZone.C:308
Foam::polyBoundaryMesh::whichPatch
label whichPatch(const label faceIndex) const
Return patch index for a given face label.
Definition: polyBoundaryMesh.C:703
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::faceZone::updateMesh
virtual void updateMesh(const mapPolyMesh &)
Update for changes in topology.
Definition: faceZone.C:419
Foam::operator<<
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:130
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::IOstream::check
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
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::primitiveMesh::nInternalFaces
label nInternalFaces() const
Definition: primitiveMeshI.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
faceZoneMesh.H
Foam::faceZoneMesh.
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::faceZone::whichFace
label whichFace(const label globalCellID) const
Helper function to re-direct to zone::localID(...)
Definition: faceZone.C:314
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::List::size
label size() const
Return the number of elements in the UList.
Foam::mapPolyMesh::reverseFaceMap
const labelList & reverseFaceMap() const
Reverse face map.
Definition: mapPolyMesh.H:497
Foam::List::setSize
void setSize(const label)
Reset size of List.
Foam::faceZone::checkParallelSync
virtual bool checkParallelSync(const bool report=false) const
Check whether all procs have faces synchronised. Return.
Definition: faceZone.C:455
Foam::faceZone::masterCells
const labelList & masterCells() const
Return labels of master cells (cells next to the master face.
Definition: faceZone.C:331
Foam::syncTools::swapBoundaryFaceList
static void swapBoundaryFaceList(const polyMesh &mesh, UList< T > &l)
Swap coupled boundary face values.
Definition: syncTools.H:430
Foam::faceZone::operator()
const primitiveFacePatch & operator()() const
Return reference to primitive patch.
Definition: faceZone.C:320
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::faceZone::calcCellLayers
void calcCellLayers() const
Calculate master and slave face layer.
Definition: faceZone.C:99
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
Foam::faceZone::calcFaceZonePatch
void calcFaceZonePatch() const
Build primitive patch.
Definition: faceZone.C:48
Foam::token::END_BLOCK
@ END_BLOCK
Definition: token.H:105
faceZone.H
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::faceZone::writeDict
virtual void writeDict(Ostream &) const
Write dictionary.
Definition: faceZone.C:560
Foam::primitiveMesh::isInternalFace
bool isInternalFace(const label faceIndex) const
Return true if given face label is internal to the mesh.
Definition: primitiveMeshI.H:70
Foam::UList
A 1D vector of objects of type <T>, where the size of the vector is known and can be used for subscri...
Definition: HashTable.H:60
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::zone::localID
label localID(const label globalID) const
Map storing the local index for every global index. Used to find.
Definition: zone.C:170
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::token::BEGIN_BLOCK
@ BEGIN_BLOCK
Definition: token.H:104
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:53
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::faceZone::flipMap
const boolList & flipMap() const
Return face flip map.
Definition: faceZone.H:249
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88
Foam::zone::checkDefinition
virtual bool checkDefinition(const bool report=false) const =0
Check zone definition. Return true if in error.