Chi-Tech
chi_mesh Namespace Reference

Namespaces

namespace  ff_interpolation
 
namespace  mesh_cutting
 
namespace  sweep_management
 

Data Structures

class  BooleanLogicalVolume
 
class  Cell
 
class  CellFace
 
struct  Edge
 
class  ExtruderMeshGenerator
 
struct  ExtrusionLayer
 
struct  Face
 
struct  FFICellIntersection
 
struct  FFIFaceEdgeIntersection
 
struct  FieldFunctionContext
 
class  FieldFunctionInterpolation
 
class  FieldFunctionInterpolationLine
 
class  FieldFunctionInterpolationPoint
 
class  FieldFunctionInterpolationSlice
 
class  FieldFunctionInterpolationVolume
 
class  FromFileMeshGenerator
 
class  GlobalCellHandler
 
class  GridFaceHistogram
 
class  LocalCellHandler
 
class  LogicalVolume
 
class  LogicalVolumeInterface
 
struct  Matrix3x3
 
class  MeshContinuum
 
class  MeshGenerator
 
class  MeshHandler
 
class  MeshModifier
 
class  OrthogonalMeshGenerator
 
struct  PolyFace
 
class  RayTracer
 
struct  RayTracerOutputInformation
 
class  RCCLogicalVolume
 
class  RPPLogicalVolume
 
class  SnapToPlaneMeshModifier
 
class  SphereLogicalVolume
 
class  SplitFileMeshGenerator
 
class  SurfaceMesh
 
class  SurfaceMesher
 
class  SurfaceMesherPassthrough
 
class  SurfaceMesherPredefined
 
class  SurfaceMeshLogicalVolume
 
struct  TensorRank2Dim3
 
class  UnpartitionedMesh
 
struct  Vector3
 
class  VertexHandler
 
class  VolumeMesher
 
class  VolumeMesherExtruder
 
class  VolumeMesherPredefinedUnpartitioned
 

Typedefs

typedef std::shared_ptr< MeshHandlerMeshHandlerPtr
 
typedef std::shared_ptr< SurfaceMeshSurfaceMeshPtr
 
typedef FieldFunctionInterpolation FFInterp
 
typedef std::shared_ptr< FFInterpFFInterpPtr
 
typedef std::shared_ptr< UnpartitionedMeshUnpartitionedMeshPtr
 
typedef UnpartitionedMeshPtr UnpartMeshPtr
 
typedef Vector3 Normal
 
typedef Vector3 Vertex
 
typedef std::shared_ptr< MeshContinuumMeshContinuumPtr
 
typedef std::shared_ptr< const MeshContinuumMeshContinuumConstPtr
 
typedef vtkSmartPointer< vtkUnstructuredGrid > vtkUGridPtr
 
typedef std::pair< vtkUGridPtr, std::string > vtkUGridPtrAndName
 

Enumerations

enum class  CellType {
  GHOST = 0 , SLAB = 1 , TRIANGLE = 4 , QUADRILATERAL = 5 ,
  POLYGON = 6 , TETRAHEDRON = 7 , HEXAHEDRON = 8 , WEDGE = 9 ,
  PYRAMID = 10 , POLYHEDRON = 20 , POINT = 99
}
 
enum  MeshAttributes : int {
  NONE = 0 , DIMENSION_1 = (1 << 0) , DIMENSION_2 = (1 << 1) , DIMENSION_3 = (1 << 2) ,
  ORTHOGONAL = (1 << 3) , EXTRUDED = (1 << 4) , UNSTRUCTURED = (1 << 5)
}
 
enum class  SurfaceMesherType { Predefined = 1 }
 
enum class  VolumeMesherType { EXTRUDER = 4 , UNPARTITIONED = 6 }
 
