Loading...
Searching...
No Matches
raySearchEngine.H
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
26Class
27 Foam::VF::raySearchEngine
28
29Description
30 Base class for ray search engines
31
32 Participating patches must be in the \c viewFactorWall group, i.e. using the
33 \c inGroups entry of the "<case>/polyMesh/boundary" file.
34
35 \verbatim
36 myPatch
37 {
38 type wall;
39 inGroups 2(wall viewFactorWall);
40 ...
41 }
42 \endverbatim
43
44 Face agglomeration can be employed, created using the \c faceAgglomerate
45 utility. The file name to be read can be user-defined:
46
47 \verbatim
48 // Name of agglomeration file; default = finalAgglom
49 agglom finalAgglom;
50 \endverbatim
51
52SourceFiles
53 raySearchEngine.C
54
55\*---------------------------------------------------------------------------*/
56
57#ifndef Foam_vf_raySearchEngine_H
58#define Foam_vf_raySearchEngine_H
59
60#include "cartesianCS.H"
61#include "mapDistribute.H"
62#include "singleCellFvMesh.H"
64
65// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
66
67namespace Foam
68{
69
70namespace VF
71{
73/*---------------------------------------------------------------------------*\
74 Class raySearchEngine Declaration
75\*---------------------------------------------------------------------------*/
76
78{
79protected:
80
81 // Protected Data
82
83 //- Reference to the mesh
84 const fvMesh& mesh_;
85
86 //- Parallel map
88
89 //- Compact to global addressing
92 //- Global numbering
94
95 //- Name of patch group to identify participating patches
97
98 //- List of participating patch IDs
100
101 //- Patch areas
103
104 //- Agglomeration flag
105 bool agglomerate_;
107 //- Agglomerated mesh representation
109
110 //- Number of original faces
111 label nFace_;
112
113 //- Number of coarse faces
114 label nCoarseFace_;
115
116 //- List of all face centres per processor
118
119 //- List of all face areas per processor
122 //- List of all face agglomeration index per processor
124
125
126 // Protected Member Functions
127
128 static void check(const labelList& nVisibleFaceFaces);
129
130 static label closestPointIndex
132 const point& p0,
133 const List<point>& pts
134 );
135
136 //- Create patch geometry based on the original mesh
137 void createGeometry();
138
139 //- Create parallel addressing - map, compact-to-global
140 void createParallelAddressing(labelList& rayEndFace) const;
142 //- Create Cartesian co-ordinate system
144 (
145 const point& origin,
146 const vector& dir
147 ) const;
148
149 //- Create patch geometry based on the agglomerated mesh
150 void createAgglomeration(const IOobject& io);
152 //- Create a set of points describing a hemisphere
153 // Note: origin is (0 0 0)
154 tmp<pointField> createHemiPoints(const label nRayPerFace) const;
155
156
157public:
158
159 static const label maxDynListLength;
160
161 //- Run-time type information
162 TypeName("raySearchEngine");
163
164 //- Selection table
166 (
169 mesh,
170 (
171 const fvMesh& mesh,
173 ),
174 (mesh, dict)
175 );
176
177 //- Selector
179 (
180 const fvMesh& mesh,
182 );
183
184
185 // Generated Methods
186
187 //- No copy construct
189
190 //- No copy assignment
191 void operator=(const raySearchEngine&) = delete;
192
194 //- Constructor
195 raySearchEngine(const fvMesh& mesh, const dictionary& dict);
196
197 //- Destructor
198 virtual ~raySearchEngine() = default;
199
200
201 // Public Member Functions
202
203 // Access
204
205 //- Reference to the mesh
206 inline const fvMesh& mesh() const noexcept;
207
208 //- Parallel map
209 inline const mapDistribute& map() const;
210
211 //- Compact to global addressing
212 inline const labelList& compactToGlobal() const noexcept;
213
214 //- Global numbering
215 inline const globalIndex& globalNumbering() const noexcept;
216
217 //- List of participating patch IDs
218 inline const labelList& patchIDs() const noexcept;
219
220 //- Patch areas
221 inline const scalarList& patchAreas() const noexcept;
222
223 //- Number of participating faces
224 inline label nParticipatingFaces() const;
225
226 //- List of all face centres per processor
227 inline const List<pointField>& allCf() const noexcept;
228
229 //- List of all face areas per processor
230 inline const List<vectorField>& allSf() const noexcept;
231
232 //- List of all face agglomeration index per processor
233 inline const List<labelField>& allAgg() const noexcept;
234
236 // Main calculation functions
237
238 //- Shoot rays; returns lists of ray start and end faces
239 virtual void shootRays
240 (
241 labelList& rayStartFaceOut,
242 labelList& rayEndFaceOut
243 ) const = 0;
244
245 //- Correct
246 virtual void correct(labelListList& visibleFaceFaces) const;
247
248 //- Create compact addressing
250 (
251 const mapDistribute& map,
252 pointField& compactCf,
253 vectorField& compactSf,
254 List<List<vector>>& compactFineSf,
255 List<List<point>>& compactFineCf,
256 DynamicList<List<point>>& compactPoints,
257 DynamicList<label>& compactPatchId
258 ) const;
259
260 //- Interpolate field
261 template<class Type>
262 void interpolate
263 (
265 const List<List<Type>>& values
266 ) const;
267};
268
269
270// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
271
272} // End namespace VF
273} // End namespace Foam
274
275// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
276
277#include "raySearchEngineI.H"
278
279// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
280
281#ifdef NoRepository
283#endif
284
285// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286
287#endif
288
289// ************************************************************************* //
Info<< nl;Info<< "Write faMesh in vtk format:"<< nl;{ vtk::uindirectPatchWriter writer(aMesh.patch(), fileName(aMesh.time().globalPath()/vtkBaseFileName));writer.writeGeometry();globalIndex procAddr(aMesh.nFaces());labelList cellIDs;if(UPstream::master()) { cellIDs.resize(procAddr.totalSize());for(const labelRange &range :procAddr.ranges()) { auto slice=cellIDs.slice(range);slice=identity(range);} } writer.beginCellData(4);writer.writeProcIDs();writer.write("cellID", cellIDs);writer.write("area", aMesh.S().field());writer.write("normal", aMesh.faceAreaNormals());writer.beginPointData(1);writer.write("normal", aMesh.pointAreaNormals());Info<< " "<< writer.output().name()<< nl;}{ vtk::lineWriter writer(aMesh.points(), aMesh.edges(), fileName(aMesh.time().globalPath()/(vtkBaseFileName+"-edges")));writer.writeGeometry();writer.beginCellData(4);writer.writeProcIDs();{ Field< scalar > fld(faMeshTools::flattenEdgeField(aMesh.magLe(), true))
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition DynamicList.H:68
Generic GeometricField class.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition IOobject.H:191
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
bool agglomerate_
Agglomeration flag.
void createGeometry()
Create patch geometry based on the original mesh.
static const label maxDynListLength
const mapDistribute & map() const
Parallel map.
label nParticipatingFaces() const
Number of participating faces.
TypeName("raySearchEngine")
Run-time type information.
globalIndex globalNumbering_
Global numbering.
const fvMesh & mesh_
Reference to the mesh.
const List< labelField > & allAgg() const noexcept
List of all face agglomeration index per processor.
label nFace_
Number of original faces.
virtual void shootRays(labelList &rayStartFaceOut, labelList &rayEndFaceOut) const =0
Shoot rays; returns lists of ray start and end faces.
static autoPtr< raySearchEngine > New(const fvMesh &mesh, const dictionary &dict)
Selector.
const globalIndex & globalNumbering() const noexcept
Global numbering.
static label closestPointIndex(const point &p0, const List< point > &pts)
List< pointField > allCf_
List of all face centres per processor.
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.
autoPtr< singleCellFvMesh > agglomMeshPtr_
Agglomerated mesh representation.
void interpolate(GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< Type > > &values) const
Interpolate field.
const labelList & compactToGlobal() const noexcept
Compact to global addressing.
void createAgglomeration(const IOobject &io)
Create patch geometry based on the agglomerated mesh.
virtual void correct(labelListList &visibleFaceFaces) const
Correct.
virtual ~raySearchEngine()=default
Destructor.
List< vectorField > allSf_
List of all face areas per processor.
raySearchEngine(const fvMesh &mesh, const dictionary &dict)
Constructor.
autoPtr< mapDistribute > mapPtr_
Parallel map.
labelList compactToGlobal_
Compact to global addressing.
static void check(const labelList &nVisibleFaceFaces)
tmp< pointField > createHemiPoints(const label nRayPerFace) const
Create a set of points describing a hemisphere.
scalarList patchAreas_
Patch areas.
const fvMesh & mesh() const noexcept
Reference to the mesh.
const List< pointField > & allCf() const noexcept
List of all face centres per processor.
declareRunTimeSelectionTable(autoPtr, raySearchEngine, mesh,(const fvMesh &mesh, const dictionary &dict),(mesh, dict))
Selection table.
label nCoarseFace_
Number of coarse faces.
const List< vectorField > & allSf() const noexcept
List of all face areas per processor.
const labelList & patchIDs() const noexcept
List of participating patch IDs.
List< labelField > allAgg_
List of all face agglomeration index per processor.
const word patchGroup_
Name of patch group to identify participating patches.
raySearchEngine(const raySearchEngine &)=delete
No copy construct.
const scalarList & patchAreas() const noexcept
Patch areas.
void operator=(const raySearchEngine &)=delete
No copy assignment.
labelList patchIDs_
List of participating patch IDs.
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.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition autoPtr.H:65
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
Mesh data needed to do the Finite Volume discretisation.
Definition fvMesh.H:85
Abstract base class with a fat-interface to all derived classes covering all possible ways in which t...
Calculates a non-overlapping list of offsets based on an input size (eg, number of cells) from differ...
Definition globalIndex.H:77
Class containing processor-to-processor mapping information.
A class for managing temporary objects.
Definition tmp.H:75
Mesh data needed to do the Finite Volume discretisation.
Definition volMesh.H:47
A class for handling words, derived from Foam::string.
Definition word.H:66
const volScalarField & p0
Definition EEqn.H:36
const auto & io
A namespace for various viewFactor model implementations.
Namespace for OpenFOAM.
List< labelList > labelListList
List of labelList.
Definition labelList.H:38
List< label > labelList
A List of labels.
Definition List.H:62
Field< vector > vectorField
Specialisation of Field<T> for vector.
Field< label > labelField
Specialisation of Field<T> for label.
Definition labelField.H:48
vector point
Point is a vector.
Definition point.H:37
const direction noexcept
Definition scalarImpl.H:265
vectorField pointField
pointField is a vectorField.
Vector< scalar > vector
Definition vector.H:57
List< scalar > scalarList
List of scalar.
Definition scalarList.H:32
Macros to ease declaration of run-time selection tables.
#define declareRunTimeSelectionTable(ptrWrapper, baseType, argNames, argList, parList)
Declare a run-time selection (variables and adder classes).
dictionary dict
const pointField & pts
#define TypeName(TypeNameString)
Declare a ClassName() with extra virtual type info.
Definition typeInfo.H:68