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 -------------------------------------------------------------------------------
10 License
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 
33 using namespace Foam::constant;
34 
35 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36 
37 namespace Foam
38 {
39 namespace VF
40 {
41  defineTypeNameAndDebug(raySearchEngine, 0);
42  defineRunTimeSelectionTable(raySearchEngine, mesh);
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 
228  Pstream::allGatherList(allCf_);
229  Pstream::allGatherList(allSf_);
230  Pstream::allGatherList(allAgg_);
231 
232  Pstream::listCombineGather(patchAreas_, plusEqOp<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 
262  Pstream::allGatherList(allCf_);
263  Pstream::allGatherList(allSf_);
264 
265  Pstream::listCombineGather(patchAreas_, plusEqOp<scalar>());
266  Pstream::broadcast(patchAreas_);
267 
268  globalNumbering_ = globalIndex(nFace_);
269 }
270 
271 
273 (
274  labelList& rayEndFace
275 ) const
276 {
277  // Construct compact numbering
278  // - return map from remote to compact indices
279  // (per processor (!= myProcNo) a map from remote index to compact index)
280  // - construct distribute map
281  // - renumber rayEndFace into compact addressing
282 
283  DebugInfo << "\nCreating map distribute" << endl;
284 
285  List<Map<label>> compactMap(Pstream::nProcs());
286  mapPtr_.reset(new mapDistribute(globalNumbering_, rayEndFace, compactMap));
287 
288  DebugInfo << "\nCreating compact-to-global addressing" << endl;
289 
290  // Invert compactMap (from processor+localface to compact) to go
291  // from compact to processor+localface (expressed as a globalIndex)
292  compactToGlobal_.resize_nocopy(mapPtr_->constructSize());
293 
294  // Local indices first
295  // Note: are not in compactMap
296  for (label i = 0; i < globalNumbering_.localSize(); ++i)
297  {
298  compactToGlobal_[i] = globalNumbering_.toGlobal(i);
299  }
300 
301  forAll(compactMap, proci)
302  {
303  const Map<label>& localToCompactMap = compactMap[proci];
304 
305  forAllConstIters(localToCompactMap, iter)
306  {
307  compactToGlobal_[*iter] =
308  globalNumbering_.toGlobal(proci, iter.key());
309  }
310  }
311 }
312 
313 
315 (
316  const point& origin,
317  const vector& dir
318 ) const
319 {
320  vector axis(Zero);
321 
322  for (direction d=0; d<3; ++d)
323  {
324  axis = dir^tensor::I.col(d);
325 
326  // Remove empty direction for 2D
327  if (mesh_.nSolutionD() == 2)
328  {
329  meshTools::constrainDirection(mesh_, mesh_.solutionD(), axis);
330  }
331 
332  if (magSqr(axis) > 0)
333  {
334  axis.normalise();
335  break;
336  }
337  }
338 
339  return coordSystem::cartesian(origin, dir, axis);
340 }
341 
342 
344 (
345  const label nRayPerFace
346 ) const
347 {
348  auto themiPts = tmp<pointField>::New(nRayPerFace);
349  auto& hemiPts = themiPts.ref();
350 
351  const label nPoints = hemiPts.size();
352 
353  if (mesh_.nSolutionD() == 3)
354  {
355  // Point in range -1 < x < 1; -1 < y < 1; 0 < z 1
356 
357  forAll(hemiPts, pointi)
358  {
359  const scalar phi = Foam::acos(1 - (pointi + 0.5)/nPoints);
360  const scalar theta =
361  mathematical::pi*(1 + Foam::sqrt(5.0))*(pointi + 0.5);
362 
363  hemiPts[pointi] =
364  vector
365  (
366  Foam::cos(theta)*Foam::sin(phi),
367  Foam::sin(theta)*Foam::sin(phi),
368  Foam::cos(phi)
369  );
370  }
371  }
372  else if (mesh_.nSolutionD() == 2)
373  {
374  // Point in range -1 < x < 1; y = 0; 0 < z < 1; _\|/_
375 
376  forAll(hemiPts, pointi)
377  {
378  const scalar theta = mathematical::pi*(pointi+0.5)/nPoints;
379  hemiPts[pointi] = vector(Foam::cos(theta), 0, Foam::sin(theta));
380  }
381  }
382 
383  return themiPts;
384 }
385 
386 
387 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
388 
390 (
391  const fvMesh& mesh,
392  const dictionary& dict
393 )
394 :
395  mesh_(mesh),
396  mapPtr_(nullptr),
397  compactToGlobal_(),
398  globalNumbering_(),
399  patchGroup_(dict.getOrDefault<word>("patchGroup", "viewFactorWall")),
400  patchIDs_(mesh_.boundaryMesh().indices(patchGroup_)),
401  patchAreas_(mesh_.boundaryMesh().nNonProcessor(), Zero),
402  agglomerate_(dict.get<bool>("agglomerate")),
403  agglomMeshPtr_(nullptr),
404  nFace_(-1),
405  nCoarseFace_(-1),
406  allCf_(UPstream::nProcs()),
407  allSf_(UPstream::nProcs()),
408  allAgg_(UPstream::nProcs())
409 {
410  Info<< "\nParticipating patches:" << endl;
411 
412  forAll(patchIDs_, i)
413  {
414  const label patchi = patchIDs_[i];
415  Info<< " " << i << ": " << mesh_.boundaryMesh()[patchi].name()
416  << endl;
417  }
418 
419  const word agglomName(dict.getOrDefault<word>("agglom", "finalAgglom"));
420 
421  IOobject agglomIO
422  (
423  agglomName,
424  mesh_.facesInstance(),
425  mesh_,
426  IOobject::MUST_READ
427  );
428 
429 
430  if (agglomerate_)
431  {
432  // Sets allCf_, allSf_, allAgg_ based on coarse mesh
433  // Sets nFace_, nCoarseFace_
434  createAgglomeration(agglomIO);
435  }
436  else
437  {
438  // Check for presence of finalAgglom - will cause problems in later
439  // calculations with viewFactor radiation model
440  if (agglomIO.typeHeaderOk<labelListIOList>())
441  {
443  << "Found agglomeration file: " << agglomIO.objectPath() << nl
444  << " This is inconsistent with the view factor calculation "
445  << "and should be removed" << nl << endl;
446  }
447 
448  // Sets allCf_, allSf_ based on fine mesh
449  // Sets nFace_; nCoarseFace_ remains unset (-1)
450  createGeometry();
451  }
452 
453  globalNumbering_ =
454  nCoarseFace_ == -1 ? globalIndex(nFace_) : globalIndex(nCoarseFace_);
455 }
456 
457 
458 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
459 
461 (
462  labelListList& visibleFaceFaces
463 ) const
464 {
465  labelList rayStartFace;
466  labelList rayEndFace;
467  shootRays(rayStartFace, rayEndFace);
468 
469  const label nFace = nParticipatingFaces();
470 
471  // Calculate number of visible faces from each local start face
472  labelList nVisibleFaceFaces(nFace, Zero);
473  for (const label facei : rayStartFace)
474  {
475  ++nVisibleFaceFaces[facei];
476  }
477 
478  check(nVisibleFaceFaces);
479 
480  createParallelAddressing(rayEndFace);
481 
482  // Set visible face-faces
483 
484  // visibleFaceFaces has:
485  // (local face, local viewed face) = compact viewed face
486  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
487 
488  visibleFaceFaces.resize_nocopy(nFace);
489  forAll(nVisibleFaceFaces, facei)
490  {
491  visibleFaceFaces[facei].resize_nocopy(nVisibleFaceFaces[facei]);
492  }
493 
494  nVisibleFaceFaces = 0;
495  forAll(rayStartFace, i)
496  {
497  const label facei = rayStartFace[i];
498  const label sloti = rayEndFace[i];
499  visibleFaceFaces[facei][nVisibleFaceFaces[facei]++] = sloti;
500  }
501 }
502 
503 
505 (
506  const mapDistribute& map,
507  pointField& compactCf,
508  vectorField& compactSf,
509  List<List<vector>>& compactFineSf,
510  List<List<point>>& compactFineCf,
511  DynamicList<List<point>>& compactPoints,
512  DynamicList<label>& compactPatchId
513 ) const
514 {
515  compactCf.resize_nocopy(map.constructSize());
516  compactSf.resize_nocopy(map.constructSize());
517  compactFineSf.resize_nocopy(map.constructSize());
518  compactFineCf.resize_nocopy(map.constructSize());
519  compactPoints.setCapacity(map.constructSize());
520  compactPatchId.setCapacity(map.constructSize());
521 
522  // Insert my local values area and centre values
523  if (agglomMeshPtr_)
524  {
525  SubList<vector>(compactSf, nCoarseFace_) = allSf_[Pstream::myProcNo()];
526  SubList<point>(compactCf, nCoarseFace_) = allCf_[Pstream::myProcNo()];
527 
528  const auto& coarseMesh = agglomMeshPtr_();
529  const auto& coarsePatches = coarseMesh.boundaryMesh();
530  const auto& coarseFaces = coarseMesh.faces();
531  const auto& coarsePoints = coarseMesh.points();
532 
533  const auto& finalAgglom = coarseMesh.patchFaceAgglomeration();
534 
535  // Insert my fine local values per coarse face
536  label sloti = 0;
537  for (const label patchi : patchIDs_)
538  {
539  const auto& fineCfp = mesh_.Cf().boundaryField()[patchi];
540  const auto& fineSfp = mesh_.Sf().boundaryField()[patchi];
541  const labelList& agglom = finalAgglom[patchi];
542 
543  if (agglom.empty()) continue;
544 
545  const label nAgglom = max(agglom) + 1;
546  const labelListList coarseToFine = invertOneToMany(nAgglom, agglom);
547  const labelList& coarsePatchFace =
548  coarseMesh.patchFaceMap()[patchi];
549 
550  const label coarseStart = coarsePatches[patchi].start();
551 
552  forAll(coarseToFine, coarsei)
553  {
554  compactPatchId.push_back(patchi);
555 
556  const vectorField localPoints
557  (
558  coarsePoints,
559  coarseFaces[coarseStart + coarsei]
560  );
561  compactPoints.push_back(localPoints);
562 
563  const label coarseFacei = coarsePatchFace[coarsei];
564  const labelList& fineFaces = coarseToFine[coarseFacei];
565 
566  List<point>& fineCf = compactFineCf[sloti];
567  fineCf.resize_nocopy(fineFaces.size());
568  fineCf = UIndirectList<point>(fineCfp, fineFaces);
569 
570  List<vector>& fineSf = compactFineSf[sloti];
571  fineSf.resize_nocopy(fineFaces.size());
572  fineSf = UIndirectList<vector>(fineSfp, fineFaces);
573 
574  ++sloti;
575  }
576  }
577  }
578  else
579  {
580  SubList<vector>(compactSf, nFace_) = allSf_[Pstream::myProcNo()];
581  SubList<point>(compactCf, nFace_) = allCf_[Pstream::myProcNo()];
582 
583  const auto& patches = mesh_.boundaryMesh();
584  const faceList& faces = mesh_.faces();
585 
586  label sloti = 0;
587 
588  for (const label patchi : patchIDs_)
589  {
590  const auto& Sfp = mesh_.Sf().boundaryField()[patchi];
591  const auto& Cfp = mesh_.Cf().boundaryField()[patchi];
592 
593  const polyPatch& pp = patches[patchi];
594 
595  forAll(pp, facei)
596  {
597  compactPatchId.push_back(patchi);
598 
599  const auto& fpts = faces[facei + pp.start()];
600  compactPoints.push_back(List<point>(mesh_.points(), fpts));
601 
602  compactFineCf[sloti] = List<point>({Cfp[facei]});
603  compactFineSf[sloti] = List<vector>({Sfp[facei]});
604  ++sloti;
605  }
606  }
607  }
608 
609 
610  // Do all swapping
611  map.distribute(compactSf);
612  map.distribute(compactCf);
613  map.distribute(compactFineCf);
614  map.distribute(compactFineSf);
615  map.distribute(compactPoints);
616  map.distribute(compactPatchId);
617 }
618 
619 
620 // ************************************************************************* //
Different types of constants.
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
Foam::surfaceFields.
static void listCombineGather(UList< T > &values, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Combines List elements.
static label closestPointIndex(const point &p0, const List< point > &pts)
coordSystem::cartesian createCoordSystem(const point &origin, const vector &dir) const
Create Cartesian co-ordinate system.
const polyBoundaryMesh & pbm
dictionary dict
dimensionedScalar acos(const dimensionedScalar &ds)
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &f1)
uint8_t direction
Definition: direction.H:46
const label maxDynListLength(viewFactorDict.getOrDefault< label >("maxDynListLength", 1000000))
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:125
void createGeometry()
Create patch geometry based on the original mesh.
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:608
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
void resize_nocopy(const label len)
Adjust allocated size of list without necessarily.
Definition: ListI.H:168
constexpr label labelMin
Definition: label.H:54
Ignore writing from objectRegistry::writeObject()
static int myProcNo(const label communicator=worldComm)
Rank of this process in the communicator (starting from masterNo()). Can be negative if the process i...
Definition: UPstream.H:1086
T returnReduce(const T &value, const BinaryOp &bop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Perform reduction on a copy, using specified binary operation.
List< labelList > labelListList
List of labelList.
Definition: labelList.H:38
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.
static void broadcast(Type &value, const label comm=UPstream::worldComm)
Broadcast content (contiguous or non-contiguous) to all communicator ranks. Does nothing in non-paral...
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
void createAgglomeration(const IOobject &io)
Create patch geometry based on the agglomerated mesh.
A Cartesian coordinate system.
Definition: cartesianCS.H:65
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
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:1077
vectorField pointField
pointField is a vectorField.
Definition: pointFieldFwd.H:38
dynamicFvMesh & mesh
dimensionedScalar cos(const dimensionedScalar &ds)
labelListList invertOneToMany(const label len, const labelUList &map)
Invert one-to-many map. Unmapped elements will be size 0.
Definition: ListOps.C:126
static const Identity< scalar > I
Definition: Identity.H:100
tmp< pointField > createHemiPoints(const label nRayPerFace) const
Create a set of points describing a hemisphere.
label nPoints
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
constexpr scalar pi(M_PI)
static void check(const labelList &nVisibleFaceFaces)
Vector< scalar > vector
Definition: vector.H:57
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
#define DebugInfo
Report an information message using Foam::Info.
dimensionedScalar sin(const dimensionedScalar &ds)
virtual void correct(labelListList &visibleFaceFaces) const
Correct.
defineRunTimeSelectionTable(reactionRateFlameArea, dictionary)
defineTypeNameAndDebug(combustionModel, 0)
static void allGatherList(UList< T > &values, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Gather data, but keep individual values separate. Uses MPI_Allgather or manual linear/tree communicat...
label nParticipatingFaces() const
Number of participating faces.
PrimitivePatch< UIndirectList< face >, const pointField & > uindirectPrimitivePatch
A PrimitivePatch with UIndirectList for the faces, const reference for the point field.
vector point
Point is a vector.
Definition: point.H:37
#define WarningInFunction
Report a warning using Foam::Warning.
raySearchEngine(const raySearchEngine &)=delete
No copy construct.
const polyBoundaryMesh & patches
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
constexpr label labelMax
Definition: label.H:55
label n
tmp< pointField > allPoints(const Triangulation &t)
Extract all points in vertex-index order.
Field< vector > vectorField
Specialisation of Field<T> for vector.
void constrainDirection(const polyMesh &mesh, const Vector< label > &dirs, vector &d)
Set the constrained components of directions/velocity to zero.
Definition: meshTools.C:680
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
A class for managing temporary objects.
Definition: HashPtrTable.H:50
const volScalarField & p0
Definition: EEqn.H:36
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
virtual void shootRays(labelList &rayStartFaceOut, labelList &rayEndFaceOut) const =0
Shoot rays; returns lists of ray start and end faces.
void createParallelAddressing(labelList &rayEndFace) const
Create parallel addressing - map, compact-to-global.
const pointField & pts
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
IOList< labelList > labelListIOList
IO for a List of labelList.