enum  VolumeMesherProperty {
  FORCE_POLYGONS = 1 , MESH_GLOBAL = 2 , PARTITION_Z = 3 , PARTITION_Y = 4 ,
  PARTITION_X = 5 , CUTS_Z = 6 , CUTS_Y = 7 , CUTS_X = 8 ,
  PARTITION_TYPE = 9 , EXTRUSION_LAYER = 10 , MATID_FROMLOGICAL = 11 , BNDRYID_FROMLOGICAL = 12 ,
  MATID_FROM_LUA_FUNCTION = 13 , BNDRYID_FROM_LUA_FUNCTION = 14
}
 

Functions

std::string CellTypeName (CellType type)
 
MeshAttributes operator| (const MeshAttributes f1, const MeshAttributes f2)
 
MeshHandlerGetCurrentHandler ()
 
size_t PushNewHandlerAndGetIndex ()
 
double ComputeLBF (std::vector< Vector3 > &points, std::vector< double > &x_cuts, std::vector< double > &y_cuts)
 
void DecomposeSurfaceMeshPxPy (const SurfaceMesh &smesh, int Px, int Py)
 
size_t CreateUnpartitioned1DOrthoMesh (std::vector< double > &vertices_1d)
 
size_t CreateUnpartitioned2DOrthoMesh (std::vector< double > &vertices_1d_x, std::vector< double > &vertices_1d_y)
 
size_t CreateUnpartitioned3DOrthoMesh (std::vector< double > &vertices_1d_x, std::vector< double > &vertices_1d_y, std::vector< double > &vertices_1d_z)
 
chi::InputParameters BooleanLogicalVolumeArgumentPair ()
 
 RegisterChiObject (chi_mesh, BooleanLogicalVolume)
 
 RegisterSyntaxBlock (chi_mesh, BooleanLogicalVolumeArgumentPair, BooleanLogicalVolumeArgumentPair)
 
 RegisterChiObject (chi_mesh, RCCLogicalVolume)
 
 RegisterChiObject (chi_mesh, RPPLogicalVolume)
 
 RegisterChiObject (chi_mesh, SphereLogicalVolume)
 
 RegisterChiObject (chi_mesh, SurfaceMeshLogicalVolume)
 
void UploadCellGeometryDiscontinuous (const chi_mesh::MeshContinuum &grid, const chi_mesh::Cell &cell, int64_t &node_counter, vtkNew< vtkPoints > &points, vtkNew< vtkUnstructuredGrid > &ugrid)
 
void UploadCellGeometryContinuous (const chi_mesh::Cell &cell, const std::vector< uint64_t > &vertex_map, vtkNew< vtkUnstructuredGrid > &ugrid)
 
void UploadFaceGeometry (const chi_mesh::CellFace &cell_face, const std::vector< uint64_t > &vertex_map, vtkNew< vtkUnstructuredGrid > &ugrid)
 
int FindHighestDimension (std::vector< vtkUGridPtrAndName > &ugrid_blocks)
 
vtkUGridPtr ConsolidateGridBlocks (std::vector< vtkUGridPtrAndName > &ugrid_blocks, const std::string &block_id_array_name="BlockID")
 
std::vector< vtkUGridPtrAndNameGetBlocksOfDesiredDimension (std::vector< vtkUGridPtrAndName > &ugrid_blocks, int desired_dimension)
 
std::vector< uint64_t > BuildBlockCellExtents (std::vector< vtkUGridPtrAndName > &ugrid_blocks, int desired_dimension)
 
void SetBlockIDArrays (std::vector< vtkUGridPtrAndName > &ugrid_blocks)
 
std::vector< int > BuildCellMaterialIDsFromField (vtkUGridPtr &ugrid, const std::string &field_name, const std::string &file_name)
 
vtkNew< vtkUnstructuredGrid > PrepareVtkUnstructuredGrid (const chi_mesh::MeshContinuum &grid, bool discontinuous=true)
 
void WritePVTUFiles (vtkNew< vtkUnstructuredGrid > &ugrid, const std::string &file_base_name)
 
 RegisterChiObject (chi_mesh, ExtruderMeshGenerator)
 
 RegisterChiObjectParametersOnly (chi_mesh, ExtrusionLayer)
 
 RegisterChiObject (chi_mesh, FromFileMeshGenerator)
 
 RegisterChiObject (chi_mesh, MeshGenerator)
 
 RegisterChiObject (chi_mesh, OrthogonalMeshGenerator)
 
 RegisterChiObject (chi_mesh, SplitFileMeshGenerator)
 
 RegisterChiObject (chi_mesh, SnapToPlaneMeshModifier)
 
