faceZoneSet.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 "faceZoneSet.H"
27 #include "mapPolyMesh.H"
28 #include "polyMesh.H"
29 
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38 
39 defineTypeNameAndDebug(faceZoneSet, 0);
40 
41 addToRunTimeSelectionTable(topoSet, faceZoneSet, word);
42 addToRunTimeSelectionTable(topoSet, faceZoneSet, size);
43 addToRunTimeSelectionTable(topoSet, faceZoneSet, set);
44 
45 
46 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
47 
49 {
50  labelList order;
51  sortedOrder(addressing_, order);
54 
58  {
60  }
61 }
62 
63 
64 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
65 
67 (
68  const polyMesh& mesh,
69  const word& name,
70  readOption r,
72 )
73 :
74  faceSet(mesh, name, 1000), // do not read faceSet
75  mesh_(mesh),
76  addressing_(0),
77  flipMap_(0)
78 {
79  const faceZoneMesh& faceZones = mesh.faceZones();
80  label zoneID = faceZones.findZoneID(name);
81 
82  if
83  (
84  (r == IOobject::MUST_READ)
86  || (r == IOobject::READ_IF_PRESENT && zoneID != -1)
87  )
88  {
89  const faceZone& fz = faceZones[zoneID];
90  addressing_ = fz;
91  flipMap_ = fz.flipMap();
92  }
93 
94  updateSet();
95 
96  check(mesh.nFaces());
97 }
98 
99 
101 (
102  const polyMesh& mesh,
103  const word& name,
104  const label size,
105  writeOption w
106 )
107 :
108  faceSet(mesh, name, size, w),
109  mesh_(mesh),
110  addressing_(0),
111  flipMap_(0)
112 {
113  updateSet();
114 }
115 
116 
118 (
119  const polyMesh& mesh,
120  const word& name,
121  const topoSet& set,
122  writeOption w
123 )
124 :
125  faceSet(mesh, name, set.size(), w),
126  mesh_(mesh),
127  addressing_(refCast<const faceZoneSet>(set).addressing()),
128  flipMap_(refCast<const faceZoneSet>(set).flipMap())
129 {
130  updateSet();
131 }
132 
133 
134 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
135 
137 {}
138 
139 
140 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
141 
142 void faceZoneSet::invert(const label maxLen)
143 {
144  // Count
145  label n = 0;
146 
147  for (label faceI = 0; faceI < maxLen; faceI++)
148  {
149  if (!found(faceI))
150  {
151  n++;
152  }
153  }
154 
155  // Fill
157  flipMap_.setSize(n);
158  n = 0;
159 
160  for (label faceI = 0; faceI < maxLen; faceI++)
161  {
162  if (!found(faceI))
163  {
164  addressing_[n] = faceI;
165  flipMap_[n] = false; //? or true?
166  n++;
167  }
168  }
169  updateSet();
170 }
171 
172 
173 void faceZoneSet::subset(const topoSet& set)
174 {
175  label nConflict = 0;
176 
177  DynamicList<label> newAddressing(addressing_.size());
178  DynamicList<bool> newFlipMap(flipMap_.size());
179 
180  Map<label> faceToIndex(addressing_.size());
181  forAll(addressing_, i)
182  {
183  faceToIndex.insert(addressing_[i], i);
184  }
185 
186  const faceZoneSet& fSet = refCast<const faceZoneSet>(set);
187 
188  forAll(fSet.addressing(), i)
189  {
190  label faceI = fSet.addressing()[i];
191 
192  Map<label>::const_iterator iter = faceToIndex.find(faceI);
193 
194  if (iter != faceToIndex.end())
195  {
196  label index = iter();
197 
198  if (fSet.flipMap()[i] != flipMap_[index])
199  {
200  nConflict++;
201  }
202  newAddressing.append(faceI);
203  newFlipMap.append(flipMap_[index]);
204  }
205  }
206 
207  if (nConflict > 0)
208  {
210  << "subset : there are " << nConflict
211  << " faces with different orientation in faceZonesSets "
212  << name() << " and " << set.name() << endl;
213  }
214 
215  addressing_.transfer(newAddressing);
216  flipMap_.transfer(newFlipMap);
217  updateSet();
218 }
219 
220 
221 void faceZoneSet::addSet(const topoSet& set)
222 {
223  label nConflict = 0;
224 
225  DynamicList<label> newAddressing(addressing_);
226  DynamicList<bool> newFlipMap(flipMap_);
227 
228  Map<label> faceToIndex(addressing_.size());
229  forAll(addressing_, i)
230  {
231  faceToIndex.insert(addressing_[i], i);
232  }
233 
234  const faceZoneSet& fSet = refCast<const faceZoneSet>(set);
235 
236  forAll(fSet.addressing(), i)
237  {
238  label faceI = fSet.addressing()[i];
239 
240  Map<label>::const_iterator iter = faceToIndex.find(faceI);
241 
242  if (iter != faceToIndex.end())
243  {
244  label index = iter();
245 
246  if (fSet.flipMap()[i] != flipMap_[index])
247  {
248  nConflict++;
249  }
250  }
251  else
252  {
253  newAddressing.append(faceI);
254  newFlipMap.append(fSet.flipMap()[i]);
255  }
256  }
257 
258  if (nConflict > 0)
259  {
261  << "addSet : there are " << nConflict
262  << " faces with different orientation in faceZonesSets "
263  << name() << " and " << set.name() << endl;
264  }
265 
266  addressing_.transfer(newAddressing);
267  flipMap_.transfer(newFlipMap);
268  updateSet();
269 }
270 
271 
273 {
274  label nConflict = 0;
275 
276  DynamicList<label> newAddressing(addressing_.size());
277  DynamicList<bool> newFlipMap(flipMap_.size());
278 
279  const faceZoneSet& fSet = refCast<const faceZoneSet>(set);
280 
281  Map<label> faceToIndex(fSet.addressing().size());
282  forAll(fSet.addressing(), i)
283  {
284  faceToIndex.insert(fSet.addressing()[i], i);
285  }
286 
287  forAll(addressing_, i)
288  {
289  label faceI = addressing_[i];
290 
291  Map<label>::const_iterator iter = faceToIndex.find(faceI);
292 
293  if (iter != faceToIndex.end())
294  {
295  label index = iter();
296 
297  if (fSet.flipMap()[index] != flipMap_[i])
298  {
299  nConflict++;
300  }
301  }
302  else
303  {
304  // Not found in fSet so add
305  newAddressing.append(faceI);
306  newFlipMap.append(fSet.flipMap()[i]);
307  }
308  }
309 
310  if (nConflict > 0)
311  {
313  << "deleteSet : there are " << nConflict
314  << " faces with different orientation in faceZonesSets "
315  << name() << " and " << set.name() << endl;
316  }
317 
318  addressing_.transfer(newAddressing);
319  flipMap_.transfer(newFlipMap);
320  updateSet();
321 }
322 
323 
325 {}
326 
327 
329 {
330  return mesh.nFaces();
331 }
332 
333 
334 //- Write using given format, version and compression
336 (
340 ) const
341 {
342  // Write shadow faceSet
343  word oldTypeName = typeName;
344  const_cast<word&>(type()) = faceSet::typeName;
345  bool ok = faceSet::writeObject(s, v, c);
346  const_cast<word&>(type()) = oldTypeName;
347 
348  // Modify faceZone
349  faceZoneMesh& faceZones = const_cast<polyMesh&>(mesh_).faceZones();
350  label zoneID = faceZones.findZoneID(name());
351 
352  if (zoneID == -1)
353  {
354  zoneID = faceZones.size();
355 
356  faceZones.setSize(zoneID+1);
357  faceZones.set
358  (
359  zoneID,
360  new faceZone
361  (
362  name(),
363  addressing_,
364  flipMap_,
365  zoneID,
366  faceZones
367  )
368  );
369  }
370  else
371  {
372  faceZones[zoneID].resetAddressing(addressing_, flipMap_);
373  }
374  faceZones.clearAddressing();
375 
376  return ok && faceZones.write();
377 }
378 
379 
381 {
382  // faceZone
383  labelList newAddressing(addressing_.size());
384  boolList newFlipMap(flipMap_.size());
385 
386  label n = 0;
387  forAll(addressing_, i)
388  {
389  label faceI = addressing_[i];
390  label newFaceI = morphMap.reverseFaceMap()[faceI];
391  if (newFaceI >= 0)
392  {
393  newAddressing[n] = newFaceI;
394  newFlipMap[n] = flipMap_[i];
395  n++;
396  }
397  }
398  newAddressing.setSize(n);
399  newFlipMap.setSize(n);
400 
401  addressing_.transfer(newAddressing);
402  flipMap_.transfer(newFlipMap);
403 
404  updateSet();
405 }
406 
407 
409 (
410  Ostream& os,
411  const primitiveMesh& mesh,
412  const label maxLen
413 ) const
414 {
415  faceSet::writeDebug(os, mesh, maxLen);
416 }
417 
418 
419 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
420 
421 } // End namespace Foam
422 
423 // ************************************************************************* //
Foam::faceZoneSet::addressing_
labelList addressing_
Definition: faceZoneSet.H:57
Foam::faceZoneSet::addressing
const labelList & addressing() const
Definition: faceZoneSet.H:107
w
volScalarField w(IOobject("w", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE), mesh, dimensionedScalar("w", dimensionSet(0, 0, 0, 0, 0, 0, 0), 0.0))
Foam::regIOobject::writeObject
virtual bool writeObject(IOstream::streamFormat, IOstream::versionNumber, IOstream::compressionType) const
Write using given format, version and compression.
Definition: regIOobjectWrite.C:37
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::HashTable< nil, word, string::hash >::resize
void resize(const label newSize)
Resize the hash table for efficiency.
Definition: HashTable.C:436
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::DynamicList< label >
faceZoneSet.H
Foam::faceZoneSet::updateSet
void updateSet()
Sort addressing and make faceSet part consistent with addressing.
Definition: faceZoneSet.C:48
Foam::addToRunTimeSelectionTable
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
Foam::IOstream::compressionType
compressionType
Enumeration for the format of data in the stream.
Definition: IOstream.H:193
mapPolyMesh.H
Foam::List::transfer
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Foam::Map< label >
Foam::faceZoneSet::updateMesh
virtual void updateMesh(const mapPolyMesh &morphMap)
Update any stored data for new labels.
Definition: faceZoneSet.C:380
Foam::faceSet
A list of face labels.
Definition: faceSet.H:48
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
polyMesh.H
Foam::faceZoneSet::flipMap
const boolList & flipMap() const
Definition: faceZoneSet.H:118
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:74
Foam::polyMesh::faceZones
const faceZoneMesh & faceZones() const
Return face zone mesh.
Definition: polyMesh.H:463
Foam::IOstream::versionNumber
Version number type.
Definition: IOstream.H:96
Foam::IOobject::MUST_READ_IF_MODIFIED
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:109
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::sortedOrder
void sortedOrder(const UList< T > &, labelList &order)
Generate the (stable) sort order for the list.
Definition: ListOpsTemplates.C:179
Foam::faceZoneSet::deleteSet
virtual void deleteSet(const topoSet &set)
Delete elements present in set.
Definition: faceZoneSet.C:272
Foam::faceSet::writeDebug
virtual void writeDebug(Ostream &os, const primitiveMesh &, const label maxLen) const
Write maxLen items with label and coordinates.
Definition: faceSet.C:165
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::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:64
Foam::faceZoneSet::invert
virtual void invert(const label maxLen)
Invert contents. (insert all members 0..maxLen-1 which were not in.
Definition: faceZoneSet.C:142
Foam::IOobject::writeOption
writeOption
Enumeration defining the write options.
Definition: IOobject.H:115
Foam::ZoneMesh< faceZone, polyMesh >
Foam::faceZoneSet
Like faceSet but updates faceZone when writing.
Definition: faceZoneSet.H:49
Foam::IOobject::name
const word & name() const
Return name.
Definition: IOobject.H:273
Foam::faceZoneSet::writeDebug
virtual void writeDebug(Ostream &os, const primitiveMesh &, const label maxLen) const
Write maxLen items with label and coordinates.
Definition: faceZoneSet.C:409
Foam::topoSet
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
Foam::faceZoneSet::faceZoneSet
faceZoneSet(const polyMesh &mesh, const word &name, readOption r=MUST_READ, writeOption w=NO_WRITE)
Construct from objectRegistry and name.
Definition: faceZoneSet.C:67
Foam::faceZoneSet::sync
virtual void sync(const polyMesh &mesh)
Sync faceZoneSet across coupled patches.
Definition: faceZoneSet.C:324
Foam::faceZoneSet::maxSize
virtual label maxSize(const polyMesh &mesh) const
Return max index+1.
Definition: faceZoneSet.C:328
Foam::HashTable::size
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::HashTable< nil, word, string::hash >::found
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:109
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::ZoneMesh::findZoneID
label findZoneID(const word &zoneName) const
Find zone index given a name.
Definition: ZoneMesh.C:348
s
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Foam::DynamicList::append
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
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::faceZoneSet::addSet
virtual void addSet(const topoSet &set)
Add elements present in set.
Definition: faceZoneSet.C:221
Foam::primitiveMesh::nFaces
label nFaces() const
Definition: primitiveMeshI.H:58
Foam::faceZoneSet::subset
virtual void subset(const topoSet &set)
Subset contents. Only elements present in both sets remain.
Definition: faceZoneSet.C:173
Foam::HashTable< nil, word, string::hash >::clearStorage
void clearStorage()
Clear the table entries and the table itself.
Definition: HashTable.C:497
Foam::faceZoneSet::writeObject
virtual bool writeObject(IOstream::streamFormat, IOstream::versionNumber, IOstream::compressionType) const
Write faceZone.
Definition: faceZoneSet.C:336
Foam::faceZoneSet::flipMap_
boolList flipMap_
Definition: faceZoneSet.H:59
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::HashSet::insert
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:116
Foam::mapPolyMesh
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
Definition: mapPolyMesh.H:158
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::UIndirectList
A List with indirect addressing.
Definition: fvMatrix.H:106
Foam::IOobject::readOption
readOption
Enumeration defining the read options.
Definition: IOobject.H:106
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::faceZoneSet::~faceZoneSet
virtual ~faceZoneSet()
Destructor.
Definition: faceZoneSet.C:136
Foam::ZoneMesh::clearAddressing
void clearAddressing()
Clear addressing.
Definition: ZoneMesh.C:394
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::IOobject::READ_IF_PRESENT
@ READ_IF_PRESENT
Definition: IOobject.H:110
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::HashSet::set
bool set(const Key &key)
Same as insert (cannot overwrite nil content)
Definition: HashSet.H:126
Foam::IOstream::streamFormat
streamFormat
Enumeration for the format of data in the stream.
Definition: IOstream.H:86
Foam::primitiveMesh
Cell-face mesh analysis engine.
Definition: primitiveMesh.H:79