Loading...
Searching...
No Matches
raySearchEngine.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) 2023-2024 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "raySearchEngine.H"
29#include "surfaceFields.H"
30#include "volFields.H"
31#include "meshTools.H"
32
33using namespace Foam::constant;
34
35// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36
37namespace Foam
38{
39namespace VF
40{
43
44 const label raySearchEngine::maxDynListLength = 1000000000;
45}
46}
47
48// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
49
51(
52 const labelList& nVisibleFaceFaces
53)
54{
55 label nRay = 0;
56 label nFaceMin = labelMax;
57 label nFaceMax = labelMin;
58 for (const label n : nVisibleFaceFaces)
59 {
60 nFaceMin = min(nFaceMin, n);
61 nFaceMax = max(nFaceMax, n);
62 nRay += n;
63 }
64
65 const label nFace = nVisibleFaceFaces.size();
66 const label nGlobalRays = returnReduce(nRay, sumOp<label>());
67
68 if (nGlobalRays == 0)
69 {
71 << "No rays identified - view factors will not be calculated"
72 << exit(FatalError);
73 }
74
75 const label globalNFacesMin = returnReduce(nFaceMin, minOp<label>());
76 const label globalNFacesMax = returnReduce(nFaceMax, maxOp<label>());
77 const label nGlobalFaces = returnReduce(nFace, sumOp<label>());
78 const scalar avgFaces = nGlobalRays/scalar(nGlobalFaces);
79
80 Info<< "\nRay summary:" << nl
81 << " Number of rays: " << nGlobalRays << nl
82 << " Number of rays-per-face (min, max, average): ("
83 << globalNFacesMin << ", "
84 << globalNFacesMax << ", "
85 << avgFaces << ")" << endl;
86}
87
88
90(
91 const point& p0,
92 const List<point>& pts
93)
94{
95 label pointi = -1;
96
97 scalar distSqr = GREAT;
98 forAll(pts, pti)
99 {
100 const scalar d2 = magSqr(pts[pti] - p0);
101 if (d2 < distSqr)
102 {
103 pointi = pti;
104 distSqr = d2;
105 }
106 }
107
108 return pointi;
109}
110
111
113{
114 Info<< "\nFace agglomeration: active" << nl
115 << " Reading file " << io.name() << endl;
116
117 // Read agglomeration map
118 const labelListIOList finalAgglom(io);
119
120 Info<< " Creating coarse mesh" << nl;
121
122 agglomMeshPtr_.reset
123 (
124 new singleCellFvMesh
125 (
126 IOobject
127 (
128 IOobject::scopedName("agglom", mesh_.name()),
129 mesh_.time().timeName(),
130 mesh_.time(),
133 ),
134 mesh_,
135 finalAgglom
136 )
137 );
138
139 const auto& coarseMesh = agglomMeshPtr_();
140
141
142 // Calculate total number of fine and coarse faces
143
144 nCoarseFace_ = 0;
145 nFace_ = 0;
146
147 const polyBoundaryMesh& finePatches = mesh_.boundaryMesh();
148 const polyBoundaryMesh& coarsePatches = coarseMesh.boundaryMesh();
149
150 for (const label patchi : patchIDs_)
151 {
152 nCoarseFace_ += coarsePatches[patchi].size();
153 nFace_ += finePatches[patchi].size();
154 }
155
156 Info<< "\nTotal number of coarse faces: "
157 << returnReduce(nCoarseFace_, sumOp<label>())
158 << endl;
159
160 Info<< "\nTotal number of fine faces: "
161 << returnReduce(nFace_, sumOp<label>())
162 << endl;
163
164 // Collect local Cf, Sf, agglom index on coarse mesh
165 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
166
167 DynamicList<point> localCf(nCoarseFace_);
168 DynamicList<vector> localSf(nCoarseFace_);
169 DynamicList<label> localAgg(nCoarseFace_);
170
171 for (const label patchi : patchIDs_)
172 {
173 const labelList& agglom = finalAgglom[patchi];
174
175 if (agglom.empty()) continue;
176
177 label nAgglom = max(agglom) + 1;
178 const labelListList coarseToFine(invertOneToMany(nAgglom, agglom));
179 const labelList& coarsePatchFace = coarseMesh.patchFaceMap()[patchi];
180
181 const pointField& coarseCf = coarseMesh.Cf().boundaryField()[patchi];
182 const vectorField& coarseSf = coarseMesh.Sf().boundaryField()[patchi];
183
184 const polyPatch& pp = finePatches[patchi];
185 patchAreas_[patchi] += sum(coarseMesh.magSf().boundaryField()[patchi]);
186
187 forAll(coarseCf, facei)
188 {
189 const label coarseFacei = coarsePatchFace[facei];
190 const label agglomi = coarseFacei + coarsePatches[patchi].start();
191
192 // Construct single coarse face
193 const labelList& fineFaces = coarseToFine[coarseFacei];
195 (
196 UIndirectList<face>(pp, fineFaces),
197 pp.points()
198 );
199
200 // Collect all points (vertices, face centres)
201 const label nFaces = cf.faceCentres().size();
202 const label nPoints = cf.localPoints().size();
203 List<point> allPoints(nFaces + nPoints);
204 SubList<point>(allPoints, nFaces) = cf.faceCentres();
205 SubList<point>(allPoints, nPoints, nFaces) = cf.localPoints();
206
207 // Coarse face centre set to closest point
208 const label pti = closestPointIndex(coarseCf[facei], allPoints);
209
210 if (pti != -1)
211 {
212 localCf.push_back(allPoints[pti]);
213 localSf.push_back(coarseSf[facei]);
214 localAgg.push_back(agglomi);
215 }
216 }
217 }
218
219 Info<< "\nAssembled coarse patch data" << endl;
220
221 // Distribute local coarse Cf and Sf for shooting rays
222 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
223
224 allCf_[Pstream::myProcNo()].transfer(localCf);
225 allSf_[Pstream::myProcNo()].transfer(localSf);
226 allAgg_[Pstream::myProcNo()].transfer(localAgg);
227
230 Pstream::allGatherList(allAgg_);
231
232 Pstream::listGather(patchAreas_, sumOp<scalar>());
233 Pstream::broadcast(patchAreas_);
234
235 globalNumbering_ = globalIndex(nCoarseFace_);
236}
237
238
240{
241 DynamicList<point> localCf(mesh_.nBoundaryFaces());
242 DynamicList<vector> localSf(mesh_.nBoundaryFaces());
243
244 const auto& pbm = mesh_.boundaryMesh();
245
246 for (const label patchi : patchIDs_)
247 {
248 localCf.push_back(pbm[patchi].faceCentres());
249 localSf.push_back(pbm[patchi].faceAreas());
250
251 patchAreas_[patchi] += sum(mesh_.magSf().boundaryField()[patchi]);
252 }
253
254 Info<< "\nAssembled patch data" << endl;
255
256 nFace_ = localCf.size();
257 nCoarseFace_ = -1;
258
259 allCf_[Pstream::myProcNo()].transfer(localCf);
260 allSf_[Pstream::myProcNo()].transfer(localSf);
261
264
265 // Pstream::listCombineGather(patchAreas_, plusEqOp<scalar>());
266 // Pstream::broadcast(patchAreas_);
267
268 // Basic type and op_sum, so can use listReduce (ie, mpiAllReduce)
269 Pstream::listReduce(patchAreas_, sumOp<scalar>());
270
271 globalNumbering_ = globalIndex(nFace_);
272}
273
274
276(
277 labelList& rayEndFace
278) const
279{
280 // Construct compact numbering
281 // - return map from remote to compact indices
282 // (per processor (!= myProcNo) a map from remote index to compact index)
283 // - construct distribute map
284 // - renumber rayEndFace into compact addressing
285
286 DebugInfo << "\nCreating map distribute" << endl;
287
288 List<Map<label>> compactMap(Pstream::nProcs());
289 mapPtr_.reset(new mapDistribute(globalNumbering_, rayEndFace, compactMap));
290
291 DebugInfo << "\nCreating compact-to-global addressing" << endl;
292
293 // Invert compactMap (from processor+localface to compact) to go
294 // from compact to processor+localface (expressed as a globalIndex)
295 compactToGlobal_.resize_nocopy(mapPtr_->constructSize());
296
297 // Local indices first
298 // Note: are not in compactMap
299 for (label i = 0; i < globalNumbering_.localSize(); ++i)
300 {
301 compactToGlobal_[i] = globalNumbering_.toGlobal(i);
302 }
303
304 forAll(compactMap, proci)
305 {
306 const Map<label>& localToCompactMap = compactMap[proci];
307
308 forAllConstIters(localToCompactMap, iter)
309 {
310 compactToGlobal_[*iter] =
311 globalNumbering_.toGlobal(proci, iter.key());
312 }
313 }
314}
315
316
317Foam::coordSystem::cartesian Foam::VF::raySearchEngine::createCoordSystem
318(
319 const point& origin,
320 const vector& dir
321) const
322{
323 vector axis(Zero);
324
325 for (direction d=0; d<3; ++d)
326 {
327 axis = dir^tensor::I.col(d);
328
329 // Remove empty direction for 2D
330 if (mesh_.nSolutionD() == 2)
331 {
332 meshTools::constrainDirection(mesh_, mesh_.solutionD(), axis);
333 }
334
335 if (magSqr(axis) > 0)
336 {
337 axis.normalise();
338 break;
339 }
340 }
341
342 return coordSystem::cartesian(origin, dir, axis);
343}
344
345
346Foam::tmp<Foam::pointField> Foam::VF::raySearchEngine::createHemiPoints
347(
348 const label nRayPerFace
349) const
350{
351 auto themiPts = tmp<pointField>::New(nRayPerFace);
352 auto& hemiPts = themiPts.ref();
353
354 const label nPoints = hemiPts.size();
355
356 if (mesh_.nSolutionD() == 3)
357 {
358 // Point in range -1 < x < 1; -1 < y < 1; 0 < z 1
359
360 forAll(hemiPts, pointi)
361 {
362 const scalar phi = Foam::acos(1 - (pointi + 0.5)/nPoints);
363 const scalar theta =
364 mathematical::pi*(1 + Foam::sqrt(5.0))*(pointi + 0.5);
365
366 hemiPts[pointi] =
367 vector
368 (
369 Foam::cos(theta)*Foam::sin(phi),
370 Foam::sin(theta)*Foam::sin(phi),
372 );
373 }
374 }
375 else if (mesh_.nSolutionD() == 2)
376 {
377 // Point in range -1 < x < 1; y = 0; 0 < z < 1; _\|/_
378
379 forAll(hemiPts, pointi)
380 {
381 const scalar theta = mathematical::pi*(pointi+0.5)/nPoints;
382 hemiPts[pointi] = vector(Foam::cos(theta), 0, Foam::sin(theta));
383 }
384 }
385
386 return themiPts;
387}
388
389
390// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
391
393(
394 const fvMesh& mesh,
395 const dictionary& dict
396)
397:
398 mesh_(mesh),
399 mapPtr_(nullptr),
400 compactToGlobal_(),
401 globalNumbering_(),
402 patchGroup_(dict.getOrDefault<word>("patchGroup", "viewFactorWall")),
403 patchIDs_(mesh_.boundaryMesh().indices(patchGroup_)),
404 patchAreas_(mesh_.boundaryMesh().nNonProcessor(), Zero),
405 agglomerate_(dict.get<bool>("agglomerate")),
406 agglomMeshPtr_(nullptr),
407 nFace_(-1),
408 nCoarseFace_(-1),
409 allCf_(UPstream::nProcs()),
410 allSf_(UPstream::nProcs()),
411 allAgg_(UPstream::nProcs())
412{
413 Info<< "\nParticipating patches:" << endl;
414
415 forAll(patchIDs_, i)
416 {
417 const label patchi = patchIDs_[i];
418 Info<< " " << i << ": " << mesh_.boundaryMesh()[patchi].name()
419 << endl;
420 }
421
422 const word agglomName(dict.getOrDefault<word>("agglom", "finalAgglom"));
423
424 IOobject agglomIO
425 (
426 agglomName,
427 mesh_.facesInstance(),
428 mesh_,
430 );
431
432
433 if (agglomerate_)
434 {
435 // Sets allCf_, allSf_, allAgg_ based on coarse mesh
436 // Sets nFace_, nCoarseFace_
437 createAgglomeration(agglomIO);
438 }
439 else
440 {
441 // Check for presence of finalAgglom - will cause problems in later
442 // calculations with viewFactor radiation model
443 if (agglomIO.typeHeaderOk<labelListIOList>())
444 {
446 << "Found agglomeration file: " << agglomIO.objectPath() << nl
447 << " This is inconsistent with the view factor calculation "
448 << "and should be removed" << nl << endl;
449 }
450
451 // Sets allCf_, allSf_ based on fine mesh
452 // Sets nFace_; nCoarseFace_ remains unset (-1)
453 createGeometry();
454 }
455
456 globalNumbering_ =
457 nCoarseFace_ == -1 ? globalIndex(nFace_) : globalIndex(nCoarseFace_);
458}
459
460
461// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
462
464(
465 labelListList& visibleFaceFaces
466) const
467{
468 labelList rayStartFace;
469 labelList rayEndFace;
470 shootRays(rayStartFace, rayEndFace);
471
472 const label nFace = nParticipatingFaces();
473
474 // Calculate number of visible faces from each local start face
475 labelList nVisibleFaceFaces(nFace, Zero);
476 for (const label facei : rayStartFace)
477 {
478 ++nVisibleFaceFaces[facei];
479 }
480
481 check(nVisibleFaceFaces);
482
483 createParallelAddressing(rayEndFace);
484
485 // Set visible face-faces
486
487 // visibleFaceFaces has:
488 // (local face, local viewed face) = compact viewed face
489 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
490
491 visibleFaceFaces.resize_nocopy(nFace);
492 forAll(nVisibleFaceFaces, facei)
493 {
494 visibleFaceFaces[facei].resize_nocopy(nVisibleFaceFaces[facei]);
495 }
496
497 nVisibleFaceFaces = 0;
498 forAll(rayStartFace, i)
499 {
500 const label facei = rayStartFace[i];
501 const label sloti = rayEndFace[i];
502 visibleFaceFaces[facei][nVisibleFaceFaces[facei]++] = sloti;
503 }
504}
505
506
508(
509 const mapDistribute& map,
510 pointField& compactCf,
511 vectorField& compactSf,
512 List<List<vector>>& compactFineSf,
513 List<List<point>>& compactFineCf,
514 DynamicList<List<point>>& compactPoints,
515 DynamicList<label>& compactPatchId
516) const
517{
518 compactCf.resize_nocopy(map.constructSize());
519 compactSf.resize_nocopy(map.constructSize());
520 compactFineSf.resize_nocopy(map.constructSize());
521 compactFineCf.resize_nocopy(map.constructSize());
522 compactPoints.setCapacity(map.constructSize());
523 compactPatchId.setCapacity(map.constructSize());
524
525 // Insert my local values area and centre values
526 if (agglomMeshPtr_)
527 {
528 SubList<vector>(compactSf, nCoarseFace_) = allSf_[Pstream::myProcNo()];
529 SubList<point>(compactCf, nCoarseFace_) = allCf_[Pstream::myProcNo()];
530
531 const auto& coarseMesh = agglomMeshPtr_();
532 const auto& coarsePatches = coarseMesh.boundaryMesh();
533 const auto& coarseFaces = coarseMesh.faces();
534 const auto& coarsePoints = coarseMesh.points();
535
536 const auto& finalAgglom = coarseMesh.patchFaceAgglomeration();
537
538 // Insert my fine local values per coarse face
539 label sloti = 0;
540 for (const label patchi : patchIDs_)
541 {
542 const auto& fineCfp = mesh_.Cf().boundaryField()[patchi];
543 const auto& fineSfp = mesh_.Sf().boundaryField()[patchi];
544 const labelList& agglom = finalAgglom[patchi];
545
546 if (agglom.empty()) continue;
547
548 const label nAgglom = max(agglom) + 1;
549 const labelListList coarseToFine = invertOneToMany(nAgglom, agglom);
550 const labelList& coarsePatchFace =
551 coarseMesh.patchFaceMap()[patchi];
552
553 const label coarseStart = coarsePatches[patchi].start();
554
555 forAll(coarseToFine, coarsei)
556 {
557 compactPatchId.push_back(patchi);
558
559 const vectorField localPoints
560 (
561 coarsePoints,
562 coarseFaces[coarseStart + coarsei]
563 );
564 compactPoints.push_back(localPoints);
565
566 const label coarseFacei = coarsePatchFace[coarsei];
567 const labelList& fineFaces = coarseToFine[coarseFacei];
568
569 List<point>& fineCf = compactFineCf[sloti];
570 fineCf.resize_nocopy(fineFaces.size());
571 fineCf = UIndirectList<point>(fineCfp, fineFaces);
572
573 List<vector>& fineSf = compactFineSf[sloti];
574 fineSf.resize_nocopy(fineFaces.size());
575 fineSf = UIndirectList<vector>(fineSfp, fineFaces);
576
577 ++sloti;
578 }
579 }
580 }
581 else
582 {
583 SubList<vector>(compactSf, nFace_) = allSf_[Pstream::myProcNo()];
584 SubList<point>(compactCf, nFace_) = allCf_[Pstream::myProcNo()];
585
586 const auto& patches = mesh_.boundaryMesh();
587 const faceList& faces = mesh_.faces();
588
589 label sloti = 0;
590
591 for (const label patchi : patchIDs_)
592 {
593 const auto& Sfp = mesh_.Sf().boundaryField()[patchi];
594 const auto& Cfp = mesh_.Cf().boundaryField()[patchi];
595
596 const polyPatch& pp = patches[patchi];
597
598 forAll(pp, facei)
599 {
600 compactPatchId.push_back(patchi);
601
602 const auto& fpts = faces[facei + pp.start()];
603 compactPoints.push_back(List<point>(mesh_.points(), fpts));
604
605 compactFineCf[sloti] = List<point>({Cfp[facei]});
606 compactFineSf[sloti] = List<vector>({Sfp[facei]});
607 ++sloti;
608 }
609 }
610 }
611
612
613 // Do all swapping
614 map.distribute(compactSf);
615 map.distribute(compactCf);
616 map.distribute(compactFineCf);
617 map.distribute(compactFineSf);
618 map.distribute(compactPoints);
619 map.distribute(compactPatchId);
620}
621
622
623// ************************************************************************* //
label n
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
const polyBoundaryMesh & pbm
void push_back(const T &val)
Copy append an element to the end of this list.
void setCapacity(const label len)
Alter the size of the underlying storage.
@ NO_READ
Nothing to be read.
@ MUST_READ
Reading required.
@ NO_WRITE
Ignore writing from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition IOobjectI.H:50
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition ListI.H:171
void push_back(const T &val)
Append an element at the end of the list.
Definition ListI.H:221
static void listGather(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather (reduce) list elements, applying bop to each list element.
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual communication.
static void listReduce(UList< T > &values, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Reduce list elements (list must be equal size on all ranks), applying bop to each list element.
static const Tensor I
Definition Tensor.H:81
bool empty() const noexcept
True if List is empty (ie, size() is zero).
Definition UList.H:701
void size(const label n)
Older name for setAddressableSize.
Definition UList.H:118
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Negative if the process is not a...
Definition UPstream.H:1706
static label nProcs(const label communicator=worldComm)
Number of ranks in parallel run (for given communicator). It is 1 for serial run.
Definition UPstream.H:1697
@ broadcast
broadcast [MPI]
Definition UPstream.H:189
label size() const noexcept
The number of entries in the list.
Definition UPtrListI.H:106
Base class for ray search engines.
void createGeometry()
Create patch geometry based on the original mesh.
static const label maxDynListLength
label nParticipatingFaces() const
Number of participating faces.
virtual void shootRays(labelList &rayStartFaceOut, labelList &rayEndFaceOut) const =0
Shoot rays; returns lists of ray start and end faces.
static label closestPointIndex(const point &p0, const List< point > &pts)
void createParallelAddressing(labelList &rayEndFace) const
Create parallel addressing - map, compact-to-global.
coordSystem::cartesian createCoordSystem(const point &origin, const vector &dir) const
Create Cartesian co-ordinate system.
void createAgglomeration(const IOobject &io)
Create patch geometry based on the agglomerated mesh.
virtual void correct(labelListList &visibleFaceFaces) const
Correct.
static void check(const labelList &nVisibleFaceFaces)
tmp< pointField > createHemiPoints(const label nRayPerFace) const
Create a set of points describing a hemisphere.
raySearchEngine(const raySearchEngine &)=delete
No copy construct.
void compactAddressing(const mapDistribute &map, pointField &compactCf, vectorField &compactSf, List< List< vector > > &compactFineSf, List< List< point > > &compactFineCf, DynamicList< List< point > > &compactPoints, DynamicList< label > &compactPatchId) const
Create compact addressing.
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
label constructSize() const noexcept
Constructed data size.
void distribute(List< T > &fld, const bool dummyTransform=true, const int tag=UPstream::msgType()) const
Distribute List data using default commsType, default flip/negate operator.
const faceList::subList faces() const
Return mesh faces for the entire boundary.
label start() const noexcept
The start label of boundary faces in the polyMesh face list.
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
const volScalarField & p0
Definition EEqn.H:36
const polyBoundaryMesh & patches
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
const auto & io
label nPoints
#define DebugInfo
Report an information message using Foam::Info.
#define WarningInFunction
Report a warning using Foam::Warning.
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
A namespace for various viewFactor model implementations.
Different types of constants.
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition meshTools.C:680
Namespace for OpenFOAM.
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
constexpr label labelMin
Definition label.H:54
dimensionedScalar sin(const dimensionedScalar &ds)
messageStream Info
Information stream (stdout output on master, null elsewhere).
List< face > faceList
List of faces.
Definition faceListFwd.H:41
constexpr label labelMax
Definition label.H:55
T returnReduce(const T &value, BinaryOp bop, const int tag=UPstream::msgType(), const int communicator=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
dimensionedScalar sqrt(const dimensionedScalar &ds)
IOList< labelList > labelListIOList
IO for a List of labelList.
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
Field< vector > vectorField
Specialisation of Field<T> for vector.
uint8_t direction
Definition direction.H:49
vector point
Point is a vector.
Definition point.H:37
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1, const label comm)
static constexpr const zero Zero
Global zero (0).
Definition zero.H:127
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vectorField pointField
pointField is a vectorField.
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition ListOps.C:125
Vector< scalar > vector
Definition vector.H:57
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
dimensionedScalar acos(const dimensionedScalar &ds)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
#define defineRunTimeSelectionTable(baseType, argNames)
Define run-time selection table.
dictionary dict
const pointField & pts
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
Definition stdFoam.H:235
Foam::surfaceFields.