DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
DNDS::Geom::UnstructuredMesh Struct Reference

#include <Mesh.hpp>

Inheritance diagram for DNDS::Geom::UnstructuredMesh:
[legend]
Collaboration diagram for DNDS::Geom::UnstructuredMesh:
[legend]

Classes

struct  ElevationInfo
 
struct  HDF5OutSetting
 
struct  WallDistOptions
 

Public Types

template<DeviceBackend B>
using t_deviceView = UnstructuredMeshDeviceView< B >
 

Public Member Functions

auto device_array_list_primary ()
 
auto device_array_list_N2CB ()
 
auto device_array_list_facial ()
 
auto device_array_list_C2F ()
 
 UnstructuredMesh (const DNDS::MPIInfo &n_mpi, int n_dim)
 
int getDim () const
 
auto getCell2NodeIndexPbiRow (index iCell)
 
auto getBnd2NodeIndexPbiRow (index iBnd)
 
template<class TPair >
void EnsureGhostMapping (TPair &pair)
 Ensure a pair has a ghost mapping on its transformer.
 
template<class TPair >
index IndexGlobal2Local (TPair &pair, DNDS::index iGlobal)
 Global-to-local conversion using father+son ghost mapping.
 
template<class TPair >
index IndexLocal2Global (TPair &pair, DNDS::index iLocal)
 Local-to-global conversion using father+son ghost mapping.
 
template<class TPair >
index IndexLocal2Global_NoSon (TPair &pair, index iLocal)
 Local-to-global conversion using only the father's global mapping (no son / ghost layer).
 
template<class TPair >
index IndexGlobal2Local_NoSon (TPair &pair, index iGlobal)
 Global-to-local conversion using only the father's global mapping (no son / ghost layer).
 
index NodeIndexGlobal2Local (DNDS::index i)
 
index NodeIndexLocal2Global (DNDS::index i)
 
index NodeIndexLocal2Global_NoSon (index i)
 
index NodeIndexGlobal2Local_NoSon (index i)
 
index CellIndexGlobal2Local (DNDS::index i)
 
index CellIndexLocal2Global (DNDS::index i)
 
index CellIndexLocal2Global_NoSon (index i)
 
index CellIndexGlobal2Local_NoSon (index i)
 
index BndIndexGlobal2Local (DNDS::index i)
 
index BndIndexLocal2Global (DNDS::index i)
 
index BndIndexLocal2Global_NoSon (index i)
 
index BndIndexGlobal2Local_NoSon (index i)
 
index FaceIndexGlobal2Local (DNDS::index i)
 
index FaceIndexLocal2Global (DNDS::index i)
 
void SetPeriodicGeometry (const tPoint &translation1, const tPoint &rotationCenter1, const tPoint &eulerAngles1, const tPoint &translation2, const tPoint &rotationCenter2, const tPoint &eulerAngles2, const tPoint &translation3, const tPoint &rotationCenter3, const tPoint &eulerAngles3)
 
void RecoverNode2CellAndNode2Bnd ()
 only requires father part of cell2node, bnd2node and coords generates node2cell and node2bnd (father part)
 
void RecoverNode2CellAndNode2BndLegacy ()
 
void RecoverCell2CellAndBnd2Cell ()
 needs to use RecoverNode2CellAndNode2Bnd before doing this. Requires node2cell.father and builds a version of its son.
 
void RecoverCell2CellAndBnd2CellLegacy ()
 
void BuildGhostPrimary (int nGhostLayers=1)
 building ghost (son) from primary (currently only cell2cell)
 
void BuildGhostPrimaryLegacy ()
 
void AdjGlobal2LocalPrimary ()
 
void AdjLocal2GlobalPrimary ()
 
void AdjGlobal2LocalPrimaryForBnd ()
 
void AdjLocal2GlobalPrimaryForBnd ()
 
void AdjGlobal2LocalFacial ()
 
void AdjLocal2GlobalFacial ()
 
void AdjGlobal2LocalC2F ()
 
void AdjLocal2GlobalC2F ()
 
void BuildGhostN2CB ()
 
void AdjGlobal2LocalN2CB ()
 
void AdjLocal2GlobalN2CB ()
 
void AssertOnN2CB ()
 
void BuildCell2CellFace ()
 
void AdjGlobal2LocalC2CFace ()
 
void AdjLocal2GlobalC2CFace ()
 
void InterpolateFace ()
 
void BuildGhostFace ()
 
void MatchFaceBoundary ()
 
void InterpolateFaceLegacy ()
 
void AssertOnFaces ()
 
void ConstructBndMesh (UnstructuredMesh &bMesh)
 
void fillRegistry (MeshConnectivity &dag) const
 Populate a MeshConnectivity registry from this mesh's currently-built adjacencies.
 
void fillRegistry (MeshConnectivity &dag, const std::unordered_set< AdjKind, AdjKindHash > &skip) const
 
ReorderRegistry buildReorderRegistry (const std::unordered_set< EntityKind > &destroyKinds={})
 Build a ReorderRegistry containing all mesh members.
 
void ReorderEntities (const ReorderInput &input)
 Reorder entities using the general framework.
 
ReorderPlan buildReorderPlan (const ReorderInput &input)
 Build a ReorderPlan without applying it.
 
ReorderInput resolveFollows (const ReorderInput &input, const ReorderRegistry &reg)
 Augment a ReorderInput with default follows and compute follow maps.
 
tLocalMatStruct GetCell2CellFaceVLocal (bool onLocalPartition=false)
 
void ObtainLocalFactFillOrdering (Direct::SerialSymLUStructure &symLU, Direct::DirectPrecControl control)
 
void ObtainSymmetricSymbolicFactorization (Direct::SerialSymLUStructure &symLU, Direct::DirectPrecControl control) const
 
void ReorderLocalCells (int nParts=1, int nPartsInner=1)
 
void ReorderLocalCellsLegacy (int nParts=1, int nPartsInner=1)
 Legacy implementation preserved for reference/fallback.
 
int NLocalParts () const
 
index LocalPartStart (int iPart) const
 
index LocalPartEnd (int iPart) const
 
index NumNode () const
 
index NumCell () const
 
index NumFace () const
 
index NumBnd () const
 
index NumNodeGhost () const
 
index NumCellGhost () const
 
index NumFaceGhost () const
 
index NumBndGhost () const
 
index NumNodeProc () const
 
index NumCellProc () const
 
index NumFaceProc () const
 
index NumBndProc () const
 
index NumCellGlobal ()
 
index NumNodeGlobal ()
 
index NumFaceGlobal ()
 
index NumBndGlobal ()
 
Elem::Element GetCellElement (index iC)
 
Elem::Element GetFaceElement (index iF)
 
Elem::Element GetBndElement (index iB)
 
t_index GetCellZone (index iC)
 
t_index GetFaceZone (index iF)
 
t_index GetBndZone (index iB)
 
MPIInfogetMPI ()
 
void BuildO2FromO1Elevation (UnstructuredMesh &meshO1)
 
void ElevatedNodesGetBoundarySmooth (const std::function< bool(t_index)> &FiFBndIdNeedSmooth)
 
void ElevatedNodesSolveInternalSmooth ()
 
void ElevatedNodesSolveInternalSmoothV1Old ()
 
void ElevatedNodesSolveInternalSmoothV1 ()
 
void ElevatedNodesSolveInternalSmoothV2 ()
 
void BuildBisectO1FormO2 (UnstructuredMesh &meshO2)
 
bool IsO1 ()
 
bool IsO2 ()
 
template<class tC2n >
void _detail_GetCoords (const tC2n &c2n, tSmallCoords &cs)
 directly load coords; gets faulty if isPeriodic!
 
template<class tC2n , class tCoordExt >
void _detail_GetCoords (const tC2n &c2n, tSmallCoords &cs, tCoordExt &coo)
 directly load coords; gets faulty if isPeriodic!
 
template<class tC2n , class tC2nPbi >
void _detail_GetCoordsOnElem (const tC2n &c2n, const tC2nPbi &c2nPbi, tSmallCoords &cs)
 specially for periodicity
 
template<class tC2n , class tC2nPbi , class tCoordExt >
void _detail_GetCoordsOnElem (const tC2n &c2n, const tC2nPbi &c2nPbi, tSmallCoords &cs, tCoordExt &coo)
 specially for periodicity
 
void GetCoordsOnCell (index iCell, tSmallCoords &cs)
 
void GetCoordsOnCell (index iCell, tSmallCoords &cs, tCoordPair &coo)
 
void GetCoordsOnFace (index iFace, tSmallCoords &cs)
 
tPoint GetCoordNodeOnCell (index iCell, rowsize ic2n)
 
tPoint GetCoordNodeOnFace (index iFace, rowsize if2n)
 
tPoint GetCoordWallDistOnCell (index iCell, rowsize ic2n)
 
tPoint GetCoordWallDistOnFace (index iFace, rowsize if2n)
 
bool CellIsFaceBack (index iCell, index iFace) const
 
index CellFaceOther (index iCell, index iFace) const
 
template<class FA , class FB , class F0 = std::function<void(void)>>
auto CellOtherCellPeriodicHandle (index iFace, rowsize if2c, FA &&fA, FB &&fB, F0 &&f0=[]() {})
 fA executes when if2c points to the donor side; fB the main side
 
void WriteSerialize (Serializer::SerializerBaseSSP serializerP, const std::string &name)
 
void ReadSerialize (Serializer::SerializerBaseSSP serializerP, const std::string &name)
 
void ReadSerializeAndDistribute (Serializer::SerializerBaseSSP serializerP, const std::string &name, const PartitionOptions &partitionOptions)
 Reads mesh from an H5 serializer using even-split distribution, then repartitions via ParMetis for locality.
 
template<class TFTrans >
void TransformCoords (TFTrans &&FTrans)
 
void RecreatePeriodicNodes ()
 
void BuildVTKConnectivity ()
 
void PrintParallelVTKHDFDataArray (std::string fname, std::string seriesName, int arraySiz, int vecArraySiz, int arraySizPoint, int vecArraySizPoint, const tFGetName &names, const tFGetData &data, const tFGetName &vectorNames, const tFGetVecData &vectorData, const tFGetName &namesPoint, const tFGetData &dataPoint, const tFGetName &vectorNamesPoint, const tFGetVecData &vectorDataPoint, double t, MPI_Comm commDup)
 
void SetHDF5OutSetting (size_t chunkSiz, int deflateLevel, bool coll_on_data, bool coll_on_meta)
 
void PrintMeshCGNS (std::string fname, const t_FBCID_2_Name &fbcid2name, const std::vector< std::string > &allNames)
 
void BuildNodeWallDist (const std::function< bool(Geom::t_index)> &fBndIsWall, WallDistOptions options=WallDistOptions{})
 
template<class F >
void op_on_device_arrays (F &&f)
 
template<typename F >
void for_each_device_member (F &&f)
 
template<DeviceBackend B>
auto deviceView ()
 
index getArrayBytes ()
 
- Public Member Functions inherited from DNDS::DeviceTransferable< UnstructuredMesh >
void to_device (DeviceBackend B)
 Mirror every registered ArrayPair to the given device backend.
 
void to_host ()
 Pull every registered pair back to host memory.
 
DeviceBackend device ()
 Consistent device backend across all registered pairs.
 
index getDeviceArrayBytes ()
 Total footprint of every registered father+son pair in bytes.
 

Static Public Member Functions

template<class TAdj , class TFn >
static void ConvertAdjEntries (TAdj &adj, index nRows, TFn &&fn)
 Apply a conversion function to every entry of an adjacency array's first nRows rows.
 
template<class TAdj , class TFn >
static void ConvertAdjEntriesOMP (TAdj &adj, index nRows, TFn &&fn)
 OpenMP-parallelized variant of ConvertAdjEntries.
 
template<class TPair , class TFn >
static void PermuteRows (TPair &pair, index nRows, TFn &&old2new)
 Permute the father rows of an ArrayPair according to a mapping function.
 

Public Attributes

MPI_int mRank {0}
 
DNDS::MPIInfo mpi
 
int dim
 
bool isPeriodic {false}
 
MeshAdjState adjPrimaryState {Adj_Unknown}
 
MeshAdjState adjFacialState {Adj_Unknown}
 
MeshAdjState adjC2FState {Adj_Unknown}
 
MeshAdjState adjN2CBState {Adj_Unknown}
 
MeshAdjState adjC2CFaceState {Adj_Unknown}
 
Periodicity periodicInfo
 
index nNodeO1 {-1}
 
MeshElevationState elevState = Elevation_Untouched
 
tCoordPair coords
 reader
 
AdjPairTracked< tAdjPaircell2node
 
AdjPairTracked< tAdjPairbnd2node
 
AdjPairTracked< tAdj2Pairbnd2cell
 
AdjPairTracked< tAdjPaircell2cell
 
tElemInfoArrayPair cellElemInfo
 
tElemInfoArrayPair bndElemInfo
 
tAdj1Pair cell2cellOrig
 
tAdj1Pair node2nodeOrig
 
tAdj1Pair bnd2bndOrig
 
tPbiPair cell2nodePbi
 periodic only, after reader
 
tPbiPair bnd2nodePbi
 
AdjPairTracked< tAdjPairnode2cell
 inverse relations
 
AdjPairTracked< tAdjPairnode2bnd
 
AdjPairTracked< tAdjPaircell2face
 interpolated
 
AdjPairTracked< tAdjPairface2node
 
AdjPairTracked< tAdj2Pairface2cell
 
tElemInfoArrayPair faceElemInfo
 
AdjPairTracked< tAdj1Pairbnd2face
 
AdjPairTracked< tAdj1Pairface2bnd
 
std::vector< indexbnd2faceV
 
std::unordered_map< index, indexface2bndM
 
tPbiPair face2nodePbi
 periodic only, after interpolated
 
AdjPairTracked< tAdjPaircell2cellFace
 constructed on demand
 
std::vector< indexnode2parentNode
 parent built
 
std::vector< indexnode2bndNode
 
std::vector< indexcell2parentCell
 
std::vector< indexvtkCell2nodeOffsets
 for parallel out
 
std::vector< uint8_tvtkCellType
 
std::vector< indexvtkCell2node
 
index vtkNodeOffset {-1}
 
index vtkCellOffset {-1}
 
index vtkCell2NodeGlobalSiz {-1}
 
index vtkNCellGlobal {-1}
 
index vtkNNodeGlobal {-1}
 
tAdjPair cell2nodePeriodicRecreated
 
tCoordPair coordsPeriodicRecreated
 
std::vector< indexnodeRecreated2nodeLocal
 
struct DNDS::Geom::UnstructuredMesh::HDF5OutSetting hdf5OutSetting
 
tCoordPair coordsElevDisp
 only elevation
 
index nTotalMoved {-1}
 
struct DNDS::Geom::UnstructuredMesh::ElevationInfo elevationInfo
 
tLocalMatStruct cell2cellFaceVLocalParts
 for cell local factorization
 
std::vector< indexlocalPartitionStarts
 
tCoordPair nodeWallDist
 wall dist:
 

Detailed Description

Definition at line 37 of file Mesh.hpp.

