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