bool CheckPlaneLineIntersect (const chi_mesh::Normal &plane_normal, const chi_mesh::Vector3 &plane_point, const chi_mesh::Vector3 &line_point_0, const chi_mesh::Vector3 &line_point_1, chi_mesh::Vector3 &intersection_point, std::pair< double, double > *weights=nullptr)
 
bool CheckLineIntersectStrip (const chi_mesh::Vector3 &strip_point0, const chi_mesh::Vector3 &strip_point1, const chi_mesh::Vector3 &strip_normal, const chi_mesh::Vector3 &line_point0, const chi_mesh::Vector3 &line_point1, chi_mesh::Vector3 &intersection_point, double *distance_to_intersection=nullptr)
 
bool CheckLineIntersectTriangle2 (const chi_mesh::Vector3 &tri_point0, const chi_mesh::Vector3 &tri_point1, const chi_mesh::Vector3 &tri_point2, const chi_mesh::Vector3 &ray_posi, const chi_mesh::Vector3 &ray_dir, chi_mesh::Vector3 &intersection_point, double *distance_to_intersection=nullptr)
 
bool CheckPointInTriangle (const chi_mesh::Vector3 &v0, const chi_mesh::Vector3 &v1, const chi_mesh::Vector3 &v2, const chi_mesh::Normal &n, const chi_mesh::Vector3 &point)
 
bool CheckPlaneTetIntersect (const chi_mesh::Normal &plane_normal, const chi_mesh::Vector3 &plane_point, const std::vector< chi_mesh::Vector3 > &tet_points)
 
void PopulateRaySegmentLengths (const chi_mesh::MeshContinuum &grid, const Cell &cell, const chi_mesh::Vector3 &line_point0, const chi_mesh::Vector3 &line_point1, const chi_mesh::Vector3 &omega, std::vector< double > &segment_lengths)
 

Detailed Description

Namespace for all meshing features

Meshes in ChiTech follow the concept of Regions. In any given region the boundaries are a collection of either line-meshes (2D) or surface-meshes (3D).

Typedef Documentation

◆ FFInterp

Definition at line 21 of file chi_runtime.h.

◆ FFInterpPtr

typedef std::shared_ptr<FFInterp> chi_mesh::FFInterpPtr

Definition at line 22 of file chi_runtime.h.

◆ MeshContinuumConstPtr

typedef std::shared_ptr<const MeshContinuum> chi_mesh::MeshContinuumConstPtr

Definition at line 46 of file chi_mesh.h.

◆ MeshContinuumPtr

typedef std::shared_ptr< MeshContinuum > chi_mesh::MeshContinuumPtr

Definition at line 45 of file chi_mesh.h.

◆ MeshHandlerPtr

typedef std::shared_ptr<MeshHandler> chi_mesh::MeshHandlerPtr

Definition at line 15 of file chi_runtime.h.

◆ Normal

Definition at line 19 of file chi_mesh.h.

◆ SurfaceMeshPtr

typedef std::shared_ptr<SurfaceMesh> chi_mesh::SurfaceMeshPtr

Definition at line 18 of file chi_runtime.h.

◆ UnpartitionedMeshPtr

Definition at line 25 of file chi_runtime.h.

◆ UnpartMeshPtr

Definition at line 26 of file chi_runtime.h.

◆ Vertex

Definition at line 20 of file chi_mesh.h.

◆ vtkUGridPtr

typedef vtkSmartPointer<vtkUnstructuredGrid> chi_mesh::vtkUGridPtr

Definition at line 36 of file chi_grid_vtk_utils.h.

◆ vtkUGridPtrAndName

typedef std::pair<vtkUGridPtr, std::string> chi_mesh::vtkUGridPtrAndName

Definition at line 37 of file chi_grid_vtk_utils.h.

Enumeration Type Documentation

◆ CellType

enum class chi_mesh::CellType
strong
Enumerator
GHOST 
SLAB 
TRIANGLE 
QUADRILATERAL 
POLYGON 
TETRAHEDRON 
HEXAHEDRON 
WEDGE 
PYRAMID 
POLYHEDRON 
POINT 

Definition at line 11 of file cell.h.

◆ MeshAttributes

Enumerator
NONE 
DIMENSION_1 
DIMENSION_2 
DIMENSION_3 
ORTHOGONAL 
EXTRUDED 
UNSTRUCTURED 

Definition at line 69 of file chi_mesh.h.

◆ SurfaceMesherType

enum class chi_mesh::SurfaceMesherType
strong
Enumerator
Predefined 

Definition at line 8 of file surfacemesher.h.

◆ VolumeMesherProperty

Enumerator
FORCE_POLYGONS 
MESH_GLOBAL 
PARTITION_Z 
PARTITION_Y 
PARTITION_X 
CUTS_Z 
CUTS_Y 
CUTS_X 
PARTITION_TYPE 
EXTRUSION_LAYER 
MATID_FROMLOGICAL 
BNDRYID_FROMLOGICAL 
MATID_FROM_LUA_FUNCTION 
BNDRYID_FROM_LUA_FUNCTION 

Definition at line 19 of file chi_volumemesher.h.

◆ VolumeMesherType

enum class chi_mesh::VolumeMesherType
strong
Enumerator
EXTRUDER 
UNPARTITIONED 

Definition at line 14 of file chi_volumemesher.h.

Function Documentation

◆ BooleanLogicalVolumeArgumentPair()

chi::InputParameters chi_mesh::BooleanLogicalVolumeArgumentPair ( )

Definition at line 58 of file BooleanLogicalVolume.cc.

◆ BuildBlockCellExtents()

std::vector< uint64_t > chi_mesh::BuildBlockCellExtents ( std::vector< vtkUGridPtrAndName > &  ugrid_blocks,
int  desired_dimension 
)

Given several unstructured grid blocks, each denoting a material id, this function sets material ids accordingly.

Definition at line 158 of file chi_grid_vtk_utils_01_reading.cc.

◆ BuildCellMaterialIDsFromField()

std::vector< int > chi_mesh::BuildCellMaterialIDsFromField ( vtkUGridPtr ugrid,
const std::string &  field_name,
const std::string &  file_name 
)

Retrieves material-ids from a field.

Definition at line 209 of file chi_grid_vtk_utils_01_reading.cc.

◆ CellTypeName()

std::string chi_mesh::CellTypeName ( CellType  type)

Provides the text name associated with a cell type.

Definition at line 12 of file cell.cc.

◆ CheckLineIntersectStrip()

bool chi_mesh::CheckLineIntersectStrip ( const chi_mesh::Vector3 strip_point0,
const chi_mesh::Vector3 strip_point1,
const chi_mesh::Vector3 strip_normal,
const chi_mesh::Vector3 line_point0,
const chi_mesh::Vector3 line_point1,
chi_mesh::Vector3 intersection_point,
double *  distance_to_intersection = nullptr 
)

Given a strip defined by two points (v0,v1) and a normal, n, (meaning infinite in the direction defined by (v1-v0).cross(n), this function determines if a line, defined from p0 to p1, intersects it. If it does then true is returned and intersection_point contains the point of intersection. If it does not then false is returned and intersection_point remains unchanged.

Definition at line 77 of file raytrace_utils.cc.

◆ CheckLineIntersectTriangle2()

bool chi_mesh::CheckLineIntersectTriangle2 ( const chi_mesh::Vector3 tri_point0,
const chi_mesh::Vector3 tri_point1,
const chi_mesh::Vector3 tri_point2,
const chi_mesh::Vector3 ray_posi,
const chi_mesh::Vector3 ray_dir,
chi_mesh::Vector3 intersection_point,
double *  distance_to_intersection = nullptr 
)

Given a triangle defined by three points, computes whether a line intersects this triangle.

Definition at line 123 of file raytrace_utils.cc.

◆ CheckPlaneLineIntersect()

bool chi_mesh::CheckPlaneLineIntersect ( const chi_mesh::Normal plane_normal,
const chi_mesh::Vector3 plane_point,
const chi_mesh::Vector3 line_point_0,
const chi_mesh::Vector3 line_point_1,
chi_mesh::Vector3 intersection_point,
std::pair< double, double > *  weights = nullptr 
)

Computes the intersection of a line with a plane.

The first step of this algorithm is to compute v0 and v1. These are vectors from the plane's reference point to each of the line-points, respectively. We then take the dot-products of these vectors with the plane normal. We then say that the vectors have a positive sense if the dot-product is positive and a negative sense if the dot-product is negative. If the senses are not equal then the line intersects the plane.

Since the face normal is a normalized vector the dot-product of v0 or v1 will give the projection of the relevant vector along the normal to the plane. We can use this projection to compute a weight associated with each vector. This also then allows us to compute the intersection point.

Parameters
plane_normalThe normal associated with the plane
plane_pointThe reference point for the plane
line_point_0The line's initial point
line_point_1The line's destination point
intersection_pointThe point to be populated with the intersection point
weightsThe weights associated with this intersection
Returns
Returns true if the line intersects the plane and false otherwise.
Author
Jan

Definition at line 35 of file raytrace_utils.cc.

◆ CheckPlaneTetIntersect()

bool chi_mesh::CheckPlaneTetIntersect ( const chi_mesh::Normal plane_normal,
const chi_mesh::Vector3 plane_point,
const std::vector< chi_mesh::Vector3 > &  tet_points 
)

This functions checks the intersection of a plane with a tetrahedron. The equation of a plane is nx(x-x0) + ny(y-y0) + nz(z-z0) = 0 Where the plane normal is (nx,ny,nz) and the plane point is (x0,y0,z0). If we form a dot product between the normal and a vector (x-x0,y-y0,z-z0) then sign of the result gives the sense to the surface. Therefore, if we encounter differing senses then the plane is indeed intersecting.

Definition at line 225 of file raytrace_utils.cc.

◆ CheckPointInTriangle()

bool chi_mesh::CheckPointInTriangle ( const chi_mesh::Vector3 v0,
const chi_mesh::Vector3 v1,
const chi_mesh::Vector3 v2,
const chi_mesh::Normal n,
const chi_mesh::Vector3 point 
)

Check whether a point lies in a triangle.

Definition at line 186 of file raytrace_utils.cc.

◆ ComputeLBF()

double chi_mesh::ComputeLBF ( std::vector< Vector3 > &  points,
std::vector< double > &  x_cuts,
std::vector< double > &  y_cuts 
)

Makes a centroid based load balance factor calculation.

Author
Jan

Definition at line 16 of file decompose_pxpy.cc.

◆ ConsolidateGridBlocks()

chi_mesh::vtkUGridPtr chi_mesh::ConsolidateGridBlocks ( std::vector< vtkUGridPtrAndName > &  ugrid_blocks,
const std::string &  block_id_array_name = "BlockID" 
)

Consolidates all blocks containing cells with the desired dimension. Thereafter it removes duplicate vertices.

Lambda to right pad an entry.

Definition at line 37 of file chi_grid_vtk_utils_01_reading.cc.

◆ CreateUnpartitioned1DOrthoMesh()

size_t chi_mesh::CreateUnpartitioned1DOrthoMesh ( std::vector< double > &  vertices)

Creates a 1D slab mesh from a set of vertices.

Definition at line 14 of file mesh_orthomacros_02_unpartitioned.cc.

◆ CreateUnpartitioned2DOrthoMesh()

size_t chi_mesh::CreateUnpartitioned2DOrthoMesh ( std::vector< double > &  vertices_1d_x,
std::vector< double > &  vertices_1d_y 
)

Creates a 2D orthogonal mesh from a set of vertices in x and y. The vertices along a dimension merely represents the divisions. They are not the complete vertices defining a cell. For example:

std::vector<chi_mesh::Vertex> vertices_x = {0.0,1.0,2.0};
std::vector<chi_mesh::Vertex> vertices_y = {0.0,1.0,2.0};
size_t CreateUnpartitioned2DOrthoMesh(std::vector< double > &vertices_1d_x, std::vector< double > &vertices_1d_y)

This code will create a 2x2 mesh with $ \vec{x} \in [0,2]^2 $.

Definition at line 96 of file mesh_orthomacros_02_unpartitioned.cc.

◆ CreateUnpartitioned3DOrthoMesh()

size_t chi_mesh::CreateUnpartitioned3DOrthoMesh ( std::vector< double > &  vertices_1d_x,
std::vector< double > &  vertices_1d_y,
std::vector< double > &  vertices_1d_z 
)

Creates a 3D orthogonal mesh from a set of vertices in x,y,z. The vertices along a dimension merely represents the divisions. They are not the complete vertices defining a cell. For example:

std::vector<chi_mesh::Vertex> vertices_x = {0.0,1.0,2.0};
std::vector<chi_mesh::Vertex> vertices_y = {0.0,1.0,2.0};
std::vector<chi_mesh::Vertex> vertices_z = {0.0,1.0,2.0};
chi_mesh::CreateUnpartitioned3DOrthoMesh(vertices_x,vertices_y,vertices_z);
size_t CreateUnpartitioned3DOrthoMesh(std::vector< double > &vertices_1d_x, std::vector< double > &vertices_1d_y, std::vector< double > &vertices_1d_z)

This code will create a 2x2 mesh with $ \vec{x} \in [0,2]^2 $.

Definition at line 210 of file mesh_orthomacros_02_unpartitioned.cc.

◆ DecomposeSurfaceMeshPxPy()

void chi_mesh::DecomposeSurfaceMeshPxPy ( const SurfaceMesh smesh,
int  px,
int  py 
)

Decomposes a 2D surface mesh using the centroids in a Px-Py fashion.

Definition at line 73 of file decompose_pxpy.cc.

◆ FindHighestDimension()

int chi_mesh::FindHighestDimension ( std::vector< vtkUGridPtrAndName > &  ugrid_blocks)

Finds the highest dimension across all the grid blocks. This is useful when a vtk-read mesh contains multiple blocks. Some of which are boundary faces.

Definition at line 16 of file chi_grid_vtk_utils_01_reading.cc.

◆ GetBlocksOfDesiredDimension()

std::vector< chi_mesh::vtkUGridPtrAndName > chi_mesh::GetBlocksOfDesiredDimension ( std::vector< vtkUGridPtrAndName > &  ugrid_blocks,
int  desired_dimension 
)

Provides a map of the different grids that have the requested dimension.

Definition at line 136 of file chi_grid_vtk_utils_01_reading.cc.

◆ GetCurrentHandler()

chi_mesh::MeshHandler & chi_mesh::GetCurrentHandler ( )

Obtains a reference to the current mesh handler from the global stack.

If the stack is empty this routine will through std::logic_error.

Author
Jan

Definition at line 13 of file chi_mesh_meshhandler_utils.cc.

◆ operator|()

MeshAttributes chi_mesh::operator| ( const MeshAttributes  f1,
const MeshAttributes  f2 
)
inline

Definition at line 80 of file chi_mesh.h.

◆ PopulateRaySegmentLengths()

void chi_mesh::PopulateRaySegmentLengths ( const chi_mesh::MeshContinuum grid,
const Cell cell,
const chi_mesh::Vector3 line_point0,
const chi_mesh::Vector3 line_point1,
const chi_mesh::Vector3 omega,
std::vector< double > &  segment_lengths 
)

Populates segment lengths along a ray. Sorted along the direction.

Definition at line 252 of file raytrace_utils.cc.

◆ PrepareVtkUnstructuredGrid()

vtkNew< vtkUnstructuredGrid > chi_mesh::PrepareVtkUnstructuredGrid ( const chi_mesh::MeshContinuum grid,
bool  discontinuous = true 
)

Uploads vertices and cells to an unstructured grid. This routine also uploads cell material ids (sub-domain ids) and partition ids.

Definition at line 21 of file chi_grid_vtk_utils_04_writing.cc.

◆ PushNewHandlerAndGetIndex()

size_t chi_mesh::PushNewHandlerAndGetIndex ( )

Adds a new mesh handler to the stack, sets it as the current handler and returns a handle to it.

Definition at line 25 of file chi_mesh_meshhandler_utils.cc.

◆ RegisterChiObject() [1/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
BooleanLogicalVolume   
)

◆ RegisterChiObject() [2/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
ExtruderMeshGenerator   
)

◆ RegisterChiObject() [3/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
FromFileMeshGenerator   
)

◆ RegisterChiObject() [4/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
MeshGenerator   
)

◆ RegisterChiObject() [5/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
OrthogonalMeshGenerator   
)

◆ RegisterChiObject() [6/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
RCCLogicalVolume   
)

◆ RegisterChiObject() [7/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
RPPLogicalVolume   
)

◆ RegisterChiObject() [8/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
SnapToPlaneMeshModifier   
)

◆ RegisterChiObject() [9/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
SphereLogicalVolume   
)

◆ RegisterChiObject() [10/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
SplitFileMeshGenerator   
)

◆ RegisterChiObject() [11/11]

chi_mesh::RegisterChiObject ( chi_mesh  ,
SurfaceMeshLogicalVolume   
)

◆ RegisterChiObjectParametersOnly()

chi_mesh::RegisterChiObjectParametersOnly ( chi_mesh  ,
ExtrusionLayer   
)

◆ RegisterSyntaxBlock()

chi_mesh::RegisterSyntaxBlock ( chi_mesh  ,
BooleanLogicalVolumeArgumentPair  ,
BooleanLogicalVolumeArgumentPair   
)

◆ SetBlockIDArrays()

void chi_mesh::SetBlockIDArrays ( std::vector< vtkUGridPtrAndName > &  ugrid_blocks)

Given several unstructured grid blocks, each denoting a material id, this function creates a VTK cell-data array called "BlockID" that holds this information.

Definition at line 183 of file chi_grid_vtk_utils_01_reading.cc.

◆ UploadCellGeometryContinuous()

void chi_mesh::UploadCellGeometryContinuous ( const chi_mesh::Cell cell,
const std::vector< uint64_t > &  vertex_map,
vtkNew< vtkUnstructuredGrid > &  ugrid 
)

Uploads vertices and cells to an unstructured grid.

Definition at line 98 of file chi_grid_vtk_utils_00_celluploading.cc.

◆ UploadCellGeometryDiscontinuous()

void chi_mesh::UploadCellGeometryDiscontinuous ( const chi_mesh::MeshContinuum grid,
const chi_mesh::Cell cell,
int64_t &  node_counter,
vtkNew< vtkPoints > &  points,
vtkNew< vtkUnstructuredGrid > &  ugrid 
)

Uploads vertices and cells to an unstructured grid.

Definition at line 10 of file chi_grid_vtk_utils_00_celluploading.cc.

◆ UploadFaceGeometry()

void chi_mesh::UploadFaceGeometry ( const chi_mesh::CellFace cell_face,
const std::vector< uint64_t > &  vertex_map,
vtkNew< vtkUnstructuredGrid > &  ugrid 
)

Uploads vertices and cells to an unstructured grid.

Definition at line 190 of file chi_grid_vtk_utils_00_celluploading.cc.

◆ WritePVTUFiles()

void chi_mesh::WritePVTUFiles ( vtkNew< vtkUnstructuredGrid > &  ugrid,
const std::string &  file_base_name 
)

Writes an unstructured grid to files (.pvtu and .vtu).

Definition at line 77 of file chi_grid_vtk_utils_04_writing.cc.