raySearchEngine Class Referenceabstract

Base class for ray search engines. More...

Inheritance diagram for raySearchEngine:
Collaboration diagram for raySearchEngine:

Public Member Functions

 TypeName ("raySearchEngine")
 Run-time type information. More...
 
 declareRunTimeSelectionTable (autoPtr, raySearchEngine, mesh,(const fvMesh &mesh, const dictionary &dict),(mesh, dict))
 Selection table. More...
 
 raySearchEngine (const raySearchEngine &)=delete
 No copy construct. More...
 
void operator= (const raySearchEngine &)=delete
 No copy assignment. More...
 
 raySearchEngine (const fvMesh &mesh, const dictionary &dict)
 Constructor. More...
 
virtual ~raySearchEngine ()=default
 Destructor. More...
 
const fvMeshmesh () const noexcept
 Reference to the mesh. More...
 
const mapDistributemap () const
 Parallel map. More...
 
const labelListcompactToGlobal () const noexcept
 Compact to global addressing. More...
 
const globalIndexglobalNumbering () const noexcept
 Global numbering. More...
 
const labelListpatchIDs () const noexcept
 List of participating patch IDs. More...
 
const scalarListpatchAreas () const noexcept
 Patch areas. More...
 
label nParticipatingFaces () const
 Number of participating faces. More...
 
const List< pointField > & allCf () const noexcept
 List of all face centres per processor. More...
 
const List< vectorField > & allSf () const noexcept
 List of all face areas per processor. More...
 
const List< labelField > & allAgg () const noexcept
 List of all face agglomeration index per processor. More...
 
virtual void shootRays (labelList &rayStartFaceOut, labelList &rayEndFaceOut) const =0
 Shoot rays; returns lists of ray start and end faces. More...
 
virtual void correct (labelListList &visibleFaceFaces) const
 Correct. More...
 
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. More...
 
template<class Type >
void interpolate (GeometricField< Type, fvPatchField, volMesh > &fld, const List< List< Type >> &values) const
 Interpolate field. More...
 

Static Public Member Functions

static autoPtr< raySearchEngineNew (const fvMesh &mesh, const dictionary &dict)
 Selector. More...
 

Static Public Attributes

static const label maxDynListLength
 

Protected Member Functions

void createGeometry ()
 Create patch geometry based on the original mesh. More...
 
void createParallelAddressing (labelList &rayEndFace) const
 Create parallel addressing - map, compact-to-global. More...
 
coordSystem::cartesian createCoordSystem (const point &origin, const vector &dir) const
 Create Cartesian co-ordinate system. More...
 
void createAgglomeration (const IOobject &io)
 Create patch geometry based on the agglomerated mesh. More...
 
tmp< pointFieldcreateHemiPoints (const label nRayPerFace) const
 Create a set of points describing a hemisphere. More...
 

Static Protected Member Functions

static void check (const labelList &nVisibleFaceFaces)
 
static label closestPointIndex (const point &p0, const List< point > &pts)
 

Protected Attributes

const fvMeshmesh_
 Reference to the mesh. More...
 
autoPtr< mapDistributemapPtr_
 Parallel map. More...
 
labelList compactToGlobal_
 Compact to global addressing. More...
 
globalIndex globalNumbering_
 Global numbering. More...
 
const word patchGroup_
 Name of patch group to identify participating patches. More...
 
labelList patchIDs_
 List of participating patch IDs. More...
 
scalarList patchAreas_
 Patch areas. More...
 
bool agglomerate_
 Agglomeration flag. More...
 
autoPtr< singleCellFvMeshagglomMeshPtr_
 Agglomerated mesh representation. More...
 
label nFace_
 Number of original faces. More...
 
label nCoarseFace_
 Number of coarse faces. More...
 
List< pointFieldallCf_
 List of all face centres per processor. More...
 
List< vectorFieldallSf_
 List of all face areas per processor. More...
 
List< labelFieldallAgg_
 List of all face agglomeration index per processor. More...
 

Detailed Description

Base class for ray search engines.

Participating patches must be in the viewFactorWall group, i.e. using the inGroups entry of the "\<case\>/polyMesh/boundary" file.

myPatch
{
    type            wall;
    inGroups        2(wall viewFactorWall);
    ...
}

Face agglomeration can be employed, created using the faceAgglomerate utility. The file name to be read can be user-defined:

// Name of agglomeration file; default = finalAgglom
agglom      finalAgglom;
Source files

Definition at line 72 of file raySearchEngine.H.

Constructor & Destructor Documentation

◆ raySearchEngine() [1/2]

raySearchEngine ( const raySearchEngine )
delete

No copy construct.

◆ raySearchEngine() [2/2]

raySearchEngine ( const fvMesh mesh,
const dictionary dict 
)

Constructor.

◆ ~raySearchEngine()

virtual ~raySearchEngine ( )
virtualdefault

Destructor.

Member Function Documentation

◆ check()

static void check ( const labelList nVisibleFaceFaces)
staticprotected

◆ closestPointIndex()

static label closestPointIndex ( const point p0,
const List< point > &  pts 
)
staticprotected

◆ createGeometry()

void createGeometry ( )
protected

Create patch geometry based on the original mesh.

◆ createParallelAddressing()

void createParallelAddressing ( labelList rayEndFace) const
protected

Create parallel addressing - map, compact-to-global.

◆ createCoordSystem()

coordSystem::cartesian createCoordSystem ( const point origin,
const vector dir 
) const
protected

