searchableSurfaces.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-2023 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 "searchableSurfaces.H"
28 #include "ListOps.H"
29 #include "Time.H"
30 #include "DynamicField.H"
31 #include "PatchTools.H"
32 #include "triSurfaceMesh.H"
33 #include "vtkWritePolyData.H"
34 
35 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
40 }
41 
42 
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
44 
45 bool Foam::searchableSurfaces::connected
46 (
47  const triSurface& s,
48  const label edgeI,
49  const pointIndexHit& hit
50 )
51 {
52  const triFace& localFace = s.localFaces()[hit.index()];
53  const edge& e = s.edges()[edgeI];
54 
55  forAll(localFace, i)
56  {
57  if (e.otherVertex(localFace[i]) != -1)
58  {
59  return true;
60  }
61  }
62 
63  return false;
64 }
65 
66 
67 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
68 
70 :
72  regionNames_(size),
73  allSurfaces_(identityMap(size))
74 {}
75 
76 
78 (
79  const IOobject& io,
80  const dictionary& topDict,
81  const bool singleRegionName
82 )
83 :
84  PtrList<searchableSurface>(topDict.size()),
85  names_(topDict.size()),
86  regionNames_(topDict.size()),
87  allSurfaces_(identityMap(topDict.size()))
88 {
89  label surfI = 0;
90  forAllConstIter(dictionary, topDict, iter)
91  {
92  const word& key = iter().keyword();
93 
94  if (!topDict.isDict(key))
95  {
97  << "Found non-dictionary entry " << iter()
98  << " in top-level dictionary " << topDict
99  << exit(FatalError);
100  }
101 
102  const dictionary& dict = topDict.subDict(key);
103 
104  names_[surfI] = key;
105  dict.readIfPresent("name", names_[surfI]);
106 
107  // Make IOobject with correct name
108  autoPtr<IOobject> namedIO(io.clone());
109  // Note: we would like to e.g. register triSurface 'sphere.stl' as
110  // 'sphere'. Unfortunately
111  // no support for having object read from different location than
112  // their object name. Maybe have stlTriSurfaceMesh which appends .stl
113  // when reading/writing?
114  namedIO().rename(key); // names_[surfI]
115 
116  // Create and hook surface
117  set
118  (
119  surfI,
121  (
122  dict.lookup("type"),
123  namedIO(),
124  dict
125  )
126  );
127  const searchableSurface& s = operator[](surfI);
128 
129  // Construct default region names by prepending surface name
130  // to region name.
131  const wordList& localNames = s.regions();
132 
133  wordList& rNames = regionNames_[surfI];
134  rNames.setSize(localNames.size());
135 
136  if (singleRegionName && localNames.size() == 1)
137  {
138  rNames[0] = names_[surfI];
139  }
140  else
141  {
142  forAll(localNames, regionI)
143  {
144  rNames[regionI] = names_[surfI] + '_' + localNames[regionI];
145  }
146  }
147 
148  // See if dictionary provides any global region names.
149  if (dict.found("regions"))
150  {
151  const dictionary& regionsDict = dict.subDict("regions");
152 
153  forAllConstIter(dictionary, regionsDict, iter)
154  {
155  const word& key = iter().keyword();
156 
157  if (regionsDict.isDict(key))
158  {
159  // Get the dictionary for region iter.keyword()
160  const dictionary& regionDict = regionsDict.subDict(key);
161 
162  label index = findIndex(localNames, key);
163 
164  if (index == -1)
165  {
167  << "Unknown region name " << key
168  << " for surface " << s.name() << endl
169  << "Valid region names are " << localNames
170  << exit(FatalError);
171  }
172 
173  rNames[index] = word(regionDict.lookup("name"));
174  }
175  }
176  }
177 
178  surfI++;
179  }
180 
181  // Trim (not really necessary since we don't allow non-dictionary entries)
183  names_.setSize(surfI);
184  regionNames_.setSize(surfI);
185  allSurfaces_.setSize(surfI);
186 }
187 
188 
189 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
190 
192 (
193  const word& wantedName
194 ) const
195 {
196  return findIndex(names_, wantedName);
197 }
198 
199 
201 (
202  const word& surfaceName,
203  const word& regionName
204 ) const
205 {
206  label surfaceIndex = findSurfaceID(surfaceName);
207 
208  return findIndex(this->operator[](surfaceIndex).regions(), regionName);
209 }
210 
211 
212 // Find any intersection
214 (
215  const pointField& start,
216  const pointField& end,
217  labelList& hitSurfaces,
218  List<pointIndexHit>& hitInfo
219 ) const
220 {
222  (
223  *this,
224  allSurfaces_,
225  start,
226  end,
227  hitSurfaces,
228  hitInfo
229  );
230 }
231 
232 
234 (
235  const pointField& start,
236  const pointField& end,
237  labelListList& hitSurfaces,
238  List<List<pointIndexHit>>& hitInfo
239 ) const
240 {
242  (
243  *this,
244  allSurfaces_,
245  start,
246  end,
247  hitSurfaces,
248  hitInfo
249  );
250 }
251 
252 
253 //Find intersections of edge nearest to both endpoints.
255 (
256  const pointField& start,
257  const pointField& end,
258  labelList& surface1,
259  List<pointIndexHit>& hit1,
260  labelList& surface2,
261  List<pointIndexHit>& hit2
262 ) const
263 {
265  (
266  *this,
267  allSurfaces_,
268  start,
269  end,
270  surface1,
271  hit1,
272  surface2,
273  hit2
274  );
275 }
276 
277 
279 (
280  const pointField& samples,
281  const scalarField& nearestDistSqr,
282  labelList& nearestSurfaces,
283  List<pointIndexHit>& nearestInfo
284 ) const
285 {
287  (
288  *this,
289  allSurfaces_,
290  samples,
291  nearestDistSqr,
292  nearestSurfaces,
293  nearestInfo
294  );
295 }
296 
297 
299 (
300  const pointField& samples,
301  const scalarField& nearestDistSqr,
302  const labelList& regionIndices,
303  labelList& nearestSurfaces,
304  List<pointIndexHit>& nearestInfo
305 ) const
306 {
308  (
309  *this,
310  allSurfaces_,
311  samples,
312  nearestDistSqr,
313  regionIndices,
314  nearestSurfaces,
315  nearestInfo
316  );
317 }
318 
319 
321 {
323  (
324  *this,
325  allSurfaces_
326  );
327 }
328 
329 
330 bool Foam::searchableSurfaces::checkClosed(const bool report) const
331 {
332  if (report)
333  {
334  Info<< "Checking for closedness." << endl;
335  }
336 
337  bool hasError = false;
338 
339  forAll(*this, surfI)
340  {
341  if (!operator[](surfI).hasVolumeType())
342  {
343  hasError = true;
344 
345  if (report)
346  {
347  Info<< " " << names()[surfI]
348  << " : not closed" << endl;
349  }
350 
351  if (isA<triSurface>(operator[](surfI)))
352  {
353  const triSurface& s = dynamic_cast<const triSurface&>
354  (
355  operator[](surfI)
356  );
357  const labelListList& edgeFaces = s.edgeFaces();
358 
359  label nSingleEdges = 0;
360  forAll(edgeFaces, edgeI)
361  {
362  if (edgeFaces[edgeI].size() == 1)
363  {
364  nSingleEdges++;
365  }
366  }
367 
368  label nMultEdges = 0;
369  forAll(edgeFaces, edgeI)
370  {
371  if (edgeFaces[edgeI].size() > 2)
372  {
373  nMultEdges++;
374  }
375  }
376 
377  if (report && (nSingleEdges != 0 || nMultEdges != 0))
378  {
379  Info<< " connected to one face : "
380  << nSingleEdges << nl
381  << " connected to >2 faces : "
382  << nMultEdges << endl;
383  }
384  }
385  }
386  }
387 
388  if (report)
389  {
390  Info<< endl;
391  }
392 
393  return returnReduce(hasError, orOp<bool>());
394 }
395 
396 
398 {
399  if (report)
400  {
401  Info<< "Checking for normal orientation." << endl;
402  }
403 
404  bool hasError = false;
405 
406  forAll(*this, surfI)
407  {
408  if (isA<triSurface>(operator[](surfI)))
409  {
410  const triSurface& s = dynamic_cast<const triSurface&>
411  (
412  operator[](surfI)
413  );
414 
415  labelHashSet borderEdge(s.size()/1000);
416  PatchTools::checkOrientation(s, false, &borderEdge);
417  // Colour all faces into zones using borderEdge
418  labelList normalZone;
419  label nZones = PatchTools::markZones(s, borderEdge, normalZone);
420  if (nZones > 1)
421  {
422  hasError = true;
423 
424  if (report)
425  {
426  Info<< " " << names()[surfI]
427  << " : has multiple orientation zones ("
428  << nZones << ")" << endl;
429  }
430  }
431  }
432  }
433 
434  if (report)
435  {
436  Info<< endl;
437  }
438 
439  return returnReduce(hasError, orOp<bool>());
440 }
441 
442 
444 (
445  const scalar maxRatio,
446  const bool report
447 ) const
448 {
449  if (report)
450  {
451  Info<< "Checking for size." << endl;
452  }
453 
454  bool hasError = false;
455 
456  forAll(*this, i)
457  {
458  const boundBox& bb = operator[](i).bounds();
459 
460  for (label j = i+1; j < size(); j++)
461  {
462  scalar ratio = bb.mag()/operator[](j).bounds().mag();
463 
464  if (ratio > maxRatio || ratio < 1.0/maxRatio)
465  {
466  hasError = true;
467 
468  if (report)
469  {
470  Info<< " " << names()[i]
471  << " bounds differ from " << names()[j]
472  << " by more than a factor 100:" << nl
473  << " bounding box : " << bb << nl
474  << " bounding box : " << operator[](j).bounds()
475  << endl;
476  }
477  break;
478  }
479  }
480  }
481 
482  if (report)
483  {
484  Info<< endl;
485  }
486 
487  return returnReduce(hasError, orOp<bool>());
488 }
489 
490 
492 (
493  const scalar tolerance,
494  const bool report
495 ) const
496 {
497  if (report)
498  {
499  Info<< "Checking for intersection." << endl;
500  }
501 
502  // cpuTime timer;
503 
504  bool hasError = false;
505 
506  forAll(*this, i)
507  {
508  if (isA<triSurfaceMesh>(operator[](i)))
509  {
510  const triSurfaceMesh& s0 = dynamic_cast<const triSurfaceMesh&>
511  (
512  operator[](i)
513  );
514  const edgeList& edges0 = s0.edges();
515  const pointField& localPoints0 = s0.localPoints();
516 
517  // Collect all the edge vectors
518  pointField start(edges0.size());
519  pointField end(edges0.size());
520  forAll(edges0, edgeI)
521  {
522  const edge& e = edges0[edgeI];
523  start[edgeI] = localPoints0[e[0]];
524  end[edgeI] = localPoints0[e[1]];
525  }
526 
527  // Test all other surfaces for intersection
528  forAll(*this, j)
529  {
530  List<pointIndexHit> hits;
531 
532  if (i == j)
533  {
534  // Slightly shorten the edges to prevent finding lots of
535  // intersections. The fast triangle intersection routine
536  // has problems with rays passing through a point of the
537  // triangle.
538  vectorField d
539  (
540  max(tolerance, 10*s0.tolerance())
541  *(end-start)
542  );
543  start += d;
544  end -= d;
545  }
546 
547  operator[](j).findLineAny(start, end, hits);
548 
549  label nHits = 0;
550  DynamicField<point> intersections(edges0.size()/100);
551  DynamicField<scalar> intersectionEdge(intersections.capacity());
552 
553  forAll(hits, edgeI)
554  {
555  if
556  (
557  hits[edgeI].hit()
558  && (i != j || !connected(s0, edgeI, hits[edgeI]))
559  )
560  {
561  intersections.append(hits[edgeI].hitPoint());
562  intersectionEdge.append(1.0*edgeI);
563  nHits++;
564  }
565  }
566 
567  if (nHits > 0)
568  {
569  if (report)
570  {
571  Info<< " " << names()[i]
572  << " intersects " << names()[j]
573  << " at " << nHits
574  << " locations."
575  << endl;
576 
577  const fileName fName
578  (
579  names()[i] + '_' + names()[j] + "_edgeIndex.vtk"
580  );
581 
582  Info<< " Writing intersection locations to "
583  << fName << endl;
584 
586  (
587  fName,
588  names()[i] + '_' + names()[j],
589  false,
590  intersections,
591  identityMap(intersections.size()),
592  edgeList(),
593  faceList(),
594  "edgeIndex",
595  true,
596  intersectionEdge
597  );
598  }
599 
600  hasError = true;
601  break;
602  }
603  }
604  }
605  }
606 
607  if (report)
608  {
609  Info<< endl;
610  }
611 
612  return returnReduce(hasError, orOp<bool>());
613 }
614 
615 
617 (
618  const scalar minQuality,
619  const bool report
620 ) const
621 {
622  if (report)
623  {
624  Info<< "Checking for triangle quality." << endl;
625  }
626 
627  bool hasError = false;
628 
629  forAll(*this, surfI)
630  {
631  if (isA<triSurface>(operator[](surfI)))
632  {
633  const triSurface& s = dynamic_cast<const triSurface&>
634  (
635  operator[](surfI)
636  );
637 
638  label nBadTris = 0;
639  forAll(s, facei)
640  {
641  const labelledTri& f = s[facei];
642 
643  scalar q = triPointRef
644  (
645  s.points()[f[0]],
646  s.points()[f[1]],
647  s.points()[f[2]]
648  ).quality();
649 
650  if (q < minQuality)
651  {
652  nBadTris++;
653  }
654  }
655 
656  if (nBadTris > 0)
657  {
658  hasError = true;
659 
660  if (report)
661  {
662  Info<< " " << names()[surfI]
663  << " : has " << nBadTris << " bad quality triangles "
664  << " (quality < " << minQuality << ")" << endl;
665  }
666  }
667  }
668  }
669 
670  if (report)
671  {
672  Info<< endl;
673  }
674 
675  return returnReduce(hasError, orOp<bool>());
676 
677 }
678 
679 
680 // Check all surfaces. Return number of failures.
682 (
683  const bool report
684 ) const
685 {
686  label noFailedChecks = 0;
687 
688  if (checkClosed(report))
689  {
690  noFailedChecks++;
691  }
692 
693  if (checkNormalOrientation(report))
694  {
695  noFailedChecks++;
696  }
697  return noFailedChecks;
698 }
699 
700 
702 (
703  const scalar maxRatio,
704  const scalar tol,
705  const scalar minQuality,
706  const bool report
707 ) const
708 {
709  label noFailedChecks = 0;
710 
711  if (maxRatio > 0 && checkSizes(maxRatio, report))
712  {
713  noFailedChecks++;
714  }
715 
716  if (checkIntersection(tol, report))
717  {
718  noFailedChecks++;
719  }
720 
721  if (checkQuality(minQuality, report))
722  {
723  noFailedChecks++;
724  }
725 
726  return noFailedChecks;
727 }
728 
729 
731 (
732  const List<wordList>& patchTypes,
733  Ostream& os
734 ) const
735 {
736  Info<< "Statistics." << endl;
737  forAll(*this, surfI)
738  {
739  Info<< " " << names()[surfI] << ':' << endl;
740 
741  const searchableSurface& s = operator[](surfI);
742 
743  Info<< " type : " << s.type() << nl
744  << " size : " << s.globalSize() << nl;
745  if (isA<triSurfaceMesh>(s))
746  {
747  const triSurfaceMesh& ts = dynamic_cast<const triSurfaceMesh&>(s);
748  Info<< " edges : " << ts.nEdges() << nl
749  << " points : " << ts.points()().size() << nl;
750  }
751  Info<< " bounds : " << s.bounds() << nl
752  << " closed : " << Switch(s.hasVolumeType()) << endl;
753 
754  if (patchTypes.size() && patchTypes[surfI].size() >= 1)
755  {
756  wordList unique(HashSet<word>(patchTypes[surfI]).sortedToc());
757  Info<< " patches : ";
758  forAll(unique, i)
759  {
760  Info<< unique[i];
761  if (i < unique.size()-1)
762  {
763  Info<< ',';
764  }
765  }
766  Info<< endl;
767  }
768  }
769  Info<< endl;
770 }
771 
772 
773 // * * * * * * * * * * * * * * Member Operators * * * * * * * * * * * * * * //
774 
776 (
777  const word& surfName
778 ) const
779 {
780  const label surfI = findSurfaceID(surfName);
781 
782  if (surfI < 0)
783  {
785  << "Surface named " << surfName << " not found." << nl
786  << "Available surface names: " << names_ << endl
787  << abort(FatalError);
788  }
789 
790  return operator[](surfI);
791 }
792 
793 
795 (
796  const word& surfName
797 )
798 {
799  const label surfI = findSurfaceID(surfName);
800 
801  if (surfI < 0)
802  {
804  << "Surface named " << surfName << " not found." << nl
805  << "Available surface names: " << names_ << endl
806  << abort(FatalError);
807  }
808 
809  return operator[](surfI);
810 }
811 
812 
813 // ************************************************************************* //
Various functions to operate on Lists.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Dynamically sized Field.
Definition: DynamicField.H:72
DynamicField< T, SizeInc, SizeMult, SizeDiv > & append(const T &)
Append an element at the end of the list.
label capacity() const
Size of the underlying storage.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
autoPtr< IOobject > clone() const
Clone.
Definition: IOobject.H:283
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
void setSize(const label)
Reset size of List.
Definition: List.C:281
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
static bool checkOrientation(const PrimitivePatch< FaceList, PointField > &, const bool report=false, labelHashSet *marked=0)
Check for orientation issues.
static label markZones(const PrimitivePatch< FaceList, PointField > &, const BoolListType &borderEdge, labelList &faceZone)
Size and fills faceZone with zone of face.
label nEdges() const
Return number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
bool set(const label) const
Is element set.
Definition: PtrListI.H:62
void setSize(const label)
Reset size of PtrList. If extending the PtrList, new entries are.
Definition: PtrList.C:131
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
friend Ostream & operator(Ostream &, const UPtrList< T > &)
Write UPtrList to Ostream.
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
A bounding box defined in terms of the points at its extremities.
Definition: boundBox.H:59
scalar mag() const
The magnitude of the bounding box span.
Definition: boundBoxI.H:96
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
bool isDict(const word &) const
Check if entry is a sub-dictionary.
Definition: dictionary.C:797
const dictionary & subDict(const word &) const
Find and return a sub-dictionary.
Definition: dictionary.C:843
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A class for handling file names.
Definition: fileName.H:82
Triangle with additional region number.
Definition: labelledTri.H:60
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
static autoPtr< searchableSurface > New(const word &surfaceType, const IOobject &io, const dictionary &dict)
Return a reference to the selected searchableSurface.
static boundBox bounds(const PtrList< searchableSurface > &allSurfaces, const labelList &surfacesToTest)
Find the boundBox of the selected surfaces.
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.
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.
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.
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.
Container for searchableSurfaces.
bool checkQuality(const scalar minQuality, const bool report) const
Check triangle quality.
label checkGeometry(const scalar maxRatio, const scalar tolerance, const scalar minQuality, const bool report) const
All geometric checks. Return number of failed checks.
boundBox bounds() const
Calculate bounding box.
label findSurfaceRegionID(const word &surfaceName, const word &regionName) const
bool checkClosed(const bool report) const
Are all surfaces closed and manifold.
bool checkSizes(const scalar maxRatio, const bool report) const
Are all bounding boxes of similar size.
const searchableSurface & operator[](const word &) const
Return const reference to searchableSurface by name.
void findAllIntersections(const pointField &start, const pointField &end, labelListList &surfaces, List< List< pointIndexHit >> &) const
Find all intersections in order from start to end. Returns for.
label checkTopology(const bool report) const
All topological checks. Return number of failed checks.
void findAnyIntersection(const pointField &start, const pointField &end, labelList &surfaces, List< pointIndexHit > &) const
Find any intersection. Return hit point information and.
void writeStats(const List< wordList > &, Ostream &) const
Write some stats.
bool checkNormalOrientation(const bool report) const
Are all (triangulated) surfaces consistent normal orientation.
bool checkIntersection(const scalar tol, const bool report) const
Do surfaces self-intersect or intersect others.
searchableSurfaces(const label)
Construct with length specified. Fill later.
void findNearest(const pointField &, const scalarField &nearestDistSqr, labelList &surfaces, List< pointIndexHit > &) const
Find nearest. Return -1 (and a miss()) or surface and nearest.
label findSurfaceID(const word &name) const
Find index of surface. Return -1 if not found.
void findNearestIntersection(const pointField &start, const pointField &end, labelList &surface1, List< pointIndexHit > &hit1, labelList &surface2, List< pointIndexHit > &hit2) const
A surface geometry formed of discrete facets, e.g. triangles and/or quadrilaterals,...
virtual tmp< pointField > points() const
Get the points that define the surface.
scalar tolerance() const
Return tolerance to use in searches.
Triangulated surface description with patch information.
Definition: triSurface.H:69
A class for handling words, derived from string.
Definition: word.H:62
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void write(const fileName &file, const word &title, const bool binary, const PointField &points, const VertexList &vertices, const LineList &lines, const FaceList &faces, const wordList &fieldNames, const boolList &fieldIsPointValues, const UPtrList< const Field< label >> &fieldLabelValues #define FieldTypeValuesConstArg(Type, nullArg))
Write VTK polygonal data to a file. Takes a PtrList of fields of labels and.
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
const doubleScalar e
Definition: doubleScalar.H:106
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
PointIndexHit< point > pointIndexHit
Definition: pointIndexHit.H:42
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
errorManip< error > abort(error &err)
Definition: errorManip.H:131
const word & regionName(const solver &region)
Definition: solver.H:209
messageStream Info
T returnReduce(const T &Value, const BinaryOp &bop, const int tag=Pstream::msgType(), const label comm=UPstream::worldComm)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
triangle< point, const point & > triPointRef
Definition: triPointRef.H:44
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
error FatalError
labelList identityMap(const label len)
Create identity map (map[i] == i) of given length.
Definition: ListOps.C:104
List< face > faceList
Definition: faceListFwd.H:41
static const char nl
Definition: Ostream.H:266
wordList patchTypes(nPatches)
face triFace(3)
labelList f(nPoints)
dictionary dict
scalarField samples(nIntervals, 0)