Member Typedef Documentation

◆ t_deviceView

Constructor & Destructor Documentation

◆ UnstructuredMesh()

DNDS::Geom::UnstructuredMesh::UnstructuredMesh ( const DNDS::MPIInfo n_mpi,
int  n_dim 
)
inline

Definition at line 179 of file Mesh.hpp.

Member Function Documentation

◆ _detail_GetCoords() [1/2]

template<class tC2n >
void DNDS::Geom::UnstructuredMesh::_detail_GetCoords ( const tC2n c2n,
tSmallCoords cs 
)
inline

directly load coords; gets faulty if isPeriodic!

Definition at line 715 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _detail_GetCoords() [2/2]

template<class tC2n , class tCoordExt >
void DNDS::Geom::UnstructuredMesh::_detail_GetCoords ( const tC2n c2n,
tSmallCoords cs,
tCoordExt coo 
)
inline

directly load coords; gets faulty if isPeriodic!

Definition at line 734 of file Mesh.hpp.

Here is the call graph for this function:

◆ _detail_GetCoordsOnElem() [1/2]

template<class tC2n , class tC2nPbi >
void DNDS::Geom::UnstructuredMesh::_detail_GetCoordsOnElem ( const tC2n c2n,
const tC2nPbi c2nPbi,
tSmallCoords cs 
)
inline

specially for periodicity

Definition at line 753 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ _detail_GetCoordsOnElem() [2/2]

void DNDS::Geom::UnstructuredMesh::_detail_GetCoordsOnElem ( const tC2n c2n,
const tC2nPbi c2nPbi,
tSmallCoords cs,
tCoordExt coo 
)
inline

specially for periodicity

Definition at line 772 of file Mesh.hpp.

Here is the call graph for this function:

◆ AdjGlobal2LocalC2CFace()

void DNDS::Geom::UnstructuredMesh::AdjGlobal2LocalC2CFace ( )

Definition at line 1080 of file Mesh.cpp.

◆ AdjGlobal2LocalC2F()

void DNDS::Geom::UnstructuredMesh::AdjGlobal2LocalC2F ( )

Definition at line 924 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjGlobal2LocalFacial()

void DNDS::Geom::UnstructuredMesh::AdjGlobal2LocalFacial ( )

Definition at line 873 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjGlobal2LocalN2CB()

void DNDS::Geom::UnstructuredMesh::AdjGlobal2LocalN2CB ( )

Definition at line 974 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjGlobal2LocalPrimary()

void DNDS::Geom::UnstructuredMesh::AdjGlobal2LocalPrimary ( )

Definition at line 805 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjGlobal2LocalPrimaryForBnd()

void DNDS::Geom::UnstructuredMesh::AdjGlobal2LocalPrimaryForBnd ( )

Definition at line 848 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjLocal2GlobalC2CFace()

void DNDS::Geom::UnstructuredMesh::AdjLocal2GlobalC2CFace ( )

Definition at line 1066 of file Mesh.cpp.

◆ AdjLocal2GlobalC2F()

void DNDS::Geom::UnstructuredMesh::AdjLocal2GlobalC2F ( )

Definition at line 910 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjLocal2GlobalFacial()

void DNDS::Geom::UnstructuredMesh::AdjLocal2GlobalFacial ( )

Definition at line 893 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjLocal2GlobalN2CB()

void DNDS::Geom::UnstructuredMesh::AdjLocal2GlobalN2CB ( )

Definition at line 960 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjLocal2GlobalPrimary()

void DNDS::Geom::UnstructuredMesh::AdjLocal2GlobalPrimary ( )

Definition at line 829 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AdjLocal2GlobalPrimaryForBnd()

void DNDS::Geom::UnstructuredMesh::AdjLocal2GlobalPrimaryForBnd ( )

Definition at line 861 of file Mesh.cpp.

Here is the caller graph for this function:

◆ AssertOnFaces()

void DNDS::Geom::UnstructuredMesh::AssertOnFaces ( )

Definition at line 1377 of file Mesh.cpp.

Here is the call graph for this function:

◆ AssertOnN2CB()

void DNDS::Geom::UnstructuredMesh::AssertOnN2CB ( )

Definition at line 988 of file Mesh.cpp.

Here is the call graph for this function:

◆ BndIndexGlobal2Local()

index DNDS::Geom::UnstructuredMesh::BndIndexGlobal2Local ( DNDS::index  i)
inline

Definition at line 336 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ BndIndexGlobal2Local_NoSon()

index DNDS::Geom::UnstructuredMesh::BndIndexGlobal2Local_NoSon ( index  i)
inline

Definition at line 339 of file Mesh.hpp.

Here is the call graph for this function:

◆ BndIndexLocal2Global()

index DNDS::Geom::UnstructuredMesh::BndIndexLocal2Global ( DNDS::index  i)
inline

Definition at line 337 of file Mesh.hpp.

Here is the call graph for this function:

◆ BndIndexLocal2Global_NoSon()

index DNDS::Geom::UnstructuredMesh::BndIndexLocal2Global_NoSon ( index  i)
inline

Definition at line 338 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildBisectO1FormO2()

void DNDS::Geom::UnstructuredMesh::BuildBisectO1FormO2 ( UnstructuredMesh meshO2)

cell2cellOrig uses new global size later, cell2cellOrigV discarded

cell2cellOrig uses new global size now

bnd2bndOrig uses new global size now

Definition at line 476 of file Mesh_Elevation.cpp.

Here is the call graph for this function:

◆ BuildCell2CellFace()

void DNDS::Geom::UnstructuredMesh::BuildCell2CellFace ( )

Definition at line 1033 of file Mesh.cpp.

Here is the call graph for this function:

◆ BuildGhostFace()

void DNDS::Geom::UnstructuredMesh::BuildGhostFace ( )

Definition at line 1265 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildGhostN2CB()

void DNDS::Geom::UnstructuredMesh::BuildGhostN2CB ( )

Definition at line 938 of file Mesh.cpp.

◆ BuildGhostPrimary()

void DNDS::Geom::UnstructuredMesh::BuildGhostPrimary ( int  nGhostLayers = 1)

building ghost (son) from primary (currently only cell2cell)

the face and bnd parts are currently only local (no comm available) only builds comm data of cell and node cells: current-father and cell2cell neighbor (face or node neighbor) nodes: needed by all cells faces/bnds: needed by all father cells

Definition at line 666 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildGhostPrimaryLegacy()

void DNDS::Geom::UnstructuredMesh::BuildGhostPrimaryLegacy ( )

Definition at line 458 of file Mesh_Legacy.cpp.

Here is the call graph for this function:

◆ BuildNodeWallDist()

void DNDS::Geom::UnstructuredMesh::BuildNodeWallDist ( const std::function< bool(Geom::t_index)> &  fBndIsWall,
WallDistOptions  options = WallDistOptions{} 
)

TODO

Definition at line 12 of file Mesh_WallDist.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildO2FromO1Elevation()

void DNDS::Geom::UnstructuredMesh::BuildO2FromO1Elevation ( UnstructuredMesh meshO1)

now cell2cell, bnd2cell are lost

jumped because cell2cell, bnd2cell is to be recovered mandatorily

cannot use sorted

cannot use sorted

Definition at line 22 of file Mesh_Elevation.cpp.

Here is the call graph for this function:

◆ buildReorderPlan()

ReorderPlan DNDS::Geom::UnstructuredMesh::buildReorderPlan ( const ReorderInput input)

Build a ReorderPlan without applying it.

Useful for external code to obtain the plan and apply it to its own arrays after the mesh reorder.

Definition at line 482 of file Mesh_Reorder.cpp.

◆ buildReorderRegistry()

ReorderRegistry DNDS::Geom::UnstructuredMesh::buildReorderRegistry ( const std::unordered_set< EntityKind > &  destroyKinds = {})

Build a ReorderRegistry containing all mesh members.

Registers all built adj arrays (as type-erased callbacks) and all companion arrays (coords, cellElemInfo, pbi, etc.). Skips adjacencies involving destroyKinds.

External code may extend the returned registry with its own arrays before passing to ReorderPlan::build.

Definition at line 227 of file Mesh_Reorder.cpp.

Here is the call graph for this function:

◆ BuildVTKConnectivity()

void DNDS::Geom::UnstructuredMesh::BuildVTKConnectivity ( )

Definition at line 639 of file Mesh_Serial_BuildCell2Cell.cpp.

Here is the call graph for this function:

◆ CellFaceOther()

index DNDS::Geom::UnstructuredMesh::CellFaceOther ( index  iCell,
index  iFace 
) const
inline

Definition at line 845 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CellIndexGlobal2Local()

index DNDS::Geom::UnstructuredMesh::CellIndexGlobal2Local ( DNDS::index  i)
inline

Definition at line 328 of file Mesh.hpp.

Here is the call graph for this function:

◆ CellIndexGlobal2Local_NoSon()

index DNDS::Geom::UnstructuredMesh::CellIndexGlobal2Local_NoSon ( index  i)
inline

Definition at line 331 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CellIndexLocal2Global()

index DNDS::Geom::UnstructuredMesh::CellIndexLocal2Global ( DNDS::index  i)
inline

Definition at line 329 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CellIndexLocal2Global_NoSon()

index DNDS::Geom::UnstructuredMesh::CellIndexLocal2Global_NoSon ( index  i)
inline

Definition at line 330 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CellIsFaceBack()

bool DNDS::Geom::UnstructuredMesh::CellIsFaceBack ( index  iCell,
index  iFace 
) const
inline

Definition at line 839 of file Mesh.hpp.

Here is the caller graph for this function:

◆ CellOtherCellPeriodicHandle()

template<class FA , class FB , class F0 = std::function<void(void)>>
auto DNDS::Geom::UnstructuredMesh::CellOtherCellPeriodicHandle ( index  iFace,
rowsize  if2c,
FA &&  fA,
FB &&  fB,
F0 &&  f0 = []() {} 
)
inline

fA executes when if2c points to the donor side; fB the main side

Template Parameters
FA
FB
F0
Parameters
iFace
if2c
fA
fB
f0
Returns
auto

Definition at line 866 of file Mesh.hpp.

◆ ConstructBndMesh()

void DNDS::Geom::UnstructuredMesh::ConstructBndMesh ( UnstructuredMesh bMesh)

Definition at line 1560 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ConvertAdjEntries()

template<class TAdj , class TFn >
static void DNDS::Geom::UnstructuredMesh::ConvertAdjEntries ( TAdj adj,
index  nRows,
TFn &&  fn 
)
inlinestatic

Apply a conversion function to every entry of an adjacency array's first nRows rows.

Replaces the recurring nested-loop pattern:

for (index i = 0; i < nRows; i++)
for (rowsize j = 0; j < adj.RowSize(i); j++)
adj(i, j) = fn(adj(i, j));
int32_t rowsize
Row-width / per-row element-count type (signed 32-bit).
Definition Defines.hpp:114
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:112
UnstructuredMeshDeviceView< B > t_deviceView
Definition Mesh.hpp:1036
auto adj
Template Parameters
TAdjAny type with RowSize(index) and operator()(index, rowsize).
TFnCallable index(index).

Definition at line 366 of file Mesh.hpp.

◆ ConvertAdjEntriesOMP()

template<class TAdj , class TFn >
static void DNDS::Geom::UnstructuredMesh::ConvertAdjEntriesOMP ( TAdj adj,
index  nRows,
TFn &&  fn 
)
inlinestatic

OpenMP-parallelized variant of ConvertAdjEntries.

Same semantics as ConvertAdjEntries but with #pragma omp parallel for over the outer row loop. Use when the adjacency array is large and the conversion function is thread-safe.

Definition at line 381 of file Mesh.hpp.

◆ device_array_list_C2F()

auto DNDS::Geom::UnstructuredMesh::device_array_list_C2F ( )
inline

Definition at line 119 of file Mesh.hpp.

Here is the caller graph for this function:

◆ device_array_list_facial()

auto DNDS::Geom::UnstructuredMesh::device_array_list_facial ( )
inline

Definition at line 109 of file Mesh.hpp.

Here is the caller graph for this function:

◆ device_array_list_N2CB()

auto DNDS::Geom::UnstructuredMesh::device_array_list_N2CB ( )
inline

Definition at line 89 of file Mesh.hpp.

Here is the caller graph for this function:

◆ device_array_list_primary()

auto DNDS::Geom::UnstructuredMesh::device_array_list_primary ( )
inline

Definition at line 71 of file Mesh.hpp.

Here is the caller graph for this function:

◆ deviceView()

template<DeviceBackend B>
auto DNDS::Geom::UnstructuredMesh::deviceView ( )
inline

Definition at line 1039 of file Mesh.hpp.

◆ ElevatedNodesGetBoundarySmooth()

void DNDS::Geom::UnstructuredMesh::ElevatedNodesGetBoundarySmooth ( const std::function< bool(t_index)> &  FiFBndIdNeedSmooth)

angle is here

Definition at line 699 of file Mesh_Elevation.cpp.

Here is the call graph for this function:

◆ ElevatedNodesSolveInternalSmooth()

void DNDS::Geom::UnstructuredMesh::ElevatedNodesSolveInternalSmooth ( )

Definition at line 203 of file Mesh_Elevation_SmoothSolver.cpp.

Here is the call graph for this function:

◆ ElevatedNodesSolveInternalSmoothV1()

void DNDS::Geom::UnstructuredMesh::ElevatedNodesSolveInternalSmoothV1 ( )

Definition at line 487 of file Mesh_Elevation_SmoothSolver.cpp.

Here is the call graph for this function:

◆ ElevatedNodesSolveInternalSmoothV1Old()

void DNDS::Geom::UnstructuredMesh::ElevatedNodesSolveInternalSmoothV1Old ( )

Definition at line 293 of file Mesh_Elevation_SmoothSolver.cpp.

Here is the call graph for this function:

◆ ElevatedNodesSolveInternalSmoothV2()

void DNDS::Geom::UnstructuredMesh::ElevatedNodesSolveInternalSmoothV2 ( )

Definition at line 89 of file Mesh_Elevation_SmoothV2.cpp.

Here is the call graph for this function:

◆ EnsureGhostMapping()

template<class TPair >
void DNDS::Geom::UnstructuredMesh::EnsureGhostMapping ( TPair &  pair)
inline

Ensure a pair has a ghost mapping on its transformer.

If trans.pLGhostMapping is already set, does nothing. Otherwise, creates a father-only ghost mapping (empty ghost set) so that IndexGlobal2Local / IndexLocal2Global work even before ghost layers are built.

Precondition
trans.pLGlobalMapping must be set (call createFatherGlobalMapping first).
Warning
Collective — calls MPI_Alltoall internally when creating the mapping.

Definition at line 230 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ FaceIndexGlobal2Local()

index DNDS::Geom::UnstructuredMesh::FaceIndexGlobal2Local ( DNDS::index  i)
inline

Definition at line 344 of file Mesh.hpp.

Here is the call graph for this function:

◆ FaceIndexLocal2Global()

index DNDS::Geom::UnstructuredMesh::FaceIndexLocal2Global ( DNDS::index  i)
inline

Definition at line 345 of file Mesh.hpp.

Here is the call graph for this function:

◆ fillRegistry() [1/2]

void DNDS::Geom::UnstructuredMesh::fillRegistry ( MeshConnectivity dag) const

Populate a MeshConnectivity registry from this mesh's currently-built adjacencies.

Registers all adjacency arrays whose father is non-null. For each entity kind that appears as a source (.from) of any registered adjacency, finds a pLGlobalMapping from any adj array for that entity kind.

Precondition
All entity kinds that have registered adjacencies must have at least one adj array with a valid pLGlobalMapping on its father. Throws if this is not satisfied.
Parameters
dagMeshConnectivity to populate (meshDim is set).

Definition at line 1743 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ fillRegistry() [2/2]

void DNDS::Geom::UnstructuredMesh::fillRegistry ( MeshConnectivity dag,
const std::unordered_set< AdjKind, AdjKindHash > &  skip 
) const

Definition at line 1749 of file Mesh.cpp.

Here is the call graph for this function:

◆ for_each_device_member()

template<typename F >
void DNDS::Geom::UnstructuredMesh::for_each_device_member ( F &&  f)
inline

Definition at line 1027 of file Mesh.hpp.

Here is the call graph for this function:

◆ getArrayBytes()

index DNDS::Geom::UnstructuredMesh::getArrayBytes ( )
inline

Definition at line 1044 of file Mesh.hpp.

Here is the call graph for this function:

◆ getBnd2NodeIndexPbiRow()

auto DNDS::Geom::UnstructuredMesh::getBnd2NodeIndexPbiRow ( index  iBnd)
inline

Definition at line 196 of file Mesh.hpp.

◆ GetBndElement()

Elem::Element DNDS::Geom::UnstructuredMesh::GetBndElement ( index  iB)
inline

Definition at line 691 of file Mesh.hpp.

◆ GetBndZone()

t_index DNDS::Geom::UnstructuredMesh::GetBndZone ( index  iB)
inline

Definition at line 695 of file Mesh.hpp.

Here is the caller graph for this function:

◆ GetCell2CellFaceVLocal()

tLocalMatStruct DNDS::Geom::UnstructuredMesh::GetCell2CellFaceVLocal ( bool  onLocalPartition = false)
Returns
cell2cell for local mesh, which do not contain the diagonal part; should be a diag-less symmetric adjacency matrix

must be local not ghost ptrs

Definition at line 182 of file Mesh_Serial_Partition.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ getCell2NodeIndexPbiRow()

auto DNDS::Geom::UnstructuredMesh::getCell2NodeIndexPbiRow ( index  iCell)
inline

Definition at line 184 of file Mesh.hpp.

Here is the caller graph for this function:

◆ GetCellElement()

Elem::Element DNDS::Geom::UnstructuredMesh::GetCellElement ( index  iC)
inline

Definition at line 689 of file Mesh.hpp.

Here is the caller graph for this function:

◆ GetCellZone()

t_index DNDS::Geom::UnstructuredMesh::GetCellZone ( index  iC)
inline

Definition at line 693 of file Mesh.hpp.

◆ GetCoordNodeOnCell()

tPoint DNDS::Geom::UnstructuredMesh::GetCoordNodeOnCell ( index  iCell,
rowsize  ic2n 
)
inline

Definition at line 811 of file Mesh.hpp.

Here is the call graph for this function:

◆ GetCoordNodeOnFace()

tPoint DNDS::Geom::UnstructuredMesh::GetCoordNodeOnFace ( index  iFace,
rowsize  if2n 
)
inline

Definition at line 818 of file Mesh.hpp.

Here is the call graph for this function:

◆ GetCoordsOnCell() [1/2]

void DNDS::Geom::UnstructuredMesh::GetCoordsOnCell ( index  iCell,
tSmallCoords cs 
)
inline

Definition at line 787 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCoordsOnCell() [2/2]

void DNDS::Geom::UnstructuredMesh::GetCoordsOnCell ( index  iCell,
tSmallCoords cs,
tCoordPair coo 
)
inline

Definition at line 795 of file Mesh.hpp.

Here is the call graph for this function:

◆ GetCoordsOnFace()

void DNDS::Geom::UnstructuredMesh::GetCoordsOnFace ( index  iFace,
tSmallCoords cs 
)
inline

Definition at line 803 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetCoordWallDistOnCell()

tPoint DNDS::Geom::UnstructuredMesh::GetCoordWallDistOnCell ( index  iCell,
rowsize  ic2n 
)
inline

Definition at line 825 of file Mesh.hpp.

Here is the call graph for this function:

◆ GetCoordWallDistOnFace()

tPoint DNDS::Geom::UnstructuredMesh::GetCoordWallDistOnFace ( index  iFace,
rowsize  if2n 
)
inline

Definition at line 832 of file Mesh.hpp.

Here is the call graph for this function:

◆ getDim()

int DNDS::Geom::UnstructuredMesh::getDim ( ) const
inline

Definition at line 182 of file Mesh.hpp.

◆ GetFaceElement()

Elem::Element DNDS::Geom::UnstructuredMesh::GetFaceElement ( index  iF)
inline

Definition at line 690 of file Mesh.hpp.

Here is the caller graph for this function:

◆ GetFaceZone()

t_index DNDS::Geom::UnstructuredMesh::GetFaceZone ( index  iF)
inline

Definition at line 694 of file Mesh.hpp.

◆ getMPI()

MPIInfo & DNDS::Geom::UnstructuredMesh::getMPI ( )
inline

Definition at line 697 of file Mesh.hpp.

Here is the caller graph for this function:

◆ IndexGlobal2Local()

template<class TPair >
index DNDS::Geom::UnstructuredMesh::IndexGlobal2Local ( TPair &  pair,
DNDS::index  iGlobal 
)
inline

Global-to-local conversion using father+son ghost mapping.

Returns
local index, or (-1 - iGlobal) when not found in the pair. UnInitIndex passes through unchanged.

Definition at line 247 of file Mesh.hpp.

Here is the caller graph for this function:

◆ IndexGlobal2Local_NoSon()

template<class TPair >
index DNDS::Geom::UnstructuredMesh::IndexGlobal2Local_NoSon ( TPair &  pair,
index  iGlobal 
)
inline

Global-to-local conversion using only the father's global mapping (no son / ghost layer).

Returns
local index if the global maps to this rank, otherwise (-1 - iGlobal). UnInitIndex passes through.

Definition at line 304 of file Mesh.hpp.

Here is the caller graph for this function:

◆ IndexLocal2Global()

template<class TPair >
index DNDS::Geom::UnstructuredMesh::IndexLocal2Global ( TPair &  pair,
DNDS::index  iLocal 
)
inline

Local-to-global conversion using father+son ghost mapping.

Returns
global index, or the decoded global when the local is a negative "not-found" encoding. UnInitIndex passes through.

Definition at line 267 of file Mesh.hpp.

Here is the caller graph for this function:

◆ IndexLocal2Global_NoSon()

template<class TPair >
index DNDS::Geom::UnstructuredMesh::IndexLocal2Global_NoSon ( TPair &  pair,
index  iLocal 
)
inline

Local-to-global conversion using only the father's global mapping (no son / ghost layer).

Returns
global index. Asserts if iLocal exceeds father size. Negative "not-found" encoding and UnInitIndex pass through.

Definition at line 285 of file Mesh.hpp.

Here is the caller graph for this function:

◆ InterpolateFace()

void DNDS::Geom::UnstructuredMesh::InterpolateFace ( )

Definition at line 1094 of file Mesh.cpp.

Here is the call graph for this function:

◆ InterpolateFaceLegacy()

void DNDS::Geom::UnstructuredMesh::InterpolateFaceLegacy ( )

sync call!!!

Definition at line 610 of file Mesh_Legacy.cpp.

Here is the call graph for this function:

◆ IsO1()

bool DNDS::Geom::UnstructuredMesh::IsO1 ( )

Definition at line 291 of file Mesh.cpp.

Here is the call graph for this function:

◆ IsO2()

bool DNDS::Geom::UnstructuredMesh::IsO2 ( )

Definition at line 319 of file Mesh.cpp.

Here is the call graph for this function:

◆ LocalPartEnd()

index DNDS::Geom::UnstructuredMesh::LocalPartEnd ( int  iPart) const
inline

Definition at line 599 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ LocalPartStart()

index DNDS::Geom::UnstructuredMesh::LocalPartStart ( int  iPart) const
inline

Definition at line 598 of file Mesh.hpp.

Here is the caller graph for this function:

◆ MatchFaceBoundary()

void DNDS::Geom::UnstructuredMesh::MatchFaceBoundary ( )

Definition at line 1340 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ NLocalParts()

int DNDS::Geom::UnstructuredMesh::NLocalParts ( ) const
inline

Definition at line 597 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NodeIndexGlobal2Local()

index DNDS::Geom::UnstructuredMesh::NodeIndexGlobal2Local ( DNDS::index  i)
inline

Definition at line 320 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ NodeIndexGlobal2Local_NoSon()

index DNDS::Geom::UnstructuredMesh::NodeIndexGlobal2Local_NoSon ( index  i)
inline

Definition at line 323 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ NodeIndexLocal2Global()

index DNDS::Geom::UnstructuredMesh::NodeIndexLocal2Global ( DNDS::index  i)
inline

Definition at line 321 of file Mesh.hpp.

Here is the call graph for this function:

◆ NodeIndexLocal2Global_NoSon()

index DNDS::Geom::UnstructuredMesh::NodeIndexLocal2Global_NoSon ( index  i)
inline

Definition at line 322 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ NumBnd()

index DNDS::Geom::UnstructuredMesh::NumBnd ( ) const
inline

Definition at line 616 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumBndGhost()

index DNDS::Geom::UnstructuredMesh::NumBndGhost ( ) const
inline

Definition at line 637 of file Mesh.hpp.

◆ NumBndGlobal()

index DNDS::Geom::UnstructuredMesh::NumBndGlobal ( )
inline
Warning
must collectively call

Definition at line 683 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumBndProc()

index DNDS::Geom::UnstructuredMesh::NumBndProc ( ) const
inline

Definition at line 658 of file Mesh.hpp.

◆ NumCell()

index DNDS::Geom::UnstructuredMesh::NumCell ( ) const
inline

Definition at line 606 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumCellGhost()

index DNDS::Geom::UnstructuredMesh::NumCellGhost ( ) const
inline

Definition at line 627 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumCellGlobal()

index DNDS::Geom::UnstructuredMesh::NumCellGlobal ( )
inline
Warning
must collectively call

Definition at line 665 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumCellProc()

index DNDS::Geom::UnstructuredMesh::NumCellProc ( ) const
inline

Definition at line 648 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumFace()

index DNDS::Geom::UnstructuredMesh::NumFace ( ) const
inline

Definition at line 611 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumFaceGhost()

index DNDS::Geom::UnstructuredMesh::NumFaceGhost ( ) const
inline

Definition at line 632 of file Mesh.hpp.

◆ NumFaceGlobal()

index DNDS::Geom::UnstructuredMesh::NumFaceGlobal ( )
inline
Warning
must collectively call

Definition at line 677 of file Mesh.hpp.

◆ NumFaceProc()

index DNDS::Geom::UnstructuredMesh::NumFaceProc ( ) const
inline

Definition at line 653 of file Mesh.hpp.

◆ NumNode()

index DNDS::Geom::UnstructuredMesh::NumNode ( ) const
inline

Definition at line 601 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumNodeGhost()

index DNDS::Geom::UnstructuredMesh::NumNodeGhost ( ) const
inline

Definition at line 622 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumNodeGlobal()

index DNDS::Geom::UnstructuredMesh::NumNodeGlobal ( )
inline
Warning
must collectively call

Definition at line 671 of file Mesh.hpp.

Here is the caller graph for this function:

◆ NumNodeProc()

index DNDS::Geom::UnstructuredMesh::NumNodeProc ( ) const
inline

Definition at line 643 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ObtainLocalFactFillOrdering()

void DNDS::Geom::UnstructuredMesh::ObtainLocalFactFillOrdering ( Direct::SerialSymLUStructure symLU,
Direct::DirectPrecControl  control 
)

Definition at line 318 of file Mesh_Serial_Partition.cpp.

Here is the call graph for this function:

◆ ObtainSymmetricSymbolicFactorization()

void DNDS::Geom::UnstructuredMesh::ObtainSymmetricSymbolicFactorization ( Direct::SerialSymLUStructure symLU,
Direct::DirectPrecControl  control 
) const

Definition at line 1734 of file Mesh.cpp.

Here is the call graph for this function:

◆ op_on_device_arrays()

template<class F >
void DNDS::Geom::UnstructuredMesh::op_on_device_arrays ( F &&  f)
inline

Definition at line 1015 of file Mesh.hpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ PermuteRows()

template<class TPair , class TFn >
static void DNDS::Geom::UnstructuredMesh::PermuteRows ( TPair &  pair,
index  nRows,
TFn &&  old2new 
)
inlinestatic

Permute the father rows of an ArrayPair according to a mapping function.

For each old row index i in [0, nRows), the row is moved to position old2new(i). Works with both CSR (variable-row-size) and fixed-row-size array pairs: CSR pairs are Decompressed before and Compressed after the permutation; fixed-size pairs are copied directly.

Template Parameters
TPairAn ArrayPair-like type with father/son, IsCSR().
TFnCallable index(index) returning the new row index.

Definition at line 409 of file Mesh.hpp.

Here is the caller graph for this function:

◆ PrintMeshCGNS()

void DNDS::Geom::UnstructuredMesh::PrintMeshCGNS ( std::string  fname,
const t_FBCID_2_Name fbcid2name,
const std::vector< std::string > &  allNames 
)

Definition at line 1551 of file Mesh_Plts.cpp.

Here is the call graph for this function:

◆ PrintParallelVTKHDFDataArray()

void DNDS::Geom::UnstructuredMesh::PrintParallelVTKHDFDataArray ( std::string  fname,
std::string  seriesName,
int  arraySiz,
int  vecArraySiz,
int  arraySizPoint,
int  vecArraySizPoint,
const tFGetName names,
const tFGetData data,
const tFGetName vectorNames,
const tFGetVecData vectorData,
const tFGetName namesPoint,
const tFGetData dataPoint,
const tFGetName vectorNamesPoint,
const tFGetVecData vectorDataPoint,
double  t,
MPI_Comm  commDup 
)

Definition at line 1335 of file Mesh_Plts.cpp.

Here is the call graph for this function:

◆ ReadSerialize()

void DNDS::Geom::UnstructuredMesh::ReadSerialize ( Serializer::SerializerBaseSSP  serializerP,
const std::string &  name 
)

Definition at line 1466 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReadSerializeAndDistribute()

void DNDS::Geom::UnstructuredMesh::ReadSerializeAndDistribute ( Serializer::SerializerBaseSSP  serializerP,
const std::string &  name,
const PartitionOptions partitionOptions 
)

Reads mesh from an H5 serializer using even-split distribution, then repartitions via ParMetis for locality.

This method enables reading a mesh written with any number of MPI ranks into any number of MPI ranks, without requiring the original partition.

Flow:

  1. Even-split read of all primary arrays (coords, cell2node, cellElemInfo, bnd2node, bndElemInfo, cell2nodePbi, bnd2nodePbi, origIndex arrays).
  2. Build distributed node2cell and cell2cell (node-neighbor) using existing RecoverNode2CellAndNode2Bnd() + RecoverCell2CellAndBnd2Cell().
  3. Filter cell2cell to facial neighbors only (shared O1 vertices >= dim).
  4. Call ParMETIS_V3_PartKway on the distributed facial graph.
  5. Redistribute cells, nodes, and boundaries to the new partition via ArrayTransformer push-based transfer.
  6. Set adjPrimaryState = Adj_PointToGlobal.

After return, the caller should proceed with the standard rebuild sequence: RecoverNode2CellAndNode2Bnd -> RecoverCell2CellAndBnd2Cell -> BuildGhostPrimary -> AdjGlobal2LocalPrimary -> ...

Parameters
serializerPH5 serializer (must be collective / non-per-rank)
namePath name in the serializer hierarchy
partitionOptionsParMetis partitioning options

Definition at line 48 of file Mesh_ReadSerializeDistributed.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ RecoverCell2CellAndBnd2Cell()

void DNDS::Geom::UnstructuredMesh::RecoverCell2CellAndBnd2Cell ( )

needs to use RecoverNode2CellAndNode2Bnd before doing this. Requires node2cell.father and builds a version of its son.

Definition at line 422 of file Mesh.cpp.

Here is the call graph for this function:

◆ RecoverCell2CellAndBnd2CellLegacy()

void DNDS::Geom::UnstructuredMesh::RecoverCell2CellAndBnd2CellLegacy ( )

warning, this is not actual final official trans, just needed temporarily

warning, this is not actual final official trans, just needed temporarily

warning, this is not actual final official trans, just needed temporarily

Definition at line 195 of file Mesh_Legacy.cpp.

Here is the call graph for this function:

◆ RecoverNode2CellAndNode2Bnd()

void DNDS::Geom::UnstructuredMesh::RecoverNode2CellAndNode2Bnd ( )

only requires father part of cell2node, bnd2node and coords generates node2cell and node2bnd (father part)

Definition at line 378 of file Mesh.cpp.

Here is the call graph for this function:

◆ RecoverNode2CellAndNode2BndLegacy()

void DNDS::Geom::UnstructuredMesh::RecoverNode2CellAndNode2BndLegacy ( )

Definition at line 138 of file Mesh_Legacy.cpp.

Here is the call graph for this function:

◆ RecreatePeriodicNodes()

void DNDS::Geom::UnstructuredMesh::RecreatePeriodicNodes ( )

Definition at line 523 of file Mesh_Serial_BuildCell2Cell.cpp.

Here is the call graph for this function:

◆ ReorderEntities()

void DNDS::Geom::UnstructuredMesh::ReorderEntities ( const ReorderInput input)

Reorder entities using the general framework.

Builds a ReorderRegistry, computes follow maps, builds a ReorderPlan, applies it, rebuilds global mappings, and updates idx states.

Precondition
All adjacencies in Adj_PointToGlobal state.
Postcondition
All (non-destroyed) adjacencies in Adj_PointToGlobal. Ghost mappings stale (caller must rebuild ghosts). Global mappings fresh on reordered entities.
Warning
Collective.

Definition at line 493 of file Mesh_Reorder.cpp.

◆ ReorderLocalCells()

void DNDS::Geom::UnstructuredMesh::ReorderLocalCells ( int  nParts = 1,
int  nPartsInner = 1 
)
Warning
RecreatePeriodicNodes and BuildVTKConnectivity results are invalid after this;
bnd mesh's cell2parentCell is invalid after this

Definition at line 662 of file Mesh_Reorder.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ReorderLocalCellsLegacy()

void DNDS::Geom::UnstructuredMesh::ReorderLocalCellsLegacy ( int  nParts = 1,
int  nPartsInner = 1 
)

Legacy implementation preserved for reference/fallback.

Definition at line 1846 of file Mesh.cpp.

Here is the call graph for this function:

◆ resolveFollows()

ReorderInput DNDS::Geom::UnstructuredMesh::resolveFollows ( const ReorderInput input,
const ReorderRegistry reg 
)

Augment a ReorderInput with default follows and compute follow maps.

When input.follows is empty and Cell is explicitly reordered, default follows (Node/Bnd follow Cell) are added automatically. Returns a finalised ReorderInput with all follows resolved to explicit maps (ready for ReorderPlan::build).

Parameters
inputOriginal input (may have empty follows).
regRegistry with global mappings and adj data.
Returns
ReorderInput with all entity kinds as explicit maps.

Definition at line 391 of file Mesh_Reorder.cpp.

◆ SetHDF5OutSetting()

void DNDS::Geom::UnstructuredMesh::SetHDF5OutSetting ( size_t  chunkSiz,
int  deflateLevel,
bool  coll_on_data,
bool  coll_on_meta 
)
inline

Definition at line 976 of file Mesh.hpp.

◆ SetPeriodicGeometry()

void DNDS::Geom::UnstructuredMesh::SetPeriodicGeometry ( const tPoint translation1,
const tPoint rotationCenter1,
const tPoint eulerAngles1,
const tPoint translation2,
const tPoint rotationCenter2,
const tPoint eulerAngles2,
const tPoint translation3,
const tPoint rotationCenter3,
const tPoint eulerAngles3 
)

Definition at line 347 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ TransformCoords()

template<class TFTrans >
void DNDS::Geom::UnstructuredMesh::TransformCoords ( TFTrans &&  FTrans)
inline

Definition at line 953 of file Mesh.hpp.

Here is the call graph for this function:

◆ WriteSerialize()

void DNDS::Geom::UnstructuredMesh::WriteSerialize ( Serializer::SerializerBaseSSP  serializerP,
const std::string &  name 
)

Definition at line 1428 of file Mesh.cpp.

Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ adjC2CFaceState

MeshAdjState DNDS::Geom::UnstructuredMesh::adjC2CFaceState {Adj_Unknown}

Definition at line 52 of file Mesh.hpp.

◆ adjC2FState

MeshAdjState DNDS::Geom::UnstructuredMesh::adjC2FState {Adj_Unknown}

Definition at line 48 of file Mesh.hpp.

◆ adjFacialState

MeshAdjState DNDS::Geom::UnstructuredMesh::adjFacialState {Adj_Unknown}

Definition at line 46 of file Mesh.hpp.

◆ adjN2CBState

MeshAdjState DNDS::Geom::UnstructuredMesh::adjN2CBState {Adj_Unknown}

Definition at line 50 of file Mesh.hpp.

◆ adjPrimaryState

MeshAdjState DNDS::Geom::UnstructuredMesh::adjPrimaryState {Adj_Unknown}

Definition at line 44 of file Mesh.hpp.

◆ bnd2bndOrig

tAdj1Pair DNDS::Geom::UnstructuredMesh::bnd2bndOrig

Definition at line 66 of file Mesh.hpp.

◆ bnd2cell

AdjPairTracked<tAdj2Pair> DNDS::Geom::UnstructuredMesh::bnd2cell

Definition at line 60 of file Mesh.hpp.

◆ bnd2face

AdjPairTracked<tAdj1Pair> DNDS::Geom::UnstructuredMesh::bnd2face

Definition at line 101 of file Mesh.hpp.

◆ bnd2faceV

std::vector<index> DNDS::Geom::UnstructuredMesh::bnd2faceV

Definition at line 104 of file Mesh.hpp.

◆ bnd2node

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::bnd2node

Definition at line 59 of file Mesh.hpp.

◆ bnd2nodePbi

tPbiPair DNDS::Geom::UnstructuredMesh::bnd2nodePbi

Definition at line 69 of file Mesh.hpp.

◆ bndElemInfo

tElemInfoArrayPair DNDS::Geom::UnstructuredMesh::bndElemInfo

Definition at line 63 of file Mesh.hpp.

◆ cell2cell

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::cell2cell

Definition at line 61 of file Mesh.hpp.

◆ cell2cellFace

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::cell2cellFace

constructed on demand

Definition at line 127 of file Mesh.hpp.

◆ cell2cellFaceVLocalParts

tLocalMatStruct DNDS::Geom::UnstructuredMesh::cell2cellFaceVLocalParts

for cell local factorization

Definition at line 172 of file Mesh.hpp.

◆ cell2cellOrig

tAdj1Pair DNDS::Geom::UnstructuredMesh::cell2cellOrig

Definition at line 64 of file Mesh.hpp.

◆ cell2face

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::cell2face

interpolated

Definition at line 97 of file Mesh.hpp.

◆ cell2node

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::cell2node

Definition at line 58 of file Mesh.hpp.

◆ cell2nodePbi

tPbiPair DNDS::Geom::UnstructuredMesh::cell2nodePbi

periodic only, after reader

Definition at line 68 of file Mesh.hpp.

◆ cell2nodePeriodicRecreated

tAdjPair DNDS::Geom::UnstructuredMesh::cell2nodePeriodicRecreated

Definition at line 143 of file Mesh.hpp.

◆ cell2parentCell

std::vector<index> DNDS::Geom::UnstructuredMesh::cell2parentCell

Definition at line 132 of file Mesh.hpp.

◆ cellElemInfo

tElemInfoArrayPair DNDS::Geom::UnstructuredMesh::cellElemInfo

Definition at line 62 of file Mesh.hpp.

◆ coords

tCoordPair DNDS::Geom::UnstructuredMesh::coords

reader

Definition at line 57 of file Mesh.hpp.

◆ coordsElevDisp

tCoordPair DNDS::Geom::UnstructuredMesh::coordsElevDisp

only elevation

Definition at line 156 of file Mesh.hpp.

◆ coordsPeriodicRecreated

tCoordPair DNDS::Geom::UnstructuredMesh::coordsPeriodicRecreated

Definition at line 144 of file Mesh.hpp.

◆ dim

int DNDS::Geom::UnstructuredMesh::dim

Definition at line 41 of file Mesh.hpp.

◆ elevationInfo

struct DNDS::Geom::UnstructuredMesh::ElevationInfo DNDS::Geom::UnstructuredMesh::elevationInfo

◆ elevState

MeshElevationState DNDS::Geom::UnstructuredMesh::elevState = Elevation_Untouched

Definition at line 55 of file Mesh.hpp.

◆ face2bnd

AdjPairTracked<tAdj1Pair> DNDS::Geom::UnstructuredMesh::face2bnd

Definition at line 102 of file Mesh.hpp.

◆ face2bndM

std::unordered_map<index, index> DNDS::Geom::UnstructuredMesh::face2bndM

Definition at line 105 of file Mesh.hpp.

◆ face2cell

AdjPairTracked<tAdj2Pair> DNDS::Geom::UnstructuredMesh::face2cell

Definition at line 99 of file Mesh.hpp.

◆ face2node

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::face2node

Definition at line 98 of file Mesh.hpp.

◆ face2nodePbi

tPbiPair DNDS::Geom::UnstructuredMesh::face2nodePbi

periodic only, after interpolated

Definition at line 107 of file Mesh.hpp.

◆ faceElemInfo

tElemInfoArrayPair DNDS::Geom::UnstructuredMesh::faceElemInfo

Definition at line 100 of file Mesh.hpp.

◆ hdf5OutSetting

struct DNDS::Geom::UnstructuredMesh::HDF5OutSetting DNDS::Geom::UnstructuredMesh::hdf5OutSetting

◆ isPeriodic

bool DNDS::Geom::UnstructuredMesh::isPeriodic {false}

Definition at line 42 of file Mesh.hpp.

◆ localPartitionStarts

std::vector<index> DNDS::Geom::UnstructuredMesh::localPartitionStarts

Definition at line 174 of file Mesh.hpp.

◆ mpi

DNDS::MPIInfo DNDS::Geom::UnstructuredMesh::mpi

Definition at line 40 of file Mesh.hpp.

◆ mRank

MPI_int DNDS::Geom::UnstructuredMesh::mRank {0}

Definition at line 39 of file Mesh.hpp.

◆ nNodeO1

index DNDS::Geom::UnstructuredMesh::nNodeO1 {-1}

Definition at line 54 of file Mesh.hpp.

◆ node2bnd

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::node2bnd

Definition at line 87 of file Mesh.hpp.

◆ node2bndNode

std::vector<index> DNDS::Geom::UnstructuredMesh::node2bndNode

Definition at line 131 of file Mesh.hpp.

◆ node2cell

AdjPairTracked<tAdjPair> DNDS::Geom::UnstructuredMesh::node2cell

inverse relations

Definition at line 86 of file Mesh.hpp.

◆ node2nodeOrig

tAdj1Pair DNDS::Geom::UnstructuredMesh::node2nodeOrig

Definition at line 65 of file Mesh.hpp.

◆ node2parentNode

std::vector<index> DNDS::Geom::UnstructuredMesh::node2parentNode

parent built

Definition at line 130 of file Mesh.hpp.

◆ nodeRecreated2nodeLocal

std::vector<index> DNDS::Geom::UnstructuredMesh::nodeRecreated2nodeLocal

Definition at line 145 of file Mesh.hpp.

◆ nodeWallDist

tCoordPair DNDS::Geom::UnstructuredMesh::nodeWallDist

wall dist:

Definition at line 177 of file Mesh.hpp.

◆ nTotalMoved

index DNDS::Geom::UnstructuredMesh::nTotalMoved {-1}

Definition at line 157 of file Mesh.hpp.

◆ periodicInfo

Periodicity DNDS::Geom::UnstructuredMesh::periodicInfo

Definition at line 53 of file Mesh.hpp.

◆ vtkCell2node

std::vector<index> DNDS::Geom::UnstructuredMesh::vtkCell2node

Definition at line 137 of file Mesh.hpp.

◆ vtkCell2NodeGlobalSiz

index DNDS::Geom::UnstructuredMesh::vtkCell2NodeGlobalSiz {-1}

Definition at line 140 of file Mesh.hpp.

◆ vtkCell2nodeOffsets

std::vector<index> DNDS::Geom::UnstructuredMesh::vtkCell2nodeOffsets

for parallel out

Definition at line 135 of file Mesh.hpp.

◆ vtkCellOffset

index DNDS::Geom::UnstructuredMesh::vtkCellOffset {-1}

Definition at line 139 of file Mesh.hpp.

◆ vtkCellType

std::vector<uint8_t> DNDS::Geom::UnstructuredMesh::vtkCellType

Definition at line 136 of file Mesh.hpp.

◆ vtkNCellGlobal

index DNDS::Geom::UnstructuredMesh::vtkNCellGlobal {-1}

Definition at line 141 of file Mesh.hpp.

◆ vtkNNodeGlobal

index DNDS::Geom::UnstructuredMesh::vtkNNodeGlobal {-1}

Definition at line 142 of file Mesh.hpp.

◆ vtkNodeOffset

index DNDS::Geom::UnstructuredMesh::vtkNodeOffset {-1}

Definition at line 138 of file Mesh.hpp.


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