Create Cartesian co-ordinate system.

◆ createAgglomeration()

void createAgglomeration ( const IOobject io)
protected

Create patch geometry based on the agglomerated mesh.

◆ createHemiPoints()

tmp<pointField> createHemiPoints ( const label  nRayPerFace) const
protected

Create a set of points describing a hemisphere.

Note: origin is (0 0 0)

◆ TypeName()

TypeName ( "raySearchEngine"  )

Run-time type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
raySearchEngine  ,
mesh  ,
(const fvMesh &mesh, const dictionary &dict ,
(mesh, dict  
)

Selection table.

◆ New()

static autoPtr<raySearchEngine> New ( const fvMesh mesh,
const dictionary dict 
)
static

Selector.

◆ operator=()

void operator= ( const raySearchEngine )
delete

No copy assignment.

◆ mesh()

const Foam::fvMesh & mesh ( ) const
inlinenoexcept

Reference to the mesh.

Definition at line 21 of file raySearchEngineI.H.

References raySearchEngine::mesh_.

◆ map()

const Foam::mapDistribute & map ( ) const
inline

Parallel map.

Definition at line 27 of file raySearchEngineI.H.

References Foam::abort(), Foam::FatalError, and FatalErrorInFunction.

Here is the call graph for this function:

◆ compactToGlobal()

const Foam::labelList & compactToGlobal ( ) const
inlinenoexcept

Compact to global addressing.

Definition at line 41 of file raySearchEngineI.H.

◆ globalNumbering()

const Foam::globalIndex & globalNumbering ( ) const
inlinenoexcept

Global numbering.

Definition at line 48 of file raySearchEngineI.H.

◆ patchIDs()

const Foam::labelList & patchIDs ( ) const
inlinenoexcept

List of participating patch IDs.

Definition at line 54 of file raySearchEngineI.H.

◆ patchAreas()

const Foam::scalarList & patchAreas ( ) const
inlinenoexcept

Patch areas.

Definition at line 60 of file raySearchEngineI.H.

◆ nParticipatingFaces()

Foam::label nParticipatingFaces ( ) const
inline

Number of participating faces.

Definition at line 66 of file raySearchEngineI.H.

◆ allCf()

const Foam::List< Foam::pointField > & allCf ( ) const
inlinenoexcept

List of all face centres per processor.

Definition at line 74 of file raySearchEngineI.H.

◆ allSf()

const Foam::List< Foam::vectorField > & allSf ( ) const
inlinenoexcept

List of all face areas per processor.

Definition at line 81 of file raySearchEngineI.H.

◆ allAgg()

const Foam::List< Foam::labelField > & allAgg ( ) const
inlinenoexcept

List of all face agglomeration index per processor.

Definition at line 88 of file raySearchEngineI.H.

◆ shootRays()

virtual void shootRays ( labelList rayStartFaceOut,
labelList rayEndFaceOut 
) const
pure virtual

Shoot rays; returns lists of ray start and end faces.

Implemented in voxel.

◆ correct()

virtual void correct ( labelListList visibleFaceFaces) const
virtual

Correct.

◆ compactAddressing()

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.

◆ interpolate()

void interpolate ( GeometricField< Type, fvPatchField, volMesh > &  fld,
const List< List< Type >> &  values 
) const

Interpolate field.

Member Data Documentation

◆ mesh_

const fvMesh& mesh_
protected

Reference to the mesh.

Definition at line 81 of file raySearchEngine.H.

Referenced by raySearchEngine::mesh().

◆ mapPtr_

autoPtr<mapDistribute> mapPtr_
mutableprotected

Parallel map.

Definition at line 86 of file raySearchEngine.H.

◆ compactToGlobal_

labelList compactToGlobal_
mutableprotected

Compact to global addressing.

Definition at line 91 of file raySearchEngine.H.

◆ globalNumbering_

globalIndex globalNumbering_
protected

Global numbering.

Definition at line 96 of file raySearchEngine.H.

◆ patchGroup_

const word patchGroup_
protected

Name of patch group to identify participating patches.

Definition at line 101 of file raySearchEngine.H.

◆ patchIDs_

labelList patchIDs_
protected

List of participating patch IDs.

Definition at line 106 of file raySearchEngine.H.

◆ patchAreas_

scalarList patchAreas_
protected

Patch areas.

Definition at line 111 of file raySearchEngine.H.

◆ agglomerate_

bool agglomerate_
protected

Agglomeration flag.

Definition at line 116 of file raySearchEngine.H.

◆ agglomMeshPtr_

autoPtr<singleCellFvMesh> agglomMeshPtr_
protected

Agglomerated mesh representation.

Definition at line 121 of file raySearchEngine.H.

◆ nFace_

label nFace_
protected

Number of original faces.

Definition at line 126 of file raySearchEngine.H.

◆ nCoarseFace_

label nCoarseFace_
protected

Number of coarse faces.

Definition at line 131 of file raySearchEngine.H.

◆ allCf_

List<pointField> allCf_
protected

List of all face centres per processor.

Definition at line 136 of file raySearchEngine.H.

◆ allSf_

List<vectorField> allSf_
protected

List of all face areas per processor.

Definition at line 141 of file raySearchEngine.H.

◆ allAgg_

List<labelField> allAgg_
protected

List of all face agglomeration index per processor.

Definition at line 146 of file raySearchEngine.H.

◆ maxDynListLength

const label maxDynListLength
static

Definition at line 193 of file raySearchEngine.H.


The documentation for this class was generated from the following files: