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.
|
#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 RefElement * | getReferenceElement () 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::TreeNode & | currentOctant () |
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::TreeNode * | getNodalSplitterNodes () 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 > ¢ered, 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::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(*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... | |
Contains all the information needed to build neighbourhood information for the balanced octree based on treeSearch.
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
[in] | in | complete sorted 2:1 balanced octree to generate mesh |
[in] | k_s | how many neighbours to check on each direction (used =1) |
[in] | pOrder | order of an element. |
[in] | commActive | MPI active communicator |
[in] | comm | MPI communicator (global) |
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
[in] | in | complete sorted 2:1 balanced octree to generate mesh |
[in] | k_s | how many neighbours to check on each direction (used =1) |
[in] | pOrder | order of an element. |
[in] | comm | MPI communicator (global) |
[in] | grainSz | prefered grain sz. (this parameter is used to perform automatic comm expansion and shrinking) |
[in] | ld_tol | load imbalance tolerance for comm expansion and shrinking |
[in] | sf_k | splitter fix _k value. (Needed by SFC_partitioinng for large p>=10,000) |
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.
[in] | in | : vector that we need to apply the stencil on. |
[in] | centered | Stencil that we need to apply on vector in, this is the centered stencil. |
[in] | backward | backward version of the centered stencil.(This is used for boundary elements. ) |
[in] | forward | foward version of the centered stencil. (This is used for boundary elements. ) |
[out] | out | output vector that after applying the stencil. |
|
inline |
Performs child to parent injection.
[in] | children | function values. (all the function values of children ordered according to the SFC ordering) |
[in] | isHanging | array of boolean variables specifying each node is hanging ot not. |
[out] | injected | values back to parent. |
|
inline |
performs the child to parent contribution (only from a single child).
[in] | in | child function values |
[in] | cnum | morton ID of the current child. |
[in] | dim | dim of the interpolation. (dim=1 for edge interpolation , dim=2 for face interpolation, dim=3 for octant to child interpolation.) |
[out] | out | interpolated values. (child to parent contribution) |
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.
[in] | vec | child var vector (nPe) |
[in] | elementID | element ID of the current element (or child octant) |
[out] | out | add the contributions to the current vector accordingly. |
Usage: This is needed when performing matrix-free matvec for FEM method.
|
inline |
Returns the index of m_uiE2NMapping (CG or DG) for a specified Edge.
[in] | elementID | element ID of the edge, that the face belongs to . |
[in] | mortonIndex | morton ID of the corner node in the cordinate change in the order of x y z. |
[out] | index | returns the indecies of the requested face. |
void ot::Mesh::createUnZippedVector | ( | std::vector< T > & | uvec | ) | const |
: creates (memory allocation) for the unzipped version of the vector.
[in] |
T* ot::Mesh::createUnZippedVector | ( | ) | const |
: creates (memory allocation) for the unzipped version of the vector.
[in] |
void ot::Mesh::createUnZippedVector | ( | std::vector< T > & | uvec, |
const T | initValue | ||
) | const |
allocate memory for variable array based on the adaptive mesh.
[in] | uvec | allocate memory for uvec (unzipped version) |
[in] | initValue | initialize the vector to the given value. |
T* ot::Mesh::createUnZippedVector | ( | const T | initValue | ) | const |
allocate memory for variable array based on the adaptive mesh.
[in] | uvec | allocate memory for uvec (unzipped version) |
[in] | initValue | initialize the vector to the given value. |
T* ot::Mesh::createVector | ( | const T | initValue | ) | const |
allocate memory for variable array based on the adaptive mesh.
[in] | initValue | initialize the vector to the given value. |
T* ot::Mesh::createVector | ( | std::function< T(T, T, T)> | func | ) | const |
allocate memory for variable array based on the adaptive mesh.
[in] | initValue | initialize the vector to the given value. |
void ot::Mesh::createVector | ( | std::vector< T > & | vec, |
const T | initValue | ||
) | const |
allocate memory for variable array based on the adaptive mesh.
[in] | vec | allocate memory for vec |
[in] | initValue | initialize the vector to the given value. |
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.
[in] | vec | mesh varaible vector |
[in] | func | function to initialize with. |
|
inline |
Returns the index of m_uiE2NMapping (CG or DG) for a specified Edge.
[in] | elementID | element ID of the edge, that the face belongs to . |
[in] | face1 | : one of the face, that an edge belongs to |
[in] | face2 | second face that an edge belongs to . Hence the edge in cosideration will be intersection of the two faces, face1 and face2. |
[in] | isInternal | return only the internal nodes if true hence the index size would be ((m_uiElementOrder-1)*(m_uiElementOrder-1)) |
[out] | index | returns the indecies of the requested face. |
|
inline |
Returns the all the internal node indices of an element.
[in] | elementID | element ID of the edge, that the face belongs to . |
[in] | isInternal | return only the internal nodes if true hence the index size would be ((m_uiElementOrder-1)*(m_uiElementOrder-1)) |
[out] | index | returns the indecies of the requested face. |
|
inline |
Returns the index of m_uiE2NMapping (CG or DG) for a specified face.
[in] | elementID | element ID of the face, that the face belongs to . |
[in] | face | : Face that you need the indexing, |
[in] | isInternal | return only the internal nodes if true hence the index size would be ((m_uiElementOrder-1)*(m_uiElementOrder-1)) |
[out] | index | returns the indecies of the requested face. |
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.
[in] | blkId | : block id. |
[in] | eleId | : element id of the block |
[in] | dir | face 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] | child | element id for the block boudary faces, (elements containing inside the block) |
[out] | fid | reference pointer to the child array for the finer elements. |
[out] | cid | reference pointer to the child array for the coarser elements. (if it was refined. (cnumbers reference to the coarser elements)) |
|
inline |
returns const e2n mapping instance
|
inline |
returns const e2n mapping instance (debuging purposes only)
void ot::Mesh::getElementalEdgeNeighbors | ( | const unsigned int | eID, |
const unsigned int | dir, | ||
unsigned int * | lookup | ||
) | const |
returns the edge neighours in specidied direction.
[in] | eID | Element ID. |
[in] | dir | direction of the edge. |
[out] | lookUp | result if the result is a same level or lower level than element ID. (only 4 neighbor) |
void ot::Mesh::getElementalFaceNeighbors | ( | const unsigned int | eID, |
const unsigned int | dir, | ||
unsigned int * | lookup | ||
) | const |
returns the face neighours in specidied direction.
[in] | eID | Element ID. |
[in] | dir | direction of the face. |
[out] | lookUp | result if the result is a same level or lower level than element ID. (only 2 neighbor) |
void ot::Mesh::getElementalVertexNeighbors | ( | const unsigned int | eID, |
const unsigned int | dir, | ||
unsigned int * | lookup | ||
) | const |
returns the vertex neighours in specidied direction.
[in] | eID | Element ID. |
[in] | dir | direction of the vertex. |
[out] | lookUp | result if the result is a same level or lower level than element ID. (only 8 neighbor) |
void ot::Mesh::getElementCoordinates | ( | unsigned int | eleID, |
double * | coords | ||
) | const |
computes the elementCoordinates (based on the nodal placement)
[in] | eleID | : element ID |
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.
[in] | vec | variable vector that we want to get the nodal values. |
[in] | elementID | element ID that we need to get the nodal values. |
[out] | nodalValues | nodal values of the specified element ID |
void ot::Mesh::getElementQMat | ( | unsigned int | currentId, |
double *& | qMat, | ||
bool | isAllocated = true |
||
) | const |
: Compute the elemental interpolation matrix.
[in] | currentId | Element ID. |
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.
[in] | eleID | element ID |
[in] | in | inpute vector |
[out] | out | output vector values are in the order of the x,y,z size : 4*NodesPerElement |
[out] | coords | get the corresponding coordinates size: 4*NodesPerElement*m_uiDim; |
[out] | neighID | face neighbor octant IDs, |
[in] | face | face direction in {OCT_DIR_LEFT,OCT_IDR_RIGHT,OCT_DIR_DOWN, OCT_DIR_UP,OCT_DIR_BACK,OCT_DIR_FRONT} |
[out] | level | the level of the neighbour octant with respect to the current octant. returns the number of face neighbours 1/4 for 3D. |
|
inline |
Returns the number of local elements in the grid. (local to the current considering processor)
|
inline |
Returns the number of pre-ghost elements
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.
T | type of the vector. |
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. |
void ot::Mesh::ghostExchangeRecvSync | ( | T * | vec, |
T * | recvNodeBuffer, | ||
MPI_Request * | recv_reqs, | ||
MPI_Status * | recv_sts | ||
) |
Perform the wait on the recv requests.
[in] | vec | adaptive mesh vector containing the values. |
[in] | recv_reqs | m_uiRecvProcList.size() recv request |
[in] | recv_sts | m_uiRecvProcList.size() recv status |
|
inline |
Perform the wait on the recv requests.
[in] | vec | adaptive mesh vector containing the values. |
[in] | send_reqs | m_uiSendProcList.size() send request |
[in] | send_sts | m_uiSendProcList.size() send status |
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.
[in] | vec | adaptive mesh vector contaiting the values. |
void ot::Mesh::interGridTransfer | ( | std::vector< T > & | vec, |
const ot::Mesh * | pMesh | ||
) |
transfer a variable vector form old grid to new grid.
[in] | vec | variable vector needs to be transfered. |
[out] | vec | transfered varaible vector |
[in] | pMesh | Mesh that we need to transfer the old varaible. |
void ot::Mesh::interGridTransfer | ( | T *& | vec, |
const ot::Mesh * | pMesh | ||
) |
transfer a variable vector form old grid to new grid.
[in] | vec | variable vector needs to be transfered. |
[out] | vec | transfered varaible vector |
[in] | pMesh | Mesh that we need to transfer the old varaible. |
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.
T | : data type of the vectors |
unzip | : unizp representation of the vectors. |
vec | : zip representation of the transfered vectors. |
pMesh | : new mesh. |
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.
[in] | elementId | element ID of the octant. |
[in] | edgeId | edge id |
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.
[in] | elementId | element ID of the octant. |
[in] | faceId | face id |
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.
[in] | eleID | element ID |
[in] | ix | i-index of the node. |
[in] | jy | j-index of the node. |
[in] | kz | k-index of the node. |
|
inline |
: Returns true if the specified node (e,i,j,k) is local.
[in] | eleID | element ID |
[in] | ix | i-index of the node. |
[in] | jy | j-index of the node. |
[in] | kz | k-index of the node. |
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).
[in] | vec | sequence of varaibles to check |
[in] | varIds | variable ids to check. (var ids to index the vec, vec[i] is a T* pointintg to one of the variable in vec. ) |
[in] | numVars | number of variables to check |
[in] | tol | wavelet tolerance Returns true if specified variable violates the specified wavelet toerlance. |
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 |
||
) |
[in] | vec | sequence of varaibles to check |
[in] | varIds | variable ids to check. (var ids to index the vec, vec[i] is a T* pointintg to one of the variable in vec. ) |
[in] | numVars | number of variables to check |
[in] | tol | wavelet tolerance Returns true if specified variable violates the specified wavelet toerlance. |
|
inline |
Performs all parent to child interpolations for the m_uiEl_i element in order to apply the stencil.
[in] | parent | function values |
[in] | cnum | child number to interpolate. |
[in] | dim | dim of the interpolation. (dim=1 for edge interpolation , dim=2 for face interpolation, dim=3 for octant to child interpolation.) |
[out] | out | interpolated values. |
void ot::Mesh::performGhostExchange | ( | std::vector< T > & | vec | ) |
Perform the ghost exchange for the vector vec.
[in] | vec | adaptive mesh vector contaiting the values. |
void ot::Mesh::performGhostExchange | ( | T * | vec | ) |
Perform the ghost exchange for the vector vec.
[in] | vec | adaptive mesh vector contaiting the values. |
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.
[in] | refinedIDs | element IDs need to be refined. (computed by isReMesh function) |
[in] | coarsenIDs | element IDs need to be coarsened. (computed by isReMesh function) |
[in] | ld_tol | tolerance value used for flexible partitioning |
[in] | sfK | spliiter fix parameter (need to specify larger value when run in super large scale) |
[in] | getWeight | function pointer which returns a uint weight values for an given octant |
void ot::Mesh::setOctreeRefineFlags | ( | unsigned int * | flags, |
unsigned int | sz | ||
) |
set refinement flags for the octree. This is non const function
[in] | flags | indicating to refine/coarsen or no change |
[in] | sz | size of the array flags, sz should be equivalent to number of local elements. |
void ot::Mesh::unzip | ( | const T * | zippedVec, |
T * | unzippedVec | ||
) |
Creates the decomposition of adaptive octree variables into blocklist variables that we computed.
[in] | zippedVec | adaptive representation of the variable array. (created by createVec function) |
[out] | unzippedVec | decomposed representation of the adaptive array. |
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.
[in] | zippedVec | adaptive representation of the variable array. (created by createVec function) |
[out] | unzippedVec | decomposed representation of the adaptive array. |
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.
[in] | vec | variable vector that needs to be written as a vtk file. |
[in] | fprefix | prefix of the output vtk file name. |
void ot::Mesh::zip | ( | const T * | unzippedVec, |
T * | zippedVec | ||
) |
Performs the compression frrom regular block grid varable list to adaptive representation.
[in] | unzippedVec | decomposed version of the adaptive array |
[out] | compressed | version of the unzippedVec. |