sampledMeshedSurface.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-2016 OpenFOAM Foundation
9 Copyright (C) 2016-2021 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
30#include "meshSearch.H"
31#include "Tuple2.H"
32#include "globalIndex.H"
33#include "treeDataCell.H"
34#include "treeDataFace.H"
35#include "meshTools.H"
37
38// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
39
40const Foam::Enum
41<
43>
44Foam::sampledMeshedSurface::samplingSourceNames_
45({
46 { samplingSource::cells, "cells" },
47 { samplingSource::insideCells, "insideCells" },
48 { samplingSource::boundaryFaces, "boundaryFaces" },
49});
50
51
52namespace Foam
53{
55 // Use shorter name only
57 (
60 word,
62 );
63 // Compatibility name (1912)
65 (
68 word,
69 sampledTriSurfaceMesh
70 );
71
72} // End namespace Foam
73
74
75// * * * * * * * * * * * * * * * Local Functions * * * * * * * * * * * * * * //
76
77namespace Foam
78{
79
80// The IOobject for reading
81inline static IOobject selectReadIO(const word& name, const Time& runTime)
82{
83 return IOobject
84 (
85 name,
86 runTime.constant(), // instance
87 "triSurface", // local
88 runTime, // registry
91 false // no register
92 );
93}
94
95} // End namespace Foam
96
97
98// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
99
100void Foam::sampledMeshedSurface::setZoneMap()
101{
102 // Populate zoneIds_ based on surfZone information
103
104 const meshedSurface& s = static_cast<const meshedSurface&>(*this);
105
106 const auto& zones = s.surfZones();
107
108 zoneIds_.resize(s.size());
109
110 // Trivial case
111 if (zoneIds_.empty() || zones.size() <= 1)
112 {
113 zoneIds_ = 0;
114 return;
115 }
116
117
118 label beg = 0;
119
120 forAll(zones, zonei)
121 {
122 const label len = min(zones[zonei].size(), zoneIds_.size() - beg);
123
124 // Assign sub-zone Ids
125 SubList<label>(zoneIds_, len, beg) = zonei;
126
127 beg += len;
128 }
129
130 // Anything remaining? Should not happen
131 {
132 const label len = (zoneIds_.size() - beg);
133
134 if (len > 0)
135 {
136 SubList<label>(zoneIds_, len, beg) = max(0, zones.size()-1);
137 }
138 }
139}
140
141
142// * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
143
144bool Foam::sampledMeshedSurface::update(const meshSearch& meshSearcher)
145{
146 // Global numbering for cells/faces
147 // - only used to uniquely identify local elements
148 globalIndex globalCells(onBoundary() ? mesh().nFaces() : mesh().nCells());
149
150 // Find the cells the triangles of the surface are in.
151 // Does approximation by looking at the face centres only
152 const pointField& fc = surface_.faceCentres();
153
154 // sqr(distance), global index
155 typedef Tuple2<scalar, label> nearInfo;
156 List<nearInfo> nearest(fc.size(), nearInfo(Foam::sqr(GREAT), labelMax));
157
158 if (sampleSource_ == samplingSource::cells)
159 {
160 // Search for nearest cell
161
162 const indexedOctree<treeDataCell>& cellTree = meshSearcher.cellTree();
163
164 forAll(fc, facei)
165 {
166 const point& pt = fc[facei];
167 auto& near = nearest[facei];
168
169 pointIndexHit info = cellTree.findNearest(pt, Foam::sqr(GREAT));
170
171 if (info.hit())
172 {
173 near.first() = magSqr(info.hitPoint()-pt);
174 near.second() = globalCells.toGlobal(info.index());
175 }
176 }
177 }
178 else if (sampleSource_ == samplingSource::insideCells)
179 {
180 // Search for cell containing point
181
182 const auto& cellTree = meshSearcher.cellTree();
183
184 forAll(fc, facei)
185 {
186 const point& pt = fc[facei];
187 auto& near = nearest[facei];
188
189 if (cellTree.bb().contains(pt))
190 {
191 const label index = cellTree.findInside(pt);
192 if (index != -1)
193 {
194 near.first() = 0;
195 near.second() = globalCells.toGlobal(index);
196 }
197 }
198 }
199 }
200 else // samplingSource::boundaryFaces
201 {
202 // Search for nearest boundary face
203 // on all non-coupled boundary faces
204
205 const auto& bndTree = meshSearcher.nonCoupledBoundaryTree();
206
207 forAll(fc, facei)
208 {
209 const point& pt = fc[facei];
210 auto& near = nearest[facei];
211
212 pointIndexHit info = bndTree.findNearest(pt, sqr(GREAT));
213
214 if (info.hit())
215 {
216 near.first() = magSqr(info.hitPoint()-pt);
217 near.second() =
218 globalCells.toGlobal
219 (
220 bndTree.shapes().faceLabels()[info.index()]
221 );
222 }
223 }
224 }
225
226
227 // See which processor has the nearest. Mark and subset
228 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
229
230 Pstream::listCombineAllGather(nearest, minFirstEqOp<scalar>{});
231
232 labelList cellOrFaceLabels(fc.size(), -1);
233
234 bitSet facesToSubset(fc.size());
235
236 forAll(nearest, facei)
237 {
238 const auto& near = nearest[facei];
239
240 const label index = near.second();
241
242 if (index == labelMax)
243 {
244 // Not found on any processor, or simply too far.
245 // How to map?
246 }
247 else if (globalCells.isLocal(index))
248 {
249 facesToSubset.set(facei);
250
251 // Mark as special if too far away
252 cellOrFaceLabels[facei] =
253 (
254 (near.first() < maxDistanceSqr_)
255 ? globalCells.toLocal(index)
256 : -1
257 );
258 }
259 }
260
261
262 if (debug)
263 {
264 Pout<< "Local out of faces:" << cellOrFaceLabels.size()
265 << " keeping:" << facesToSubset.count() << endl;
266 }
267
268
269 // Subset the surface
270 meshedSurface& s = static_cast<meshedSurface&>(*this);
271
272 labelList pointMap;
274
275 s = surface_.subsetMesh(facesToSubset, pointMap, faceMap);
276
277
278 // Ensure zoneIds_ are indeed correct
279 setZoneMap();
280
281
282 // Subset cellOrFaceLabels (for compact faces)
283 cellOrFaceLabels = labelUIndList(cellOrFaceLabels, faceMap)();
284
285
286 // Collect the samplePoints and sampleElements
287 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
288
290 {
291 // With point interpolation
292
293 samplePoints_ = points();
294 sampleElements_.resize(pointMap.size(), -1);
295
296 // Store any face per point (without using pointFaces())
297 labelList pointToFace(std::move(pointMap));
298
299 forAll(s, facei)
300 {
301 const face& f = s[facei];
302
303 for (const label labi : f)
304 {
305 pointToFace[labi] = facei;
306 }
307 }
308
309
310 if (sampleSource_ == samplingSource::cells)
311 {
312 // sampleElements_ : per surface point, the cell
313 // samplePoints_ : per surface point, a location inside the cell
314
315 forAll(samplePoints_, pointi)
316 {
317 // Use _copy_ of point
318 const point pt = samplePoints_[pointi];
319
320 const label celli = cellOrFaceLabels[pointToFace[pointi]];
321
322 sampleElements_[pointi] = celli;
323
324 if
325 (
326 celli >= 0
327 && !mesh().pointInCell(pt, celli, meshSearcher.decompMode())
328 )
329 {
330 // Point not inside cell
331 // Find nearest point on faces of cell
332
333 scalar minDistSqr = VGREAT;
334
335 for (const label facei : mesh().cells()[celli])
336 {
337 const face& f = mesh().faces()[facei];
338
339 pointHit info =
340 f.nearestPoint
341 (
342 pt,
343 mesh().points()
344 );
345
346 if (info.distance() < minDistSqr)
347 {
348 minDistSqr = info.distance();
349 samplePoints_[pointi] = info.rawPoint();
350 }
351 }
352 }
353 }
354 }
355 else if (sampleSource_ == samplingSource::insideCells)
356 {
357 // sampleElements_ : per surface point the cell
358 // samplePoints_ : per surface point a location inside the cell
359
360 forAll(samplePoints_, pointi)
361 {
362 const label celli = cellOrFaceLabels[pointToFace[pointi]];
363
364 sampleElements_[pointi] = celli;
365 }
366 }
367 else // samplingSource::boundaryFaces
368 {
369 // sampleElements_ : per surface point, the boundary face containing
370 // the location
371 // samplePoints_ : per surface point, a location on the boundary
372
373 forAll(samplePoints_, pointi)
374 {
375 const point& pt = samplePoints_[pointi];
376
377 const label facei = cellOrFaceLabels[pointToFace[pointi]];
378
379 sampleElements_[pointi] = facei;
380
381 if (facei >= 0)
382 {
383 samplePoints_[pointi] =
384 mesh().faces()[facei].nearestPoint
385 (
386 pt,
387 mesh().points()
388 ).rawPoint();
389 }
390 }
391 }
392 }
393 else
394 {
395 // if sampleSource_ == cells:
396 // sampleElements_ : per surface face, the cell
397 // samplePoints_ : n/a
398 // if sampleSource_ == insideCells:
399 // sampleElements_ : -1 or per surface face, the cell
400 // samplePoints_ : n/a
401 // if sampleSource_ == boundaryFaces:
402 // sampleElements_ : per surface face, the boundary face
403 // samplePoints_ : n/a
404
405 sampleElements_.transfer(cellOrFaceLabels);
406 samplePoints_.clear();
407 }
408
409
410 if (debug)
411 {
412 this->clearOut();
413 OFstream str(mesh().time().path()/"surfaceToMesh.obj");
414 Info<< "Dumping correspondence from local surface (points or faces)"
415 << " to mesh (cells or faces) to " << str.name() << endl;
416
417 const vectorField& centres =
418 (
419 onBoundary()
420 ? mesh().faceCentres()
421 : mesh().cellCentres()
422 );
423
425 {
426 label vertI = 0;
427 forAll(samplePoints_, pointi)
428 {
429 meshTools::writeOBJ(str, points()[pointi]);
430 ++vertI;
431
432 meshTools::writeOBJ(str, samplePoints_[pointi]);
433 ++vertI;
434
435 label elemi = sampleElements_[pointi];
436 meshTools::writeOBJ(str, centres[elemi]);
437 ++vertI;
438 str << "l " << vertI-2 << ' ' << vertI-1 << ' ' << vertI << nl;
439 }
440 }
441 else
442 {
443 label vertI = 0;
444 forAll(sampleElements_, triI)
445 {
446 meshTools::writeOBJ(str, faceCentres()[triI]);
447 ++vertI;
448
449 label elemi = sampleElements_[triI];
450 meshTools::writeOBJ(str, centres[elemi]);
451 ++vertI;
452 str << "l " << vertI-1 << ' ' << vertI << nl;
453 }
454 }
455 }
456
457 needsUpdate_ = false;
458 return true;
459}
460
461
462// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
463
465(
466 const word& name,
467 const polyMesh& mesh,
468 const word& surfaceName,
469 const samplingSource sampleSource
470)
471:
474 surfaceName_(surfaceName),
475 surface_
476 (
477 selectReadIO(surfaceName, mesh.time()),
478 dictionary::null
479 ),
480 sampleSource_(sampleSource),
481 needsUpdate_(true),
482 keepIds_(true),
483 zoneIds_(),
484 sampleElements_(),
485 samplePoints_(),
486 maxDistanceSqr_(Foam::sqr(GREAT)),
487 defaultValues_()
488{}
489
490
492(
493 const word& name,
494 const polyMesh& mesh,
495 const dictionary& dict
496)
497:
500 surfaceName_
501 (
502 meshedSurface::findFile
503 (
504 selectReadIO(dict.get<word>("surface"), mesh.time()),
505 dict
506 ).name()
507 ),
508 surface_
509 (
510 selectReadIO(dict.get<word>("surface"), mesh.time()),
511 dict
512 ),
513 sampleSource_(samplingSourceNames_.get("source", dict)),
514 needsUpdate_(true),
515 keepIds_(dict.getOrDefault("keepIds", true)),
516 zoneIds_(),
517 sampleElements_(),
518 samplePoints_(),
519 maxDistanceSqr_(Foam::sqr(GREAT)),
520 defaultValues_(dict.subOrEmptyDict("defaultValue"))
521{
522 if (dict.readIfPresent("maxDistance", maxDistanceSqr_))
523 {
524 // Info<< "Limit to maxDistance " << maxDistanceSqr_ << nl;
525 // if (defaultValues_.empty())
526 // {
527 // Info<< "defaultValues = zero" << nl;
528 // }
529 // else
530 // {
531 // defaultValues_.writeEntry(Info);
532 // }
533
534 maxDistanceSqr_ = Foam::sqr(maxDistanceSqr_);
535 }
536
537 wordRes includePatches;
538 dict.readIfPresent("patches", includePatches);
539 includePatches.uniq();
540
541 // Could also shift this to the reader itself,
542 // but not yet necessary.
543
544 if (!includePatches.empty())
545 {
546 Info<< "Subsetting surface " << surfaceName_
547 << " to patches: " << flatOutput(includePatches) << nl;
548
549 const surfZoneList& zones = surface_.surfZones();
550
551 const labelList zoneIndices
552 (
554 (
555 zones,
556 includePatches,
557 wordRes(),
559 )
560 );
561
562 // Faces to subset
563 bitSet includeMap(surface_.size());
564
565 for (const label zonei : zoneIndices)
566 {
567 const surfZone& zn = zones[zonei];
568 includeMap.set(zn.range());
569 }
570
571 if (includeMap.none())
572 {
574 << "Patch selection results in an empty surface"
575 << " - ignoring" << nl;
576 }
577 else if (!includeMap.all())
578 {
579 meshedSurface subSurf(surface_.subsetMesh(includeMap));
580
581 if (subSurf.empty())
582 {
584 << "Bad surface subset (empty)"
585 << " - skip and hope for the best" << nl;
586 }
587 else
588 {
589 // Replace
590 surface_.transfer(subSurf);
591 }
592 }
593 }
594}
595
596
597// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
598
600{
601 return needsUpdate_;
602}
603
604
606{
607 // already marked as expired
608 if (needsUpdate_)
609 {
610 return false;
611 }
612
615 zoneIds_.clear();
616
617 sampleElements_.clear();
618 samplePoints_.clear();
619
620 needsUpdate_ = true;
621 return true;
622}
623
624
626{
627 if (!needsUpdate_)
628 {
629 return false;
630 }
631
632 // Calculate surface and mesh overlap bounding box
633 treeBoundBox bb(surface_.points(), surface_.meshPoints());
634
635 // Check for overlap with (global!) mesh bb
636 const bool intersect = bb.intersect(mesh().bounds());
637
638 if (!intersect)
639 {
640 // Surface and mesh do not overlap at all. Guarantee a valid
641 // bounding box so we don't get any 'invalid bounding box' errors.
642
644 << "Surface " << surfaceName_
645 << " does not overlap bounding box of mesh " << mesh().bounds()
646 << endl;
647
648 bb = treeBoundBox(mesh().bounds());
649 const vector span(bb.span());
650
651 bb.min() += (0.5-1e-6)*span;
652 bb.max() -= (0.5-1e-6)*span;
653 }
654 else
655 {
656 // Extend a bit
657 const vector span(bb.span());
658 bb.min() -= 0.5*span;
659 bb.max() += 0.5*span;
660
661 bb.inflate(1e-6);
662 }
663
664 // Mesh search engine, no triangulation of faces.
665 meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
666
667 return update(meshSearcher);
668}
669
670
672{
673 if (!needsUpdate_)
674 {
675 return false;
676 }
677
678 // Mesh search engine on subset, no triangulation of faces.
679 meshSearch meshSearcher(mesh(), bb, polyMesh::FACE_PLANES);
680
681 return update(meshSearcher);
682}
683
684
686(
687 const interpolation<scalar>& sampler
688) const
689{
690 return sampleOnFaces(sampler);
691}
692
693
695(
696 const interpolation<vector>& sampler
697) const
698{
699 return sampleOnFaces(sampler);
700}
701
702
704(
705 const interpolation<sphericalTensor>& sampler
706) const
707{
708 return sampleOnFaces(sampler);
709}
710
711
713(
714 const interpolation<symmTensor>& sampler
715) const
716{
717 return sampleOnFaces(sampler);
718}
719
720
722(
723 const interpolation<tensor>& sampler
724) const
725{
726 return sampleOnFaces(sampler);
727}
728
729
731(
732 const interpolation<scalar>& interpolator
733) const
734{
735 return sampleOnPoints(interpolator);
736}
737
738
740(
741 const interpolation<vector>& interpolator
742) const
743{
744 return sampleOnPoints(interpolator);
745}
746
748(
749 const interpolation<sphericalTensor>& interpolator
750) const
751{
752 return sampleOnPoints(interpolator);
753}
754
755
757(
758 const interpolation<symmTensor>& interpolator
759) const
760{
761 return sampleOnPoints(interpolator);
762}
763
764
766(
767 const interpolation<tensor>& interpolator
768) const
769{
770 return sampleOnPoints(interpolator);
771}
772
773
775{
776 os << "meshedSurface: " << name() << " :"
777 << " surface:" << surfaceName_;
778
779 if (level)
780 {
781 os << " faces:" << faces().size()
782 << " points:" << points().size()
783 << " zoneids:" << zoneIds().size();
784 }
785}
786
787
788// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addNamedToRunTimeSelectionTable(baseType, thisType, argNames, lookupName)
Add to construction table with 'lookupName' as the key.
Enum is a wrapper around a list of names/values that represent particular enumeration (or int) values...
Definition: Enum.H:61
Minimal example by using system/controlDict.functions:
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
void transfer(List< T > &list)
Definition: List.C:447
void resize(const label len)
Adjust allocated size of list.
Definition: ListI.H:139
MeshedSurface subsetMesh(const UList< bool > &include, labelList &pointMap, labelList &faceMap) const
Return a new surface subsetted on the selected faces.
const surfZoneList & surfZones() const
Const access to the surface zones.
void transfer(pointField &pointLst, List< Face > &faceLst)
Transfer the components.
label size() const
The surface size is the number of faces.
virtual void clear()
Clear all storage.
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:62
static void listCombineAllGather(const List< commsStruct > &comms, List< T > &values, const CombineOp &cop, const int tag, const label comm)
After completion all processors have the same data.
const word & constant() const
Return constant name.
Definition: TimePathsI.H:96
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:80
bool empty() const noexcept
True if the UList is empty (ie, size() is zero)
Definition: UListI.H:427
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
A bitSet stores bits (elements with only two states) in packed internal format and supports a variety...
Definition: bitSet.H:66
void set(const bitSet &bitset)
Set specified bits from another bitset.
Definition: bitSetI.H:590
bool none() const
True if no bits in this bitset are set.
Definition: bitSetI.H:488
bool all() const
True if all bits in this bitset are set or if the set is empty.
Definition: bitSetI.H:461
const point & min() const
Minimum describing the bounding box.
Definition: boundBoxI.H:91
const point & max() const
Maximum describing the bounding box.
Definition: boundBoxI.H:97
void inflate(const scalar s)
Inflate box by factor*mag(span) in all dimensions.
Definition: boundBox.C:175
bool intersect(const boundBox &bb)
Intersection (union) with the second box.
Definition: boundBox.C:191
vector span() const
The bounding box span (from minimum to maximum)
Definition: boundBoxI.H:127
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Abstract base class for volume field interpolation.
Definition: interpolation.H:60
Various (local, not parallel) searches on polyMesh; uses (demand driven) octree to search.
Definition: meshSearch.H:61
scalar print()
Print to screen.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:81
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1108
const boundBox & bounds() const
Return mesh bounding box.
Definition: polyMesh.H:462
const vectorField & faceCentres() const
const vectorField & cellCentres() const
A sampledSurface from a meshed surface. It samples on the points/faces of the meshed surface.
virtual bool expire()
Mark the surface as needing an update.
virtual bool needsUpdate() const
Does the surface need an update?
virtual bool update()
Update the surface as required.
samplingSource
Types of sampling regions.
An abstract class for surfaces with sampling.
bool isPointData() const noexcept
Using interpolation to surface points.
virtual void clearGeom() const
Additional cleanup when clearing the geometry.
bool interpolate() const noexcept
Same as isPointData()
A surface zone on a MeshedSurface.
Definition: surfZone.H:59
labelRange range() const
The start/size range of this zone.
Definition: surfZone.H:152
A class for managing temporary objects.
Definition: tmp.H:65
Standard boundBox with extra functionality for use in octree.
Definition: treeBoundBox.H:89
A List of wordRe with additional matching capabilities.
Definition: wordRes.H:54
static wordRes uniq(const UList< wordRe > &input)
Return a wordRes with duplicate entries filtered out.
Definition: wordRes.C:32
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
fileName path(UMean.rootPath()/UMean.caseName()/"graphs"/UMean.instance())
mesh update()
dynamicFvMesh & mesh
engineTime & runTime
OBJstream os(runTime.globalPath()/outputName)
const labelList nFaces(UPstream::listGatherValues< label >(aMesh.nFaces()))
const pointField & points
const cellShapeList & cells
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))
#define WarningInFunction
Report a warning using Foam::Warning.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
labelList findMatching(const StringListType &input, const wordRes::filter &pred, AccessOp aop=identityOp())
Return ids for items with matching names.
Namespace for OpenFOAM.
Pair< int > faceMap(const label facePi, const face &faceP, const label faceNi, const face &faceN)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
List< label > labelList
A List of labels.
Definition: List.H:66
dimensionedSymmTensor sqr(const dimensionedVector &dv)
PointIndexHit< point > pointIndexHit
A PointIndexHit for 3D points.
Definition: pointIndexHit.H:46
messageStream Info
Information stream (stdout output on master, null elsewhere)
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:44
vector point
Point is a vector.
Definition: point.H:43
constexpr label labelMax
Definition: label.H:61
static IOobject selectReadIO(const word &name, const Time &runTime)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
MeshedSurface< face > meshedSurface
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
FlatOutput::OutputAdaptor< Container, Delimiters > flatOutput(const Container &obj, Delimiters delim)
Global flatOutput() function with specified output delimiters.
Definition: FlatOutput.H:215
PointHit< point > pointHit
A PointIndexHit for 3D points.
Definition: pointHit.H:44
Field< vector > vectorField
Specialisation of Field<T> for vector.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
UIndirectList< label > labelUIndList
UIndirectList of labels.
Definition: IndirectList.H:68
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
dictionary dict
volScalarField & e
Definition: createFields.H:11
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333
Extract name (as a word) from an object, typically using its name() method.
Definition: word.H:238