Loading...
Searching...
No Matches
searchableSurfaceCollection.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) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2016-2022 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
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
31#include "Time.H"
32#include "ListOps.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
40 (
43 dict
44 );
46 (
49 dict,
50 collection
51 );
52}
53
54
55// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
56
57void Foam::searchableSurfaceCollection::findNearest
58(
59 const pointField& samples,
60 scalarField& minDistSqr,
61 List<pointIndexHit>& nearestInfo,
62 labelList& nearestSurf
63) const
64{
65 // Initialise
66 nearestInfo.setSize(samples.size());
67 nearestInfo = pointIndexHit();
68 nearestSurf.setSize(samples.size());
69 nearestSurf = -1;
70
71 List<pointIndexHit> hitInfo(samples.size());
72
73 const scalarField localMinDistSqr(samples.size(), GREAT);
74
75 forAll(subGeom_, surfI)
76 {
77 subGeom_[surfI].findNearest
78 (
79 cmptDivide // Transform then divide
80 (
81 transform_[surfI].localPosition(samples),
82 scale_[surfI]
83 ),
84 localMinDistSqr,
85 hitInfo
86 );
87
88 forAll(hitInfo, pointi)
89 {
90 if (hitInfo[pointi].hit())
91 {
92 // Rework back into global coordinate sys. Multiply then
93 // transform
94 point globalPt = transform_[surfI].globalPosition
95 (
97 (
98 hitInfo[pointi].point(),
99 scale_[surfI]
100 )
101 );
102
103 scalar distSqr = globalPt.distSqr(samples[pointi]);
104
105 if (distSqr < minDistSqr[pointi])
106 {
107 minDistSqr[pointi] = distSqr;
108 nearestInfo[pointi].hitPoint(globalPt);
109 nearestInfo[pointi].setIndex
110 (
111 hitInfo[pointi].index()
112 + indexOffset_[surfI]
113 );
114 nearestSurf[pointi] = surfI;
115 }
116 }
117 }
118 }
119}
120
121
122// Sort hits into per-surface bins. Misses are rejected. Maintains map back
123// to position
124void Foam::searchableSurfaceCollection::sortHits
125(
126 const List<pointIndexHit>& info,
127 List<List<pointIndexHit>>& surfInfo,
128 labelListList& infoMap
129) const
130{
131 // Count hits per surface.
132 labelList nHits(subGeom_.size(), Zero);
133
134 forAll(info, pointi)
135 {
136 if (info[pointi].hit())
137 {
138 label index = info[pointi].index();
139 label surfI = findLower(indexOffset_, index+1);
140 nHits[surfI]++;
141 }
142 }
143
144 // Per surface the hit
145 surfInfo.setSize(subGeom_.size());
146 // Per surface the original position
147 infoMap.setSize(subGeom_.size());
148
149 forAll(surfInfo, surfI)
150 {
151 surfInfo[surfI].setSize(nHits[surfI]);
152 infoMap[surfI].setSize(nHits[surfI]);
153 }
154 nHits = 0;
155
156 forAll(info, pointi)
157 {
158 if (info[pointi].hit())
159 {
160 label index = info[pointi].index();
161 label surfI = findLower(indexOffset_, index+1);
162
163 // Store for correct surface and adapt indices back to local
164 // ones
165 label localI = nHits[surfI]++;
166 surfInfo[surfI][localI] = pointIndexHit
167 (
168 info[pointi],
169 index-indexOffset_[surfI]
170 );
171 infoMap[surfI][localI] = pointi;
173 }
174}
175
176
177// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
178
179Foam::searchableSurfaceCollection::searchableSurfaceCollection
180(
181 const IOobject& io,
182 const dictionary& dict
183)
184:
186 instance_(dict.size()),
187 scale_(dict.size()),
188 transform_(dict.size()),
189 subGeom_(dict.size()),
190 mergeSubRegions_(dict.get<bool>("mergeSubRegions")),
191 indexOffset_(dict.size()+1)
192{
193 Info<< "SearchableCollection : " << name() << endl;
194
195 label surfI = 0;
196 label startIndex = 0;
197 for (const entry& dEntry : dict)
198 {
199 if (dEntry.isDict())
200 {
201 instance_[surfI] = dEntry.keyword();
202
203 const dictionary& sDict = dEntry.dict();
204
205 sDict.readEntry("scale", scale_[surfI]);
206
207 // Mandatory 'transform' sub-dictionary
208 transform_.set
209 (
210 surfI,
211 new coordSystem::cartesian(sDict, "transform")
212 );
213
214 const word subGeomName(sDict.get<word>("surface"));
215 //Pout<< "Trying to find " << subGeomName << endl;
216
217 searchableSurface& s =
218 io.db().lookupObjectRef<searchableSurface>(subGeomName);
219
220 // I don't know yet how to handle the globalSize combined with
221 // regionOffset. Would cause non-consecutive indices locally
222 // if all indices offset by globalSize() of the local region...
223 if (s.size() != s.globalSize())
224 {
226 << "Cannot use a distributed surface in a collection."
227 << exit(FatalError);
228 }
229
230 subGeom_.set(surfI, &s);
231
232 indexOffset_[surfI] = startIndex;
233 startIndex += subGeom_[surfI].size();
234
235 Info<< " instance : " << instance_[surfI] << endl;
236 Info<< " surface : " << s.name() << endl;
237 Info<< " scale : " << scale_[surfI] << endl;
238 Info<< " transform: " << transform_[surfI] << endl;
239
240 surfI++;
241 }
242 }
243 indexOffset_[surfI] = startIndex;
244
245 instance_.setSize(surfI);
246 scale_.setSize(surfI);
247 transform_.setSize(surfI);
248 subGeom_.setSize(surfI);
249 indexOffset_.setSize(surfI+1);
250
251 // Bounds is the overall bounds - prepare for min/max ops
252 bounds().reset();
253
254 forAll(subGeom_, surfI)
255 {
256 const boundBox& surfBb = subGeom_[surfI].bounds();
257
258 // Transform back to global coordinate sys.
259 const point surfBbMin = transform_[surfI].globalPosition
260 (
262 (
263 surfBb.min(),
264 scale_[surfI]
265 )
266 );
267 const point surfBbMax = transform_[surfI].globalPosition
268 (
270 (
271 surfBb.max(),
272 scale_[surfI]
273 )
274 );
275
276 bounds().min() = min(bounds().min(), surfBbMin);
277 bounds().max() = max(bounds().max(), surfBbMax);
278 }
279}
280
281
282// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
285{}
286
287
288// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
289
291{
292 if (regions_.empty())
293 {
294 regionOffset_.setSize(subGeom_.size());
295
296 DynamicList<word> allRegions;
297 forAll(subGeom_, surfI)
298 {
299 regionOffset_[surfI] = allRegions.size();
300
301 if (mergeSubRegions_)
302 {
303 // Single name regardless how many regions subsurface has
304 allRegions.append(instance_[surfI] + "_" + Foam::name(surfI));
305 }
306 else
307 {
308 const wordList& subRegions = subGeom_[surfI].regions();
309
310 for (const word& regionName : subRegions)
311 {
312 allRegions.append(instance_[surfI] + "_" + regionName);
313 }
314 }
316 regions_.transfer(allRegions);
317 }
318 return regions_;
319}
320
321
323{
324 return indexOffset_.last();
325}
326
327
330{
331 auto tctrs = tmp<pointField>::New(size());
332 auto& ctrs = tctrs.ref();
333
334 // Append individual coordinates
335 label coordI = 0;
336
337 forAll(subGeom_, surfI)
338 {
339 const pointField subCoords(subGeom_[surfI].coordinates());
340
341 forAll(subCoords, i)
342 {
343 ctrs[coordI++] = transform_[surfI].globalPosition
344 (
346 (
347 subCoords[i],
348 scale_[surfI]
349 )
350 );
352 }
353
354 return tctrs;
355}
356
357
359(
360 pointField& centres,
361 scalarField& radiusSqr
362) const
363{
364 centres.setSize(size());
365 radiusSqr.setSize(centres.size());
366
367 // Append individual coordinates
368 label coordI = 0;
369
370 forAll(subGeom_, surfI)
371 {
372 scalar maxScale = cmptMax(scale_[surfI]);
373
374 pointField subCentres;
375 scalarField subRadiusSqr;
376 subGeom_[surfI].boundingSpheres(subCentres, subRadiusSqr);
377
378 forAll(subCentres, i)
379 {
380 centres[coordI] = transform_[surfI].globalPosition
381 (
383 (
384 subCentres[i],
385 scale_[surfI]
386 )
387 );
388 radiusSqr[coordI] = maxScale*subRadiusSqr[i];
389 coordI++;
390 }
391 }
392}
393
394
395Foam::tmp<Foam::pointField>
397{
398 // Get overall size
399 label nPoints = 0;
400
401 forAll(subGeom_, surfI)
402 {
403 nPoints += subGeom_[surfI].points()().size();
404 }
405
406 auto tpts = tmp<pointField>::New(nPoints);
407 auto& pts = tpts.ref();
408
409 // Append individual coordinates
410 nPoints = 0;
411
412 forAll(subGeom_, surfI)
413 {
414 const pointField subCoords(subGeom_[surfI].points());
415
416 forAll(subCoords, i)
417 {
418 pts[nPoints++] = transform_[surfI].globalPosition
419 (
421 (
422 subCoords[i],
423 scale_[surfI]
424 )
425 );
427 }
428
429 return tpts;
430}
431
432
433void Foam::searchableSurfaceCollection::findNearest
434(
435 const pointField& samples,
436 const scalarField& nearestDistSqr,
437 List<pointIndexHit>& nearestInfo
438) const
439{
440 // How to scale distance?
441 scalarField minDistSqr(nearestDistSqr);
442
443 labelList nearestSurf;
444 findNearest
445 (
446 samples,
447 minDistSqr,
448 nearestInfo,
449 nearestSurf
450 );
451}
452
453
455(
456 const pointField& start,
457 const pointField& end,
459) const
460{
461 info.setSize(start.size());
462 info = pointIndexHit();
463
464 // Current nearest (to start) intersection
465 pointField nearest(end);
466
467 List<pointIndexHit> hitInfo(start.size());
468
469 forAll(subGeom_, surfI)
470 {
471 // Starting point
473 (
474 transform_[surfI].localPosition
475 (
476 start
477 ),
478 scale_[surfI]
479 );
480
481 // Current best end point
483 (
484 transform_[surfI].localPosition
485 (
486 nearest
487 ),
488 scale_[surfI]
489 );
490
491 subGeom_[surfI].findLine(e0, e1, hitInfo);
492
493 forAll(hitInfo, pointi)
494 {
495 if (hitInfo[pointi].hit())
496 {
497 // Transform back to global coordinate sys.
498 nearest[pointi] = transform_[surfI].globalPosition
499 (
501 (
502 hitInfo[pointi].point(),
503 scale_[surfI]
504 )
505 );
506 info[pointi] = hitInfo[pointi];
507 info[pointi].point() = nearest[pointi];
508 info[pointi].setIndex
509 (
510 hitInfo[pointi].index()
511 + indexOffset_[surfI]
512 );
513 }
514 }
515 }
516
517
518 // Debug check
519 if (false)
520 {
521 forAll(info, pointi)
522 {
523 if (info[pointi].hit())
524 {
525 vector n(end[pointi] - start[pointi]);
526 scalar magN = mag(n);
527
528 if (magN > SMALL)
529 {
530 n /= mag(n);
531
532 scalar s = ((info[pointi].point()-start[pointi])&n);
533
534 if (s < 0 || s > 1)
535 {
537 << "point:" << info[pointi]
538 << " s:" << s
539 << " outside vector "
540 << " start:" << start[pointi]
541 << " end:" << end[pointi]
542 << abort(FatalError);
543 }
545 }
546 }
547 }
548}
549
550
552(
553 const pointField& start,
554 const pointField& end,
556) const
557{
558 // To be done ...
559 findLine(start, end, info);
560}
561
562
564(
565 const pointField& start,
566 const pointField& end,
568) const
569{
570 // To be done. Assume for now only one intersection.
571 List<pointIndexHit> nearestInfo;
572 findLine(start, end, nearestInfo);
573
574 info.setSize(start.size());
575 forAll(info, pointi)
576 {
577 if (nearestInfo[pointi].hit())
578 {
579 info[pointi].setSize(1);
580 info[pointi][0] = nearestInfo[pointi];
581 }
582 else
584 info[pointi].clear();
585 }
586 }
587}
588
589
591(
592 const List<pointIndexHit>& info,
593 labelList& region
594) const
595{
596 if (subGeom_.size() == 0)
597 {}
598 else if (subGeom_.size() == 1)
599 {
600 if (mergeSubRegions_)
601 {
602 region.setSize(info.size());
603 region = regionOffset_[0];
604 }
605 else
606 {
607 subGeom_[0].getRegion(info, region);
608 }
609 }
610 else
611 {
612 // Multiple surfaces. Sort by surface.
613
614 // Per surface the hit
615 List<List<pointIndexHit>> surfInfo;
616 // Per surface the original position
617 List<List<label>> infoMap;
618 sortHits(info, surfInfo, infoMap);
619
620 region.setSize(info.size());
621 region = -1;
622
623 // Do region tests
624
625 if (mergeSubRegions_)
626 {
627 // Actually no need for surfInfo. Just take region for surface.
628 forAll(infoMap, surfI)
629 {
630 const labelList& map = infoMap[surfI];
631 forAll(map, i)
632 {
633 region[map[i]] = regionOffset_[surfI];
634 }
635 }
636 }
637 else
638 {
639 forAll(infoMap, surfI)
640 {
641 labelList surfRegion;
642 subGeom_[surfI].getRegion(surfInfo[surfI], surfRegion);
643
644 const labelList& map = infoMap[surfI];
645 forAll(map, i)
646 {
647 region[map[i]] = regionOffset_[surfI] + surfRegion[i];
649 }
650 }
651 }
652}
653
654
656(
657 const List<pointIndexHit>& info,
658 vectorField& normal
659) const
660{
661 if (subGeom_.size() == 0)
662 {}
663 else if (subGeom_.size() == 1)
664 {
665 subGeom_[0].getNormal(info, normal);
666 }
667 else
668 {
669 // Multiple surfaces. Sort by surface.
670
671 // Per surface the hit
672 List<List<pointIndexHit>> surfInfo;
673 // Per surface the original position
674 List<List<label>> infoMap;
675 sortHits(info, surfInfo, infoMap);
676
677 normal.setSize(info.size());
678
679 // Do region tests
680 forAll(surfInfo, surfI)
681 {
682 vectorField surfNormal;
683 subGeom_[surfI].getNormal(surfInfo[surfI], surfNormal);
684
685 // Transform back to global coordinate sys.
686 surfNormal = transform_[surfI].globalVector(surfNormal);
687
688 const labelList& map = infoMap[surfI];
689 forAll(map, i)
690 {
691 normal[map[i]] = surfNormal[i];
692 }
693 }
694 }
695}
696
697
699(
700 const pointField& points,
701 List<volumeType>& volType
702) const
705 << "Volume type not supported for collection."
706 << exit(FatalError);
707}
708
709
711(
712 const List<treeBoundBox>& bbs,
713 const bool keepNonLocal,
715 autoPtr<mapDistribute>& pointMap
716)
717{
718 forAll(subGeom_, surfI)
719 {
720 // Note:Transform the bounding boxes? Something like
721 // pointField bbPoints =
722 // cmptDivide
723 // (
724 // transform_[surfI].localPosition(bbs[i].points()),
725 // scale_[surfI]
726 // );
727 // treeBoundBox newBb(bbPoints);
728
729 // Note: what to do with faceMap, pointMap from multiple surfaces?
730 subGeom_[surfI].distribute
731 (
732 bbs,
733 keepNonLocal,
735 pointMap
736 );
737 }
738}
739
740
742{
743 forAll(subGeom_, surfI)
744 {
745 subGeom_[surfI].setField
746 (
747 static_cast<const labelList&>
748 (
749 SubList<label>
750 (
751 values,
752 subGeom_[surfI].size(),
753 indexOffset_[surfI]
755 )
756 );
757 }
758}
759
760
762(
763 const List<pointIndexHit>& info,
764 labelList& values
765) const
766{
767 if (subGeom_.size() == 0)
768 {}
769 else if (subGeom_.size() == 1)
770 {
771 subGeom_[0].getField(info, values);
772 }
773 else
774 {
775 // Multiple surfaces. Sort by surface.
776
777 // Per surface the hit
778 List<List<pointIndexHit>> surfInfo;
779 // Per surface the original position
780 List<List<label>> infoMap;
781 sortHits(info, surfInfo, infoMap);
782
783 // Do surface tests
784 forAll(surfInfo, surfI)
785 {
786 labelList surfValues;
787 subGeom_[surfI].getField(surfInfo[surfI], surfValues);
788
789 if (surfValues.size())
790 {
791 // Size values only when we have a surface that supports it.
792 values.setSize(info.size());
793
794 const labelList& map = infoMap[surfI];
795 forAll(map, i)
796 {
797 values[map[i]] = surfValues[i];
798 }
799 }
800 }
801 }
802}
803
804
805// ************************************************************************* //
Various functions to operate on Lists.
label n
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
void append(const T &val)
Copy append an element to the end of this list.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
const word & name() const noexcept
Return the object name.
Definition IOobjectI.H:205
InfoProxy< IOobject > info() const noexcept
Return info proxy, for printing information to a stream.
Definition IOobject.H:1041
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition List.H:72
void transfer(List< T > &list)
Transfer the contents of the argument List into this list and annul the argument list.
Definition List.C:347
void setSize(label n)
Alias for resize().
Definition List.H:536
A non-owning sub-view of a List (allocated or unallocated storage).
Definition SubList.H:61
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
A bounding box defined in terms of min/max extrema points.
Definition boundBox.H:71
const point & max() const noexcept
Maximum describing the bounding box.
Definition boundBoxI.H:168
const point & min() const noexcept
Minimum describing the bounding box.
Definition boundBoxI.H:162
void reset()
Reset to an inverted box.
Definition boundBoxI.H:295
A Cartesian coordinate system.
Definition cartesianCS.H:68
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
T get(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a T. FatalIOError if not found, or if the number of tokens is incorrect.
bool readEntry(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX, IOobjectOption::readOption readOpt=IOobjectOption::MUST_READ) const
Find entry and assign to T val. FatalIOError if it is found and the number of tokens is incorrect,...
A keyword and a list of tokens is an 'entry'.
Definition entry.H:66
Set of transformed searchableSurfaces. Does not do boolean operations so when meshing might find part...
virtual void findLine(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Find first intersection on segment from start to end.
virtual void distribute(const List< treeBoundBox > &, const bool keepNonLocal, autoPtr< mapDistribute > &faceMap, autoPtr< mapDistribute > &pointMap)
Set bounds of surface. Bounds currently set as list of.
virtual void findLineAll(const pointField &start, const pointField &end, List< List< pointIndexHit > > &) const
Get all intersections in order from start to end.
virtual void findLineAny(const pointField &start, const pointField &end, List< pointIndexHit > &) const
Return any intersection on segment from start to end.
virtual label size() const
Range of local indices that can be returned.
virtual void boundingSpheres(pointField &centres, scalarField &radiusSqr) const
Get bounding spheres (centre and radius squared), one per element.
virtual void getVolumeType(const pointField &points, List< volumeType > &volType) const
Determine type (inside/outside/mixed) for point.
virtual void setField(const labelList &values)
WIP. Store element-wise field.
virtual void getField(const List< pointIndexHit > &, labelList &) const
WIP. From a set of hits (points and.
virtual void getNormal(const List< pointIndexHit > &, vectorField &normal) const
From a set of points and indices get the normal.
virtual void getRegion(const List< pointIndexHit > &, labelList &region) const
From a set of points and indices get the region.
virtual const wordList & regions() const
Names of regions.
virtual tmp< pointField > coordinates() const
Get representative set of element coordinates.
virtual tmp< pointField > points() const
Get the points that define the surface.
Base class of (analytical or triangulated) surface. Encapsulates all the search routines....
virtual const boundBox & bounds() const
Return const reference to boundBox.
A class for managing temporary objects.
Definition tmp.H:75
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition tmp.H:215
A class for handling words, derived from Foam::string.
Definition word.H:66
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
PtrList< coordinateSystem > coordinates(solidRegions.size())
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
const pointField & points
label nPoints
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))
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition HashOps.H:164
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
const char * end
Definition SVGTools.H:223
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
List< word > wordList
List of word.
Definition fileName.H:60
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
dimensioned< Type > cmptDivide(const dimensioned< Type > &, const dimensioned< Type > &)
messageStream Info
Information stream (stdout output on master, null elsewhere).
label findLower(const ListType &input, const T &val, const label start, const ComparePredicate &comp)
Binary search to find the index of the last element in a sorted list that is less than value.
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
errorManip< error > abort(error &err)
Definition errorManip.H:139
Field< vector > vectorField
Specialisation of Field<T> for vector.
vector point
Point is a vector.
Definition point.H:37
void cmptMax(FieldField< Field, typename FieldField< Field, Type >::cmptType > &cf, const FieldField< Field, Type > &f)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
dimensioned< Type > cmptMultiply(const dimensioned< Type > &, const dimensioned< Type > &)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition exprTraits.C:127
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
PointIndexHit< point > pointIndexHit
A PointIndexHit with a 3D point.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict
const pointField & pts
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
scalarField samples(nIntervals, Zero)