triSurface.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2024 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 #include "plane.H"
35 #include "tensor2D.H"
36 #include "symmTensor2D.H"
37 #include "transform.H"
38 #include "OSspecific.H"
39 #include "scalarMatrices.H"
40 
41 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
42 
43 namespace Foam
44 {
46 }
47 
48 
49 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
50 
52 {
53  fileName foamName(d.caseName() + ".ftr");
54 
55  // Search back through the time directories list to find the time
56  // closest to and lower than current time
57 
58  instantList ts = d.times();
59  label i;
60 
61  for (i=ts.size()-1; i>=0; i--)
62  {
63  if (ts[i].value() <= d.userTimeValue())
64  {
65  break;
66  }
67  }
68 
69  // Noting that the current directory has already been searched
70  // for mesh data, start searching from the previously stored time directory
71 
72  if (i>=0)
73  {
74  for (label j=i; j>=0; j--)
75  {
76  if (isFile(d.path()/ts[j].name()/typeName/foamName))
77  {
78  if (debug)
79  {
80  Pout<< " triSurface::triSurfInstance(const Time& d)"
81  << "reading " << foamName
82  << " from " << ts[j].name()/typeName
83  << endl;
84  }
85 
86  return ts[j].name();
87  }
88  }
89  }
90 
91  if (debug)
92  {
93  Pout<< " triSurface::triSurfInstance(const Time& d)"
94  << "reading " << foamName
95  << " from constant/" << endl;
96  }
97  return d.constant();
98 }
99 
100 
101 Foam::List<Foam::labelledTri> Foam::triSurface::convertToTri
102 (
103  const faceList& faces,
104  const label defaultRegion
105 )
106 {
107  List<labelledTri> triFaces(faces.size());
108 
109  forAll(triFaces, facei)
110  {
111  const face& f = faces[facei];
112 
113  if (f.size() != 3)
114  {
116  << "Face at position " << facei
117  << " does not have three vertices:" << f
118  << abort(FatalError);
119  }
120 
121  labelledTri& tri = triFaces[facei];
122 
123  tri[0] = f[0];
124  tri[1] = f[1];
125  tri[2] = f[2];
126  tri.region() = defaultRegion;
127  }
128 
129  return triFaces;
130 }
131 
132 
133 Foam::List<Foam::labelledTri> Foam::triSurface::convertToTri
134 (
135  const triFaceList& faces,
136  const label defaultRegion
137 )
138 {
139  List<labelledTri> triFaces(faces.size());
140 
141  forAll(triFaces, facei)
142  {
143  const triFace& f = faces[facei];
144 
145  labelledTri& tri = triFaces[facei];
146 
147  tri[0] = f[0];
148  tri[1] = f[1];
149  tri[2] = f[2];
150  tri.region() = defaultRegion;
151  }
152 
153  return triFaces;
154 }
155 
156 
157 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
158 
159 void Foam::triSurface::printTriangle
160 (
161  Ostream& os,
162  const string& pre,
163  const labelledTri& f,
164  const pointField& points
165 )
166 {
167  os
168  << pre.c_str() << "vertex numbers:"
169  << f[0] << ' ' << f[1] << ' ' << f[2] << endl
170  << pre.c_str() << "vertex coords :"
171  << points[f[0]] << ' ' << points[f[1]] << ' ' << points[f[2]]
172  << pre.c_str() << "region :" << f.region() << endl
173  << endl;
174 }
175 
176 
177 Foam::string Foam::triSurface::getLineNoComment(IFstream& is)
178 {
179  string line;
180  do
181  {
182  is.getLine(line);
183  }
184  while ((line.empty() || line[0] == '#') && is.good());
185 
186  return line;
187 }
188 
189 
190 Foam::scalar Foam::triSurface::pointNormalWeight
191 (
192  const triFace& f,
193  const label pi,
194  const vector& fa,
195  const pointField& points
196 ) const
197 {
198  const label index = findIndex(f, pi);
199 
200  if (index == -1)
201  {
203  << "Point not in face" << abort(FatalError);
204  }
205 
206  const vector e1 = points[f[index]] - points[f[f.fcIndex(index)]];
207  const vector e2 = points[f[index]] - points[f[f.rcIndex(index)]];
208 
209  return mag(fa)/(magSqr(e1)*magSqr(e2) + vSmall);
210 }
211 
212 
213 Foam::tmp<Foam::vectorField> Foam::triSurface::weightedPointNormals() const
214 {
215  // Weighted average of normals of faces attached to the vertex
216  // Weight = fA / (mag(e0)^2 * mag(e1)^2);
217 
218  tmp<vectorField> tpointNormals
219  (
220  new vectorField(nPoints(), Zero)
221  );
222  vectorField& pointNormals = tpointNormals.ref();
223 
224  const pointField& points = this->points();
225  const labelListList& pointFaces = this->pointFaces();
226  const labelList& meshPoints = this->meshPoints();
227 
228  forAll(pointFaces, pi)
229  {
230  const labelList& pFaces = pointFaces[pi];
231 
232  forAll(pFaces, fi)
233  {
234  const label facei = pFaces[fi];
235  const triFace& f = operator[](facei);
236 
237  const vector fa = f.area(points);
238  const scalar weight =
239  pointNormalWeight(f, meshPoints[pi], fa, points);
240 
241  pointNormals[pi] += weight*fa;
242  }
243 
244  pointNormals[pi] /= mag(pointNormals[pi]) + vSmall;
245  }
246 
247  return tpointNormals;
248 }
249 
250 
251 Foam::tmp<Foam::triadField> Foam::triSurface::pointCoordSys
252 (
253  const vectorField& pointNormals
254 ) const
255 {
256  tmp<triadField> tpointCoordSys(new triadField(nPoints()));
257  triadField& pointCoordSys = tpointCoordSys.ref();
258 
259  const pointField& points = this->points();
260  const labelList& meshPoints = this->meshPoints();
261 
262  forAll(meshPoints, pi)
263  {
264  const point& pt = points[meshPoints[pi]];
265  const vector& normal = pointNormals[pi];
266 
267  if (mag(normal) < small)
268  {
269  pointCoordSys[pi] = triad::unset;
270  continue;
271  }
272 
273  plane p(pt, normal);
274 
275  // Pick a point in plane
276  vector dir1 = pt - p.aPoint();
277  dir1 /= mag(dir1);
278 
279  vector dir2 = dir1 ^ normal;
280  dir2 /= mag(dir2);
281 
282  pointCoordSys[pi] = triad(dir1, dir2, normal);
283  }
284 
285  return pointCoordSys;
286 }
287 
288 
289 bool Foam::triSurface::read(Istream& is)
290 {
291  is >> patches_ >> storedPoints() >> storedFaces();
292 
293  return true;
294 }
295 
296 
297 bool Foam::triSurface::read
298 (
299  const fileName& name,
300  const word& ext,
301  const bool check
302 )
303 {
304  if (check && !exists(name))
305  {
307  << "Cannot read " << name << exit(FatalError);
308  }
309 
310  if (ext == "gz")
311  {
312  fileName unzipName = name.lessExt();
313 
314  // Do not check for existence. Let IFstream do the unzipping.
315  return read(unzipName, unzipName.ext(), false);
316  }
317  else if (ext == "ftr")
318  {
319  return read(IFstream(name)());
320  }
321  else if (ext == "stl")
322  {
323  return readSTL(name);
324  }
325  else if (ext == "stlb")
326  {
327  return readSTLBINARY(name);
328  }
329  else if (ext == "gts")
330  {
331  return readGTS(name);
332  }
333  else if (ext == "obj")
334  {
335  return readOBJ(name);
336  }
337  else if (ext == "off")
338  {
339  return readOFF(name);
340  }
341  else if (ext == "tri")
342  {
343  return readTRI(name);
344  }
345  else if (ext == "ac")
346  {
347  return readAC(name);
348  }
349  else if (ext == "nas")
350  {
351  return readNAS(name);
352  }
353  else if (ext == "vtk")
354  {
355  return readVTK(name);
356  }
357  else
358  {
360  << "unknown file extension " << ext
361  << ". Supported extensions are '.ftr', '.stl', '.stlb', '.gts'"
362  << ", '.obj', '.ac', '.off', '.nas', '.tri' and '.vtk'"
363  << exit(FatalError);
364 
365  return false;
366  }
367 }
368 
369 
370 void Foam::triSurface::write
371 (
372  const fileName& name,
373  const word& ext,
374  const bool sort
375 ) const
376 {
377  if (ext == "ftr")
378  {
379  return write(OFstream(name)());
380  }
381  else if (ext == "stl")
382  {
383  return writeSTLASCII(sort, OFstream(name)());
384  }
385  else if (ext == "stlb")
386  {
387  ofstream outFile(name.c_str(), std::ios::binary);
388 
389  writeSTLBINARY(outFile);
390  }
391  else if (ext == "gts")
392  {
393  return writeGTS(sort, OFstream(name)());
394  }
395  else if (ext == "obj")
396  {
397  writeOBJ(sort, OFstream(name)());
398  }
399  else if (ext == "off")
400  {
401  writeOFF(sort, OFstream(name)());
402  }
403  else if (ext == "vtk")
404  {
405  writeVTK(sort, OFstream(name)());
406  }
407  else if (ext == "tri")
408  {
409  writeTRI(sort, OFstream(name)());
410  }
411  else if (ext == "ac")
412  {
413  writeAC(OFstream(name)());
414  }
415  else if (ext == "smesh")
416  {
417  writeSMESH(sort, OFstream(name)());
418  }
419  else
420  {
422  << "unknown file extension " << ext
423  << " for file " << name
424  << ". Supported extensions are '.ftr', '.stl', '.stlb', "
425  << "'.gts', '.obj', '.vtk'"
426  << ", '.off', '.dx', '.smesh', '.ac' and '.tri'"
427  << exit(FatalError);
428  }
429 }
430 
431 
432 Foam::surfacePatchList Foam::triSurface::calcPatches(labelList& faceMap) const
433 {
434  // Sort according to region numbers of labelledTri
435  SortableList<label> sortedRegion(size());
436 
437  forAll(sortedRegion, facei)
438  {
439  sortedRegion[facei] = operator[](facei).region();
440  }
441  sortedRegion.sort();
442 
443  faceMap = sortedRegion.indices();
444 
445  // Extend regions
446  label maxRegion = patches_.size()-1; // for non-compacted regions
447 
448  if (faceMap.size())
449  {
450  maxRegion = max
451  (
452  maxRegion,
453  operator[](faceMap.last()).region()
454  );
455  }
456 
457  // Get new region list
458  surfacePatchList newPatches(maxRegion + 1);
459 
460  // Fill patch sizes
461  forAll(*this, facei)
462  {
463  label region = operator[](facei).region();
464 
465  newPatches[region].size()++;
466  }
467 
468 
469  // Fill rest of patch info
470 
471  label startFacei = 0;
472  forAll(newPatches, newPatchi)
473  {
474  surfacePatch& newPatch = newPatches[newPatchi];
475 
476  newPatch.index() = newPatchi;
477 
478  label oldPatchi = newPatchi;
479 
480  // start of patch
481  newPatch.start() = startFacei;
482 
483 
484  // Take over any information from existing patches
485  if ((oldPatchi < patches_.size()) && (patches_[oldPatchi].name() != ""))
486  {
487  newPatch.name() = patches_[oldPatchi].name();
488  }
489  else
490  {
491  newPatch.name() = word("patch") + name(newPatchi);
492  }
493 
494  if
495  (
496  (oldPatchi < patches_.size())
497  && (patches_[oldPatchi].geometricType() != "")
498  )
499  {
500  newPatch.geometricType() = patches_[oldPatchi].geometricType();
501  }
502  else
503  {
504  newPatch.geometricType() = "empty";
505  }
506 
507  startFacei += newPatch.size();
508  }
509 
510  return newPatches;
511 }
512 
513 
514 void Foam::triSurface::setDefaultPatches()
515 {
517 
518  // Get names, types and sizes
519  surfacePatchList newPatches(calcPatches(faceMap));
520 
521  // Take over names and types (but not size)
522  patches_.setSize(newPatches.size());
523 
524  forAll(newPatches, patchi)
525  {
526  patches_[patchi].index() = patchi;
527  patches_[patchi].name() = newPatches[patchi].name();
528  patches_[patchi].geometricType() = newPatches[patchi].geometricType();
529  }
530 }
531 
532 
533 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
534 
536 :
538  patches_(0),
539  sortedEdgeFacesPtr_(nullptr),
540  edgeOwnerPtr_(nullptr)
541 {}
542 
543 
545 (
546  const List<labelledTri>& triangles,
548  const pointField& points
549 )
550 :
551  ParentType(triangles, points),
552  patches_(patches),
553  sortedEdgeFacesPtr_(nullptr),
554  edgeOwnerPtr_(nullptr)
555 {}
556 
557 
559 (
560  List<labelledTri>& triangles,
563  const bool reuse
564 )
565 :
566  ParentType(triangles, points, reuse),
567  patches_(patches),
568  sortedEdgeFacesPtr_(nullptr),
569  edgeOwnerPtr_(nullptr)
570 {}
571 
572 
574 (
575  List<labelledTri>&& triangles,
578 )
579 :
580  ParentType(move(triangles), move(points)),
581  patches_(patches),
582  sortedEdgeFacesPtr_(nullptr),
583  edgeOwnerPtr_(nullptr)
584 {}
585 
586 
588 (
589  const List<labelledTri>& triangles,
590  const pointField& points
591 )
592 :
593  ParentType(triangles, points),
594  patches_(),
595  sortedEdgeFacesPtr_(nullptr),
596  edgeOwnerPtr_(nullptr)
597 {
598  setDefaultPatches();
599 }
600 
601 
603 (
604  const triFaceList& triangles,
605  const pointField& points
606 )
607 :
608  ParentType(convertToTri(triangles, 0), points),
609  patches_(),
610  sortedEdgeFacesPtr_(nullptr),
611  edgeOwnerPtr_(nullptr)
612 {
613  setDefaultPatches();
614 }
615 
616 
618 :
620  patches_(),
621  sortedEdgeFacesPtr_(nullptr),
622  edgeOwnerPtr_(nullptr)
623 {
624  word ext = name.ext();
625 
626  read(name, ext);
627 
628  setDefaultPatches();
629 }
630 
631 
633 :
635  patches_(),
636  sortedEdgeFacesPtr_(nullptr),
637  edgeOwnerPtr_(nullptr)
638 {
639  read(is);
640 
641  setDefaultPatches();
642 }
643 
644 
646 :
648  patches_(),
649  sortedEdgeFacesPtr_(nullptr),
650  edgeOwnerPtr_(nullptr)
651 {
652  fileName foamFile(d.caseName() + ".ftr");
653 
654  fileName foamPath(d.path()/triSurfInstance(d)/typeName/foamFile);
655 
656  IFstream foamStream(foamPath);
657 
658  read(foamStream);
659 
660  setDefaultPatches();
661 }
662 
663 
665 :
666  ParentType(ts, ts.points()),
667  patches_(ts.patches()),
668  sortedEdgeFacesPtr_(nullptr),
669  edgeOwnerPtr_(nullptr)
670 {}
671 
672 
674 :
675  ParentType(move(ts), move(ts.points())),
676  patches_(move(ts.patches())),
677  sortedEdgeFacesPtr_(nullptr),
678  edgeOwnerPtr_(nullptr)
679 {}
680 
681 
682 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
683 
685 {
686  clearOut();
687 }
688 
689 
690 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
691 
693 {
694  ParentType::clearTopology();
695  deleteDemandDrivenData(sortedEdgeFacesPtr_);
696  deleteDemandDrivenData(edgeOwnerPtr_);
697 }
698 
699 
701 {
702  ParentType::clearPatchMeshAddr();
703 }
704 
705 
707 {
708  ParentType::clearOut();
709 
710  clearTopology();
711  clearPatchMeshAddr();
712 }
713 
714 
716 {
717  if (!sortedEdgeFacesPtr_)
718  {
719  calcSortedEdgeFaces();
720  }
721 
722  return *sortedEdgeFacesPtr_;
723 }
724 
725 
727 {
728  if (!edgeOwnerPtr_)
729  {
730  calcEdgeOwner();
731  }
732 
733  return *edgeOwnerPtr_;
734 }
735 
736 
738 {
739  // Remove all geometry dependent data
740  deleteDemandDrivenData(sortedEdgeFacesPtr_);
741 
742  // Adapt for new point position
743  ParentType::clearGeom();
744 
745  // Copy new points
746  storedPoints() = newPoints;
747 }
748 
749 
750 void Foam::triSurface::scalePoints(const scalar scaleFactor)
751 {
752  // avoid bad scaling
753  if (scaleFactor > 0 && scaleFactor != 1.0)
754  {
755  // Remove all geometry dependent data
756  clearTopology();
757 
758  // Adapt for new point position
759  ParentType::clearGeom();
760 
761  storedPoints() *= scaleFactor;
762  }
763 }
764 
765 
766 void Foam::triSurface::checkTriangles(const bool verbose)
767 {
768  // Simple check on indices ok.
769  const label maxPointi = points().size() - 1;
770 
771  forAll(*this, facei)
772  {
773  const triSurface::FaceType& f = (*this)[facei];
774 
775  forAll(f, fp)
776  {
777  if (f[fp] < 0 || f[fp] > maxPointi)
778  {
780  << "triangle " << f
781  << " uses point indices outside point range 0.."
782  << maxPointi
783  << exit(FatalError);
784  }
785  }
786  }
787 
788  // Two phase process
789  // 1. mark invalid faces
790  // 2. pack
791  // Done to keep numbering constant in phase 1
792 
793  // List of valid triangles
794  boolList valid(size(), true);
795  bool hasInvalid = false;
796 
797  forAll(*this, facei)
798  {
799  const labelledTri& f = (*this)[facei];
800 
801  if ((f[0] == f[1]) || (f[0] == f[2]) || (f[1] == f[2]))
802  {
803  // 'degenerate' triangle check
804  valid[facei] = false;
805  hasInvalid = true;
806 
807  if (verbose)
808  {
810  << "triangle " << facei
811  << " does not have three unique vertices:\n";
812  printTriangle(Warning, " ", f, points());
813  }
814  }
815  else
816  {
817  // duplicate triangle check
818  const labelList& fEdges = faceEdges()[facei];
819 
820  // Check if faceNeighbours use same points as this face.
821  // Note: discards normal information - sides of baffle are merged.
822 
823  forAll(fEdges, fp)
824  {
825  const labelList& eFaces = edgeFaces()[fEdges[fp]];
826 
827  forAll(eFaces, i)
828  {
829  label neighbour = eFaces[i];
830 
831  if (neighbour > facei)
832  {
833  // lower numbered faces already checked
834  const labelledTri& n = (*this)[neighbour];
835 
836  if
837  (
838  ((f[0] == n[0]) || (f[0] == n[1]) || (f[0] == n[2]))
839  && ((f[1] == n[0]) || (f[1] == n[1]) || (f[1] == n[2]))
840  && ((f[2] == n[0]) || (f[2] == n[1]) || (f[2] == n[2]))
841  )
842  {
843  valid[facei] = false;
844  hasInvalid = true;
845 
846  if (verbose)
847  {
849  << "triangles share the same vertices:\n"
850  << " face 1 :" << facei << endl;
851  printTriangle(Warning, " ", f, points());
852 
853  Warning
854  << endl
855  << " face 2 :"
856  << neighbour << endl;
857  printTriangle(Warning, " ", n, points());
858  }
859 
860  break;
861  }
862  }
863  }
864  }
865  }
866  }
867 
868  if (hasInvalid)
869  {
870  // Pack
871  label newFacei = 0;
872  forAll(*this, facei)
873  {
874  if (valid[facei])
875  {
876  const labelledTri& f = (*this)[facei];
877  (*this)[newFacei++] = f;
878  }
879  }
880 
881  if (verbose)
882  {
884  << "Removing " << size() - newFacei
885  << " illegal faces." << endl;
886  }
887  (*this).setSize(newFacei);
888 
889  // Topology can change because of renumbering
890  clearOut();
891  }
892 }
893 
894 
895 void Foam::triSurface::checkEdges(const bool verbose)
896 {
897  const labelListList& eFaces = edgeFaces();
898 
899  forAll(eFaces, edgeI)
900  {
901  const labelList& myFaces = eFaces[edgeI];
902 
903  if (myFaces.empty())
904  {
906  << "Edge " << edgeI << " with vertices " << edges()[edgeI]
907  << " has no edgeFaces"
908  << exit(FatalError);
909  }
910  else if (myFaces.size() > 2 && verbose)
911  {
913  << "Edge " << edgeI << " with vertices " << edges()[edgeI]
914  << " has more than 2 faces connected to it : " << myFaces
915  << endl;
916  }
917  }
918 }
919 
920 
921 void Foam::triSurface::cleanup(const bool verbose)
922 {
923  // Merge points (already done for STL, TRI)
924  stitchTriangles(small, verbose);
925 
926  // Merging points might have changed geometric factors
927  clearOut();
928 
929  checkTriangles(verbose);
930 
931  checkEdges(verbose);
932 }
933 
934 
936 (
937  const boolList& borderEdge,
938  const label facei,
939  const label currentZone,
941 ) const
942 {
943  // List of faces whose faceZone has been set.
944  labelList changedFaces(1, facei);
945 
946  while (true)
947  {
948  // Pick up neighbours of changedFaces
949  DynamicList<label> newChangedFaces(2*changedFaces.size());
950 
951  forAll(changedFaces, i)
952  {
953  label facei = changedFaces[i];
954 
955  const labelList& fEdges = faceEdges()[facei];
956 
957  forAll(fEdges, i)
958  {
959  label edgeI = fEdges[i];
960 
961  if (!borderEdge[edgeI])
962  {
963  const labelList& eFaces = edgeFaces()[edgeI];
964 
965  forAll(eFaces, j)
966  {
967  label nbrFacei = eFaces[j];
968 
969  if (faceZone[nbrFacei] == -1)
970  {
971  faceZone[nbrFacei] = currentZone;
972  newChangedFaces.append(nbrFacei);
973  }
974  else if (faceZone[nbrFacei] != currentZone)
975  {
977  << "Zones " << faceZone[nbrFacei]
978  << " at face " << nbrFacei
979  << " connects to zone " << currentZone
980  << " at face " << facei
981  << abort(FatalError);
982  }
983  }
984  }
985  }
986  }
987 
988  if (newChangedFaces.empty())
989  {
990  break;
991  }
992 
993  changedFaces.transfer(newChangedFaces);
994  }
995 }
996 
997 
999 (
1000  const boolList& borderEdge,
1002 ) const
1003 {
1004  faceZone.setSize(size());
1005  faceZone = -1;
1006 
1007  if (borderEdge.size() != nEdges())
1008  {
1010  << "borderEdge boolList not same size as number of edges" << endl
1011  << "borderEdge:" << borderEdge.size() << endl
1012  << "nEdges :" << nEdges()
1013  << exit(FatalError);
1014  }
1015 
1016  label zoneI = 0;
1017 
1018  for (label startFacei = 0;; zoneI++)
1019  {
1020  // Find first non-coloured face
1021  for (; startFacei < size(); startFacei++)
1022  {
1023  if (faceZone[startFacei] == -1)
1024  {
1025  break;
1026  }
1027  }
1028 
1029  if (startFacei >= size())
1030  {
1031  break;
1032  }
1033 
1034  faceZone[startFacei] = zoneI;
1035 
1036  markZone(borderEdge, startFacei, zoneI, faceZone);
1037  }
1038 
1039  return zoneI;
1040 }
1041 
1042 
1044 (
1045  const boolList& include,
1046  labelList& pointMap,
1048 ) const
1049 {
1050  const List<labelledTri>& locFaces = localFaces();
1051 
1052  label facei = 0;
1053  label pointi = 0;
1054 
1055  faceMap.setSize(locFaces.size());
1056 
1057  pointMap.setSize(nPoints());
1058 
1059  boolList pointHad(nPoints(), false);
1060 
1061  forAll(include, oldFacei)
1062  {
1063  if (include[oldFacei])
1064  {
1065  // Store new faces compact
1066  faceMap[facei++] = oldFacei;
1067 
1068  // Renumber labels for face
1069  const triSurface::FaceType& f = locFaces[oldFacei];
1070 
1071  forAll(f, fp)
1072  {
1073  label labI = f[fp];
1074  if (!pointHad[labI])
1075  {
1076  pointHad[labI] = true;
1077  pointMap[pointi++] = labI;
1078  }
1079  }
1080  }
1081  }
1082 
1083  // Trim
1084  faceMap.setSize(facei);
1085  pointMap.setSize(pointi);
1086 }
1087 
1088 
1090 (
1091  const boolList& include,
1092  labelList& pointMap,
1094 ) const
1095 {
1096  const pointField& locPoints = localPoints();
1097  const List<labelledTri>& locFaces = localFaces();
1098 
1099  // Fill pointMap, faceMap
1100  subsetMeshMap(include, pointMap, faceMap);
1101 
1102 
1103  // Create compact coordinate list and forward mapping array
1104  pointField newPoints(pointMap.size());
1105  labelList oldToNew(locPoints.size());
1106  forAll(pointMap, pointi)
1107  {
1108  newPoints[pointi] = locPoints[pointMap[pointi]];
1109  oldToNew[pointMap[pointi]] = pointi;
1110  }
1111 
1112  // Renumber triangle node labels and compact
1113  List<labelledTri> newTriangles(faceMap.size());
1114 
1115  forAll(faceMap, facei)
1116  {
1117  // Get old vertex labels
1118  const labelledTri& tri = locFaces[faceMap[facei]];
1119 
1120  newTriangles[facei][0] = oldToNew[tri[0]];
1121  newTriangles[facei][1] = oldToNew[tri[1]];
1122  newTriangles[facei][2] = oldToNew[tri[2]];
1123  newTriangles[facei].region() = tri.region();
1124  }
1125 
1126  // Construct subsurface
1127  return triSurface(newTriangles, patches(), newPoints, true);
1128 }
1129 
1130 
1132 {
1133  faceList faces(size());
1134 
1135  forAll(*this, facei)
1136  {
1137  faces[facei] = operator[](facei).triFaceFace();
1138  }
1139 
1140  return faces;
1141 }
1142 
1143 
1145 {
1146  // Curvature calculation is an implementation of the algorithm from:
1147  // "Estimating Curvatures and their Derivatives on Triangle Meshes"
1148  // by S. Rusinkiewicz
1149 
1150  const pointField& points = this->points();
1151  const labelList& meshPoints = this->meshPoints();
1152  const Map<label>& meshPointMap = this->meshPointMap();
1153 
1154  const vectorField pointNormals
1155  (
1156  this->weightedPointNormals()
1157  );
1158 
1159  const triadField pointCoordSys
1160  (
1161  this->pointCoordSys(pointNormals)
1162  );
1163 
1164  tmp<scalarField> tcurvaturePointField(new scalarField(nPoints(), 0.0));
1165  scalarField& curvaturePointField = tcurvaturePointField.ref();
1166 
1167  List<symmTensor2D> pointFundamentalTensors(nPoints(), Zero);
1168  scalarList accumulatedWeights(nPoints(), 0.0);
1169 
1170  forAll(*this, fi)
1171  {
1172  const triFace& f = operator[](fi);
1173  const edgeList fEdges = f.edges();
1174 
1175  // Calculate the edge vectors and the normal differences
1176  vectorField edgeVectors(f.size(), Zero);
1177  vectorField normalDifferences(f.size(), Zero);
1178 
1179  forAll(fEdges, feI)
1180  {
1181  const edge& e = fEdges[feI];
1182 
1183  edgeVectors[feI] = e.vec(points);
1184  normalDifferences[feI] =
1185  pointNormals[meshPointMap[e[0]]]
1186  - pointNormals[meshPointMap[e[1]]];
1187  }
1188 
1189  // Set up a local coordinate system for the face
1190  const vector& e0 = edgeVectors[0];
1191  const vector eN = f.area(points);
1192  const vector e1 = (e0 ^ eN);
1193 
1194  if (magSqr(eN) < rootVSmall)
1195  {
1196  continue;
1197  }
1198 
1199  triad faceCoordSys(e0, e1, eN);
1200  faceCoordSys.normalise();
1201 
1202  // Construct the matrix to solve
1204  scalarDiagonalMatrix Z(3, 0);
1205 
1206  // Least Squares
1207  for (label i = 0; i < 3; ++i)
1208  {
1209  const scalar x = edgeVectors[i] & faceCoordSys[0];
1210  const scalar y = edgeVectors[i] & faceCoordSys[1];
1211 
1212  T(0, 0) += sqr(x);
1213  T(1, 0) += x*y;
1214  T(1, 1) += sqr(x) + sqr(y);
1215  T(2, 1) += x*y;
1216  T(2, 2) += sqr(y);
1217 
1218  const scalar dndx = normalDifferences[i] & faceCoordSys[0];
1219  const scalar dndy = normalDifferences[i] & faceCoordSys[1];
1220 
1221  Z[0] += dndx*x;
1222  Z[1] += dndx*y + dndy*x;
1223  Z[2] += dndy*y;
1224  }
1225 
1226  // Perform Cholesky decomposition and back substitution.
1227  // Decomposed matrix is in T and solution is in Z.
1228  LUsolve(T, Z);
1229  const symmTensor2D secondFundamentalTensor(Z[0], Z[1], Z[2]);
1230 
1231  // Loop over the face points adding the contribution of the face
1232  // curvature to the points.
1233  forAll(f, fpi)
1234  {
1235  const label patchPointIndex = meshPointMap[f[fpi]];
1236 
1237  const triad& ptCoordSys = pointCoordSys[patchPointIndex];
1238 
1239  if (!ptCoordSys.set())
1240  {
1241  continue;
1242  }
1243 
1244  // Rotate faceCoordSys to ptCoordSys
1245  const tensor rotTensor = rotationTensor
1246  (
1247  ptCoordSys[2],
1248  faceCoordSys[2]
1249  );
1250  const triad rotatedFaceCoordSys = rotTensor & tensor(faceCoordSys);
1251 
1252  // Project the face curvature onto the point plane
1253 
1254  const vector2D cmp1
1255  (
1256  ptCoordSys[0] & rotatedFaceCoordSys[0],
1257  ptCoordSys[0] & rotatedFaceCoordSys[1]
1258  );
1259  const vector2D cmp2
1260  (
1261  ptCoordSys[1] & rotatedFaceCoordSys[0],
1262  ptCoordSys[1] & rotatedFaceCoordSys[1]
1263  );
1264 
1265  const tensor2D projTensor
1266  (
1267  cmp1,
1268  cmp2
1269  );
1270 
1271  const symmTensor2D projectedFundamentalTensor
1272  (
1273  projTensor.x() & (secondFundamentalTensor & projTensor.x()),
1274  projTensor.x() & (secondFundamentalTensor & projTensor.y()),
1275  projTensor.y() & (secondFundamentalTensor & projTensor.y())
1276  );
1277 
1278  // Calculate weight
1279  const scalar weight = pointNormalWeight
1280  (
1281  f,
1282  meshPoints[patchPointIndex],
1283  f.area(points),
1284  points
1285  );
1286 
1287  // Sum contribution of face to this point
1288  pointFundamentalTensors[patchPointIndex] +=
1289  weight*projectedFundamentalTensor;
1290 
1291  accumulatedWeights[patchPointIndex] += weight;
1292  }
1293  }
1294 
1295  forAll(curvaturePointField, pi)
1296  {
1297  pointFundamentalTensors[pi] /= (accumulatedWeights[pi] + small);
1298 
1299  const vector2D principalCurvatures =
1300  eigenValues(pointFundamentalTensors[pi]);
1301 
1302  // scalar curvature =
1303  // (principalCurvatures[0] + principalCurvatures[1])/2;
1304  const scalar curvature = max
1305  (
1306  mag(principalCurvatures[0]),
1307  mag(principalCurvatures[1])
1308  );
1309  // scalar curvature = principalCurvatures[0]*principalCurvatures[1];
1310 
1311  curvaturePointField[meshPoints[pi]] = curvature;
1312  }
1313 
1314  return tcurvaturePointField;
1315 }
1316 
1317 
1318 void Foam::triSurface::write
1319 (
1320  const fileName& name,
1321  const bool sortByRegion
1322 ) const
1323 {
1324  write(name, name.ext(), sortByRegion);
1325 }
1326 
1327 
1328 void Foam::triSurface::write(Ostream& os) const
1329 {
1330  os << patches() << endl;
1331 
1332  // Note: Write with global point numbering
1333  os << points() << nl
1334  << static_cast<const List<labelledTri>&>(*this) << endl;
1335 
1336  // Check state of Ostream
1337  os.check("triSurface::write(Ostream&)");
1338 }
1339 
1340 
1341 void Foam::triSurface::write(const Time& d) const
1342 {
1343  fileName foamFile(d.caseName() + ".ftr");
1344 
1345  fileName foamPath(d.path()/triSurfInstance(d)/typeName/foamFile);
1346 
1347  OFstream foamStream(foamPath);
1348 
1349  write(foamStream);
1350 }
1351 
1352 
1354 {
1355  // Unfortunately nPoints constructs meshPoints() so do compact version
1356  // ourselves.
1357  PackedBoolList pointIsUsed(points().size());
1358 
1359  label nPoints = 0;
1361 
1362  forAll(*this, facei)
1363  {
1364  const triSurface::FaceType& f = operator[](facei);
1365 
1366  forAll(f, fp)
1367  {
1368  label pointi = f[fp];
1369  if (pointIsUsed.set(pointi, 1))
1370  {
1371  bb.min() = ::Foam::min(bb.min(), points()[pointi]);
1372  bb.max() = ::Foam::max(bb.max(), points()[pointi]);
1373  nPoints++;
1374  }
1375  }
1376  }
1377 
1378  os << "Triangles : " << size() << endl
1379  << "Vertices : " << nPoints << endl
1380  << "Bounding Box : " << bb << endl;
1381 }
1382 
1383 
1384 // * * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * //
1385 
1387 {
1389  clearOut();
1390  storedPoints() = ts.points();
1391  patches_ = ts.patches();
1392 }
1393 
1394 
1396 {
1397  List<labelledTri>::operator=(move(ts));
1398  clearOut();
1399  storedPoints() = move(ts.points());
1400  patches_ = move(ts.patches());
1401 }
1402 
1403 
1404 // * * * * * * * * * * * * * * * IOstream Operators * * * * * * * * * * * * //
1405 
1407 {
1408  sm.write(os);
1409  return os;
1410 }
1411 
1412 
1413 // ************************************************************************* //
scalar y
Functions used by OpenFOAM that are specific to POSIX compliant operating systems and need to be repl...
label n
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
DynamicList< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
Definition: DynamicListI.H:296
Input from file stream.
Definition: IFstream.H:85
virtual bool check(const char *operation) const
Check IOstream status for given operation.
Definition: IOstream.C:92
An Istream is an abstract base class for all input systems (streams, files, token lists etc)....
Definition: Istream.H:60
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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
Output to file stream.
Definition: OFstream.H:86
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
A bit-packed bool list.
void set(const PackedList< 1 > &)
Set specified bits.
const Field< PointType > & points() const
Return reference to global points.
std::remove_reference< ::Foam::List< labelledTri > >::type::value_type FaceType
Templated 2D symmetric tensor derived from VectorSpace adding construction from 4 components,...
Definition: SymmTensor2D.H:56
A templated 2D square symmetric matrix of objects of <T>, where the n x n matrix dimension is known a...
Templated 2D tensor derived from VectorSpace adding construction from 4 components,...
Definition: Tensor2D.H:59
Vector2D< Cmpt > y() const
Definition: Tensor2DI.H:106
Vector2D< Cmpt > x() const
Definition: Tensor2DI.H:100
static const word & constant()
Return constant name.
Definition: TimePaths.H:122
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
instantList times() const
Search the case for valid time directories.
Definition: Time.C:651
fileName path() const
Explicitly inherit path from TimePaths to disambiguate from.
Definition: TimePaths.H:138
scalar userTimeValue() const
Return current user time value.
Definition: Time.C:830
const fileName & caseName() const
Explicitly inherit caseName from TimePaths to disambiguate from.
Definition: TimePaths.H:106
label rcIndex(const label i) const
Return the reverse circular index, i.e. the previous index.
Definition: UListI.H:65
label fcIndex(const label i) const
Return the forward circular index, i.e. the next index.
Definition: UListI.H:58
bool empty() const
Return true if the UList is empty (ie, size() is zero)
Definition: UListI.H:325
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
const point & min() const
Minimum point defining the bounding box.
Definition: boundBoxI.H:60
const point & max() const
Maximum point defining the bounding box.
Definition: boundBoxI.H:66
static const boundBox invertedBox
A very large inverted boundBox: min/max == +/- vGreat.
Definition: boundBox.H:82
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A subset of mesh faces organised as a primitive patch.
Definition: faceZone.H:59
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
A class for handling file names.
Definition: fileName.H:82
Triangle with additional region number.
Definition: labelledTri.H:60
label region() const
Return region label.
Definition: labelledTriI.H:68
A class for handling character strings derived from std::string.
Definition: string.H:79
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A triangular face using a FixedList of labels corresponding to mesh vertices.
Definition: triFace.H:71
Triangulated surface description with patch information.
Definition: triSurface.H:69
virtual void scalePoints(const scalar)
Scale points. A non-positive factor is ignored.
Definition: triSurface.C:750
static fileName triSurfInstance(const Time &)
Name of triSurface directory to use.
Definition: triSurface.C:51
triSurface()
Construct null.
Definition: triSurface.C:535
triSurface subsetMesh(const boolList &include, labelList &pointMap, labelList &faceMap) const
Return new surface. Returns pointMap, faceMap from.
Definition: triSurface.C:1090
label markZones(const boolList &borderEdge, labelList &faceZone) const
(size and) fills faceZone with zone of face. Zone is area
Definition: triSurface.C:999
tmp< scalarField > curvature() const
Return the curvature of surface at the points.
Definition: triSurface.C:1144
void operator=(const triSurface &)
Definition: triSurface.C:1386
void cleanup(const bool verbose)
Remove non-valid triangles.
Definition: triSurface.C:921
void clearTopology()
Definition: triSurface.C:692
const geometricSurfacePatchList & patches() const
Definition: triSurface.H:322
void checkTriangles(const bool verbose)
Check/remove duplicate/degenerate triangles.
Definition: triSurface.C:766
void checkEdges(const bool verbose)
Check triply (or more) connected edges.
Definition: triSurface.C:895
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:936
void writeStats(Ostream &) const
Write some statistics.
Definition: triSurface.C:1353
faceList faces() const
Return the list of triangles as a faceList.
Definition: triSurface.C:1131
const labelList & edgeOwner() const
If 2 face neighbours: label of face where ordering of edge.
Definition: triSurface.C:726
void subsetMeshMap(const boolList &include, labelList &pointMap, labelList &faceMap) const
'Create' sub mesh, including only faces for which
Definition: triSurface.C:1044
void clearPatchMeshAddr()
Definition: triSurface.C:700
const labelListList & sortedEdgeFaces() const
Return edge-face addressing sorted (for edges with more than.
Definition: triSurface.C:715
virtual ~triSurface()
Destructor.
Definition: triSurface.C:684
virtual void setPoints(const pointField &)
Move points.
Definition: triSurface.C:737
Representation of a 3D Cartesian coordinate system as a Vector of vectors.
Definition: triad.H:68
static const triad unset
Definition: triad.H:99
bool set(const direction d) const
Is the vector in the direction d set.
Definition: triadI.H:61
void normalise()
Normalise each set axis vector to have a unit magnitude.
Definition: triadI.H:111
A class for handling words, derived from string.
Definition: word.H:62
Template functions to aid in the implementation of demand driven data.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
label patchi
const pointField & points
label nPoints
const fvPatchList & patches
#define WarningInFunction
Report a warning using Foam::Warning.
bool valid(const PtrList< ModelType > &l)
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
bool isFile(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist as a file in the file system?
Definition: POSIX.C:555
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
static const zero Zero
Definition: zero.H:97
const doubleScalar e
Definition: doubleScalar.H:106
bool read(const char *, int32_t &)
Definition: int32IO.C:85
tensor rotationTensor(const vector &n1, const vector &n2)
Rotational transformation tensor from unit vector n1 to n2.
Definition: transform.H:47
List< label > labelList
A List of labels.
Definition: labelList.H:56
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Tensor< scalar > tensor
Tensor of scalars.
Definition: tensor.H:51
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
void deleteDemandDrivenData(DataType *&dataPtr)
errorManip< error > abort(error &err)
Definition: errorManip.H:131
bool exists(const fileName &, const bool checkVariants=true, const bool followLink=true)
Does the name exist (as directory or file) in the file system?
Definition: POSIX.C:520
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:42
vector point
Point is a vector.
Definition: point.H:41
Field< triad > triadField
Specialisation of Field<T> for triad.
Definition: triadField.H:49
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
List< triFace > triFaceList
list of triFaces
Definition: triFaceList.H:42
Vector< scalar > vector
A scalar version of the templated Vector.
Definition: vector.H:49
void sort(UList< T > &)
Definition: UList.C:115
List< labelList > labelListList
A List of labelList.
Definition: labelList.H:57
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedVector eigenValues(const dimensionedTensor &dt)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
Field< vector > vectorField
Specialisation of Field<T> for vector.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
void writeVTK(OFstream &os, const Type &value)
prefixOSstream Pout(cout, "Pout")
Definition: IOstreams.H:53
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
Ostream & operator<<(Ostream &os, const fvConstraints &constraints)
error FatalError
void LUsolve(scalarSquareMatrix &matrix, List< Type > &source)
Solve the matrix using LU decomposition with pivoting.
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > magSqr(const dimensioned< Type > &)
List< surfacePatch > surfacePatchList
messageStream Warning
face triFace(3)
labelList f(nPoints)
Info<< "Finished reading KIVA file"<< endl;cellShapeList cellShapes(nPoints);labelList cellZoning(nPoints, -1);const cellModel &hex=*(cellModeller::lookup("hex"));labelList hexLabels(8);label activeCells=0;labelList pointMap(nPoints);forAll(pointMap, i){ pointMap[i]=i;}for(label i=0;i< nPoints;i++){ if(f[i] > 0.0) { hexLabels[0]=i;hexLabels[1]=i1tab[i];hexLabels[2]=i3tab[i1tab[i]];hexLabels[3]=i3tab[i];hexLabels[4]=i8tab[i];hexLabels[5]=i1tab[i8tab[i]];hexLabels[6]=i3tab[i1tab[i8tab[i]]];hexLabels[7]=i3tab[i8tab[i]];cellShapes[activeCells]=cellShape(hex, hexLabels);edgeList edges=cellShapes[activeCells].edges();forAll(edges, ei) { if(edges[ei].mag(points)< small) { label start=pointMap[edges[ei].start()];while(start !=pointMap[start]) { start=pointMap[start];} label end=pointMap[edges[ei].end()];while(end !=pointMap[end]) { end=pointMap[end];} label minLabel=min(start, end);pointMap[start]=pointMap[end]=minLabel;} } cellZoning[activeCells]=idreg[i];activeCells++;}}cellShapes.setSize(activeCells);cellZoning.setSize(activeCells);forAll(cellShapes, celli){ cellShape &cs=cellShapes[celli];forAll(cs, i) { cs[i]=pointMap[cs[i]];} cs.collapse();}label bcIDs[11]={-1, 0, 2, 4, -1, 5, -1, 6, 7, 8, 9};const label nBCs=12;const word *kivaPatchTypes[nBCs]={ &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &wallPolyPatch::typeName, &symmetryPolyPatch::typeName, &wedgePolyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &polyPatch::typeName, &symmetryPolyPatch::typeName, &mergedCyclicPolyPatch::typeName};enum patchTypeNames{ PISTON, VALVE, LINER, CYLINDERHEAD, AXIS, WEDGE, INFLOW, OUTFLOW, PRESIN, PRESOUT, SYMMETRYPLANE, CYCLIC};const char *kivaPatchNames[nBCs]={ "piston", "valve", "liner", "cylinderHead", "axis", "wedge", "inflow", "outflow", "presin", "presout", "symmetryPlane", "cyclic"};List< SLList< face > > pFaces[nBCs]
Definition: readKivaGrid.H:229
volScalarField & p
3D tensor transformation operations.