conformationSurfaces.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "conformationSurfaces.H"
30 #include "conformalVoronoiMesh.H"
31 #include "triSurface.H"
33 
34 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35 
36 namespace Foam
37 {
38  defineTypeNameAndDebug(conformationSurfaces, 0);
39 }
40 
41 
42 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
43 
44 void Foam::conformationSurfaces::hasBoundedVolume
45 (
46  List<volumeType>& referenceVolumeTypes
47 ) const
48 {
49  vector sum(Zero);
50  label totalTriangles = 0;
51 
52  forAll(surfaces_, s)
53  {
54  const searchableSurface& surface(allGeometry_[surfaces_[s]]);
55 
56  if
57  (
58  surface.hasVolumeType()
59  && (
60  normalVolumeTypes_[regionOffset_[s]]
62  )
63  )
64  {
65  pointField pts(1, locationInMesh_);
66 
67  List<volumeType> vTypes
68  (
69  pts.size(),
71  );
72 
73  surface.getVolumeType(pts, vTypes);
74 
75  referenceVolumeTypes[s] = vTypes[0];
76 
77  Info<< " is " << referenceVolumeTypes[s].str()
78  << " surface " << surface.name()
79  << endl;
80  }
81 
82  if (isA<triSurface>(surface))
83  {
84  const triSurface& triSurf = refCast<const triSurface>(surface);
85 
86  const pointField& surfPts = triSurf.points();
87 
88  Info<< " Checking " << surface.name() << endl;
89 
90  label nBaffles = 0;
91 
92  Info<< " Index = " << surfaces_[s] << endl;
93  Info<< " Offset = " << regionOffset_[s] << endl;
94 
95  for (const labelledTri& f : triSurf)
96  {
97  const label patchID = f.region() + regionOffset_[s];
98 
99  // Don't include baffle surfaces in the calculation
100  if
101  (
102  normalVolumeTypes_[patchID]
104  )
105  {
106  sum += f.areaNormal(surfPts);
107  }
108  else
109  {
110  ++nBaffles;
111  }
112  }
113  Info<< " has " << nBaffles << " baffles out of "
114  << triSurf.size() << " triangles" << nl;
115 
116  totalTriangles += triSurf.size();
117  }
118  }
119 
120  Info<< " Sum of all the surface normals (if near zero, surface is"
121  << " probably closed):" << nl
122  << " Note: Does not include baffle surfaces in calculation" << nl
123  << " Sum = " << sum/(totalTriangles + SMALL) << nl
124  << " mag(Sum) = " << mag(sum)/(totalTriangles + SMALL)
125  << endl;
126 }
127 
128 
129 void Foam::conformationSurfaces::readFeatures
130 (
131  const label surfI,
132  const dictionary& featureDict,
133  const word& surfaceName,
134  label& featureIndex
135 )
136 {
137  const word featureMethod =
138  featureDict.getOrDefault<word>("featureMethod", "none");
139 
140  if (featureMethod == "extendedFeatureEdgeMesh")
141  {
142  fileName feMeshName
143  (
144  featureDict.get<fileName>("extendedFeatureEdgeMesh")
145  );
146 
147  Info<< " features: " << feMeshName << endl;
148 
149  features_.set
150  (
151  featureIndex,
152  new extendedFeatureEdgeMesh
153  (
154  IOobject
155  (
156  feMeshName,
157  runTime_.time().constant(),
158  "extendedFeatureEdgeMesh",
159  runTime_.time(),
162  )
163  )
164  );
165 
166  featureIndex++;
167  }
168  else if (featureMethod == "extractFeatures")
169  {
170  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
171 
172  Info<< " features: " << surface.name()
173  << " of type " << surface.type()
174  << ", id: " << featureIndex << endl;
175 
176  autoPtr<searchableSurfaceFeatures> ssFeatures
177  (
179  );
180 
181  if (ssFeatures().hasFeatures())
182  {
183  features_.set
184  (
185  featureIndex,
186  ssFeatures().features()
187  );
188 
189  featureIndex++;
190  }
191  else
192  {
194  << surface.name() << " of type "
195  << surface.type() << " does not have features"
196  << endl;
197  }
198  }
199  else if (featureMethod == "none")
200  {
201  // Currently nothing to do
202  }
203  else
204  {
206  << "No valid featureMethod found for surface " << surfaceName
207  << nl << "Use \"extendedFeatureEdgeMesh\" "
208  << "or \"extractFeatures\"."
209  << exit(FatalError);
210  }
211 }
212 
213 void Foam::conformationSurfaces::readFeatures
214 (
215  const dictionary& featureDict,
216  const word& surfaceName,
217  label& featureIndex
218 )
219 {
220  const word featureMethod =
221  featureDict.getOrDefault<word>("featureMethod", "none");
222 
223  if (featureMethod == "extendedFeatureEdgeMesh")
224  {
225  fileName feMeshName
226  (
227  featureDict.get<fileName>("extendedFeatureEdgeMesh")
228  );
229 
230  Info<< " features: " << feMeshName << ", id: " << featureIndex
231  << endl;
232 
233  features_.set
234  (
235  featureIndex,
236  new extendedFeatureEdgeMesh
237  (
238  IOobject
239  (
240  feMeshName,
241  runTime_.time().constant(),
242  "extendedFeatureEdgeMesh",
243  runTime_.time(),
246  )
247  )
248  );
249 
250  featureIndex++;
251  }
252  else if (featureMethod == "none")
253  {
254  // Currently nothing to do
255  }
256  else
257  {
259  << "No valid featureMethod found for surface " << surfaceName
260  << nl << "Use \"extendedFeatureEdgeMesh\" "
261  << "or \"extractFeatures\"."
262  << exit(FatalError);
263  }
264 }
265 
266 
267 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
268 
269 Foam::conformationSurfaces::conformationSurfaces
270 (
271  const Time& runTime,
272  Random& rndGen,
273  const searchableSurfaces& allGeometry,
274  const dictionary& surfaceConformationDict
275 )
276 :
277  runTime_(runTime),
278  allGeometry_(allGeometry),
279  features_(),
280  locationInMesh_(surfaceConformationDict.get<point>("locationInMesh")),
281  surfaces_(),
282  allGeometryToSurfaces_(),
283  normalVolumeTypes_(),
284  patchNames_(),
285  surfZones_(),
286  regionOffset_(),
287  patchInfo_(),
288  globalBounds_(),
289  referenceVolumeTypes_()
290 {
291  const dictionary& surfacesDict
292  (
293  surfaceConformationDict.subDict("geometryToConformTo")
294  );
295 
296  const dictionary& additionalFeaturesDict
297  (
298  surfaceConformationDict.subDict("additionalFeatures")
299  );
300 
301 
302  // Wildcard specification : loop over all surface, all regions
303  // and try to find a match.
304 
305  // Count number of surfaces.
306  label surfI = 0;
307  for (const word& geomName : allGeometry_.names())
308  {
309  if (surfacesDict.found(geomName))
310  {
311  ++surfI;
312  }
313  }
314 
315  const label nAddFeat = additionalFeaturesDict.size();
316 
317  Info<< nl << "Reading geometryToConformTo" << endl;
318 
319  allGeometryToSurfaces_.setSize(allGeometry_.size(), -1);
320 
321  normalVolumeTypes_.setSize(surfI);
322  surfaces_.setSize(surfI);
323  surfZones_.setSize(surfI);
324 
325  // Features may be attached to host surfaces or independent
326  features_.setSize(surfI + nAddFeat);
327 
328  label featureI = 0;
329 
330  regionOffset_.setSize(surfI, 0);
331 
332  PtrList<dictionary> globalPatchInfo(surfI);
333  List<Map<autoPtr<dictionary>>> regionPatchInfo(surfI);
334  List<sideVolumeType> globalVolumeTypes(surfI);
335  List<Map<sideVolumeType>> regionVolumeTypes(surfI);
336 
337  wordHashSet unmatchedKeys(surfacesDict.toc());
338 
339  surfI = 0;
340  forAll(allGeometry_.names(), geomI)
341  {
342  const word& geomName = allGeometry_.names()[geomI];
343 
344  const entry* ePtr = surfacesDict.findEntry(geomName, keyType::REGEX);
345 
346  if (ePtr)
347  {
348  const dictionary& dict = ePtr->dict();
349  unmatchedKeys.erase(ePtr->keyword());
350 
351  surfaces_[surfI] = geomI;
352 
353  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
354 
355  // Surface zones
356  if (dict.found("faceZone"))
357  {
358  surfZones_.set
359  (
360  surfI,
361  new surfaceZonesInfo
362  (
363  surface,
364  dict,
365  allGeometry_.regionNames()[surfaces_[surfI]]
366  )
367  );
368  }
369 
370  allGeometryToSurfaces_[surfaces_[surfI]] = surfI;
371 
372  Info<< nl << " " << geomName << endl;
373 
374  const wordList& regionNames =
375  allGeometry_.regionNames()[surfaces_[surfI]];
376 
377  patchNames_.append(regionNames);
378 
379  globalVolumeTypes[surfI] =
380  (
382  [
383  dict.getOrDefault<word>
384  (
385  "meshableSide",
386  "inside"
387  )
388  ]
389  );
390 
391  if (!globalVolumeTypes[surfI])
392  {
393  if (!surface.hasVolumeType())
394  {
396  << "Non-baffle surface "
397  << surface.name()
398  << " does not allow inside/outside queries."
399  << " This usually is an error." << endl;
400  }
401  }
402 
403  // Load patch info
404  if (dict.found("patchInfo"))
405  {
406  globalPatchInfo.set
407  (
408  surfI,
409  dict.subDict("patchInfo").clone()
410  );
411  }
412 
413  readFeatures
414  (
415  surfI,
416  dict,
417  geomName,
418  featureI
419  );
420 
421  const wordList& rNames = surface.regions();
422 
423  if (dict.found("regions"))
424  {
425  const dictionary& regionsDict = dict.subDict("regions");
426 
427  forAll(rNames, regionI)
428  {
429  const word& regionName = rNames[regionI];
430 
431  if (regionsDict.found(regionName))
432  {
433  Info<< " region " << regionName << endl;
434 
435  // Get the dictionary for region
436  const dictionary& regionDict = regionsDict.subDict
437  (
438  regionName
439  );
440 
441  if (regionDict.found("patchInfo"))
442  {
443  regionPatchInfo[surfI].insert
444  (
445  regionI,
446  regionDict.subDict("patchInfo").clone()
447  );
448  }
449 
450  regionVolumeTypes[surfI].insert
451  (
452  regionI,
454  [
455  regionDict.getOrDefault<word>
456  (
457  "meshableSide",
458  extendedFeatureEdgeMesh::
459  sideVolumeTypeNames_
460  [
461  globalVolumeTypes[surfI]
462  ]
463  )
464  ]
465  );
466 
467  readFeatures(regionDict, regionName, featureI);
468  }
469  }
470  }
471 
472  surfI++;
473  }
474  }
475 
476 
477  if (unmatchedKeys.size() > 0)
478  {
479  IOWarningInFunction(surfacesDict)
480  << "Not all entries in conformationSurfaces dictionary were used."
481  << " The following entries were not used : "
482  << unmatchedKeys.sortedToc()
483  << endl;
484  }
485 
486 
487  // Calculate local to global region offset
488  label nRegions = 0;
489 
490  forAll(surfaces_, surfI)
491  {
492  regionOffset_[surfI] = nRegions;
493 
494  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
495  nRegions += surface.regions().size();
496  }
497 
498  // Rework surface specific information into information per global region
499  patchInfo_.setSize(nRegions);
500  normalVolumeTypes_.setSize(nRegions);
501 
502  forAll(surfaces_, surfI)
503  {
504  const searchableSurface& surface = allGeometry_[surfaces_[surfI]];
505 
506  label nRegions = surface.regions().size();
507 
508  // Initialise to global (i.e. per surface)
509  for (label i = 0; i < nRegions; i++)
510  {
511  label globalRegionI = regionOffset_[surfI] + i;
512  normalVolumeTypes_[globalRegionI] = globalVolumeTypes[surfI];
513  if (globalPatchInfo.set(surfI))
514  {
515  patchInfo_.set
516  (
517  globalRegionI,
518  globalPatchInfo[surfI].clone()
519  );
520  }
521  }
522 
523  forAllConstIters(regionVolumeTypes[surfI], iter)
524  {
525  label globalRegionI = regionOffset_[surfI] + iter.key();
526 
527  normalVolumeTypes_[globalRegionI] =
528  regionVolumeTypes[surfI][iter.key()];
529  }
530 
531  const Map<autoPtr<dictionary>>& localInfo = regionPatchInfo[surfI];
532  forAllConstIters(localInfo, iter)
533  {
534  label globalRegionI = regionOffset_[surfI] + iter.key();
535 
536  patchInfo_.set(globalRegionI, iter()().clone());
537  }
538  }
539 
540 
541 
542  if (!additionalFeaturesDict.empty())
543  {
544  Info<< nl << "Reading additionalFeatures" << endl;
545  }
546 
547  for (const entry& dEntry : additionalFeaturesDict)
548  {
549  const word& featureName = dEntry.keyword();
550  const dictionary& featureSubDict = dEntry.dict();
551 
552  Info<< nl << " " << featureName << endl;
553 
554  readFeatures(featureSubDict, featureName, featureI);
555  }
556 
557  // Remove unnecessary space from the features list
558  features_.setSize(featureI);
559 
560  globalBounds_ = treeBoundBox
561  (
562  searchableSurfacesQueries::bounds(allGeometry_, surfaces_)
563  );
564 
565  // Extend the global bounds to stop the bound box sitting on the surfaces
566  // to be conformed to
567  //globalBounds_.inflate(rndGen, 1e-4);
568 
569  globalBounds_.grow(1e-4*globalBounds_.span());
570 
571  // Look at all surfaces at determine whether the locationInMesh point is
572  // inside or outside each, to establish a signature for the domain to be
573  // meshed.
574 
575  referenceVolumeTypes_.setSize
576  (
577  surfaces_.size(),
579  );
580 
581  Info<< endl
582  << "Testing for locationInMesh " << locationInMesh_ << endl;
583 
584  hasBoundedVolume(referenceVolumeTypes_);
585 
586  if (debug)
587  {
588  Info<< "Names = " << allGeometry_.names() << endl;
589  Info<< "Surfaces = " << surfaces_ << endl;
590  Info<< "AllGeom to Surfaces = " << allGeometryToSurfaces_ << endl;
591  Info<< "Volume types = " << normalVolumeTypes_ << endl;
592  Info<< "Patch names = " << patchNames_ << endl;
593  Info<< "Region Offset = " << regionOffset_ << endl;
594 
595  forAll(features_, fI)
596  {
597  Info<< features_[fI].name() << endl;
598  }
599  }
600 }
601 
602 
603 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
604 
605 bool Foam::conformationSurfaces::overlaps(const boundBox& bb) const
606 {
607  forAll(surfaces_, s)
608  {
609  if (allGeometry_[surfaces_[s]].overlaps(bb))
610  {
611  return true;
612  }
613  }
614 
615  return false;
616 }
617 
618 
620 (
621  const pointField& samplePts
622 ) const
623 {
624  return wellInside(samplePts, scalarField(samplePts.size(), Zero));
625 }
626 
627 
629 (
630  const point& samplePt
631 ) const
632 {
633  return wellInside(pointField(1, samplePt), scalarField(1, Zero))[0];
634 }
635 
636 
638 (
639  const pointField& samplePts
640 ) const
641 {
642  return wellOutside(samplePts, scalarField(samplePts.size(), Zero));
643 }
644 
645 
647 (
648  const point& samplePt
649 ) const
650 {
651  return wellOutside(pointField(1, samplePt), scalarField(1, Zero))[0];
652  //return !inside(samplePt);
653 }
654 
655 
657 (
658  const pointField& samplePts,
659  const scalarField& testDistSqr,
660  const bool testForInside
661 ) const
662 {
663  List<List<volumeType>> surfaceVolumeTests
664  (
665  surfaces_.size(),
666  List<volumeType>
667  (
668  samplePts.size(),
670  )
671  );
672 
673  // Get lists for the volumeTypes for each sample wrt each surface
674  forAll(surfaces_, s)
675  {
676  const searchableSurface& surface(allGeometry_[surfaces_[s]]);
677 
678  const label regionI = regionOffset_[s];
679 
680  if (normalVolumeTypes_[regionI] != extendedFeatureEdgeMesh::BOTH)
681  {
682  surface.getVolumeType(samplePts, surfaceVolumeTests[s]);
683  }
684  }
685 
686  // Compare the volumeType result for each point wrt to each surface with the
687  // reference value and if the points are inside the surface by a given
688  // distanceSquared
689 
690  // Assume that the point is wellInside until demonstrated otherwise.
691  Field<bool> insideOutsidePoint(samplePts.size(), testForInside);
692 
693  //Check if the points are inside the surface by the given distance squared
694 
695  labelList hitSurfaces;
696  List<pointIndexHit> hitInfo;
698  (
699  allGeometry_,
700  surfaces_,
701  samplePts,
702  testDistSqr,
703  hitSurfaces,
704  hitInfo
705  );
706 
707  forAll(samplePts, i)
708  {
709  const pointIndexHit& pHit = hitInfo[i];
710 
711  if (pHit.hit())
712  {
713  // If the point is within range of the surface, then it can't be
714  // well (in|out)side
715  insideOutsidePoint[i] = false;
716 
717  continue;
718  }
719 
720  forAll(surfaces_, s)
721  {
722  const label regionI = regionOffset_[s];
723 
724  if (normalVolumeTypes_[regionI] == extendedFeatureEdgeMesh::BOTH)
725  {
726  continue;
727  }
728 
729  const searchableSurface& surface(allGeometry_[surfaces_[s]]);
730 
731  if
732  (
733  !surface.hasVolumeType()
734  //&& surfaceVolumeTests[s][i] == volumeType::UNKNOWN
735  )
736  {
737  pointField sample(1, samplePts[i]);
738  scalarField nearestDistSqr(1, GREAT);
739  List<pointIndexHit> info;
740 
741  surface.findNearest(sample, nearestDistSqr, info);
742 
743  const vector hitDir =
744  normalised
745  (
746  info[0].point() - samplePts[i]
747  );
748 
749  pointIndexHit surfHit;
750  label hitSurface;
751 
752  findSurfaceNearestIntersection
753  (
754  samplePts[i],
755  info[0].point() - 1e-3*mag(hitDir)*hitDir,
756  surfHit,
757  hitSurface
758  );
759 
760  if (surfHit.hit() && hitSurface != surfaces_[s])
761  {
762  continue;
763  }
764  }
765 
766  if (surfaceVolumeTests[s][i] == volumeType::OUTSIDE)
767  {
768  if
769  (
770  normalVolumeTypes_[regionI]
772  )
773  {
774  insideOutsidePoint[i] = !testForInside;
775  break;
776  }
777  }
778  else if (surfaceVolumeTests[s][i] == volumeType::INSIDE)
779  {
780  if
781  (
782  normalVolumeTypes_[regionI]
784  )
785  {
786  insideOutsidePoint[i] = !testForInside;
787  break;
788  }
789  }
790  }
791  }
792 
793  return insideOutsidePoint;
794 }
795 
796 
798 (
799  const pointField& samplePts,
800  const scalarField& testDistSqr
801 ) const
802 {
803  return wellInOutSide(samplePts, testDistSqr, true);
804 }
805 
806 
808 (
809  const point& samplePt,
810  scalar testDistSqr
811 ) const
812 {
813  return wellInside(pointField(1, samplePt), scalarField(1, testDistSqr))[0];
814 }
815 
816 
818 (
819  const pointField& samplePts,
820  const scalarField& testDistSqr
821 ) const
822 {
823  return wellInOutSide(samplePts, testDistSqr, false);
824 }
825 
826 
828 (
829  const point& samplePt,
830  scalar testDistSqr
831 ) const
832 {
833  return wellOutside(pointField(1, samplePt), scalarField(1, testDistSqr))[0];
834 }
835 
836 
838 (
839  const point& start,
840  const point& end
841 ) const
842 {
843  labelList hitSurfaces;
844  List<pointIndexHit> hitInfo;
845 
847  (
848  allGeometry_,
849  surfaces_,
850  pointField(1, start),
851  pointField(1, end),
852  hitSurfaces,
853  hitInfo
854  );
855 
856  return hitInfo[0].hit();
857 }
858 
859 
861 (
862  const point& start,
863  const point& end,
864  pointIndexHit& surfHit,
865  label& hitSurface
866 ) const
867 {
868  labelList hitSurfaces;
869  List<pointIndexHit> hitInfo;
870 
872  (
873  allGeometry_,
874  surfaces_,
875  pointField(1, start),
876  pointField(1, end),
877  hitSurfaces,
878  hitInfo
879  );
880 
881  surfHit = hitInfo[0];
882 
883  if (surfHit.hit())
884  {
885  // hitSurfaces has returned the index of the entry in surfaces_ that was
886  // found, not the index of the surface in allGeometry_, translating this
887  // to allGeometry_
888 
889  hitSurface = surfaces_[hitSurfaces[0]];
890  }
891 }
892 
893 
895 (
896  const point& start,
897  const point& end,
898  List<pointIndexHit>& surfHit,
899  labelList& hitSurface
900 ) const
901 {
902  labelListList hitSurfaces;
903  List<List<pointIndexHit>> hitInfo;
904 
906  (
907  allGeometry_,
908  surfaces_,
909  pointField(1, start),
910  pointField(1, end),
911  hitSurfaces,
912  hitInfo
913  );
914 
915  surfHit = hitInfo[0];
916 
917  hitSurface.setSize(hitSurfaces[0].size());
918 
919  forAll(hitSurfaces[0], surfI)
920  {
921  // hitSurfaces has returned the index of the entry in surfaces_ that was
922  // found, not the index of the surface in allGeometry_, translating this
923  // to allGeometry_
924 
925  hitSurface[surfI] = surfaces_[hitSurfaces[0][surfI]];
926  }
927 }
928 
929 
931 (
932  const point& start,
933  const point& end,
934  pointIndexHit& surfHit,
935  label& hitSurface
936 ) const
937 {
938  labelList hitSurfacesStart;
939  List<pointIndexHit> hitInfoStart;
940  labelList hitSurfacesEnd;
941  List<pointIndexHit> hitInfoEnd;
942 
944  (
945  allGeometry_,
946  surfaces_,
947  pointField(1, start),
948  pointField(1, end),
949  hitSurfacesStart,
950  hitInfoStart,
951  hitSurfacesEnd,
952  hitInfoEnd
953  );
954 
955  surfHit = hitInfoStart[0];
956 
957  if (surfHit.hit())
958  {
959  // hitSurfaces has returned the index of the entry in surfaces_ that was
960  // found, not the index of the surface in allGeometry_, translating this
961  // to allGeometry_
962 
963  hitSurface = surfaces_[hitSurfacesStart[0]];
964  }
965 }
966 
967 
969 (
970  const point& sample,
971  scalar nearestDistSqr,
972  pointIndexHit& surfHit,
973  label& hitSurface
974 ) const
975 {
976  labelList hitSurfaces;
977  List<pointIndexHit> surfaceHits;
978 
980  (
981  allGeometry_,
982  surfaces_,
983  pointField(1, sample),
984  scalarField(1, nearestDistSqr),
985  hitSurfaces,
986  surfaceHits
987  );
988 
989  surfHit = surfaceHits[0];
990 
991  if (surfHit.hit())
992  {
993  // hitSurfaces has returned the index of the entry in surfaces_ that was
994  // found, not the index of the surface in allGeometry_, translating this
995  // to allGeometry_
996 
997  hitSurface = surfaces_[hitSurfaces[0]];
998  }
999 }
1000 
1001 
1003 (
1004  const pointField& samples,
1005  const scalarField& nearestDistSqr,
1006  List<pointIndexHit>& surfaceHits,
1007  labelList& hitSurfaces
1008 ) const
1009 {
1011  (
1012  allGeometry_,
1013  surfaces_,
1014  samples,
1015  nearestDistSqr,
1016  hitSurfaces,
1017  surfaceHits
1018  );
1019 
1020  forAll(surfaceHits, i)
1021  {
1022  if (surfaceHits[i].hit())
1023  {
1024  // hitSurfaces has returned the index of the entry in surfaces_ that
1025  // was found, not the index of the surface in allGeometry_,
1026  // translating this to the surface in allGeometry_.
1027 
1028  hitSurfaces[i] = surfaces_[hitSurfaces[i]];
1029  }
1030  }
1031 }
1032 
1033 
1035 (
1036  const point& sample,
1037  scalar nearestDistSqr,
1038  pointIndexHit& fpHit,
1039  label& featureHit
1040 ) const
1041 {
1042  // Work arrays
1043  scalar minDistSqr = nearestDistSqr;
1044  pointIndexHit hitInfo;
1045 
1046  forAll(features_, testI)
1047  {
1048  features_[testI].nearestFeaturePoint
1049  (
1050  sample,
1051  minDistSqr,
1052  hitInfo
1053  );
1054 
1055  if (hitInfo.hit())
1056  {
1057  minDistSqr = hitInfo.point().distSqr(sample);
1058  fpHit = hitInfo;
1059  featureHit = testI;
1060  }
1061  }
1062 }
1063 
1064 
1066 (
1067  const point& sample,
1068  scalar nearestDistSqr,
1069  pointIndexHit& edgeHit,
1070  label& featureHit
1071 ) const
1072 {
1073  pointField samples(1, sample);
1074  scalarField nearestDistsSqr(1, nearestDistSqr);
1075 
1076  List<pointIndexHit> edgeHits;
1077  labelList featuresHit;
1078 
1079  findEdgeNearest
1080  (
1081  samples,
1082  nearestDistsSqr,
1083  edgeHits,
1084  featuresHit
1085  );
1086 
1087  edgeHit = edgeHits[0];
1088  featureHit = featuresHit[0];
1089 }
1090 
1091 
1093 (
1094  const pointField& samples,
1095  const scalarField& nearestDistsSqr,
1096  List<pointIndexHit>& edgeHits,
1097  labelList& featuresHit
1098 ) const
1099 {
1100  // Initialise
1101  featuresHit.setSize(samples.size());
1102  featuresHit = -1;
1103  edgeHits.setSize(samples.size());
1104 
1105  // Work arrays
1106  scalarField minDistSqr(nearestDistsSqr);
1107  List<pointIndexHit> hitInfo(samples.size());
1108 
1109  forAll(features_, testI)
1110  {
1111  features_[testI].nearestFeatureEdge
1112  (
1113  samples,
1114  minDistSqr,
1115  hitInfo
1116  );
1117 
1118  // Update minDistSqr and arguments
1119  forAll(hitInfo, pointi)
1120  {
1121  if (hitInfo[pointi].hit())
1122  {
1123  minDistSqr[pointi] =
1124  hitInfo[pointi].point().distSqr(samples[pointi]);
1125 
1126  edgeHits[pointi] = hitInfo[pointi];
1127  featuresHit[pointi] = testI;
1128  }
1129  }
1130  }
1131 }
1132 
1133 
1135 (
1136  const point& sample,
1137  scalar nearestDistSqr,
1138  List<pointIndexHit>& edgeHits,
1139  List<label>& featuresHit
1140 ) const
1141 {
1142  // Initialise
1143  featuresHit.setSize(extendedFeatureEdgeMesh::nEdgeTypes);
1144  featuresHit = -1;
1145  edgeHits.setSize(extendedFeatureEdgeMesh::nEdgeTypes);
1146 
1147  // Work arrays
1148  scalarField minDistSqr(extendedFeatureEdgeMesh::nEdgeTypes, nearestDistSqr);
1149  List<pointIndexHit> hitInfo(extendedFeatureEdgeMesh::nEdgeTypes);
1150 
1151  forAll(features_, testI)
1152  {
1153  features_[testI].nearestFeatureEdgeByType
1154  (
1155  sample,
1156  minDistSqr,
1157  hitInfo
1158  );
1159 
1160  // Update minDistSqr and arguments
1161  forAll(hitInfo, typeI)
1162  {
1163  if (hitInfo[typeI].hit())
1164  {
1165  minDistSqr[typeI] = hitInfo[typeI].point().distSqr(sample);
1166  edgeHits[typeI] = hitInfo[typeI];
1167  featuresHit[typeI] = testI;
1168  }
1169  }
1170  }
1171 }
1172 
1173 
1175 (
1176  const point& sample,
1177  const scalar searchRadiusSqr,
1178  List<List<pointIndexHit>>& edgeHitsByFeature,
1179  List<label>& featuresHit
1180 ) const
1181 {
1182  // Initialise
1183  //featuresHit.setSize(features_.size());
1184  //featuresHit = -1;
1185  //edgeHitsByFeature.setSize(features_.size());
1186 
1187  // Work arrays
1188  List<pointIndexHit> hitInfo(extendedFeatureEdgeMesh::nEdgeTypes);
1189 
1190  forAll(features_, testI)
1191  {
1192  features_[testI].allNearestFeatureEdges
1193  (
1194  sample,
1195  searchRadiusSqr,
1196  hitInfo
1197  );
1198 
1199  bool anyHit = false;
1200  forAll(hitInfo, hitI)
1201  {
1202  if (hitInfo[hitI].hit())
1203  {
1204  anyHit = true;
1205  break;
1206  }
1207  }
1208 
1209  if (anyHit)
1210  {
1211  edgeHitsByFeature.append(hitInfo);
1212  featuresHit.append(testI);
1213  }
1214  }
1215 }
1216 
1217 
1218 void Foam::conformationSurfaces::writeFeatureObj(const fileName& prefix) const
1219 {
1220  OFstream ftStr(runTime_.time().path()/prefix + "_allFeatures.obj");
1221 
1222  Pout<< nl << "Writing all features to " << ftStr.name() << endl;
1223 
1224  label verti = 0;
1225 
1226  forAll(features_, i)
1227  {
1228  const extendedFeatureEdgeMesh& fEM(features_[i]);
1229  const pointField pts(fEM.points());
1230  const edgeList eds(fEM.edges());
1231 
1232  ftStr << "g " << fEM.name() << endl;
1233 
1234  forAll(eds, j)
1235  {
1236  const edge& e = eds[j];
1237 
1238  meshTools::writeOBJ(ftStr, pts[e[0]]); verti++;
1239  meshTools::writeOBJ(ftStr, pts[e[1]]); verti++;
1240  ftStr << "l " << verti-1 << ' ' << verti << endl;
1241  }
1242  }
1243 }
1244 
1245 
1247 (
1248  const point& ptA,
1249  const point& ptB
1250 ) const
1251 {
1252  pointIndexHit surfHit;
1253  label hitSurface;
1254 
1255  findSurfaceAnyIntersection(ptA, ptB, surfHit, hitSurface);
1256 
1257  return getPatchID(hitSurface, surfHit);
1258 }
1259 
1260 
1261 Foam::label Foam::conformationSurfaces::findPatch(const point& pt) const
1262 {
1263  pointIndexHit surfHit;
1264  label hitSurface;
1265 
1266  findSurfaceNearest(pt, sqr(GREAT), surfHit, hitSurface);
1267 
1268  return getPatchID(hitSurface, surfHit);
1269 }
1270 
1271 
1273 (
1274  const label hitSurface,
1275  const pointIndexHit& surfHit
1276 ) const
1277 {
1278  if (!surfHit.hit())
1279  {
1280  return -1;
1281  }
1282 
1283  labelList surfLocalRegion;
1284 
1285  allGeometry_[hitSurface].getRegion
1286  (
1287  List<pointIndexHit>(1, surfHit),
1288  surfLocalRegion
1289  );
1290 
1291  const label patchID =
1292  surfLocalRegion[0]
1293  + regionOffset_[allGeometryToSurfaces_[hitSurface]];
1294 
1295  return patchID;
1296 }
1297 
1298 
1301 (
1302  const label hitSurface,
1303  const pointIndexHit& surfHit
1304 ) const
1305 {
1306  const label patchID = getPatchID(hitSurface, surfHit);
1307 
1308  if (patchID == -1)
1309  {
1311  }
1312 
1313  return normalVolumeTypes_[patchID];
1314 }
1315 
1316 
1318 (
1319  const label hitSurface,
1320  const List<pointIndexHit>& surfHit,
1321  vectorField& normal
1322 ) const
1323 {
1324  allGeometry_[hitSurface].getNormal(surfHit, normal);
1325 
1326  const label patchID = regionOffset_[allGeometryToSurfaces_[hitSurface]];
1327 
1328  // Now flip sign of normal depending on mesh side
1329  if (normalVolumeTypes_[patchID] == extendedFeatureEdgeMesh::OUTSIDE)
1330  {
1331  normal *= -1;
1332  }
1333 }
1334 
1335 
1336 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Field< bool > inside(const pointField &samplePts) const
Check if points are inside surfaces to conform to.
void findEdgeNearestByType(const point &sample, scalar nearestDistSqr, List< pointIndexHit > &edgeHit, List< label > &featuresHit) const
Find the nearest point on each type of feature edge.
dictionary dict
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
Field< bool > outside(const pointField &samplePts) const
Check if points are outside surfaces to conform to.
static const Enum< sideVolumeType > sideVolumeTypeNames_
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
Field< bool > wellInside(const pointField &samplePts, const scalarField &testDistSqr) const
Check if point is inside surfaces to conform to by at least.
void append(const T &val)
Append an element at the end of the list.
Definition: List.H:521
List< edge > edgeList
List of edge.
Definition: edgeList.H:32
static void findNearest(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &)
Find nearest. Return -1 (and a miss()) or surface and nearest.
static autoPtr< searchableSurfaceFeatures > New(const searchableSurface &surface, const dictionary &dict)
Return a reference to the selected searchableSurfaceFeatures.
std::enable_if< std::is_same< bool, TypeT >::value, bool >::type set(const label i, bool val=true)
A bitSet::set() method for a list of bool.
Definition: List.H:493
dimensionedSymmTensor sqr(const dimensionedVector &dv)
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
Random rndGen
Definition: createFields.H:23
scalarField samples(nIntervals, Zero)
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
Definition: pointIndexHit.H:58
void findFeaturePointNearest(const point &sample, scalar nearestDistSqr, pointIndexHit &fpHit, label &featureHit) const
Find the nearest point on any feature edge.
extendedFeatureEdgeMesh::sideVolumeType meshableSide(const label hitSurface, const pointIndexHit &surfHit) const
Is the surface a baffle.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
wordList regionNames
static void findAllIntersections(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelListList &surfaces, List< List< pointIndexHit >> &surfaceHits)
Find all intersections in order from start to end. Returns for.
Ignore writing from objectRegistry::writeObject()
quaternion normalised(const quaternion &q)
Return the normalised (unit) quaternion of the given quaternion.
Definition: quaternionI.H:674
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:441
autoPtr< dictionary > clone() const
Construct and return clone.
Definition: dictionary.C:165
Unknown state.
Definition: volumeType.H:64
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
static void findNearestIntersection(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &surface2, List< pointIndexHit > &hit2)
Find intersections of edge nearest to both endpoints.
bool findSurfaceAnyIntersection(const point &start, const point &end) const
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
void findAllNearestEdges(const point &sample, const scalar searchRadiusSqr, List< List< pointIndexHit >> &edgeHitsByFeature, List< label > &featuresHit) const
Find the nearest points on each feature edge that is within.
void getNormal(const label hitSurface, const List< pointIndexHit > &surfHit, vectorField &normal) const
sideVolumeType
Normals point to the outside.
void writeFeatureObj(const fileName &prefix) const
Write all components of all the extendedFeatureEdgeMeshes as.
void findSurfaceNearest(const point &sample, scalar nearestDistSqr, pointIndexHit &surfHit, label &hitSurface) const
Find the nearest point to the sample and return it to the.
void findEdgeNearest(const point &sample, scalar nearestDistSqr, pointIndexHit &edgeHit, label &featureHit) const
Find the nearest point on any feature edge.
Vector< scalar > vector
Definition: vector.H:57
HashSet< word, Hash< word > > wordHashSet
A HashSet of words, uses string hasher.
Definition: HashSet.H:73
A location inside the volume.
Definition: volumeType.H:65
Field< bool > wellInOutSide(const pointField &samplePts, const scalarField &testDistSqr, bool testForInside) const
Check if point is closer to the surfaces to conform to than.
static void findAnyIntersection(const PtrList< searchableSurface > &, const labelList &surfacesToTest, const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &)
Find any intersection. Return hit point information and.
A location outside the volume.
Definition: volumeType.H:66
const wordList surface
Standard surface field types (scalar, vector, tensor, etc)
not sure when this may be used
int debug
Static debugging option.
constexpr auto end(C &c) -> decltype(c.end())
Return iterator to the end of the container c.
Definition: stdFoam.H:201
defineTypeNameAndDebug(combustionModel, 0)
Field< bool > wellOutside(const pointField &samplePts, const scalarField &testDistSqr) const
Check if point is outside surfaces to conform to by at least.
labelList f(nPoints)
bool erase(T *item)
Remove the specified element from the list and delete it.
Definition: ILList.C:101
void findSurfaceAllIntersections(const point &start, const point &end, List< pointIndexHit > &surfHit, labelList &hitSurface) const
List< word > wordList
List of word.
Definition: fileName.H:59
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
messageStream Info
Information stream (stdout output on master, null elsewhere)
void findSurfaceNearestIntersection(const point &start, const point &end, pointIndexHit &surfHit, label &hitSurface) const
Finding the nearestIntersection of the surface to start.
Field< vector > vectorField
Specialisation of Field<T> for vector.
#define IOWarningInFunction(ios)
Report an IO warning using Foam::Warning.
T getOrDefault(const word &keyword, const T &deflt, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T, or return the given default value. FatalIOError if it is found and the number of...
List< label > labelList
A List of labels.
Definition: List.H:62
label getPatchID(const label hitSurface, const pointIndexHit &surfHit) const
Get the region number of a hit surface.
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
bool overlaps(const boundBox &bb) const
Check if the supplied bound box overlaps any part of any of.
Regular expression.
Definition: keyType.H:83
static constexpr label nEdgeTypes
Number of possible feature edge types (i.e. number of slices)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static boundBox bounds(const PtrList< searchableSurface > &allSurfaces, const labelUList &surfacesToTest)
Find the boundBox of the selected surfaces.
label findPatch(const point &ptA, const point &ptB) const
Find which patch is intersected by the line from one point to.
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127