topoSet.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-2025 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 "topoSet.H"
27 #include "polyMesh.H"
28 #include "polyTopoChangeMap.H"
29 #include "Time.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
45 (
46  const word& setType,
47  const polyMesh& mesh,
48  const word& name,
49  readOption r,
50  writeOption w
51 )
52 {
53  wordConstructorTable::iterator cstrIter =
54  wordConstructorTablePtr_->find(setType);
55 
56  if (cstrIter == wordConstructorTablePtr_->end())
57  {
59  << "Unknown set type " << setType
60  << endl << endl
61  << "Valid set types : " << endl
62  << wordConstructorTablePtr_->sortedToc()
63  << exit(FatalError);
64  }
65 
66  return autoPtr<topoSet>(cstrIter()(mesh, name, r, w));
67 }
68 
69 
71 (
72  const word& setType,
73  const polyMesh& mesh,
74  const word& name,
75  const label size,
76  writeOption w
77 )
78 {
79  sizeConstructorTable::iterator cstrIter =
80  sizeConstructorTablePtr_->find(setType);
81 
82  if (cstrIter == sizeConstructorTablePtr_->end())
83  {
85  << "Unknown set type " << setType
86  << endl << endl
87  << "Valid set types : " << endl
88  << sizeConstructorTablePtr_->sortedToc()
89  << exit(FatalError);
90  }
91 
92  return autoPtr<topoSet>(cstrIter()(mesh, name, size, w));
93 }
94 
95 
97 (
98  const word& setType,
99  const polyMesh& mesh,
100  const word& name,
101  const topoSet& set,
102  writeOption w
103 )
104 {
105  setConstructorTable::iterator cstrIter =
106  setConstructorTablePtr_->find(setType);
107 
108  if (cstrIter == setConstructorTablePtr_->end())
109  {
111  << "Unknown set type " << setType
112  << endl << endl
113  << "Valid set types : " << endl
114  << setConstructorTablePtr_->sortedToc()
115  << exit(FatalError);
116  }
117 
118  return autoPtr<topoSet>(cstrIter()(mesh, name, set, w));
119 }
120 
121 
123 (
124  const polyMesh& mesh,
125  const word& name
126 )
127 {
129 }
130 
131 
132 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
133 
135 {
136  // Iterate over map to see if anything changed
137  bool changed = false;
138 
139  forAllConstIter(labelHashSet, *this, iter)
140  {
141  if ((iter.key() < 0) || (iter.key() > map.size()))
142  {
144  << "Illegal content " << iter.key() << " of set:" << name()
145  << " of type " << type() << endl
146  << "Value should be between 0 and " << map.size()-1
147  << abort(FatalError);
148  }
149 
150  const label newCelli = map[iter.key()];
151 
152  if (newCelli != iter.key())
153  {
154  changed = true;
155 
156  break;
157  }
158  }
159 
160  // Relabel (use second Map to prevent overlapping)
161  if (changed)
162  {
163  labelHashSet newSet(2*size());
164 
165  forAllConstIter(labelHashSet, *this, iter)
166  {
167  const label newCelli = map[iter.key()];
168 
169  if (newCelli >= 0)
170  {
171  newSet.insert(newCelli);
172  }
173  }
174 
175  transfer(newSet);
176  }
177 }
178 
179 
180 void Foam::topoSet::check(const label maxLabel)
181 {
182  forAllConstIter(topoSet, *this, iter)
183  {
184  if ((iter.key() < 0) || (iter.key() > maxLabel))
185  {
187  << "Illegal content " << iter.key() << " of set:" << name()
188  << " of type " << type() << endl
189  << "Value should be between 0 and " << maxLabel
190  << abort(FatalError);
191  }
192  }
193 }
194 
195 
196 void Foam::topoSet::writeDebug
197 (
198  Ostream& os,
199  const label maxElem,
201  label& elemI
202 ) const
203 {
204  label n = 0;
205 
206  for (; (iter != end()) && (n < maxElem); ++iter)
207  {
208  if ((n != 0) && ((n % 10) == 0))
209  {
210  os << endl;
211  }
212  os << iter.key() << ' ';
213 
214  n++;
215  elemI++;
216  }
217 }
218 
219 
220 void Foam::topoSet::writeDebug
221 (
222  Ostream& os,
223  const pointField& coords,
224  const label maxElem,
226  label& elemI
227 ) const
228 {
229  label n = 0;
230 
231  for (; (iter != end()) && (n < maxElem); ++iter)
232  {
233  if ((n != 0) && ((n % 3) == 0))
234  {
235  os << endl;
236  }
237  os << iter.key() << coords[iter.key()] << ' ';
238 
239  n++;
240  elemI++;
241  }
242 }
243 
244 
245 void Foam::topoSet::writeDebug
246 (
247  Ostream& os,
248  const pointField& coords,
249  const label maxLen
250 ) const
251 {
252  // Bounding box of contents.
253  boundBox bb(pointField(coords, toc()), true);
254 
255  os << "Set bounding box: min = "
256  << bb.min() << " max = " << bb.max() << " meters. " << endl << endl;
257 
258  label n = 0;
259 
260  topoSet::const_iterator iter = begin();
261 
262  if (size() <= maxLen)
263  {
264  writeDebug(os, coords, maxLen, iter, n);
265  }
266  else
267  {
268  label halfLen = maxLen/2;
269 
270  os << "Size larger than " << maxLen << ". Printing first and last "
271  << halfLen << " elements:" << endl << endl;
272 
273  writeDebug(os, coords, halfLen, iter, n);
274 
275  os << nl << " .." << nl << endl;
276 
277  for (; n < size() - halfLen; ++n)
278  {
279  ++iter;
280  }
281 
282  writeDebug(os, coords, halfLen, iter, n);
283  }
284 }
285 
286 
287 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
288 
289 Foam::topoSet::topoSet(const IOobject& obj, const word& wantedType)
290 :
291  regIOobject(obj)
292 {
293  if
294  (
297  || (
299  && headerOk()
300  )
301  )
302  {
303  if (readStream(wantedType).good())
304  {
305  readStream(wantedType) >> static_cast<labelHashSet&>(*this);
306 
307  close();
308  }
309  }
310 }
311 
312 
314 (
315  const polyMesh& mesh,
316  const word& wantedType,
317  const word& name,
318  readOption r,
319  writeOption w
320 )
321 :
323  (
324  IOobject
325  (
326  name,
327  mesh.time().findInstance
328  (
329  mesh.dbDir()/polyMesh::meshSubDir/"sets",
330  word::null,
331  r, // IOobject::MUST_READ,
332  mesh.facesInstance()
333  ),
334  polyMesh::meshSubDir/"sets",
335  mesh,
336  r,
337  w
338  )
339  )
340 {
341  if
342  (
345  || (
347  && headerOk()
348  )
349  )
350  {
351  if (readStream(wantedType).good())
352  {
353  readStream(wantedType) >> static_cast<labelHashSet&>(*this);
354 
355  close();
356  }
357  }
358 }
359 
360 
362 (
363  const polyMesh& mesh,
364  const word& name,
365  const label size,
366  writeOption w
367 )
368 :
370  (
371  IOobject
372  (
373  name,
374  mesh.time().findInstance
375  (
376  mesh.dbDir()/polyMesh::meshSubDir/"sets",
377  word::null,
378  IOobject::NO_READ,
379  mesh.facesInstance()
380  ),
381  polyMesh::meshSubDir/"sets",
382  mesh,
383  IOobject::NO_READ,
384  w
385  )
386  ),
387  labelHashSet(size)
388 {}
389 
390 
392 (
393  const polyMesh& mesh,
394  const word& name,
395  const labelHashSet& set,
396  writeOption w
397 )
398 :
400  (
401  IOobject
402  (
403  name,
404  mesh.time().findInstance
405  (
406  mesh.dbDir()/polyMesh::meshSubDir/"sets",
407  word::null,
408  IOobject::NO_READ,
409  mesh.facesInstance()
410  ),
411  polyMesh::meshSubDir/"sets",
412  mesh,
413  IOobject::NO_READ,
414  w
415  )
416  ),
417  labelHashSet(set)
418 {}
419 
420 
421 Foam::topoSet::topoSet(const IOobject& obj, const label size)
422 :
423  regIOobject(obj),
424  labelHashSet(size)
425 {}
426 
427 
429 :
430  regIOobject(obj),
431  labelHashSet(set)
432 {}
433 
434 
435 
436 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
437 
439 {}
440 
441 
442 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
443 
444 void Foam::topoSet::invert(const label maxLen)
445 {
446  // Keep copy of current set.
447  labelHashSet currentSet(*this);
448 
449  clear();
450  resize(2*(maxLen - currentSet.size()));
451 
452  for (label celli = 0; celli < maxLen; celli++)
453  {
454  if (!currentSet.found(celli))
455  {
456  insert(celli);
457  }
458  }
459 }
460 
461 
463 {
464  // Keep copy of current set.
465  labelHashSet currentSet(*this);
466 
467  clear();
468  resize(2*min(currentSet.size(), set.size()));
469 
470  forAllConstIter(labelHashSet, currentSet, iter)
471  {
472  if (set.found(iter.key()))
473  {
474  // element present in both currentSet and set.
475  insert(iter.key());
476  }
477  }
478 }
479 
480 
482 {
483  forAllConstIter(topoSet, set, iter)
484  {
485  insert(iter.key());
486  }
487 }
488 
489 
491 {
492  forAllConstIter(topoSet, set, iter)
493  {
494  erase(iter.key());
495  }
496 }
497 
498 
500 {
502 }
503 
504 
505 void Foam::topoSet::writeDebug(Ostream& os, const label maxLen) const
506 {
507  label n = 0;
508 
509  topoSet::const_iterator iter = begin();
510 
511  if (size() <= maxLen)
512  {
513  writeDebug(os, maxLen, iter, n);
514  }
515  else
516  {
517  label halfLen = maxLen/2;
518 
519  os << "Size larger than " << maxLen << ". Printing first and last "
520  << halfLen << " elements:" << endl << endl;
521 
522  writeDebug(os, halfLen, iter, n);
523 
524  os << nl << " .." << nl << endl;
525 
526  for (; n < size() - halfLen; ++n)
527  {
528  ++iter;
529  }
530 
531  writeDebug(os, halfLen, iter, n);
532  }
533 }
534 
535 
537 {
538  return (os << *this).good();
539 }
540 
541 
543 {
545 }
546 
547 
548 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
549 
551 {
553 }
554 
555 
556 // ************************************************************************* //
label n
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
void operator=(const HashSet< label, Hash< label > > &)
Assignment operator.
Definition: HashSet.C:139
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:109
HashTable< nil, label, Hash< label > >::const_iterator const_iterator
Definition: HashSet.H:67
An STL-conforming const_iterator.
Definition: HashTable.H:498
const Key & key() const
Return the Key corresponding to the iterator.
Definition: HashTableI.H:254
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
bool found(const Key &) const
Return true if hashedEntry is found in table.
Definition: HashTable.C:138
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
readOption
Enumeration defining the read options.
Definition: IOobject.H:117
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
readOption readOpt() const
Definition: IOobject.H:357
bool good() const
Definition: IOobject.H:470
const word & name() const
Return name.
Definition: IOobject.H:307
writeOption
Enumeration defining the write options.
Definition: IOobject.H:126
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
A class for handling file names.
Definition: fileName.H:82
const fileName & dbDir() const
Local directory path of this objectRegistry relative to the time.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
const fileName & facesInstance() const
Return the current instance directory for faces.
Definition: polyMesh.C:994
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:273
Class containing mesh-to-mesh mapping information after a change in polyMesh topology.
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:55
void close()
Close Istream.
bool headerOk()
Read and check header info.
Definition: regIOobject.C:453
Istream & readStream(const word &, const bool read=true)
Return Istream and check object type against that given.
General set of labels of mesh quantity (points, cells, faces).
Definition: topoSet.H:61
virtual bool writeData(Ostream &) const
Write contents.
Definition: topoSet.C:536
virtual void invert(const label maxLen)
Invert contents.
Definition: topoSet.C:444
static fileName localPath(const polyMesh &mesh, const word &name)
Name of file set will use.
Definition: topoSet.C:123
void operator=(const topoSet &)
Copy labelHashSet part only.
Definition: topoSet.C:550
virtual ~topoSet()
Destructor.
Definition: topoSet.C:438
void check(const label maxLabel)
Check validity of contents.
Definition: topoSet.C:180
void updateLabels(const labelList &map)
Update map from map. Used to update cell/face labels.
Definition: topoSet.C:134
topoSet(const IOobject &, const word &wantedType)
Construct from IOobject as explicitly passed type.
Definition: topoSet.C:289
virtual void deleteSet(const topoSet &set)
Delete elements present in set.
Definition: topoSet.C:490
virtual void sync(const polyMesh &mesh)
Sync set across coupled patches.
Definition: topoSet.C:499
virtual void addSet(const topoSet &set)
Add elements present in set.
Definition: topoSet.C:481
static autoPtr< topoSet > New(const word &setType, const polyMesh &mesh, const word &name, readOption r=MUST_READ, writeOption w=NO_WRITE)
Return a pointer to a toposet read from file.
Definition: topoSet.C:45
virtual void subset(const topoSet &set)
Subset contents. Only elements present in both sets remain.
Definition: topoSet.C:462
virtual void topoChange(const polyTopoChangeMap &map)
Update any stored data for new labels. Not implemented.
Definition: topoSet.C:542
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
tUEqn clear()
srcOptions erase("case")
void insert(const scalar, DynamicList< floatScalar > &)
Append scalar to given DynamicList.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
errorManip< error > abort(error &err)
Definition: errorManip.H:131
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
error FatalError
static const char nl
Definition: Ostream.H:267
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
triSurfaceToAgglom resize(surfacesMesh.size())