Dendro  5.01
Dendro in Greek language means tree. The Dendro library is a large scale (262K cores on ORNL's Titan) distributed memory adaptive octree framework. The main goal of Dendro is to perform large scale multiphysics simulations efficeiently in mordern supercomputers. Dendro consists of efficient parallel data structures and algorithms to perform variational ( finite element) methods and finite difference mthods on 2:1 balanced arbitary adaptive octrees which enables the users to perform simulations raning from black holes (binary black hole mergers) to blood flow in human body, where applications ranging from relativity, astrophysics to biomedical engineering.
Public Member Functions | List of all members
ot::Mesh Class Reference

#include <mesh.h>

Public Member Functions

 Mesh (std::vector< ot::TreeNode > &in, unsigned int k_s, unsigned int pOrder, unsigned int activeNpes, MPI_Comm comm, bool pBlockSetup=true, SM_TYPE smType=SM_TYPE::FDM, unsigned int grainSz=DENDRO_DEFAULT_GRAIN_SZ, double ld_tol=DENDRO_DEFAULT_LB_TOL, unsigned int sf_k=DENDRO_DEFAULT_SF_K)
 parallel mesh constructor More...
 
 Mesh (std::vector< ot::TreeNode > &in, unsigned int k_s, unsigned int pOrder, MPI_Comm comm, bool pBlockSetup=true, SM_TYPE smType=SM_TYPE::FDM, unsigned int grainSz=DENDRO_DEFAULT_GRAIN_SZ, double ld_tol=DENDRO_DEFAULT_LB_TOL, unsigned int sf_k=DENDRO_DEFAULT_SF_K)
 parallel mesh constructor More...
 
 ~Mesh ()
 destructor for mesh (releases the allocated variables in the class. )
 
void performBlocksSetup ()
 Perform the blocks initialization so that we can apply the stencil for the grid as a sequnce of finite number of regular grids. note that this should be called after performing all E2N and E2N mapping.
 
void buildF2EMap ()
 computes the face to element map. needs to be called after e2e and e2n maps has built.
 
void flagBlocks ()
 
bool isBlockSetep ()
 : returns if the block setup has performed or not
 
SM_TYPE getScatterMapType ()
 : returns if the scatter map typed set
 
unsigned int getNumLocalMeshElements () const
 
unsigned int getNumPreGhostElements () const
 
unsigned int getNumPostGhostElements () const
 Returns the number of post-ghost elements.
 
unsigned int getNumLocalMeshNodes () const
 return the number of nodes local to the mesh
 
unsigned int getNumPreMeshNodes () const
 return the number of pre ghost mesh nodes
 
unsigned int getNumPostMeshNodes () const
 return the number of post ghost mesh nodes
 
unsigned int getElementPreGhostBegin () const
 return the begin location of element pre ghost
 
unsigned int getElementPreGhostEnd () const
 return the end location of element pre ghost
 
unsigned int getElementLocalBegin () const
 return the begin location of element local
 
unsigned int getElementLocalEnd () const
 return the end location of element local
 
unsigned int getElementPostGhostBegin () const
 return the begin location of element post ghost
 
unsigned int getElementPostGhostEnd () const
 return the end location of element post ghost
 
unsigned int getNodePreGhostBegin () const
 return the begin location of pre ghost nodes
 
unsigned int getNodePreGhostEnd () const
 return the end location of pre ghost nodes
 
unsigned int getNodeLocalBegin () const
 return the location of local node begin
 
unsigned int getNodeLocalEnd () const
 return the location of local node end
 
unsigned int getNodePostGhostBegin () const
 return the location of post node begin
 
unsigned int getNodePostGhostEnd () const
 return the location of post node end
 
unsigned int getDegOfFreedom () const
 returns the dof for a partition
 
unsigned int getDegOfFreedomUnZip () const
 returns the dof for a partition
 
const std::vector< ot::TreeNode > & getAllElements () const
 returns the pointer to All elements array.
 
const std::vector< unsigned int > & getLevel1GhostElementIndices () const
 returns the Level 1 ghost element indices.
 
const std::vector< ot::TreeNode > & getSplitterElements () const
 returns the splitter elements of the mesh local elements.
 
const std::vector< ot::TreeNode > & getAllLocalNodes () const
 returns all local nodes(vertices)
 
const std::vector< unsigned int > & getE2EMapping () const
 returns const e2e mapping instance.
 
const std::vector< unsigned int > & getE2NMapping () const
 
const std::vector< unsigned int > & getDG2CGMap () const
 returns dg to cg map
 
const std::vector< unsigned int > & getCG2DGMap () const
 returns cg to dg map
 
const std::vector< unsigned int > & getE2NMapping_DG () const
 
const std::vector< ot::Block > & getLocalBlockList () const
 returns const list of local blocks (regular grids) for the consdering mesh.
 
unsigned int getNumDirections () const
 return the number of directions in the E2E mapping.
 
unsigned int getNumNodesPerElement () const
 return the number of nodes per element.
 
unsigned int getElementOrder () const
 returns the order of an element
 
MPI_Comm getMPICommunicator () const
 returns the communicator (acitve)
 
MPI_Comm getMPIGlobalCommunicator () const
 returns the global communicator
 
unsigned int getMPIRankGlobal () const
 returns the rank w.r.t. global comm
 
unsigned int getMPICommSizeGlobal () const
 returns the comm size w.r.t. global comm
 
unsigned int getMPIRank () const
 returns the rank
 
unsigned int getMPICommSize () const
 returns the comm size:
 
const RefElementgetReferenceElement () const
 returns const pointer to reference element
 
const unsigned int getSendProcListSize () const
 returns the send proc list size
 
const unsigned int getRecvProcListSize () const
 returns the recv proc list size
 
const std::vector< unsigned int > & getNodalSendCounts () const
 returns the nodal send counts
 
const std::vector< unsigned int > & getNodalSendOffsets () const
 returns the nodal send offsets
 
const std::vector< unsigned int > & getNodalRecvCounts () const
 returns the nodal recv counts
 
const std::vector< unsigned int > & getNodalRecvOffsets () const
 returns the nodal recv offsets
 
const std::vector< unsigned int > & getSendProcList () const
 returns the send proc. list
 
const std::vector< unsigned int > & getRecvProcList () const
 returns the recv proc. list
 
const std::vector< unsigned int > & getSendNodeSM () const
 return Scatter map for node send
 
const std::vector< unsigned int > & getRecvNodeSM () const
 return Scatter map for node send
 
void dg2eijk (unsigned int dg_index, unsigned int &e, unsigned int &i, unsigned int &j, unsigned int &k) const
 : Decompose the DG index to element id and it's i,j,k values.
 
unsigned int getMortonchildNum (unsigned int eleID) const
 returns the morton child number
 
bool isActive () const
 returns true if mesh is active
 
void waitAll () const
 waiting for all the mesh instances both active and inactive. This should not be called if not needed. this is a BARRIER.
 
void setOctreeRefineFlags (unsigned int *flags, unsigned int sz)
 set refinement flags for the octree. This is non const function More...
 
void getElementalFaceNeighbors (const unsigned int eID, const unsigned int dir, unsigned int *lookup) const
 returns the face neighours in specidied direction. More...
 
void getElementalEdgeNeighbors (const unsigned int eID, const unsigned int dir, unsigned int *lookup) const
 returns the edge neighours in specidied direction. More...
 
void getElementalVertexNeighbors (const unsigned int eID, const unsigned int dir, unsigned int *lookup) const
 returns the vertex neighours in specidied direction. More...
 
void getElementQMat (unsigned int currentId, double *&qMat, bool isAllocated=true) const
 : Compute the elemental interpolation matrix. More...
 
template<typename T >
void getUnzipElementalNodalValues (const T *uzipVec, unsigned int blkID, unsigned int ele, T *out, bool isPadded=true) const
 Get the elemental nodal values using unzip representation of the array. More...
 
template<ot::WaveletDA::LoopType type>
void init ()
 Initilizes the wavelet DA loop depending on the WaveletDA flags specified.
 
template<ot::WaveletDA::LoopType type>
bool nextAvailable ()
 Check whether the next element is available.
 
template<ot::WaveletDA::LoopType type>
void next ()
 Increment the counters to access the next element in the mesh.
 
const ot::TreeNodecurrentOctant ()
 Return the current element as an octant.
 
unsigned int currentIndex ()
 Returns the current element index.
 
void currentElementNeighbourIndexList (unsigned int *neighList)
 Returns the current neighbour list (Element) information. Note that for 3D it is 8 neighbours for each octant and for 2D it is 4.
 
void currentElementNodeList (unsigned int *nodeList)
 : Returns the node index list belongs to current element. NodeList size should be m_uiNpE;
 
void currentElementNodeList_DG (unsigned int *nodeList)
 Returns the node index list belogns to the currentl element in DG indexing. NodeList size should be m_uiNpE;.
 
void faceNodesIndex (unsigned int elementID, unsigned int face, std::vector< unsigned int > &index, bool isInternal) const
 Returns the index of m_uiE2NMapping (CG or DG) for a specified face. More...
 
void edgeNodeIndex (unsigned int elementID, unsigned int face1, unsigned int face2, std::vector< unsigned int > &index, bool isInternal) const
 Returns the index of m_uiE2NMapping (CG or DG) for a specified Edge. More...
 
void cornerNodeIndex (unsigned int elementID, unsigned int mortonIndex, unsigned int &index) const
 Returns the index of m_uiE2NMapping (CG or DG) for a specified Edge. More...
 
void elementNodeIndex (unsigned int elementID, std::vector< unsigned int > &index, bool isInternal) const
 Returns the all the internal node indices of an element. More...
 
bool isEdgeHanging (unsigned int elementId, unsigned int edgeId, unsigned int &cnum) const
 : Returns true or false (bases on specified edge is hanging or not.)for a given element id , and edge id. More...
 
bool isFaceHanging (unsigned int elementId, unsigned int faceId, unsigned int &cnum) const
 : Returns true or false (bases on specified edge is hanging or not.)for a given element id , and edge id. More...
 
bool isNodeHanging (unsigned int eleID, unsigned int ix, unsigned int jy, unsigned int kz) const
 : Returns true if the specified node (e,i,j,k) is hanging. More...
 
bool isNodeLocal (unsigned int eleID, unsigned int ix, unsigned int jy, unsigned int kz) const
 : Returns true if the specified node (e,i,j,k) is local. More...
 
bool isBoundaryOctant (unsigned int ele) const
 
DendroIntL getGhostExcgTotalSendNodeCount () const
 
DendroIntL getGhostExcgTotalRecvNodeCount () const
 
const ot::TreeNodegetNodalSplitterNodes () const
 get plitter nodes for each processor.
 
template<typename T >
T * createVector () const
 allocate memory for variable array based on the adaptive mesh
 
template<typename T >
void createVector (std::vector< T > &vec) const
 allocate memory for variable array based on the adaptive mesh.
 
template<typename T >
T * createVector (const T initValue) const
 allocate memory for variable array based on the adaptive mesh. More...
 
template<typename T >
T * createVector (std::function< T(T, T, T)> func) const
 allocate memory for variable array based on the adaptive mesh. More...
 
template<typename T >
void createVector (std::vector< T > &vec, const T initValue) const
 allocate memory for variable array based on the adaptive mesh. More...
 
template<typename T >
void createVector (std::vector< T > &vec, std::function< T(T, T, T)> func) const
 create and initialize local elements based on a given function. More...
 
template<typename T >
void createUnZippedVector (std::vector< T > &uvec) const
 : creates (memory allocation) for the unzipped version of the vector. More...
 
template<typename T >
T * createUnZippedVector () const
 : creates (memory allocation) for the unzipped version of the vector. More...
 
template<typename T >
void createUnZippedVector (std::vector< T > &uvec, const T initValue) const
 allocate memory for variable array based on the adaptive mesh. More...
 
template<typename T >
T * createUnZippedVector (const T initValue) const
 allocate memory for variable array based on the adaptive mesh. More...
 
void parent2ChildInterpolation (const double *in, double *out, unsigned int cnum, unsigned int dim=3) const
 Performs all parent to child interpolations for the m_uiEl_i element in order to apply the stencil. More...
 
void child2ParentInterpolation (const double *in, double *out, unsigned int cnum, unsigned int dim=3) const
 performs the child to parent contribution (only from a single child). More...
 
void child2ParentInjection (double *in, double *out, bool *isHanging) const
 Performs child to parent injection. More...
 
template<typename T >
void unzip (const T *zippedVec, T *unzippedVec)
 Creates the decomposition of adaptive octree variables into blocklist variables that we computed. More...
 
template<typename T >
void unzip_async (T *zippedVec, T *unzippedVec, MPI_Request *send_reqs, MPI_Request *recv_reqs, MPI_Status *send_sts, MPI_Status *recv_sts)
 perform the unzip operation overlap with the communication. GHOST exchage done inside. User doesn't need to perofrm ghost exchange. More...
 
template<typename T >
void zip (const T *unzippedVec, T *zippedVec)
 Performs the compression frrom regular block grid varable list to adaptive representation. More...
 
template<typename T , unsigned int length, unsigned int offsetCentered, unsigned int offsetBackward, unsigned int offsetForward>
void applyStencil (const std::vector< T > &in, std::vector< T > &out, const Stencil< T, length, offsetCentered > &centered, const Stencil< T, length, offsetBackward > &backward, const Stencil< T, length, offsetForward > &forward)
 Apply a given stencil to for provided variable array. More...
 
template<typename T >
void performGhostExchange (std::vector< T > &vec)
 Perform the ghost exchange for the vector vec. More...
 
template<typename T >
void performGhostExchange (T *vec)
 Perform the ghost exchange for the vector vec. More...
 
template<typename T >
void ghostExchangeStart (T *vec, T *sendNodeBuffer, T *recvNodeBuffer, MPI_Request *send_reqs, MPI_Request *recv_reqs)
 Perform the ghost asynchronous send for the vector vec Note: this is a non-blocking asynchronous communication. User is resposible to call the (synchronous call such as MPI_WaitAll) when overlapping the communication and computation. More...
 
template<typename T >
void ghostExchangeRecvSync (T *vec, T *recvNodeBuffer, MPI_Request *recv_reqs, MPI_Status *recv_sts)
 Perform the wait on the recv requests. More...
 
void ghostExchangeSendSync (MPI_Request *send_reqs, MPI_Status *send_sts)
 Perform the wait on the recv requests. More...
 
template<typename T >
void vectorToVTK (const std::vector< T > &vec, char *fprefix, double pTime=0.0, unsigned int nCycle=0) const
 write out function values to a vtk file. More...
 
template<typename T >
bool isReMesh (const T **vec, const unsigned int *varIds, const unsigned int numVars, double tol, double amr_coarse_fac=DENDRO_AMR_COARSEN_FAC)
 : determine whether any refinement or coarsening need for a specified set of elements. (This uses zipped version of the varibles which needs to satiesfy some constraints. Not every element is eligible for refinemenet or coarsening). More...
 
template<typename T >
bool isReMeshUnzip (const T **unzippedVec, const unsigned int *varIds, const unsigned int numVars, std::function< double(double, double, double)> wavelet_tol, double amr_coarse_fac=DENDRO_AMR_COARSEN_FAC, double coarsen_hx=DENDRO_REMESH_UNZIP_SCALE_FAC)
 
ot::MeshReMesh (unsigned int grainSz=DENDRO_DEFAULT_GRAIN_SZ, double ld_tol=DENDRO_DEFAULT_LB_TOL, unsigned int sfK=DENDRO_DEFAULT_SF_K, unsigned int(*getWeight)(const ot::TreeNode *)=NULL)
 : Remesh the mesh with the new computed elements. More...
 
template<typename T >
void interGridTransfer (std::vector< T > &vec, const ot::Mesh *pMesh)
 transfer a variable vector form old grid to new grid. More...
 
template<typename T >
void interGridTransfer (T *&vec, const ot::Mesh *pMesh)
 transfer a variable vector form old grid to new grid. More...
 
template<typename T >
void interGridTransferUnzip (T *&unzip, T *&vec, const ot::Mesh *pMesh)
 performs the intergrid transfer operation using the unzip representation. Compared to elemental intergrid transfer this uses information from neighbouring elements while the transfer occurs. More...
 
template<typename T >
void getElementNodalValues (const T *vec, T *nodalValues, unsigned int elementID) const
 : Returns the nodal values of a given element for a given variable vector. More...
 
template<typename T >
void computeElementalContribution (const T *in, T *out, unsigned int elementID) const
 : Computes the contribution of elemental nodal values to the parent elements if it is hanging. Note: internal nodes for the elements cannnot be hagging. Only the face edge nodes are possible for hanging. More...
 
void getElementCoordinates (unsigned int eleID, double *coords) const
 computes the elementCoordinates (based on the nodal placement) More...
 
template<typename T >
int getFaceNeighborValues (unsigned int eleID, const T *in, T *out, T *coords, unsigned int *neighID, unsigned int face, NeighbourLevel &level) const
 computes the face neighbor points for additional computations for a specified direction. More...
 
EType getElementType (unsigned int eleID)
 returns the type of the element.
 
int getBlkBdyParentNodeIndices (unsigned int blkId, unsigned int eleId, unsigned int dir, unsigned int *nid, unsigned int *child, unsigned int *fid, unsigned int *cid)
 compute the block boundary parent nodal locations. More...
 

Detailed Description

Contains all the information needed to build neighbourhood information for the balanced octree based on treeSearch.

Constructor & Destructor Documentation

◆ Mesh() [1/2]

ot::Mesh::Mesh ( std::vector< ot::TreeNode > &  in,
unsigned int  k_s,
unsigned int  pOrder,
unsigned int  activeNpes,
MPI_Comm  comm,
bool  pBlockSetup = true,
SM_TYPE  smType = SM_TYPE::FDM,
unsigned int  grainSz = DENDRO_DEFAULT_GRAIN_SZ,
double  ld_tol = DENDRO_DEFAULT_LB_TOL,
unsigned int  sf_k = DENDRO_DEFAULT_SF_K 
)

parallel mesh constructor

Parameters
[in]incomplete sorted 2:1 balanced octree to generate mesh
[in]k_show many neighbours to check on each direction (used =1)
[in]pOrderorder of an element.
[in]commActiveMPI active communicator
[in]commMPI communicator (global)

◆ Mesh() [2/2]

ot::Mesh::Mesh ( std::vector< ot::TreeNode > &  in,
unsigned int  k_s,
unsigned int  pOrder,
MPI_Comm  comm,
bool  pBlockSetup = true,
SM_TYPE  smType = SM_TYPE::FDM,
unsigned int  grainSz = DENDRO_DEFAULT_GRAIN_SZ,
double  ld_tol = DENDRO_DEFAULT_LB_TOL,
unsigned int  sf_k = DENDRO_DEFAULT_SF_K 
)

parallel mesh constructor

Parameters
[in]incomplete sorted 2:1 balanced octree to generate mesh
[in]k_show many neighbours to check on each direction (used =1)
[in]pOrderorder of an element.
[in]commMPI communicator (global)
[in]grainSzprefered grain sz. (this parameter is used to perform automatic comm expansion and shrinking)
[in]ld_tolload imbalance tolerance for comm expansion and shrinking
[in]sf_ksplitter fix _k value. (Needed by SFC_partitioinng for large p>=10,000)

Member Function Documentation

◆ applyStencil()

template<typename T , unsigned int length, unsigned int offsetCentered, unsigned int offsetBackward, unsigned int offsetForward>
void ot::Mesh::applyStencil ( const std::vector< T > &  in,
std::vector< T > &  out,
const Stencil< T, length, offsetCentered > &  centered,
const Stencil< T, length, offsetBackward > &  backward,
const Stencil< T, length, offsetForward > &  forward 
)

Apply a given stencil to for provided variable array.

Parameters
[in]in: vector that we need to apply the stencil on.
[in]centeredStencil that we need to apply on vector in, this is the centered stencil.
[in]backwardbackward version of the centered stencil.(This is used for boundary elements. )
[in]forwardfoward version of the centered stencil. (This is used for boundary elements. )
[out]outoutput vector that after applying the stencil.

◆ child2ParentInjection()

void ot::Mesh::child2ParentInjection ( double *  in,
double *  out,
bool *  isHanging 
) const
inline

Performs child to parent injection.

Parameters
[in]childrenfunction values. (all the function values of children ordered according to the SFC ordering)
[in]isHangingarray of boolean variables specifying each node is hanging ot not.
[out]injectedvalues back to parent.

◆ child2ParentInterpolation()

void ot::Mesh::child2ParentInterpolation ( const double *  in,
double *  out,
unsigned int  cnum,
unsigned int  dim = 3 
) const
inline

performs the child to parent contribution (only from a single child).

Parameters
[in]inchild function values
[in]cnummorton ID of the current child.
[in]dimdim of the interpolation. (dim=1 for edge interpolation , dim=2 for face interpolation, dim=3 for octant to child interpolation.)
[out]outinterpolated values. (child to parent contribution)

◆ computeElementalContribution()

template<typename T >
void ot::Mesh::computeElementalContribution ( const T *  in,
T *  out,
unsigned int  elementID 
) const

: Computes the contribution of elemental nodal values to the parent elements if it is hanging. Note: internal nodes for the elements cannnot be hagging. Only the face edge nodes are possible for hanging.

: input is the elemental nodal values.

Parameters
[in]vecchild var vector (nPe)
[in]elementIDelement ID of the current element (or child octant)
[out]outadd the contributions to the current vector accordingly.

Usage: This is needed when performing matrix-free matvec for FEM method.

◆ cornerNodeIndex()

void ot::Mesh::cornerNodeIndex ( unsigned int  elementID,
unsigned int  mortonIndex,
unsigned int &  index 
) const
inline

Returns the index of m_uiE2NMapping (CG or DG) for a specified Edge.

Parameters
[in]elementIDelement ID of the edge, that the face belongs to .
[in]mortonIndexmorton ID of the corner node in the cordinate change in the order of x y z.
[out]indexreturns the indecies of the requested face.

◆ createUnZippedVector() [1/4]

template<typename T >
void ot::Mesh::createUnZippedVector ( std::vector< T > &  uvec) const

: creates (memory allocation) for the unzipped version of the vector.

Parameters
[in]

◆ createUnZippedVector() [2/4]

template<typename T >
T* ot::Mesh::createUnZippedVector ( ) const

: creates (memory allocation) for the unzipped version of the vector.

Parameters
[in]

◆ createUnZippedVector() [3/4]

template<typename T >
void ot::Mesh::createUnZippedVector ( std::vector< T > &  uvec,
const T  initValue 
) const

allocate memory for variable array based on the adaptive mesh.

Parameters
[in]uvecallocate memory for uvec (unzipped version)
[in]initValueinitialize the vector to the given value.

◆ createUnZippedVector() [4/4]

template<typename T >
T* ot::Mesh::createUnZippedVector ( const T  initValue) const

allocate memory for variable array based on the adaptive mesh.

Parameters
[in]uvecallocate memory for uvec (unzipped version)
[in]initValueinitialize the vector to the given value.

◆ createVector() [1/4]

template<typename T >
T* ot::Mesh::createVector ( const T  initValue) const

allocate memory for variable array based on the adaptive mesh.

Parameters
[in]initValueinitialize the vector to the given value.

◆ createVector() [2/4]

template<typename T >
T* ot::Mesh::createVector ( std::function< T(T, T, T)>  func) const

allocate memory for variable array based on the adaptive mesh.

Parameters
[in]initValueinitialize the vector to the given value.

◆ createVector() [3/4]

template<typename T >
void ot::Mesh::createVector ( std::vector< T > &  vec,
const T  initValue 
) const

allocate memory for variable array based on the adaptive mesh.

Parameters
[in]vecallocate memory for vec
[in]initValueinitialize the vector to the given value.

◆ createVector() [4/4]

template<typename T >
void ot::Mesh::createVector ( std::vector< T > &  vec,
std::function< T(T, T, T)>  func 
) const

create and initialize local elements based on a given function.

Parameters
[in]vecmesh varaible vector
[in]funcfunction to initialize with.

◆ edgeNodeIndex()

void ot::Mesh::edgeNodeIndex ( unsigned int  elementID,
unsigned int  face1,
unsigned int  face2,
std::vector< unsigned int > &  index,
bool  isInternal 
) const
inline

Returns the index of m_uiE2NMapping (CG or DG) for a specified Edge.

Parameters
[in]elementIDelement ID of the edge, that the face belongs to .
[in]face1: one of the face, that an edge belongs to
[in]face2second face that an edge belongs to . Hence the edge in cosideration will be intersection of the two faces, face1 and face2.
[in]isInternalreturn only the internal nodes if true hence the index size would be ((m_uiElementOrder-1)*(m_uiElementOrder-1))
[out]indexreturns the indecies of the requested face.

◆ elementNodeIndex()

void ot::Mesh::elementNodeIndex ( unsigned int  elementID,
std::vector< unsigned int > &  index,
bool  isInternal 
) const
inline

Returns the all the internal node indices of an element.

Parameters
[in]elementIDelement ID of the edge, that the face belongs to .
[in]isInternalreturn only the internal nodes if true hence the index size would be ((m_uiElementOrder-1)*(m_uiElementOrder-1))
[out]indexreturns the indecies of the requested face.

◆ faceNodesIndex()

void ot::Mesh::faceNodesIndex ( unsigned int  elementID,
unsigned int  face,
std::vector< unsigned int > &  index,
bool  isInternal 
) const
inline

Returns the index of m_uiE2NMapping (CG or DG) for a specified face.

Parameters
[in]elementIDelement ID of the face, that the face belongs to .
[in]face: Face that you need the indexing,
[in]isInternalreturn only the internal nodes if true hence the index size would be ((m_uiElementOrder-1)*(m_uiElementOrder-1))
[out]indexreturns the indecies of the requested face.

◆ getBlkBdyParentNodeIndices()

int ot::Mesh::getBlkBdyParentNodeIndices ( unsigned int  blkId,
unsigned int  eleId,
unsigned int  dir,
unsigned int *  nid,
unsigned int *  child,
unsigned int *  fid,
unsigned int *  cid 
)

compute the block boundary parent nodal locations.

Parameters
[in]blkId: block id.
[in]eleId: element id of the block
[in]dirface direction.
[out]nid(assumed to be allocated) nodal id values of the parent which containing, block elements and the parent element in the paddinf region.
[out]childelement id for the block boudary faces, (elements containing inside the block)
[out]fidreference pointer to the child array for the finer elements.
[out]cidreference pointer to the child array for the coarser elements. (if it was refined. (cnumbers reference to the coarser elements))

◆ getE2NMapping()

const std::vector<unsigned int>& ot::Mesh::getE2NMapping ( ) const
inline

returns const e2n mapping instance

◆ getE2NMapping_DG()

const std::vector<unsigned int>& ot::Mesh::getE2NMapping_DG ( ) const
inline

returns const e2n mapping instance (debuging purposes only)

◆ getElementalEdgeNeighbors()

void ot::Mesh::getElementalEdgeNeighbors ( const unsigned int  eID,
const unsigned int  dir,
unsigned int *  lookup 
) const

returns the edge neighours in specidied direction.

Parameters
[in]eIDElement ID.
[in]dirdirection of the edge.
[out]lookUpresult if the result is a same level or lower level than element ID. (only 4 neighbor)

◆ getElementalFaceNeighbors()

void ot::Mesh::getElementalFaceNeighbors ( const unsigned int  eID,
const unsigned int  dir,
unsigned int *  lookup 
) const

returns the face neighours in specidied direction.

Parameters
[in]eIDElement ID.
[in]dirdirection of the face.
[out]lookUpresult if the result is a same level or lower level than element ID. (only 2 neighbor)

◆ getElementalVertexNeighbors()

void ot::Mesh::getElementalVertexNeighbors ( const unsigned int  eID,
const unsigned int  dir,
unsigned int *  lookup 
) const

returns the vertex neighours in specidied direction.

Parameters
[in]eIDElement ID.
[in]dirdirection of the vertex.
[out]lookUpresult if the result is a same level or lower level than element ID. (only 8 neighbor)

◆ getElementCoordinates()

void ot::Mesh::getElementCoordinates ( unsigned int  eleID,
double *  coords 
) const

computes the elementCoordinates (based on the nodal placement)

Parameters
[in]eleID: element ID

◆ getElementNodalValues()

template<typename T >
void ot::Mesh::getElementNodalValues ( const T *  vec,
T *  nodalValues,
unsigned int  elementID 
) const

: Returns the nodal values of a given element for a given variable vector.

Parameters
[in]vecvariable vector that we want to get the nodal values.
[in]elementIDelement ID that we need to get the nodal values.
[out]nodalValuesnodal values of the specified element ID

◆ getElementQMat()

void ot::Mesh::getElementQMat ( unsigned int  currentId,
double *&  qMat,
bool  isAllocated = true 
) const

: Compute the elemental interpolation matrix.

Parameters
[in]currentIdElement ID.

◆ getFaceNeighborValues()

template<typename T >
int ot::Mesh::getFaceNeighborValues ( unsigned int  eleID,
const T *  in,
T *  out,
T *  coords,
unsigned int *  neighID,
unsigned int  face,
NeighbourLevel &  level 
) const

computes the face neighbor points for additional computations for a specified direction.

Parameters
[in]eleIDelement ID
[in]ininpute vector
[out]outoutput vector values are in the order of the x,y,z size : 4*NodesPerElement
[out]coordsget the corresponding coordinates size: 4*NodesPerElement*m_uiDim;
[out]neighIDface neighbor octant IDs,
[in]faceface direction in {OCT_DIR_LEFT,OCT_IDR_RIGHT,OCT_DIR_DOWN, OCT_DIR_UP,OCT_DIR_BACK,OCT_DIR_FRONT}
[out]levelthe level of the neighbour octant with respect to the current octant. returns the number of face neighbours 1/4 for 3D.

◆ getNumLocalMeshElements()

unsigned int ot::Mesh::getNumLocalMeshElements ( ) const
inline

Returns the number of local elements in the grid. (local to the current considering processor)

◆ getNumPreGhostElements()

unsigned int ot::Mesh::getNumPreGhostElements ( ) const
inline

Returns the number of pre-ghost elements

◆ getUnzipElementalNodalValues()

template<typename T >
void ot::Mesh::getUnzipElementalNodalValues ( const T *  uzipVec,
unsigned int  blkID,
unsigned int  ele,
T *  out,
bool  isPadded = true 
) const

Get the elemental nodal values using unzip representation of the array.

Template Parameters
Ttype of the vector.
Parameters
uzipVec: unzip vector
blkID: block ID the element belongs to
ele: element ID of the block
out: output values. (allocated with corresponding padding width)
isPadded: true if the we need the elemental values with padding.

◆ ghostExchangeRecvSync()

template<typename T >
void ot::Mesh::ghostExchangeRecvSync ( T *  vec,
T *  recvNodeBuffer,
MPI_Request *  recv_reqs,
MPI_Status *  recv_sts 
)

Perform the wait on the recv requests.

Parameters
[in]vecadaptive mesh vector containing the values.
[in]recv_reqsm_uiRecvProcList.size() recv request
[in]recv_stsm_uiRecvProcList.size() recv status

◆ ghostExchangeSendSync()

void ot::Mesh::ghostExchangeSendSync ( MPI_Request *  send_reqs,
MPI_Status *  send_sts 
)
inline

Perform the wait on the recv requests.

Parameters
[in]vecadaptive mesh vector containing the values.
[in]send_reqsm_uiSendProcList.size() send request
[in]send_stsm_uiSendProcList.size() send status

◆ ghostExchangeStart()

template<typename T >
void ot::Mesh::ghostExchangeStart ( T *  vec,
T *  sendNodeBuffer,
T *  recvNodeBuffer,
MPI_Request *  send_reqs,
MPI_Request *  recv_reqs 
)

Perform the ghost asynchronous send for the vector vec Note: this is a non-blocking asynchronous communication. User is resposible to call the (synchronous call such as MPI_WaitAll) when overlapping the communication and computation.

Parameters
[in]vecadaptive mesh vector contaiting the values.

◆ interGridTransfer() [1/2]

template<typename T >
void ot::Mesh::interGridTransfer ( std::vector< T > &  vec,
const ot::Mesh pMesh 
)

transfer a variable vector form old grid to new grid.

Parameters
[in]vecvariable vector needs to be transfered.
[out]vectransfered varaible vector
[in]pMeshMesh that we need to transfer the old varaible.

◆ interGridTransfer() [2/2]

template<typename T >
void ot::Mesh::interGridTransfer ( T *&  vec,
const ot::Mesh pMesh 
)

transfer a variable vector form old grid to new grid.

Parameters
[in]vecvariable vector needs to be transfered.
[out]vectransfered varaible vector
[in]pMeshMesh that we need to transfer the old varaible.

◆ interGridTransferUnzip()

template<typename T >
void ot::Mesh::interGridTransferUnzip ( T *&  unzip,
T *&  vec,
const ot::Mesh pMesh 
)

performs the intergrid transfer operation using the unzip representation. Compared to elemental intergrid transfer this uses information from neighbouring elements while the transfer occurs.

Template Parameters
T: data type of the vectors
Parameters
unzip: unizp representation of the vectors.
vec: zip representation of the transfered vectors.
pMesh: new mesh.

◆ isEdgeHanging()

bool ot::Mesh::isEdgeHanging ( unsigned int  elementId,
unsigned int  edgeId,
unsigned int &  cnum 
) const

: Returns true or false (bases on specified edge is hanging or not.)for a given element id , and edge id.

Parameters
[in]elementIdelement ID of the octant.
[in]edgeIdedge id

◆ isFaceHanging()

bool ot::Mesh::isFaceHanging ( unsigned int  elementId,
unsigned int  faceId,
unsigned int &  cnum 
) const

: Returns true or false (bases on specified edge is hanging or not.)for a given element id , and edge id.

Parameters
[in]elementIdelement ID of the octant.
[in]faceIdface id

◆ isNodeHanging()

bool ot::Mesh::isNodeHanging ( unsigned int  eleID,
unsigned int  ix,
unsigned int  jy,
unsigned int  kz 
) const

: Returns true if the specified node (e,i,j,k) is hanging.

Parameters
[in]eleIDelement ID
[in]ixi-index of the node.
[in]jyj-index of the node.
[in]kzk-index of the node.

◆ isNodeLocal()

bool ot::Mesh::isNodeLocal ( unsigned int  eleID,
unsigned int  ix,
unsigned int  jy,
unsigned int  kz 
) const
inline

: Returns true if the specified node (e,i,j,k) is local.

Parameters
[in]eleIDelement ID
[in]ixi-index of the node.
[in]jyj-index of the node.
[in]kzk-index of the node.

◆ isReMesh()

template<typename T >
bool ot::Mesh::isReMesh ( const T **  vec,
const unsigned int *  varIds,
const unsigned int  numVars,
double  tol,
double  amr_coarse_fac = DENDRO_AMR_COARSEN_FAC 
)

: determine whether any refinement or coarsening need for a specified set of elements. (This uses zipped version of the varibles which needs to satiesfy some constraints. Not every element is eligible for refinemenet or coarsening).

Parameters
[in]vecsequence of varaibles to check
[in]varIdsvariable ids to check. (var ids to index the vec, vec[i] is a T* pointintg to one of the variable in vec. )
[in]numVarsnumber of variables to check
[in]tolwavelet tolerance Returns true if specified variable violates the specified wavelet toerlance.
Note
: this method will flag every element in the mesh with OCT_NO_CHANGE, OCT_SPLIT, OCT_COARSE.

◆ isReMeshUnzip()

template<typename T >
bool ot::Mesh::isReMeshUnzip ( const T **  unzippedVec,
const unsigned int *  varIds,
const unsigned int  numVars,
std::function< double(double, double, double)>  wavelet_tol,
double  amr_coarse_fac = DENDRO_AMR_COARSEN_FAC,
double  coarsen_hx = DENDRO_REMESH_UNZIP_SCALE_FAC 
)
Parameters
[in]vecsequence of varaibles to check
[in]varIdsvariable ids to check. (var ids to index the vec, vec[i] is a T* pointintg to one of the variable in vec. )
[in]numVarsnumber of variables to check
[in]tolwavelet tolerance Returns true if specified variable violates the specified wavelet toerlance.
Note
: this method will flag every element in the mesh with OCT_NO_CHANGE, OCT_SPLIT, OCT_COARSE.

◆ parent2ChildInterpolation()

void ot::Mesh::parent2ChildInterpolation ( const double *  in,
double *  out,
unsigned int  cnum,
unsigned int  dim = 3 
) const
inline

Performs all parent to child interpolations for the m_uiEl_i element in order to apply the stencil.

Parameters
[in]parentfunction values
[in]cnumchild number to interpolate.
[in]dimdim of the interpolation. (dim=1 for edge interpolation , dim=2 for face interpolation, dim=3 for octant to child interpolation.)
[out]outinterpolated values.

◆ performGhostExchange() [1/2]

template<typename T >
void ot::Mesh::performGhostExchange ( std::vector< T > &  vec)

Perform the ghost exchange for the vector vec.

Parameters
[in]vecadaptive mesh vector contaiting the values.

◆ performGhostExchange() [2/2]

template<typename T >
void ot::Mesh::performGhostExchange ( T *  vec)

Perform the ghost exchange for the vector vec.

Parameters
[in]vecadaptive mesh vector contaiting the values.

◆ ReMesh()

ot::Mesh * ot::Mesh::ReMesh ( unsigned int  grainSz = DENDRO_DEFAULT_GRAIN_SZ,
double  ld_tol = DENDRO_DEFAULT_LB_TOL,
unsigned int  sfK = DENDRO_DEFAULT_SF_K,
unsigned int(*)(const ot::TreeNode *)  getWeight = NULL 
)

: Remesh the mesh with the new computed elements.

Note
assumes that refinedIDs and corasenIDs are sorted. (This is automatically done by the isRemesh Fucntion)
Parameters
[in]refinedIDselement IDs need to be refined. (computed by isReMesh function)
[in]coarsenIDselement IDs need to be coarsened. (computed by isReMesh function)
[in]ld_toltolerance value used for flexible partitioning
[in]sfKspliiter fix parameter (need to specify larger value when run in super large scale)
[in]getWeightfunction pointer which returns a uint weight values for an given octant

◆ setOctreeRefineFlags()

void ot::Mesh::setOctreeRefineFlags ( unsigned int *  flags,
unsigned int  sz 
)

set refinement flags for the octree. This is non const function

Parameters
[in]flagsindicating to refine/coarsen or no change
[in]szsize of the array flags, sz should be equivalent to number of local elements.

◆ unzip()

template<typename T >
void ot::Mesh::unzip ( const T *  zippedVec,
T *  unzippedVec 
)

Creates the decomposition of adaptive octree variables into blocklist variables that we computed.

Author
Milinda Fernando
Parameters
[in]zippedVecadaptive representation of the variable array. (created by createVec function)
[out]unzippedVecdecomposed representation of the adaptive array.
Note
this routine assumes that for both arrays memory has been allocated. Routine is responsible only to fill up the unzipped entries.

◆ unzip_async()

template<typename T >
void ot::Mesh::unzip_async ( T *  zippedVec,
T *  unzippedVec,
MPI_Request *  send_reqs,
MPI_Request *  recv_reqs,
MPI_Status *  send_sts,
MPI_Status *  recv_sts 
)

perform the unzip operation overlap with the communication. GHOST exchage done inside. User doesn't need to perofrm ghost exchange.

Parameters
[in]zippedVecadaptive representation of the variable array. (created by createVec function)
[out]unzippedVecdecomposed representation of the adaptive array.
Note
this routine assumes that for both arrays memory has been allocated (No prior ghost exchange required). Routine is responsible only to fill up the unzipped entries.

◆ vectorToVTK()

template<typename T >
void ot::Mesh::vectorToVTK ( const std::vector< T > &  vec,
char *  fprefix,
double  pTime = 0.0,
unsigned int  nCycle = 0 
) const

write out function values to a vtk file.

Parameters
[in]vecvariable vector that needs to be written as a vtk file.
[in]fprefixprefix of the output vtk file name.

◆ zip()

template<typename T >
void ot::Mesh::zip ( const T *  unzippedVec,
T *  zippedVec 
)

Performs the compression frrom regular block grid varable list to adaptive representation.

Author
Milinda Fernando
Parameters
[in]unzippedVecdecomposed version of the adaptive array
[out]compressedversion of the unzippedVec.

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