The SKIRT project
advanced radiative transfer for astrophysics
Classes | Public Types | Public Member Functions | Protected Member Functions | Private Types | Private Member Functions | Private Attributes | Friends | List of all members
TetraMeshSpatialGrid Class Reference

#include <TetraMeshSpatialGrid.hpp>

Inheritance diagram for TetraMeshSpatialGrid:
Inheritance graph
[legend]

Classes

struct  Face
 
class  Tetra
 

Public Types

enum class  Policy : int {
  Uniform , CentralPeak , DustDensity , ElectronDensity ,
  GasDensity , File
}
 

Public Member Functions

int cellIndex (Position bfr) const override
 
Position centralPositionInCell (int m) const override
 
std::unique_ptr< PathSegmentGeneratorcreatePathSegmentGenerator () const override
 
double diagonal (int m) const override
 
string filename () const
 
int numCells () const override
 
int numSamples () const
 
Policy policy () const
 
Position randomPositionInCell (int m) const override
 
bool refine () const
 
double volume (int m) const override
 
void writeGridPlotFiles (const SimulationItem *probe) const override
 
- Public Member Functions inherited from BoxSpatialGrid
Box boundingBox () const override
 
int dimension () const override
 
double maxX () const
 
double maxY () const
 
double maxZ () const
 
double minX () const
 
double minY () const
 
double minZ () const
 
virtual Box boundingBox () const =0
 
virtual int cellIndex (Position bfr) const =0
 
virtual Position centralPositionInCell (int m) const =0
 
virtual std::unique_ptr< PathSegmentGeneratorcreatePathSegmentGenerator () const =0
 
virtual double diagonal (int m) const =0
 
virtual int dimension () const =0
 
virtual int numCells () const =0
 
virtual Position randomPositionInCell (int m) const =0
 
virtual double volume (int m) const =0
 
virtual void writeGridPlotFiles (const SimulationItem *probe) const
 
- Public Member Functions inherited from SimulationItem
template<class T >
T * find (bool setup=true) const
 
template<class T >
T * interface (int levels=-999999, bool setup=true) const
 
virtual string itemName () const
 
void setup ()
 
string typeAndName () const
 
- Public Member Functions inherited from Item
 Item (const Item &)=delete
 
virtual ~Item ()
 
void addChild (Item *child)
 
const vector< Item * > & children () const
 
virtual void clearItemListProperty (const PropertyDef *property)
 
void destroyChild (Item *child)
 
virtual bool getBoolProperty (const PropertyDef *property) const
 
virtual vector< double > getDoubleListProperty (const PropertyDef *property) const
 
virtual double getDoubleProperty (const PropertyDef *property) const
 
virtual string getEnumProperty (const PropertyDef *property) const
 
virtual int getIntProperty (const PropertyDef *property) const
 
virtual vector< Item * > getItemListProperty (const PropertyDef *property) const
 
virtual ItemgetItemProperty (const PropertyDef *property) const
 
virtual string getStringProperty (const PropertyDef *property) const
 
int getUtilityProperty (string name) const
 
virtual void insertIntoItemListProperty (const PropertyDef *property, int index, Item *item)
 
Itemoperator= (const Item &)=delete
 
Itemparent () const
 
virtual void removeFromItemListProperty (const PropertyDef *property, int index)
 
virtual void setBoolProperty (const PropertyDef *property, bool value)
 
virtual void setDoubleListProperty (const PropertyDef *property, vector< double > value)
 
virtual void setDoubleProperty (const PropertyDef *property, double value)
 
virtual void setEnumProperty (const PropertyDef *property, string value)
 
virtual void setIntProperty (const PropertyDef *property, int value)
 
virtual void setItemProperty (const PropertyDef *property, Item *item)
 
virtual void setStringProperty (const PropertyDef *property, string value)
 
void setUtilityProperty (string name, int value)
 
virtual string type () const
 
- Public Member Functions inherited from Box
 Box ()
 
 Box (double xmin, double ymin, double zmin, double xmax, double ymax, double zmax)
 
 Box (Vec rmin, Vec rmax)
 
void cellIndices (int &i, int &j, int &k, Vec r, int nx, int ny, int nz) const
 
Vec center () const
 
bool contains (const Box &box) const
 
bool contains (double x, double y, double z) const
 
bool contains (Vec r) const
 
double diagonal () const
 
void extend (const Box &box)
 
const Boxextent () const
 
void extent (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) const
 
Vec fracPos (double xfrac, double yfrac, double zfrac) const
 
Vec fracPos (int xd, int yd, int zd, int xn, int yn, int zn) const
 
bool intersects (const Box &box) const
 
bool intersects (Vec r, const Vec k, double &smin, double &smax) const
 
bool intersects (Vec rc, double r) const
 
Vec rmax () const
 
Vec rmin () const
 
double volume () const
 
Vec widths () const
 
double xmax () const
 
double xmin () const
 
double xwidth () const
 
double ymax () const
 
double ymin () const
 
double ywidth () const
 
double zmax () const
 
double zmin () const
 
double zwidth () const
 

Protected Member Functions

 TetraMeshSpatialGrid ()
 
void setupSelfBefore () override
 
- Protected Member Functions inherited from BoxSpatialGrid
 BoxSpatialGrid ()
 
void setupSelfBefore () override
 
- Protected Member Functions inherited from SpatialGrid
 SpatialGrid ()
 
Randomrandom () const
 
void setupSelfBefore () override
 
virtual void write_xy (SpatialGridPlotFile *outfile) const
 
virtual void write_xyz (SpatialGridPlotFile *outfile) const
 
virtual void write_xz (SpatialGridPlotFile *outfile) const
 
virtual void write_yz (SpatialGridPlotFile *outfile) const
 
- Protected Member Functions inherited from SimulationItem
 SimulationItem ()
 
virtual bool offersInterface (const std::type_info &interfaceTypeInfo) const
 
virtual void setupSelfAfter ()
 
virtual void setupSelfBefore ()
 
- Protected Member Functions inherited from Item
 Item ()
 
- Protected Member Functions inherited from Box
void setExtent (const Box &extent)
 
void setExtent (double xmin, double ymin, double zmin, double xmax, double ymax, double zmax)
 

Private Types

using BaseType = BoxSpatialGrid
 
using FourFaces = std::array< Face, 4 >
 
using FourIndices = std::array< int, 4 >
 
using ItemType = TetraMeshSpatialGrid
 

Private Member Functions

void addCorners ()
 
void buildDelaunay (tetgenio &out)
 
void buildMesh ()
 
void buildSearch ()
 
void refineDelaunay (tetgenio &in, tetgenio &out)
 
void removeInvalid ()
 
void storeTetrahedra (const tetgenio &final, bool storeVertices)
 

Private Attributes

double _eps
 
string _filename
 
Log_log
 
int _numCells
 
int _numSamples
 
int _numVertices
 
Policy _policy
 
bool _refine
 
BoxSearch _search
 
vector< Tetra_tetrahedra
 
vector< Vec_vertices
 

Friends

class ItemRegistry
 
class MySegmentGenerator
 

Detailed Description

TetraMeshSpatialGrid is a concrete subclass of SpatialGrid representing a 3D tetrahedral mesh generated from a set of vertices. The resulting grid fully covers the domain by including the 8 corners of the domain as vertices. Since a tetrahedral mesh always fills the convex hull of the vertices, the domain is always fully covered. This grid will thus always have at least 8 vertices. The grid is constructed using the open-source library TetGen version 1.6.0 (released on August 31, 2020). TetGen is an advanced C++ tetrahedral mesh generator with many features. This class uses the following key features of TetGen:

It should be noted that TetGen is a single-threaded library, but the algorithms used are generally quite fast. The refinement process is by far the most time-consuming part of the mesh generation.

The 3D Delaunay tetrahedralization often contains cells less suited for a computational grid. While the 2D Delaunay triangulation maximizes the smallest angle in each triangle, resulting in more regular cells, the 3D version lacks this property. Consequently, 3D tetrahedralizations frequently contain elongated, irregular cells. To address this, Delaunay refinement algorithms can be used. TetGen implements its own refinement process, applying local mesh operations like vertex smoothing, edge and face swapping, edge contraction, and vertex insertion. This refinement aims to optimize quality metrics such as the radius-edge ratio and the minimum dihedral angle. All refinement parameters are kept at their default values provided by TetGen. However, it is uncertain whether this refinement process will greatly improve the grids for radiative transfer applications.

The positions of the vertices used for generating the tetrahedral mesh are determined by the policy property. The available policies are:

Vertices are removed if they lie outside the simulation domain or are too close to another.

This item type is displayed only if the Boolean expression "Level2" evaluates to true after replacing the names by true or false depending on their presence.

Member Typedef Documentation

◆ FourFaces

using TetraMeshSpatialGrid::FourFaces = std::array<Face, 4>
private

Alias for a fixed array of 4 Faces.

◆ FourIndices

using TetraMeshSpatialGrid::FourIndices = std::array<int, 4>
private

Alias for a fixed array of 4 integer indices.

Member Enumeration Documentation

◆ Policy

enum class TetraMeshSpatialGrid::Policy : int
strong

The enumeration type indicating the policy for determining the positions of the vertices.

Uniform : "random from uniform distribution" .

CentralPeak : "random from distribution with a steep central peak" .

DustDensity : "random from dust density distribution" .

ElectronDensity : "random from electron density distribution" .

GasDensity : "random from gas density distribution" .

File : "loaded from text column data file" .

Constructor & Destructor Documentation

◆ TetraMeshSpatialGrid()

TetraMeshSpatialGrid::TetraMeshSpatialGrid ( )
inlineprotected

Default constructor for concrete Item subclass TetraMeshSpatialGrid : "a tetrahedral spatial grid" .

Member Function Documentation

◆ addCorners()

void TetraMeshSpatialGrid::addCorners ( )
private

This private function adds the 8 corners of the domain to the vertex list. This way the full domain will be tetrahedralized.

◆ buildDelaunay()

void TetraMeshSpatialGrid::buildDelaunay ( tetgenio &  out)
private

This private function constructs the Delaunay tetrahedralization using the vertices obtained from the policy. The output is placed inside the tetgenio reference.

◆ buildMesh()

void TetraMeshSpatialGrid::buildMesh ( )
private

This private function builds the tetrahedral mesh. It starts by constructing the Delaunay tetrahedralization and optionally refines it if the refine option is set to true.

◆ buildSearch()

void TetraMeshSpatialGrid::buildSearch ( )
private

This private function builds the search data structure for the tetrahedral mesh.

◆ cellIndex()

int TetraMeshSpatialGrid::cellIndex ( Position  bfr) const
overridevirtual

This function returns the index of the cell that contains the position r. It uses a search data structure to accelerate the search. If no cell is found to contain this position, the function returns -1.

Implements SpatialGrid.

◆ centralPositionInCell()

Position TetraMeshSpatialGrid::centralPositionInCell ( int  m) const
overridevirtual

This function returns the centroid of the tetrahedron with index m.

Implements SpatialGrid.

◆ createPathSegmentGenerator()

std::unique_ptr< PathSegmentGenerator > TetraMeshSpatialGrid::createPathSegmentGenerator ( ) const
overridevirtual

This function creates and returns ownership of a path segment generator suitable for the tetrahedral spatial grid, implemented as a private PathSegmentGenerator subclass. The algorithm for constructing the path is taken from Maria et al. (2017).

The algorithm uses Plücker coordinates to identify the exit face of a ray inside a given tetrahedron. Plücker coordinates are a set of six values that describe a directed line in 3D space:

πR={k:k×r}={UR:VR}

where r is a position along the ray, and k is the ray's direction. The Plücker product is defined as:

πRπS=URVS+USVR

This product determines the relative orientation between two rays:

πRπS={>0S goes counterclockwise around R<0S goes clockwise around R=0S intersects or is parallel to R

A ray exits a tetrahedron through a particular face if the Plücker products with all three clockwise-ordered edges of that face are negative. The algorithm in Maria et al. (2017) optimizes this by requiring only two Plücker products to be evaluated. Our implementation is described below.

In the first step, the function checks whether the start point is inside the domain. If so, the current point is simply initialized to the start point. If not, the function computes the path segment to the first intersection with one of the domain walls and moves the current point inside the domain. Next, the function determines the cell index of the tetrahedron containing the current point. If none is found, the path is terminated. Before the traversal algorithm can commence, a non-leaving face must be identified. This face acts as the entry face for the ray. Note that this face does not necessarily have to be the actual entry face. This task is handled by the findEnteringFace function of the Tetra class.

Next, the traversal algorithm begins. The entering face is labeled as face 0, with its opposing vertex labeled as vertex 0. We start by evaluating the Plücker product of the ray with the edge 10. Based on whether this product is positive or negative, the next Plücker product to evaluate is determined. If the product is negative (i.e., clockwise orientation), we evaluate the product with the edge in the clockwise direction viewed from vertex 0. The same applies for a positive product. The exit face is determined using the decision tree described in Maria et al. (2017). The distance traveled through the cell is calculated using a simple line-plane intersection:

si=n(vr)nk

where n is the outward-pointing normal of the face, v is any vertex on the exit face, and r is the current position.

Although the algorithm described in Maria et al. (2017) works even if one of the Plücker products is zero, we revert to a plane intersection algorithm in such cases. This approach is similar to the one used in the VoronoiMeshSnapshot class, where the closest intersection distance with all faces is found.

The algorithm continues until the exit face lies on the convex hull boundary. At this point, the path is terminated. If the exit face is not found, which should only rarely happen due to computational inaccuracies, the current point is advanced by a small distance, and the cell index is recalculated.

Implements SpatialGrid.

◆ diagonal()

double TetraMeshSpatialGrid::diagonal ( int  m) const
overridevirtual

This function returns the approximate diagonal of the cell with index m. For a tetrahedron, it takes the square root of the average of the squared edge lengths.

Implements SpatialGrid.

◆ filename()

TetraMeshSpatialGrid::filename ( ) const
inline

This function returns the value of the discoverable string property filename : "the name of the file containing the vertex positions" .

This property is relevant only if the Boolean expression "policyFile" evaluates to true after replacing the names by true or false depending on their presence.

◆ numCells()

int TetraMeshSpatialGrid::numCells ( ) const
overridevirtual

This function returns the number of cells in the grid.

Implements SpatialGrid.

◆ numSamples()

TetraMeshSpatialGrid::numSamples ( ) const
inline

This function returns the value of the discoverable integer property numSamples : "the number of random positions to be used as vertices" .

The minimum value for this property is "4" .

The default value for this property is given by the conditional value expression "500" .

This property is relevant only if the Boolean expression "policyUniform|policyCentralPeak|policyDustDensity|" "policyElectronDensity|policyGasDensity" evaluates to true after replacing the names by true or false depending on their presence.

◆ policy()

TetraMeshSpatialGrid::policy ( ) const
inline

This function returns the value of the discoverable Policy enumeration property policy : "the policy for determining the positions of the vertices" .

The default value for this property is given by the conditional value expression "DustDensity" .

◆ randomPositionInCell()

Position TetraMeshSpatialGrid::randomPositionInCell ( int  m) const
overridevirtual

This function returns a random location from the tetrahedron with index m.

Implements SpatialGrid.

◆ refine()

TetraMeshSpatialGrid::refine ( ) const
inline

This function returns the value of the discoverable Boolean property refine : "refine the grid by performing local mesh operations" .

The default value for this property is given by the conditional value expression "false" .

◆ refineDelaunay()

void TetraMeshSpatialGrid::refineDelaunay ( tetgenio &  in,
tetgenio &  out 
)
private

This private function refines the Delaunay tetrahedralization to improve cell quality. The refinement process is controlled by TetGen with default quality parameters. The input is the initial Delaunay tetrahedralization in the in tetgenio reference. The output, with the refined Delaunay tetrahedralization, is placed inside the out tetgenio reference.

◆ removeInvalid()

void TetraMeshSpatialGrid::removeInvalid ( )
private

This private function removes vertices that are outside the domain or too close to other vertices.

◆ setupSelfBefore()

void TetraMeshSpatialGrid::setupSelfBefore ( )
overrideprotectedvirtual

This function verifies that the attributes are correctly set, generates or retrieves vertex positions based on the configured policy, builds (and optionally refines) the tetrahedralization, and constructs the search structure to optimize the CellIndex function.

Reimplemented from BoxSpatialGrid.

◆ storeTetrahedra()

void TetraMeshSpatialGrid::storeTetrahedra ( const tetgenio &  final,
bool  storeVertices 
)
private

This private function stores the tetrahedra and vertices from the final tetgenio container into the TetraMeshSpatialGrid members. The input is the tetgenio reference with the final tetrahedralization. The storeVertices parameter indicates whether to overwrite the vertices with those from the tetgenio container. This function also logs some cell statistics after it finishes transferring the data.

◆ volume()

double TetraMeshSpatialGrid::volume ( int  m) const
overridevirtual

This function returns the volume of the cell with index m.

Implements SpatialGrid.

◆ writeGridPlotFiles()

void TetraMeshSpatialGrid::writeGridPlotFiles ( const SimulationItem probe) const
overridevirtual

This function outputs the grid plot files. It writes each tetrahedral face as a triangle to the grid plot file.

Reimplemented from SpatialGrid.


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