triSurface.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2016-2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "triSurface.H"
30 #include "demandDrivenData.H"
31 #include "Time.H"
32 #include "surfZoneList.H"
33 #include "MeshedSurface.H"
34 #include "ListOps.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40  defineTypeNameAndDebug(triSurface, 0);
41 }
42 
43 
44 // * * * * * * * * * * * * * Static Member Functions * * * * * * * * * * * * //
45 
47 {
48  fileName foamName(d.caseName() + ".ftr");
49 
50  // Search back through the time directories list to find the time
51  // closest to and lower than current time
52 
53  instantList ts = d.times();
54  label i;
55 
56  for (i=ts.size()-1; i>=0; i--)
57  {
58  if (ts[i].value() <= d.timeOutputValue())
59  {
60  break;
61  }
62  }
63 
64  // Noting that the current directory has already been searched
65  // for mesh data, start searching from the previously stored time directory
66 
67  if (i>=0)
68  {
69  for (label j=i; j>=0; j--)
70  {
71  if (isFile(d.path()/ts[j].name()/typeName/foamName))
72  {
73  if (debug)
74  {
75  Pout<< " triSurface::triSurfInstance(const Time& d)"
76  << "reading " << foamName
77  << " from " << ts[j].name()/typeName
78  << endl;
79  }
80 
81  return ts[j].name();
82  }
83  }
84  }
85 
86  if (debug)
87  {
88  Pout<< " triSurface::triSurfInstance(const Time& d)"
89  << "reading " << foamName
90  << " from constant/" << endl;
91  }
92  return d.constant();
93 }
94 
95 
96 Foam::List<Foam::labelledTri> Foam::triSurface::convertToTri
97 (
98  const faceList& faces,
99  const label defaultRegion
100 )
101 {
102  List<labelledTri> triFaces(faces.size());
103 
104  forAll(triFaces, facei)
105  {
106  const face& f = faces[facei];
107 
108  if (f.size() != 3)
109  {
111  << "Face at position " << facei
112  << " does not have three vertices:" << f
113  << abort(FatalError);
114  }
115 
116  labelledTri& tri = triFaces[facei];
117 
118  tri[0] = f[0];
119  tri[1] = f[1];
120  tri[2] = f[2];
121  tri.region() = defaultRegion;
122  }
123 
124  return triFaces;
125 }
126 
127 
128 Foam::List<Foam::labelledTri> Foam::triSurface::convertToTri
129 (
130  const triFaceList& faces,
131  const label defaultRegion
132 )
133 {
134  List<labelledTri> triFaces(faces.size());
135 
136  forAll(triFaces, facei)
137  {
138  const triFace& f = faces[facei];
139 
140  labelledTri& tri = triFaces[facei];
141 
142  tri[0] = f[0];
143  tri[1] = f[1];
144  tri[2] = f[2];
145  tri.region() = defaultRegion;
146  }
147 
148  return triFaces;
149 }
150 
151 
152 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
153 
154 // Remove non-triangles, double triangles.
155 void Foam::triSurface::checkTriangles(const bool verbose)
156 {
157  // Simple check on indices ok.
158  const label maxPointi = points().size() - 1;
159 
160  for (const triSurface::FaceType& f : *this)
161  {
162  for (const label verti : f)
163  {
164  if (verti < 0 || verti > maxPointi)
165  {
167  << "triangle " << f
168  << " uses point indices outside point range 0.."
169  << maxPointi
170  << exit(FatalError);
171  }
172  }
173  }
174 
175  // Two phase process
176  // 1. mark invalid faces
177  // 2. pack
178  // Done to keep numbering constant in phase 1
179 
180  // List of valid triangles
181  bitSet valid(size(), true);
182 
183  forAll(*this, facei)
184  {
185  const labelledTri& f = (*this)[facei];
186 
187  if ((f[0] == f[1]) || (f[0] == f[2]) || (f[1] == f[2]))
188  {
189  // 'degenerate' triangle check
190  valid.unset(facei);
191 
192  if (verbose)
193  {
195  << "triangle " << facei
196  << " does not have three unique vertices:\n";
197  printTriangle(Warning, " ", f, points());
198  }
199  }
200  else
201  {
202  // duplicate triangle check
203  const labelList& fEdges = faceEdges()[facei];
204 
205  // Check if faceNeighbours use same points as this face.
206  // Note: discards normal information - sides of baffle are merged.
207 
208  for (const label edgei : fEdges)
209  {
210  const labelList& eFaces = edgeFaces()[edgei];
211 
212  for (const label neighbour : eFaces)
213  {
214  if (neighbour > facei)
215  {
216  // lower numbered faces already checked
217  const labelledTri& n = (*this)[neighbour];
218 
219  if
220  (
221  ((f[0] == n[0]) || (f[0] == n[1]) || (f[0] == n[2]))
222  && ((f[1] == n[0]) || (f[1] == n[1]) || (f[1] == n[2]))
223  && ((f[2] == n[0]) || (f[2] == n[1]) || (f[2] == n[2]))
224  )
225  {
226  valid.unset(facei);
227 
228  if (verbose)
229  {
231  << "triangles share the same vertices:\n"
232  << " face 1 :" << facei << endl;
233  printTriangle(Warning, " ", f, points());
234 
235  Warning
236  << endl
237  << " face 2 :"
238  << neighbour << endl;
239  printTriangle(Warning, " ", n, points());
240  }
241 
242  break;
243  }
244  }
245  }
246  }
247  }
248  }
249 
250  if (!valid.all())
251  {
252  // Compact
253  label newFacei = 0;
254  for (const label facei : valid)
255  {
256  (*this)[newFacei++] = (*this)[facei];
257  }
258 
259  if (verbose)
260  {
262  << "Removing " << size() - newFacei
263  << " illegal faces." << endl;
264  }
265  (*this).setSize(newFacei);
266 
267  // Topology can change because of renumbering
268  clearOut();
269  }
270 }
271 
272 
273 // Check/fix edges with more than two triangles
274 void Foam::triSurface::checkEdges(const bool verbose)
275 {
276  const labelListList& eFaces = edgeFaces();
277 
278  forAll(eFaces, edgei)
279  {
280  const labelList& myFaces = eFaces[edgei];
281 
282  if (myFaces.empty())
283  {
285  << "Edge " << edgei << " with vertices " << edges()[edgei]
286  << " has no edgeFaces"
287  << exit(FatalError);
288  }
289  else if (myFaces.size() > 2 && verbose)
290  {
292  << "Edge " << edgei << " with vertices " << edges()[edgei]
293  << " has more than 2 faces connected to it : " << myFaces
294  << endl;
295  }
296  }
297 }
298 
299 
300 // Returns patch info. Sets faceMap to the indexing according to patch
301 // numbers. Patch numbers start at 0.
303 Foam::triSurface::calcPatches(labelList& faceMap) const
304 {
305  // Determine the sorted order:
306  // use sortedOrder directly (the intermediate list is discarded anyhow)
307 
308  labelList regions(size());
309  forAll(regions, facei)
310  {
311  regions[facei] = operator[](facei).region();
312  }
313  sortedOrder(regions, faceMap);
314  regions.clear();
315 
316  // Extend regions
317  label maxRegion = patches_.size()-1; // for non-compacted regions
318 
319  if (faceMap.size())
320  {
321  maxRegion = max
322  (
323  maxRegion,
324  operator[](faceMap.last()).region()
325  );
326  }
327 
328  // Get new region list
329  surfacePatchList newPatches(maxRegion + 1);
330 
331  // Fill patch sizes
332  forAll(*this, facei)
333  {
334  label region = operator[](facei).region();
335 
336  newPatches[region].size()++;
337  }
338 
339 
340  // Fill rest of patch info
341 
342  label startFacei = 0;
343  forAll(newPatches, newPatchi)
344  {
345  surfacePatch& newPatch = newPatches[newPatchi];
346 
347  newPatch.index() = newPatchi;
348  newPatch.start() = startFacei;
349 
350  // Take over any information from existing patches
351  if
352  (
353  newPatchi < patches_.size()
354  && !patches_[newPatchi].name().empty()
355  )
356  {
357  newPatch.name() = patches_[newPatchi].name();
358  }
359  else
360  {
361  newPatch.name() = word("patch") + Foam::name(newPatchi);
362  }
363 
364  if
365  (
366  newPatchi < patches_.size()
367  && !patches_[newPatchi].geometricType().empty()
368  )
369  {
370  newPatch.geometricType() = patches_[newPatchi].geometricType();
371  }
372  else
373  {
374  newPatch.geometricType() = geometricSurfacePatch::emptyType;
375  }
376 
377  startFacei += newPatch.size();
378  }
379 
380  return newPatches;
381 }
382 
383 
384 void Foam::triSurface::setDefaultPatches()
385 {
387 
388  // Get names, types and sizes
389  surfacePatchList newPatches(calcPatches(faceMap));
390 
391  // Take over names and types (but not size)
392  patches_.setSize(newPatches.size());
393 
394  forAll(newPatches, patchi)
395  {
396  patches_[patchi].index() = patchi;
397  patches_[patchi].name() = newPatches[patchi].name();
398  patches_[patchi].geometricType() = newPatches[patchi].geometricType();
399  }
400 }
401 
402 
403 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
404 
406 :
408  patches_(),
409  sortedEdgeFacesPtr_(nullptr),
410  edgeOwnerPtr_(nullptr)
411 {}
412 
413 
415 :
416  ParentType(surf, surf.points()),
417  patches_(surf.patches()),
418  sortedEdgeFacesPtr_(nullptr),
419  edgeOwnerPtr_(nullptr)
420 {}
421 
422 
424 :
425  triSurface()
426 {
427  transfer(surf);
428 }
429 
430 
432 (
433  const List<labelledTri>& triangles,
435  const pointField& pts
436 )
437 :
438  ParentType(triangles, pts),
439  patches_(patches),
440  sortedEdgeFacesPtr_(nullptr),
441  edgeOwnerPtr_(nullptr)
442 {}
443 
444 
446 (
447  List<labelledTri>& triangles,
449  pointField& pts,
450  const bool reuse
451 )
452 :
453  ParentType(triangles, pts, reuse),
454  patches_(patches),
455  sortedEdgeFacesPtr_(nullptr),
456  edgeOwnerPtr_(nullptr)
457 {}
458 
459 
461 (
462  const List<labelledTri>& triangles,
463  const pointField& pts
464 )
465 :
466  ParentType(triangles, pts),
467  patches_(),
468  sortedEdgeFacesPtr_(nullptr),
469  edgeOwnerPtr_(nullptr)
470 {
471  setDefaultPatches();
472 }
473 
474 
476 (
477  const triFaceList& triangles,
478  const pointField& pts
479 )
480 :
481  ParentType(convertToTri(triangles, 0), pts),
482  patches_(),
483  sortedEdgeFacesPtr_(nullptr),
484  edgeOwnerPtr_(nullptr)
485 {
486  setDefaultPatches();
487 }
488 
489 
491 (
492  const fileName& name,
493  const scalar scaleFactor
494 )
495 :
496  triSurface(name, name.ext(), scaleFactor)
497 {}
498 
499 
501 (
502  const fileName& name,
503  const word& ext,
504  const scalar scaleFactor
505 )
506 :
508  patches_(),
509  sortedEdgeFacesPtr_(nullptr),
510  edgeOwnerPtr_(nullptr)
511 {
512  read(name, ext);
513  scalePoints(scaleFactor);
514  setDefaultPatches();
515 }
516 
517 
518 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
519 
521 {
522  clearOut();
523 }
524 
525 
526 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
527 
529 {
530  ParentType::clearTopology();
531  deleteDemandDrivenData(sortedEdgeFacesPtr_);
532  deleteDemandDrivenData(edgeOwnerPtr_);
533 }
534 
535 
537 {
538  ParentType::clearPatchMeshAddr();
539 }
540 
541 
543 {
544  ParentType::clearOut();
545  clearTopology();
546  clearPatchMeshAddr();
547 }
548 
549 
551 {
552  if (this == &surf)
553  {
554  return; // Self-swap is a no-op
555  }
556 
557  clearOut();
558  surf.clearOut();
559 
560  FaceListType::swap(static_cast<FaceListType&>(surf));
561  storedPoints().swap(surf.storedPoints());
562  patches_.swap(surf.patches());
563 }
564 
565 
567 {
568  if (!sortedEdgeFacesPtr_)
569  {
570  calcSortedEdgeFaces();
571  }
572 
573  return *sortedEdgeFacesPtr_;
574 }
575 
576 
578 {
579  if (!edgeOwnerPtr_)
580  {
581  calcEdgeOwner();
582  }
583 
584  return *edgeOwnerPtr_;
585 }
586 
587 
589 {
590  // Remove all geometry dependent data
591  deleteDemandDrivenData(sortedEdgeFacesPtr_);
592 
593  // Adapt for new point positions
594  ParentType::movePoints(pts);
595 
596  // Copy new points
597  storedPoints() = pts;
598 }
599 
600 
602 {
603  // Remove all geometry dependent data
604  deleteDemandDrivenData(sortedEdgeFacesPtr_);
605 
606  // Adapt for new point positions
607  ParentType::movePoints(pts);
608 
609  // Move/swap new points
610  storedPoints().swap(pts);
611 }
612 
613 
614 void Foam::triSurface::scalePoints(const scalar scaleFactor)
615 {
616  // Avoid bad scaling
617  if (scaleFactor > SMALL && scaleFactor != 1.0)
618  {
619  // Remove all geometry dependent data
620  clearTopology();
621 
622  // Adapt for new point positions
623  ParentType::movePoints(pointField());
624 
625  storedPoints() *= scaleFactor;
626  }
627 }
628 
629 
630 // Remove non-triangles, double triangles.
631 void Foam::triSurface::cleanup(const bool verbose)
632 {
633  // Merge points (already done for STL, TRI)
634  stitchTriangles(SMALL, verbose);
635 
636  // Merging points might have changed geometric factors
637  clearOut();
638 
639  checkTriangles(verbose);
640 
641  checkEdges(verbose);
642 }
643 
644 
647 {
648  surfacePatchList patches(calcPatches(faceMap));
649 
650  surfZoneList zones(patches.size());
651  forAll(patches, patchi)
652  {
653  zones[patchi] = surfZone(patches[patchi]);
654  }
655 
656  return zones;
657 }
658 
659 
661 {
662  plainFaces.setSize(size());
663 
664  forAll(*this, facei)
665  {
666  plainFaces[facei] = this->operator[](facei);
667  }
668 }
669 
670 
671 // Finds area, starting at facei, delimited by borderEdge. Marks all visited
672 // faces (from face-edge-face walk) with currentZone.
674 (
675  const boolList& borderEdge,
676  const label facei,
677  const label currentZone,
679 ) const
680 {
681  // List of faces whose faceZone has been set.
682  labelList changedFaces(1, facei);
683 
684  while (true)
685  {
686  // Pick up neighbours of changedFaces
687  DynamicList<label> newChangedFaces(2*changedFaces.size());
688 
689  for (const label facei : changedFaces)
690  {
691  const labelList& fEdges = faceEdges()[facei];
692 
693  for (const label edgei : fEdges)
694  {
695  if (!borderEdge[edgei])
696  {
697  const labelList& eFaces = edgeFaces()[edgei];
698 
699  for (const label nbrFacei : eFaces)
700  {
701  if (faceZone[nbrFacei] == -1)
702  {
703  faceZone[nbrFacei] = currentZone;
704  newChangedFaces.append(nbrFacei);
705  }
706  else if (faceZone[nbrFacei] != currentZone)
707  {
709  << "Zones " << faceZone[nbrFacei]
710  << " at face " << nbrFacei
711  << " connects to zone " << currentZone
712  << " at face " << facei
713  << abort(FatalError);
714  }
715  }
716  }
717  }
718  }
719 
720  if (newChangedFaces.empty())
721  {
722  break;
723  }
724 
725  changedFaces.transfer(newChangedFaces);
726  }
727 }
728 
729 
730 // Finds areas delimited by borderEdge (or 'real' edges).
731 // Fills faceZone accordingly
733 (
734  const boolList& borderEdge,
736 ) const
737 {
738  faceZone.setSize(size());
739  faceZone = -1;
740 
741  if (borderEdge.size() != nEdges())
742  {
744  << "borderEdge boolList not same size as number of edges" << endl
745  << "borderEdge:" << borderEdge.size() << endl
746  << "nEdges :" << nEdges()
747  << exit(FatalError);
748  }
749 
750  label zoneI = 0;
751 
752  for (label startFacei = 0;; zoneI++)
753  {
754  // Find first non-coloured face
755  for (; startFacei < size(); startFacei++)
756  {
757  if (faceZone[startFacei] == -1)
758  {
759  break;
760  }
761  }
762 
763  if (startFacei >= size())
764  {
765  break;
766  }
767 
768  faceZone[startFacei] = zoneI;
769 
770  markZone(borderEdge, startFacei, zoneI, faceZone);
771  }
772 
773  return zoneI;
774 }
775 
776 
778 (
779  const boolList& include,
780  labelList& pointMap,
782 ) const
783 {
784  const List<labelledTri>& locFaces = localFaces();
785 
786  label facei = 0;
787  label pointi = 0;
788 
789  faceMap.setSize(locFaces.size());
790 
791  pointMap.setSize(nPoints());
792 
793  bitSet pointHad(nPoints());
794 
795  forAll(include, oldFacei)
796  {
797  if (include[oldFacei])
798  {
799  // Store new faces compact
800  faceMap[facei++] = oldFacei;
801 
802  // Renumber labels for face
803  const triSurface::FaceType& f = locFaces[oldFacei];
804 
805  for (const label verti : f)
806  {
807  if (pointHad.set(verti))
808  {
809  pointMap[pointi++] = verti;
810  }
811  }
812  }
813  }
814 
815  // Trim
816  faceMap.setSize(facei);
817  pointMap.setSize(pointi);
818 }
819 
820 
822 (
823  const boolList& include,
824  labelList& pointMap,
826 ) const
827 {
828  const pointField& locPoints = localPoints();
829  const List<labelledTri>& locFaces = localFaces();
830 
831  // Fill pointMap, faceMap
832  subsetMeshMap(include, pointMap, faceMap);
833 
834 
835  // Create compact coordinate list and forward mapping array
836  pointField newPoints(pointMap.size());
837  labelList oldToNew(locPoints.size());
838  forAll(pointMap, pointi)
839  {
840  newPoints[pointi] = locPoints[pointMap[pointi]];
841  oldToNew[pointMap[pointi]] = pointi;
842  }
843 
844  // Renumber triangle node labels and compact
845  List<labelledTri> newTriangles(faceMap.size());
846 
847  forAll(faceMap, facei)
848  {
849  // Get old vertex labels
850  const labelledTri& tri = locFaces[faceMap[facei]];
851 
852  newTriangles[facei][0] = oldToNew[tri[0]];
853  newTriangles[facei][1] = oldToNew[tri[1]];
854  newTriangles[facei][2] = oldToNew[tri[2]];
855  newTriangles[facei].region() = tri.region();
856  }
857 
858  // Construct sub-surface
859  return triSurface(newTriangles, patches(), newPoints, true);
860 }
861 
862 
864 {
865  labelList pointMap, faceMap;
866  return subsetMesh(include, pointMap, faceMap);
867 }
868 
869 
871 {
872  clearOut(); // Topology changes
873 
874  this->storedFaces().swap(faceLst);
875 }
876 
877 
879 {
880  clearOut();
881 
882  FaceListType::transfer(surf.storedFaces());
883  storedPoints().transfer(surf.storedPoints());
884  patches_.transfer(surf.patches());
885 
886  surf.clearOut();
887 }
888 
889 
891 {
892  // Transcribe zone -> patch info
893  auto patches = ListOps::create<geometricSurfacePatch>
894  (
895  surf.surfZones(),
897  );
898 
899  // Fairly ugly, but the only way to get at the data safely
900  List<labelledTri> fcs;
901  pointField pts;
902 
903  surf.swapFaces(fcs);
904  surf.swapPoints(pts);
905 
906  clearOut();
907  surf.clear();
908 
909  triSurface s(fcs, patches, pts, true);
910  swap(s);
911 }
912 
913 
914 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
915 
917 {
918  clearOut();
919 
920  FaceListType::operator=(static_cast<const FaceListType&>(surf));
921  storedPoints() = surf.points();
922  patches_ = surf.patches();
923 }
924 
925 
927 {
928  transfer(surf);
929 }
930 
931 
933 {
934  transfer(surf);
935 }
936 
937 
938 // ************************************************************************* //
Foam::expressions::patchExpr::debug
int debug
Static debugging option.
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:74
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
Foam::TimeState::timeOutputValue
scalar timeOutputValue() const
Return current time value.
Definition: TimeStateI.H:30
Foam::MeshedSurface::surfZones
const surfZoneList & surfZones() const
Const access to the surface zones.
Definition: MeshedSurface.H:370
Foam::triSurface::subsetMeshMap
void subsetMeshMap(const boolList &include, labelList &pointMap, labelList &faceMap) const
'Create' sub mesh
Definition: triSurface.C:778
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:73
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:94
Foam::PrimitivePatch::points
const Field< PointType > & points() const
Return reference to global points.
Definition: PrimitivePatch.H:300
Foam::triSurface::clearOut
void clearOut()
Definition: triSurface.C:542
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:62
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.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))
Definition: gmvOutputSpray.H:25
Foam::bitSet::all
bool all() const
True if all bits in this bitset are set or if the set is empty.
Definition: bitSetI.H:454
Foam::bitSet
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:64
Foam::triSurface::markZone
void markZone(const boolList &borderEdge, const label facei, const label currentZone, labelList &faceZone) const
Fill faceZone with currentZone for every face reachable.
Definition: triSurface.C:674
Foam::DynamicList< label >
Foam::triSurface::triFaceFaces
void triFaceFaces(List< face > &plainFaceList) const
Create a list of faces from the triFaces.
Definition: triSurface.C:660
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::bitSet::unset
bitSet & unset(const bitSet &other)
Definition: bitSetI.H:601
Foam::fileName::name
static std::string name(const std::string &str)
Return basename (part beyond last /), including its extension.
Definition: fileNameI.H:209
Foam::Warning
messageStream Warning
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::surfacePatchList
List< surfacePatch > surfacePatchList
Definition: surfacePatchList.H:46
Foam::bitSet::set
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:563
Foam::isFile
bool isFile(const fileName &name, const bool checkGzip=true, const bool followLink=true)
Does the name exist as a FILE in the file system?
Definition: MSwindows.C:658
Foam::List::append
void append(const T &val)
Append an element at the end of the list.
Definition: ListI.H:182
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:337
Foam::word::ext
word ext() const
Return file name extension (part after last .)
Definition: word.C:126
Foam::MeshedSurface::clear
virtual void clear()
Clear all storage.
Definition: MeshedSurface.C:514
Foam::triSurface::clearPatchMeshAddr
void clearPatchMeshAddr()
Definition: triSurface.C:536
Foam::Pout
prefixOSstream Pout
An Ostream wrapper for parallel output to std::cout.
triSurface.H
Foam::triSurface::patches
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
Foam::triSurface::~triSurface
virtual ~triSurface()
Destructor.
Definition: triSurface.C:520
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::triSurface::swapFaces
void swapFaces(List< labelledTri > &faceLst)
Swap the list of faces being addressed.
Definition: triSurface.C:870
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:42
n
label n
Definition: TABSMDCalcMethod2.H:31
Foam::triSurface::scalePoints
virtual void scalePoints(const scalar scaleFactor)
Scale points. A non-positive factor is ignored.
Definition: triSurface.C:614
Foam::triSurface::checkTriangles
void checkTriangles(const bool verbose)
Check/remove duplicate/degenerate triangles.
Definition: triSurface.C:155
Foam::triSurface::subsetMesh
triSurface subsetMesh(const boolList &include, labelList &pointMap, labelList &faceMap) const
Return new surface.
Definition: triSurface.C:822
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:62
Foam::triSurface::movePoints
virtual void movePoints(const pointField &pts)
Move points.
Definition: triSurface.C:588
Foam::Field< vector >
Foam::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:70
Foam::triSurface::swapPoints
virtual void swapPoints(pointField &pts)
Swap points. Similar to movePoints, but returns the old points.
Definition: triSurface.C:601
Foam::faceZone
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:65
Foam::name
word name(const complex &c)
Return string representation of complex.
Definition: complex.C:76
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:91
Foam::TimePaths::times
instantList times() const
Search the case for valid time directories.
Definition: TimePaths.C:149
Foam::triSurface::transfer
void transfer(triSurface &surf)
Alter contents by transferring (triangles, points) components.
Definition: triSurface.C:878
Foam::IOstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.C:39
Foam::triSurface::triSurface
triSurface()
Construct null.
Definition: triSurface.C:405
Foam::triSurface::storedPoints
pointField & storedPoints()
Non-const access to global points.
Definition: triSurface.H:198
Foam::triSurface::sortedZones
List< surfZone > sortedZones(labelList &faceMap) const
Sort faces according to zoneIds.
Definition: triSurface.C:646
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::MeshedSurface::swapPoints
void swapPoints(pointField &points)
Swap the stored points.
Definition: MeshedSurface.C:1219
Foam::FatalError
error FatalError
Foam::triSurface::storedFaces
List< Face > & storedFaces()
Non-const access to the faces.
Definition: triSurface.H:204
Foam::triSurface::clearTopology
void clearTopology()
Definition: triSurface.C:528
Foam::triFaceList
List< triFace > triFaceList
list of triFaces
Definition: triFaceList.H:44
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:137
Foam::triSurface::edgeOwner
const labelList & edgeOwner() const
If 2 face neighbours: label of face where ordering of edge.
Definition: triSurface.C:577
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::triSurface::cleanup
void cleanup(const bool verbose)
Remove non-valid triangles.
Definition: triSurface.C:631
Time.H
Foam::triSurface::markZones
label markZones(const boolList &borderEdge, labelList &faceZone) const
(size and) fills faceZone with zone of face. Zone is area
Definition: triSurface.C:733
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:355
Foam::Time::caseName
const fileName & caseName() const
Return case name.
Definition: TimePathsI.H:54
Foam::geometricSurfacePatch::fromIdentifier
Helper to convert identifier types as an operation.
Definition: geometricSurfacePatch.H:81
Foam::triSurface::swap
void swap(triSurface &surf)
Definition: triSurface.C:550
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:303
f
labelList f(nPoints)
Foam::List< instant >
Foam::surfZone
A surface zone on a MeshedSurface.
Definition: surfZone.H:65
points
const pointField & points
Definition: gmvOutputHeader.H:1
Foam::geometricSurfacePatch::emptyType
static const word emptyType
The name for an 'empty' type.
Definition: geometricSurfacePatch.H:78
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::labelledTri
Triangle with additional region number.
Definition: labelledTri.H:58
surfZoneList.H
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:74
ListOps.H
Various functions to operate on Lists.
Foam::triSurface::sortedEdgeFaces
const labelListList & sortedEdgeFaces() const
Return edge-face addressing sorted (for edges with more than.
Definition: triSurface.C:566
Foam::sortedOrder
labelList sortedOrder(const UList< T > &input)
Return the (stable) sort order for the list.
Foam::triSurface::operator=
void operator=(const triSurface &surf)
Copy assignment.
Definition: triSurface.C:916
Foam::MeshedSurface::swapFaces
void swapFaces(List< Face > &faces)
Swap the stored faces.
Definition: MeshedSurface.C:1210
Foam::List::setSize
void setSize(const label newSize)
Alias for resize(const label)
Definition: ListI.H:146
Foam::TimePaths::constant
const word & constant() const
Return constant name.
Definition: TimePathsI.H:88
Foam::MeshedSurface
A surface geometry mesh with zone information, not to be confused with the similarly named surfaceMes...
Definition: triSurfaceTools.H:80
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
triFace
face triFace(3)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:294
Foam::triSurface::triSurfInstance
static fileName triSurfInstance(const Time &)
Name of triSurface directory to use.
Definition: triSurface.C:46
MeshedSurface.H
Foam::PrimitivePatch< labelledTri, ::Foam::List, pointField, point >
Foam::triSurface::checkEdges
void checkEdges(const bool verbose)
Check triply (or more) connected edges.
Definition: triSurface.C:274