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 | 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 "triSurface.H"
27 #include "demandDrivenData.H"
28 #include "IFstream.H"
29 #include "OFstream.H"
30 #include "Time.H"
31 #include "boundBox.H"
32 #include "SortableList.H"
33 #include "PackedBoolList.H"
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 defineTypeNameAndDebug(triSurface, 0);
40 }
41 
42 
44 {
45  fileName foamName(d.caseName() + ".ftr");
46 
47  // Search back through the time directories list to find the time
48  // closest to and lower than current time
49 
50  instantList ts = d.times();
51  label i;
52 
53  for (i=ts.size()-1; i>=0; i--)
54  {
55  if (ts[i].value() <= d.timeOutputValue())
56  {
57  break;
58  }
59  }
60 
61  // Noting that the current directory has already been searched
62  // for mesh data, start searching from the previously stored time directory
63 
64  if (i>=0)
65  {
66  for (label j=i; j>=0; j--)
67  {
68  if (isFile(d.path()/ts[j].name()/typeName/foamName))
69  {
70  if (debug)
71  {
72  Pout<< " triSurface::triSurfInstance(const Time& d)"
73  << "reading " << foamName
74  << " from " << ts[j].name()/typeName
75  << endl;
76  }
77 
78  return ts[j].name();
79  }
80  }
81  }
82 
83  if (debug)
84  {
85  Pout<< " triSurface::triSurfInstance(const Time& d)"
86  << "reading " << foamName
87  << " from constant/" << endl;
88  }
89  return d.constant();
90 }
91 
92 
94 (
95  const faceList& faces,
96  const label defaultRegion
97 )
98 {
99  List<labelledTri> triFaces(faces.size());
100 
101  forAll(triFaces, faceI)
102  {
103  const face& f = faces[faceI];
104 
105  if (f.size() != 3)
106  {
108  << "Face at position " << faceI
109  << " does not have three vertices:" << f
110  << abort(FatalError);
111  }
112 
113  labelledTri& tri = triFaces[faceI];
114 
115  tri[0] = f[0];
116  tri[1] = f[1];
117  tri[2] = f[2];
118  tri.region() = defaultRegion;
119  }
120 
121  return triFaces;
122 }
123 
124 
126 (
127  const triFaceList& faces,
128  const label defaultRegion
129 )
130 {
131  List<labelledTri> triFaces(faces.size());
132 
133  forAll(triFaces, faceI)
134  {
135  const triFace& f = faces[faceI];
136 
137  labelledTri& tri = triFaces[faceI];
138 
139  tri[0] = f[0];
140  tri[1] = f[1];
141  tri[2] = f[2];
142  tri.region() = defaultRegion;
143  }
144 
145  return triFaces;
146 }
147 
148 
149 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
150 
152 (
153  Ostream& os,
154  const string& pre,
155  const labelledTri& f,
156  const pointField& points
157 )
158 {
159  os
160  << pre.c_str() << "vertex numbers:"
161  << f[0] << ' ' << f[1] << ' ' << f[2] << endl
162  << pre.c_str() << "vertex coords :"
163  << points[f[0]] << ' ' << points[f[1]] << ' ' << points[f[2]]
164  << pre.c_str() << "region :" << f.region() << endl
165  << endl;
166 }
167 
168 
170 {
171  string line;
172  do
173  {
174  is.getLine(line);
175  }
176  while ((line.empty() || line[0] == '#') && is.good());
177 
178  return line;
179 }
180 
181 
182 // Remove non-triangles, double triangles.
183 void Foam::triSurface::checkTriangles(const bool verbose)
184 {
185  // Simple check on indices ok.
186  const label maxPointI = points().size() - 1;
187 
188  forAll(*this, faceI)
189  {
190  const triSurface::FaceType& f = (*this)[faceI];
191 
192  forAll(f, fp)
193  {
194  if (f[fp] < 0 || f[fp] > maxPointI)
195  {
197  << "triangle " << f
198  << " uses point indices outside point range 0.."
199  << maxPointI
200  << exit(FatalError);
201  }
202  }
203  }
204 
205  // Two phase process
206  // 1. mark invalid faces
207  // 2. pack
208  // Done to keep numbering constant in phase 1
209 
210  // List of valid triangles
211  boolList valid(size(), true);
212  bool hasInvalid = false;
213 
214  forAll(*this, faceI)
215  {
216  const labelledTri& f = (*this)[faceI];
217 
218  if ((f[0] == f[1]) || (f[0] == f[2]) || (f[1] == f[2]))
219  {
220  // 'degenerate' triangle check
221  valid[faceI] = false;
222  hasInvalid = true;
223 
224  if (verbose)
225  {
227  << "triangle " << faceI
228  << " does not have three unique vertices:\n";
229  printTriangle(Warning, " ", f, points());
230  }
231  }
232  else
233  {
234  // duplicate triangle check
235  const labelList& fEdges = faceEdges()[faceI];
236 
237  // Check if faceNeighbours use same points as this face.
238  // Note: discards normal information - sides of baffle are merged.
239 
240  forAll(fEdges, fp)
241  {
242  const labelList& eFaces = edgeFaces()[fEdges[fp]];
243 
244  forAll(eFaces, i)
245  {
246  label neighbour = eFaces[i];
247 
248  if (neighbour > faceI)
249  {
250  // lower numbered faces already checked
251  const labelledTri& n = (*this)[neighbour];
252 
253  if
254  (
255  ((f[0] == n[0]) || (f[0] == n[1]) || (f[0] == n[2]))
256  && ((f[1] == n[0]) || (f[1] == n[1]) || (f[1] == n[2]))
257  && ((f[2] == n[0]) || (f[2] == n[1]) || (f[2] == n[2]))
258  )
259  {
260  valid[faceI] = false;
261  hasInvalid = true;
262 
263  if (verbose)
264  {
266  << "triangles share the same vertices:\n"
267  << " face 1 :" << faceI << endl;
268  printTriangle(Warning, " ", f, points());
269 
270  Warning
271  << endl
272  << " face 2 :"
273  << neighbour << endl;
274  printTriangle(Warning, " ", n, points());
275  }
276 
277  break;
278  }
279  }
280  }
281  }
282  }
283  }
284 
285  if (hasInvalid)
286  {
287  // Pack
288  label newFaceI = 0;
289  forAll(*this, faceI)
290  {
291  if (valid[faceI])
292  {
293  const labelledTri& f = (*this)[faceI];
294  (*this)[newFaceI++] = f;
295  }
296  }
297 
298  if (verbose)
299  {
301  << "Removing " << size() - newFaceI
302  << " illegal faces." << endl;
303  }
304  (*this).setSize(newFaceI);
305 
306  // Topology can change because of renumbering
307  clearOut();
308  }
309 }
310 
311 
312 // Check/fix edges with more than two triangles
313 void Foam::triSurface::checkEdges(const bool verbose)
314 {
315  const labelListList& eFaces = edgeFaces();
316 
317  forAll(eFaces, edgeI)
318  {
319  const labelList& myFaces = eFaces[edgeI];
320 
321  if (myFaces.empty())
322  {
324  << "Edge " << edgeI << " with vertices " << edges()[edgeI]
325  << " has no edgeFaces"
326  << exit(FatalError);
327  }
328  else if (myFaces.size() > 2 && verbose)
329  {
331  << "Edge " << edgeI << " with vertices " << edges()[edgeI]
332  << " has more than 2 faces connected to it : " << myFaces
333  << endl;
334  }
335  }
336 }
337 
338 
339 // Read triangles, points from Istream
340 bool Foam::triSurface::read(Istream& is)
341 {
342  is >> patches_ >> storedPoints() >> storedFaces();
343 
344  return true;
345 }
346 
347 
348 // Read from file in given format
350 (
351  const fileName& name,
352  const word& ext,
353  const bool check
354 )
355 {
356  if (check && !exists(name))
357  {
359  << "Cannnot read " << name << exit(FatalError);
360  }
361 
362  if (ext == "gz")
363  {
364  fileName unzipName = name.lessExt();
365 
366  // Do not check for existence. Let IFstream do the unzipping.
367  return read(unzipName, unzipName.ext(), false);
368  }
369  else if (ext == "ftr")
370  {
371  return read(IFstream(name)());
372  }
373  else if (ext == "stl")
374  {
375  return readSTL(name);
376  }
377  else if (ext == "stlb")
378  {
379  return readSTLBINARY(name);
380  }
381  else if (ext == "gts")
382  {
383  return readGTS(name);
384  }
385  else if (ext == "obj")
386  {
387  return readOBJ(name);
388  }
389  else if (ext == "off")
390  {
391  return readOFF(name);
392  }
393  else if (ext == "tri")
394  {
395  return readTRI(name);
396  }
397  else if (ext == "ac")
398  {
399  return readAC(name);
400  }
401  else if (ext == "nas")
402  {
403  return readNAS(name);
404  }
405  else if (ext == "vtk")
406  {
407  return readVTK(name);
408  }
409  else
410  {
412  << "unknown file extension " << ext
413  << ". Supported extensions are '.ftr', '.stl', '.stlb', '.gts'"
414  << ", '.obj', '.ac', '.off', '.nas', '.tri' and '.vtk'"
415  << exit(FatalError);
416 
417  return false;
418  }
419 }
420 
421 
422 // Write to file in given format
424 (
425  const fileName& name,
426  const word& ext,
427  const bool sort
428 ) const
429 {
430  if (ext == "ftr")
431  {
432  return write(OFstream(name)());
433  }
434  else if (ext == "stl")
435  {
436  return writeSTLASCII(sort, OFstream(name)());
437  }
438  else if (ext == "stlb")
439  {
440  ofstream outFile(name.c_str(), std::ios::binary);
441 
442  writeSTLBINARY(outFile);
443  }
444  else if (ext == "gts")
445  {
446  return writeGTS(sort, OFstream(name)());
447  }
448  else if (ext == "obj")
449  {
450  writeOBJ(sort, OFstream(name)());
451  }
452  else if (ext == "off")
453  {
454  writeOFF(sort, OFstream(name)());
455  }
456  else if (ext == "vtk")
457  {
458  writeVTK(sort, OFstream(name)());
459  }
460  else if (ext == "tri")
461  {
462  writeTRI(sort, OFstream(name)());
463  }
464  else if (ext == "dx")
465  {
466  writeDX(sort, OFstream(name)());
467  }
468  else if (ext == "ac")
469  {
470  writeAC(OFstream(name)());
471  }
472  else if (ext == "smesh")
473  {
474  writeSMESH(sort, OFstream(name)());
475  }
476  else
477  {
479  << "unknown file extension " << ext
480  << " for file " << name
481  << ". Supported extensions are '.ftr', '.stl', '.stlb', "
482  << "'.gts', '.obj', '.vtk'"
483  << ", '.off', '.dx', '.smesh', '.ac' and '.tri'"
484  << exit(FatalError);
485  }
486 }
487 
488 
489 // Returns patch info. Sets faceMap to the indexing according to patch
490 // numbers. Patch numbers start at 0.
492 {
493  // Sort according to region numbers of labelledTri
494  SortableList<label> sortedRegion(size());
495 
496  forAll(sortedRegion, faceI)
497  {
498  sortedRegion[faceI] = operator[](faceI).region();
499  }
500  sortedRegion.sort();
501 
502  faceMap = sortedRegion.indices();
503 
504  // Extend regions
505  label maxRegion = patches_.size()-1; // for non-compacted regions
506 
507  if (faceMap.size())
508  {
509  maxRegion = max
510  (
511  maxRegion,
512  operator[](faceMap.last()).region()
513  );
514  }
515 
516  // Get new region list
517  surfacePatchList newPatches(maxRegion + 1);
518 
519  // Fill patch sizes
520  forAll(*this, faceI)
521  {
522  label region = operator[](faceI).region();
523 
524  newPatches[region].size()++;
525  }
526 
527 
528  // Fill rest of patch info
529 
530  label startFaceI = 0;
531  forAll(newPatches, newPatchI)
532  {
533  surfacePatch& newPatch = newPatches[newPatchI];
534 
535  newPatch.index() = newPatchI;
536 
537  label oldPatchI = newPatchI;
538 
539  // start of patch
540  newPatch.start() = startFaceI;
541 
542 
543  // Take over any information from existing patches
544  if ((oldPatchI < patches_.size()) && (patches_[oldPatchI].name() != ""))
545  {
546  newPatch.name() = patches_[oldPatchI].name();
547  }
548  else
549  {
550  newPatch.name() = word("patch") + name(newPatchI);
551  }
552 
553  if
554  (
555  (oldPatchI < patches_.size())
556  && (patches_[oldPatchI].geometricType() != "")
557  )
558  {
559  newPatch.geometricType() = patches_[oldPatchI].geometricType();
560  }
561  else
562  {
563  newPatch.geometricType() = "empty";
564  }
565 
566  startFaceI += newPatch.size();
567  }
568 
569  return newPatches;
570 }
571 
572 
574 {
576 
577  // Get names, types and sizes
578  surfacePatchList newPatches(calcPatches(faceMap));
579 
580  // Take over names and types (but not size)
581  patches_.setSize(newPatches.size());
582 
583  forAll(newPatches, patchI)
584  {
585  patches_[patchI].index() = patchI;
586  patches_[patchI].name() = newPatches[patchI].name();
587  patches_[patchI].geometricType() = newPatches[patchI].geometricType();
588  }
589 }
590 
591 
592 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
593 
595 :
596  ParentType(List<Face>(), pointField()),
597  patches_(0),
598  sortedEdgeFacesPtr_(NULL),
599  edgeOwnerPtr_(NULL)
600 {}
601 
602 
603 
605 (
606  const List<labelledTri>& triangles,
608  const pointField& points
609 )
610 :
611  ParentType(triangles, points),
612  patches_(patches),
613  sortedEdgeFacesPtr_(NULL),
614  edgeOwnerPtr_(NULL)
615 {}
616 
617 
619 (
620  List<labelledTri>& triangles,
623  const bool reUse
624 )
625 :
626  ParentType(triangles, points, reUse),
627  patches_(patches),
628  sortedEdgeFacesPtr_(NULL),
629  edgeOwnerPtr_(NULL)
630 {}
631 
632 
634 (
635  const Xfer<List<labelledTri> >& triangles,
637  const Xfer<List<point> >& points
638 )
639 :
640  ParentType(triangles, points),
641  patches_(patches),
642  sortedEdgeFacesPtr_(NULL),
643  edgeOwnerPtr_(NULL)
644 {}
645 
646 
648 (
649  const List<labelledTri>& triangles,
650  const pointField& points
651 )
652 :
653  ParentType(triangles, points),
654  patches_(),
655  sortedEdgeFacesPtr_(NULL),
656  edgeOwnerPtr_(NULL)
657 {
659 }
660 
661 
663 (
664  const triFaceList& triangles,
665  const pointField& points
666 )
667 :
668  ParentType(convertToTri(triangles, 0), points),
669  patches_(),
670  sortedEdgeFacesPtr_(NULL),
671  edgeOwnerPtr_(NULL)
672 {
674 }
675 
676 
677 Foam::triSurface::triSurface(const fileName& name)
678 :
679  ParentType(List<Face>(), pointField()),
680  patches_(),
681  sortedEdgeFacesPtr_(NULL),
682  edgeOwnerPtr_(NULL)
683 {
684  word ext = name.ext();
685 
686  read(name, ext);
687 
689 }
690 
691 
692 Foam::triSurface::triSurface(Istream& is)
693 :
694  ParentType(List<Face>(), pointField()),
695  patches_(),
696  sortedEdgeFacesPtr_(NULL),
697  edgeOwnerPtr_(NULL)
698 {
699  read(is);
700 
702 }
703 
704 
705 Foam::triSurface::triSurface(const Time& d)
706 :
707  ParentType(List<Face>(), pointField()),
708  patches_(),
709  sortedEdgeFacesPtr_(NULL),
710  edgeOwnerPtr_(NULL)
711 {
712  fileName foamFile(d.caseName() + ".ftr");
713 
714  fileName foamPath(d.path()/triSurfInstance(d)/typeName/foamFile);
715 
716  IFstream foamStream(foamPath);
717 
718  read(foamStream);
719 
721 }
722 
723 
724 Foam::triSurface::triSurface(const triSurface& ts)
725 :
726  ParentType(ts, ts.points()),
727  patches_(ts.patches()),
728  sortedEdgeFacesPtr_(NULL),
729  edgeOwnerPtr_(NULL)
730 {}
731 
732 
733 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
734 
736 {
737  clearOut();
738 }
739 
740 
741 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
742 
744 {
745  ParentType::clearTopology();
746  deleteDemandDrivenData(sortedEdgeFacesPtr_);
747  deleteDemandDrivenData(edgeOwnerPtr_);
748 }
749 
750 
752 {
753  ParentType::clearPatchMeshAddr();
754 }
755 
756 
758 {
759  ParentType::clearOut();
760 
761  clearTopology();
762  clearPatchMeshAddr();
763 }
764 
765 
767 {
768  if (!sortedEdgeFacesPtr_)
769  {
770  calcSortedEdgeFaces();
771  }
772 
773  return *sortedEdgeFacesPtr_;
774 }
775 
776 
778 {
779  if (!edgeOwnerPtr_)
780  {
781  calcEdgeOwner();
782  }
783 
784  return *edgeOwnerPtr_;
785 }
786 
787 
788 void Foam::triSurface::movePoints(const pointField& newPoints)
789 {
790  // Remove all geometry dependent data
791  deleteDemandDrivenData(sortedEdgeFacesPtr_);
792 
793  // Adapt for new point position
794  ParentType::movePoints(newPoints);
795 
796  // Copy new points
797  storedPoints() = newPoints;
798 }
799 
800 
801 void Foam::triSurface::scalePoints(const scalar scaleFactor)
802 {
803  // avoid bad scaling
804  if (scaleFactor > 0 && scaleFactor != 1.0)
805  {
806  // Remove all geometry dependent data
807  clearTopology();
808 
809  // Adapt for new point position
810  ParentType::movePoints(pointField());
811 
812  storedPoints() *= scaleFactor;
813  }
814 }
815 
816 
817 // Remove non-triangles, double triangles.
818 void Foam::triSurface::cleanup(const bool verbose)
819 {
820  // Merge points (already done for STL, TRI)
821  stitchTriangles(SMALL, verbose);
822 
823  // Merging points might have changed geometric factors
824  clearOut();
825 
826  checkTriangles(verbose);
827 
828  checkEdges(verbose);
829 }
830 
831 
832 // Finds area, starting at faceI, delimited by borderEdge. Marks all visited
833 // faces (from face-edge-face walk) with currentZone.
835 (
836  const boolList& borderEdge,
837  const label faceI,
838  const label currentZone,
839  labelList& faceZone
840 ) const
841 {
842  // List of faces whose faceZone has been set.
843  labelList changedFaces(1, faceI);
844 
845  while (true)
846  {
847  // Pick up neighbours of changedFaces
848  DynamicList<label> newChangedFaces(2*changedFaces.size());
849 
850  forAll(changedFaces, i)
851  {
852  label faceI = changedFaces[i];
853 
854  const labelList& fEdges = faceEdges()[faceI];
855 
856  forAll(fEdges, i)
857  {
858  label edgeI = fEdges[i];
859 
860  if (!borderEdge[edgeI])
861  {
862  const labelList& eFaces = edgeFaces()[edgeI];
863 
864  forAll(eFaces, j)
865  {
866  label nbrFaceI = eFaces[j];
867 
868  if (faceZone[nbrFaceI] == -1)
869  {
870  faceZone[nbrFaceI] = currentZone;
871  newChangedFaces.append(nbrFaceI);
872  }
873  else if (faceZone[nbrFaceI] != currentZone)
874  {
876  << "Zones " << faceZone[nbrFaceI]
877  << " at face " << nbrFaceI
878  << " connects to zone " << currentZone
879  << " at face " << faceI
880  << abort(FatalError);
881  }
882  }
883  }
884  }
885  }
886 
887  if (newChangedFaces.empty())
888  {
889  break;
890  }
891 
892  changedFaces.transfer(newChangedFaces);
893  }
894 }
895 
896 
897 // Finds areas delimited by borderEdge (or 'real' edges).
898 // Fills faceZone accordingly
900 (
901  const boolList& borderEdge,
902  labelList& faceZone
903 ) const
904 {
905  faceZone.setSize(size());
906  faceZone = -1;
907 
908  if (borderEdge.size() != nEdges())
909  {
911  << "borderEdge boolList not same size as number of edges" << endl
912  << "borderEdge:" << borderEdge.size() << endl
913  << "nEdges :" << nEdges()
914  << exit(FatalError);
915  }
916 
917  label zoneI = 0;
918 
919  for (label startFaceI = 0;; zoneI++)
920  {
921  // Find first non-coloured face
922  for (; startFaceI < size(); startFaceI++)
923  {
924  if (faceZone[startFaceI] == -1)
925  {
926  break;
927  }
928  }
929 
930  if (startFaceI >= size())
931  {
932  break;
933  }
934 
935  faceZone[startFaceI] = zoneI;
936 
937  markZone(borderEdge, startFaceI, zoneI, faceZone);
938  }
939 
940  return zoneI;
941 }
942 
943 
945 (
946  const boolList& include,
947  labelList& pointMap,
949 ) const
950 {
951  const List<labelledTri>& locFaces = localFaces();
952 
953  label faceI = 0;
954  label pointI = 0;
955 
956  faceMap.setSize(locFaces.size());
957 
958  pointMap.setSize(nPoints());
959 
960  boolList pointHad(nPoints(), false);
961 
962  forAll(include, oldFaceI)
963  {
964  if (include[oldFaceI])
965  {
966  // Store new faces compact
967  faceMap[faceI++] = oldFaceI;
968 
969  // Renumber labels for face
970  const triSurface::FaceType& f = locFaces[oldFaceI];
971 
972  forAll(f, fp)
973  {
974  label labI = f[fp];
975  if (!pointHad[labI])
976  {
977  pointHad[labI] = true;
978  pointMap[pointI++] = labI;
979  }
980  }
981  }
982  }
983 
984  // Trim
985  faceMap.setSize(faceI);
986  pointMap.setSize(pointI);
987 }
988 
989 
991 (
992  const boolList& include,
993  labelList& pointMap,
995 ) const
996 {
997  const pointField& locPoints = localPoints();
998  const List<labelledTri>& locFaces = localFaces();
999 
1000  // Fill pointMap, faceMap
1001  subsetMeshMap(include, pointMap, faceMap);
1002 
1003 
1004  // Create compact coordinate list and forward mapping array
1005  pointField newPoints(pointMap.size());
1006  labelList oldToNew(locPoints.size());
1007  forAll(pointMap, pointi)
1008  {
1009  newPoints[pointi] = locPoints[pointMap[pointi]];
1010  oldToNew[pointMap[pointi]] = pointi;
1011  }
1012 
1013  // Renumber triangle node labels and compact
1014  List<labelledTri> newTriangles(faceMap.size());
1015 
1016  forAll(faceMap, facei)
1017  {
1018  // Get old vertex labels
1019  const labelledTri& tri = locFaces[faceMap[facei]];
1020 
1021  newTriangles[facei][0] = oldToNew[tri[0]];
1022  newTriangles[facei][1] = oldToNew[tri[1]];
1023  newTriangles[facei][2] = oldToNew[tri[2]];
1024  newTriangles[facei].region() = tri.region();
1025  }
1026 
1027  // Construct subsurface
1028  return triSurface(newTriangles, patches(), newPoints, true);
1029 }
1030 
1031 
1033 (
1034  const fileName& name,
1035  const bool sortByRegion
1036 ) const
1037 {
1038  write(name, name.ext(), sortByRegion);
1039 }
1040 
1041 
1042 void Foam::triSurface::write(Ostream& os) const
1043 {
1044  os << patches() << endl;
1045 
1046  //Note: Write with global point numbering
1047  os << points() << nl
1048  << static_cast<const List<labelledTri>&>(*this) << endl;
1049 
1050  // Check state of Ostream
1051  os.check("triSurface::write(Ostream&)");
1052 }
1053 
1054 
1055 void Foam::triSurface::write(const Time& d) const
1056 {
1057  fileName foamFile(d.caseName() + ".ftr");
1058 
1059  fileName foamPath(d.path()/triSurfInstance(d)/typeName/foamFile);
1060 
1061  OFstream foamStream(foamPath);
1062 
1063  write(foamStream);
1064 }
1065 
1066 
1067 void Foam::triSurface::writeStats(Ostream& os) const
1068 {
1069  // Unfortunately nPoints constructs meshPoints() so do compact version
1070  // ourselves.
1071  PackedBoolList pointIsUsed(points().size());
1072 
1073  label nPoints = 0;
1074  boundBox bb = boundBox::invertedBox;
1075 
1076  forAll(*this, faceI)
1077  {
1078  const triSurface::FaceType& f = operator[](faceI);
1079 
1080  forAll(f, fp)
1081  {
1082  label pointI = f[fp];
1083  if (pointIsUsed.set(pointI, 1))
1084  {
1085  bb.min() = ::Foam::min(bb.min(), points()[pointI]);
1086  bb.max() = ::Foam::max(bb.max(), points()[pointI]);
1087  nPoints++;
1088  }
1089  }
1090  }
1091 
1092  os << "Triangles : " << size() << endl
1093  << "Vertices : " << nPoints << endl
1094  << "Bounding Box : " << bb << endl;
1095 }
1096 
1097 
1098 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1099 
1100 void Foam::triSurface::operator=(const triSurface& ts)
1101 {
1103  clearOut();
1104  storedPoints() = ts.points();
1105  patches_ = ts.patches();
1106 }
1107 
1108 
1109 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1110 
1111 Foam::Ostream& Foam::operator<<(Ostream& os, const triSurface& sm)
1112 {
1113  sm.write(os);
1114  return os;
1115 }
1116 
1117 
1118 // ************************************************************************* //
Foam::pointField
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
Foam::triSurface::convertToTri
static List< labelledTri > convertToTri(const faceList &, const label defaultRegion=0)
Convert faces to labelledTri. All get same region.
Definition: triSurface.C:95
Foam::triSurface::subsetMeshMap
void subsetMeshMap(const boolList &include, labelList &pointMap, labelList &faceMap) const
'Create' sub mesh, including only faces for which
Definition: triSurface.C:955
Foam::faceMap
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
Definition: blockMeshMergeFast.C:90
Foam::triSurface::clearOut
void clearOut()
Clear all data.
Definition: triSurface.C:757
Foam::triSurface::writeStats
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1089
Foam::fileName
A class for handling file names.
Definition: fileName.H:69
Foam::exists
bool exists(const fileName &, const bool checkGzip=true)
Does the name exist (as DIRECTORY or FILE) in the file system?
Definition: POSIX.C:608
Foam::labelList
List< label > labelList
A List of labels.
Definition: labelList.H:56
Foam::triSurface::checkTriangles
void checkTriangles(const bool verbose=false)
Check/remove duplicate/degenerate triangles.
Definition: triSurface.C:186
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
writeOBJ
void writeOBJ(Ostream &os, label &vertI, const tetPoints &tet)
Definition: Test-tetTetOverlap.C:38
Foam::triSurface::read
bool read(Istream &)
Read in Foam format.
Definition: triSurface.C:349
demandDrivenData.H
Template functions to aid in the implementation of demand driven data.
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::boundBox::invertedBox
static const boundBox invertedBox
A very large inverted boundBox: min/max == +/- VGREAT.
Definition: boundBox.H:79
Foam::Warning
messageStream Warning
Foam::geometricSurfacePatchList
List< geometricSurfacePatch > geometricSurfacePatchList
Definition: geometricSurfacePatchList.H:44
triSurface.H
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::string
A class for handling character strings derived from std::string.
Definition: string.H:74
Foam::triSurface::clearPatchMeshAddr
void clearPatchMeshAddr()
Clear patch addressing.
Definition: triSurface.C:751
Foam::List::operator=
void operator=(const UList< T > &)
Assignment from UList operator. Takes linear time.
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:42
Foam::Xfer
A simple container for copying or transferring objects of type <T>.
Definition: Xfer.H:85
Foam::triSurface::~triSurface
~triSurface()
Definition: triSurface.C:735
Foam::writeVTK
void writeVTK(OFstream &os, const Type &value)
OFstream.H
Foam::triSurface::movePoints
virtual void movePoints(const pointField &)
Move points.
Definition: triSurface.C:789
nPoints
label nPoints
Definition: gmvOutputHeader.H:2
Foam::deleteDemandDrivenData
void deleteDemandDrivenData(DataPtr &dataPtr)
Definition: demandDrivenData.H:40
List::size
int size() const
Definition: ListI.H:83
n
label n
Definition: TABSMDCalcMethod2.H:31
SortableList.H
Foam::triSurface::cleanup
void cleanup(const bool verbose=false)
Remove non-valid triangles.
Definition: triSurface.C:820
Foam::triSurface::subsetMesh
triSurface subsetMesh(const boolList &include, labelList &pointMap, labelList &faceMap) const
Return new surface. Returns pointMap, faceMap from.
Definition: triSurface.C:1013
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::triSurface
Triangulated surface description with patch information.
Definition: triSurface.H:57
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::nl
static const char nl
Definition: Ostream.H:260
PackedBoolList.H
Foam::OSstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: OSstream.H:88
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:837
Foam::triSurface::operator=
void operator=(const triSurface &)
Definition: triSurface.C:1122
Foam::triSurface::triSurface
triSurface()
Construct null.
Definition: triSurface.C:608
Foam::PrimitivePatch< labelledTri, List, pointField, point >::FaceType
labelledTri FaceType
Definition: PrimitivePatchTemplate.H:98
IFstream.H
Foam::operator<<
Ostream & operator<<(Ostream &, const edgeMesh &)
Definition: edgeMeshIO.C:130
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::triSurface::clearTopology
void clearTopology()
Clear topology.
Definition: triSurface.C:743
Foam::triFaceList
List< triFace > triFaceList
list of triFaces
Definition: triFaceList.H:42
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::abort
errorManip< error > abort(error &err)
Definition: errorManip.H:131
Foam::isFile
bool isFile(const fileName &, const bool checkGzip=true)
Does the name exist as a FILE in the file system?
Definition: POSIX.C:622
Foam::boolList
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
Foam::triSurface::edgeOwner
const labelList & edgeOwner() const
If 2 face neighbours: label of face where ordering of edge.
Definition: triSurface.C:777
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
Foam::triSurface::calcPatches
surfacePatchList calcPatches(labelList &faceMap) const
Sort faces according to region. Returns patch list.
Definition: triSurface.C:501
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::triSurface::setDefaultPatches
void setDefaultPatches()
Sets default values for patches.
Definition: triSurface.C:587
boundBox.H
Foam::labelListList
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
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:906
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:318
Foam::triSurface::checkEdges
void checkEdges(const bool verbose=false)
Check triply (or more) connected edges.
Definition: triSurface.C:320
Foam::triSurface::scalePoints
virtual void scalePoints(const scalar &)
Scale points. A non-positive factor is ignored.
Definition: triSurface.C:803
Foam::Pout
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
Definition: faceListFwd.H:43
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
points
const pointField & points
Definition: gmvOutputHeader.H:1
patches
patches[0]
Definition: createSingleCellMesh.H:36
List
Definition: Test.C:19
Foam::triSurface::sortedEdgeFaces
const labelListList & sortedEdgeFaces() const
Return edge-face addressing sorted (for edges with more than.
Definition: triSurface.C:766
Foam::sort
void sort(UList< T > &)
Definition: UList.C:107
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::triSurface::write
void write(const fileName &, const word &ext, const bool sort) const
Generic write routine. Chooses writer based on extension.
Definition: triSurface.C:433
Foam::surfacePatchList
List< surfacePatch > surfacePatchList
Definition: surfacePatchList.H:44
Foam::defineTypeNameAndDebug
defineTypeNameAndDebug(combustionModel, 0)
Foam::IOstream::name
virtual const fileName & name() const
Return the name of the stream.
Definition: IOstream.H:297
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
triFace
face triFace(3)
WarningInFunction
#define WarningInFunction
Report a warning using Foam::Warning.
Definition: messageStream.H:259
Foam::triSurface::getLineNoComment
static string getLineNoComment(IFstream &)
Read non-comment line.
Definition: triSurface.C:172
write
Tcoeff write()
Foam::triSurface::triSurfInstance
static fileName triSurfInstance(const Time &)
Name of triSurface directory to use.
Definition: triSurface.C:43
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Foam::triSurface::printTriangle
static void printTriangle(Ostream &, const Foam::string &pre, const labelledTri &, const pointField &)
Helper function to print triangle info.
Definition: triSurface.C:156
Foam::PrimitivePatch
A list of faces which address into the list of points.
Definition: PrimitivePatchTemplate.H:88