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 | |
std::vector< ot::TreeNode > | getLocalOctants () |
: returns the vector containing only local elements | |
const std::vector< ot::TreeNode > & | getBlocks () |
: get octree to block decomposition blocks | |
DA (std::vector< ot::TreeNode > &balOct, MPI_Comm comm, unsigned int order, unsigned int grainSz=100, double sfc_tol=0.3, SM_TYPE smType=SM_TYPE::FEM_CG) | |
: Constructor for the DA data structures More... | |
DA (ot::Mesh *pMesh) | |
: Create a DA from a specified mesh More... | |
template<typename T > | |
DA (std::function< void(T, T, T, T *)>func, unsigned int dofSz, MPI_Comm comm, unsigned int order, double interp_tol, unsigned int grainSz=100, double sfc_tol=0.3, SM_TYPE smType=SM_TYPE::FEM_CG) | |
DA (const Point *pts, unsigned int numPts, Point pt_min, Point pt_max, MPI_Comm comm, unsigned int order, unsigned int grainSz=100, double sfc_tol=0.3, SM_TYPE smType=SM_TYPE::FEM_CG) | |
Construct a new DA object from specified point set. More... | |
~DA () | |
deconstructor for the DA class. | |
unsigned int | getLocalNodalSz () const |
returns the local nodal size | |
unsigned int | getPreNodalSz () const |
returns the pre ghost nodal size | |
unsigned int | getPostNodalSz () const |
returns the post nodal size | |
unsigned int | getTotalNodalSz () const |
returns the total nodal size (this includes the ghosted region as well.) | |
unsigned int | getLocalElemSz () const |
returns the local elemental size | |
unsigned int | getTotalElemSz () const |
returns the local elemental size (includes the ghost elements as well) | |
bool | isActive () const |
see if the current DA is active | |
unsigned int | getNumNodesPerElement () const |
get number of nodes per element | |
unsigned int | getElementOrder () const |
get element order | |
MPI_Comm | getGlobalComm () const |
: returns the global MPI communicator | |
const std::vector< DendroIntL > | getNodeLocalToGlobalMap () const |
: returns node local to node global map | |
const ot::Mesh * | getMesh () const |
returns the mesh | |
const std::vector< unsigned int > | getOctFlags () const |
: Returns octDA flags | |
MPI_Comm | getCommActive () const |
: returns the active MPI sub com of the global communicator | |
unsigned int | getNpesAll () const |
: global mpi com. size | |
unsigned int | getNpesActive () const |
: number of processors active | |
unsigned int | getRankAll () const |
: rank with respect to the global comm. | |
unsigned int | getRankActive () const |
: rank w.r.t active comm. | |
bool | isBoundaryOctant (unsigned int eleID) const |
: returns true if specified eleID is a boundary element false if the eleID is local and not a boundary element. More... | |
void | getElementalCoords (unsigned int eleID, double *coords) const |
computes the elementCoordinates (based on the nodal placement) More... | |
const RefElement * | getReferenceElement () const |
get a constant pointer for the reference element | |
unsigned int | getElementSize () const |
: get the number of local elements. | |
unsigned int | getPreGhostElementSize () const |
: get the number of pre-ghost element size | |
unsigned int | getPostGhostElementSize () const |
: get the number of post-ghost element size | |
unsigned int | getPreAndPostGhostNodeSize () const |
: get the number of pre and post ghost elements | |
unsigned int | getNumLayer1GhostEleSz () const |
number of layer 1 ghost elements. | |
unsigned int | getGhostedElementSize () const |
: get the number of pre and post ghost elements | |
unsigned int | getMaxDepth () const |
: get the max depth of the octree | |
unsigned int | getDimension () const |
: get the dimensionality of the octree | |
ot::TreeNode | getMinTreeNode () const |
get min local node (nodal) | |
ot::TreeNode | getMaxTreeNode () const |
get max localnode (nodal) | |
void | computeTreeNodeOwnerProc (const ot::TreeNode *pNodes, unsigned int n, int *ownerranks) |
template<typename T > | |
int | createVector (T *&local, bool isElemental=false, bool isGhosted=false, unsigned int dof=1) const |
Creates a ODA vector. More... | |
template<typename T > | |
int | createVector (std::vector< T > &local, bool isElemental=false, bool isGhosted=false, unsigned int dof=1) const |
Creates a ODA vector std::vector<T> More... | |
template<typename T > | |
void | destroyVector (T *&local) const |
deallocates the memory allocated for a vector More... | |
template<typename T > | |
void | destroyVector (std::vector< T > &local) const |
template<ot::DA_FLAGS::LoopType type> | |
void | init () |
initialize the loop counters. | |
unsigned int | curr () |
get the current elmenent in the iteration | |
template<ot::DA_FLAGS::LoopType type> | |
unsigned int | end () |
Returns if the current element has reached end or not. | |
template<ot::DA_FLAGS::LoopType type> | |
void | next () |
: increment it to the next element. | |
template<typename T > | |
void | readFromGhostBegin (T *vec, unsigned int dof=1) |
Initiate the ghost nodal value exchange. | |
template<typename T > | |
void | readFromGhostEnd (T *vec, unsigned int dof=1) |
Sync the ghost element exchange. | |
template<typename T > | |
void | writeToGhostsBegin (T *vec, unsigned int dof=1) |
Initiate accumilation across ghost elements. | |
template<typename T > | |
void | writeToGhostsEnd (T *vec, DA_FLAGS::WriteMode mode, unsigned int dof=1) |
Sync accumilation across ghost elements. More... | |
template<typename T > | |
void | nodalVecToGhostedNodal (const T *in, T *&out, bool isAllocated=false, unsigned int dof=1) const |
convert nodal local vector with ghosted buffer regions. More... | |
template<typename T > | |
void | ghostedNodalToNodalVec (const T *gVec, T *&local, bool isAllocated=false, unsigned int dof=1) const |
convert ghosted nodal vector to local vector (without ghosting) More... | |
template<typename T > | |
void | getElementNodalValues (const T *in, T *eleVecOut, unsigned int eleID, unsigned int dof=1) const |
computes the element nodal values using interpolation if needed. More... | |
template<typename T > | |
void | eleVecToVecAccumilation (T *out, const T *eleVecIn, unsigned int eleID, unsigned int dof=1) const |
computes the elemental vec to global vec accumilation More... | |
void | getOctreeBoundaryNodeIndices (std::vector< unsigned int > &bdyIndex, std::vector< double > &coords, bool isGhosted=false) |
Computes the octree writable boundary nodes. More... | |
int | getNodeIndices (DendroIntL *nodeIdx, unsigned int ele, bool isGhosted) const |
Returns the local indices to the nodes of the current element. More... | |
int | getGlobalNodeIndices (DendroIntL *nodeIdx, unsigned int ele) const |
Returns the global node indices of a given element,. More... | |
template<typename T > | |
void | setVectorByFunction (T *local, std::function< void(T, T, T, T *)>func, bool isElemental=false, bool isGhosted=false, unsigned int dof=1) const |
initialize a variable vector to a function depends on spatial coords. More... | |
template<typename T > | |
void | setVectorByScalar (T *local, const T *value, bool isElemental=false, bool isGhosted=false, unsigned int dof=1) const |
initialize a variable vector to a function depends on spatial coords. More... | |
template<typename T > | |
void | vecTopvtu (T *local, const char *fPrefix, char **nodalVarNames=NULL, bool isElemental=false, bool isGhosted=false, unsigned int dof=1) |
write the vec to pvtu file More... | |
unsigned int | getLevel (unsigned int ele) const |
return the level of current octant More... | |
ot::TreeNode | getOctant (unsigned int ele) const |
return the TreeNode of the current octnat. More... | |
template<typename T > | |
T * | getVecPointerToDof (T *in, unsigned int dofInex, bool isElemental=false, bool isGhosted=false) const |
returns a pointer to a dof index, More... | |
template<typename T > | |
void | copyVectors (T *dest, const T *source, bool isElemental=false, bool isGhosted=false, unsigned int dof=1) const |
copy vecotor to sorce to destination, assumes the same number of dof. More... | |
template<typename T > | |
void | copyVector (T *dest, const T *source, bool isElemental=false, bool isGhosted=false) const |
more premitive copy, from source pointer to the dest pointer More... | |
ot::DA * | remesh (const DA_FLAGS::Refine *flags, unsigned int sz, unsigned int grainSz=100, double ld_bal=0.3, unsigned int sfK=2, unsigned int(*getWeight)(const ot::TreeNode *)=NULL) const |
: Performs remesh based on the DA_FLAGS::Refine, which specifies no change, refine or coarsen. More... | |
template<typename T > | |
void | intergridTransfer (const T *varIn, T *&varOut, const ot::DA *newDA, bool isElemental=false, bool isGhosted=false, unsigned int dof=1) |
performs grid transfer operations after the remesh. More... | |
template<typename T > | |
int | getFaceNeighborValues (unsigned int eleID, const T *in, T *out, T *coords, unsigned int *neighID, unsigned int face, NeighbourLevel &level, unsigned int dof) const |
computes the face neighbor points for additional computations for a specified direction. More... | |
unsigned int | getMortonChildNum (unsigned int eleID) const |
computes the child number of the given octant More... | |
ot::DA::DA | ( | std::vector< ot::TreeNode > & | balOct, |
MPI_Comm | comm, | ||
unsigned int | order, | ||
unsigned int | grainSz = 100 , |
||
double | sfc_tol = 0.3 , |
||
SM_TYPE | smType = SM_TYPE::FEM_CG |
||
) |
ot::DA::DA | ( | std::function< void(T, T, T, T *)> | func, |
unsigned int | dofSz, | ||
MPI_Comm | comm, | ||
unsigned int | order, | ||
double | interp_tol, | ||
unsigned int | grainSz = 100 , |
||
double | sfc_tol = 0.3 , |
||
SM_TYPE | smType = SM_TYPE::FEM_CG |
||
) |
Construct a DA from a function
[in] | func | function to capture |
[in] | dofSz | size of the degrees of freedoms |
[in] | comm | MPI communicator |
[in] | order | element order |
[in] | interp_tol | interpolation tolerance for func approximation/ capturing |
[in] | grainSz | number of elements per core. |
[in] | sfc_tol | flexible partitioning tolerance. |
[in] | SM_TYPE | scatter map type |
ot::DA::DA | ( | const Point * | pts, |
unsigned int | numPts, | ||
Point | pt_min, | ||
Point | pt_max, | ||
MPI_Comm | comm, | ||
unsigned int | order, | ||
unsigned int | grainSz = 100 , |
||
double | sfc_tol = 0.3 , |
||
SM_TYPE | smType = SM_TYPE::FEM_CG |
||
) |
Construct a new DA object from specified point set.
[in] | pts | : input points |
numPts | : number of points | |
pt_min | : global min of point | |
pt_max | : global max of point | |
comm | : MPI communicator | |
order | : element order | |
grainSz | : number of elements per core | |
sfc_tol | : SFC partition tolerance | |
smType | : scatter map type |
void ot::DA::computeTreeNodeOwnerProc | ( | const ot::TreeNode * | pNodes, |
unsigned int | n, | ||
int * | ownerranks | ||
) |
[in] | pNodes | Nodes that needs to be searched across. |
[in] | n | number of nodes. |
void ot::DA::copyVector | ( | T * | dest, |
const T * | source, | ||
bool | isElemental = false , |
||
bool | isGhosted = false |
||
) | const |
more premitive copy, from source pointer to the dest pointer
[in/out] | dest: destination pointer | |
[in] | source | source pointer |
[in] | isElemental | true if this is an elemental vector/ false otherwise |
[in] | isGhosted | true if this is a ghosted vector |
void ot::DA::copyVectors | ( | T * | dest, |
const T * | source, | ||
bool | isElemental = false , |
||
bool | isGhosted = false , |
||
unsigned int | dof = 1 |
||
) | const |
copy vecotor to sorce to destination, assumes the same number of dof.
[in/out] | dest: destination pointer | |
[in] | source | source pointer |
[in] | isElemental | true if this is an elemental vector/ false otherwise |
[in] | isGhosted | true if this is a ghosted vector |
[in] | dof | degrees of freedoms |
int ot::DA::createVector | ( | T *& | local, |
bool | isElemental = false , |
||
bool | isGhosted = false , |
||
unsigned int | dof = 1 |
||
) | const |
Creates a ODA vector.
[in] | local | : VecType pointer |
[in] | isElemental | True if creating a elemental vector (cell data vector) false for a nodal vector |
[in] | isGhosted | True will allocate ghost nodal values as well, false will only allocate memory for local nodes. |
[in] | dof | degrees of freedoms |
int ot::DA::createVector | ( | std::vector< T > & | local, |
bool | isElemental = false , |
||
bool | isGhosted = false , |
||
unsigned int | dof = 1 |
||
) | const |
Creates a ODA vector std::vector<T>
[in] | local | : VecType pointer |
[in] | isElemental | True if creating a elemental vector (cell data vector) false for a nodal vector |
[in] | isGhosted | True will allocate ghost nodal values as well, false will only allocate memory for local nodes. |
[in] | dof | degrees of freedoms |
void ot::DA::destroyVector | ( | T *& | local | ) | const |
deallocates the memory allocated for a vector
void ot::DA::eleVecToVecAccumilation | ( | T * | out, |
const T * | eleVecIn, | ||
unsigned int | eleID, | ||
unsigned int | dof = 1 |
||
) | const |
computes the elemental vec to global vec accumilation
[out] | out | output (need to be ghosted ) |
[in] | eleVecIn | eleVec values. |
[in] | eleID | element ID |
[in] | dof | degrees of freedoms |
void ot::DA::getElementalCoords | ( | unsigned int | eleID, |
double * | coords | ||
) | const |
computes the elementCoordinates (based on the nodal placement)
[in] | eleID | : element ID |
void ot::DA::getElementNodalValues | ( | const T * | in, |
T * | eleVecOut, | ||
unsigned int | eleID, | ||
unsigned int | dof = 1 |
||
) | const |
computes the element nodal values using interpolation if needed.
[in] | in | input vector (need to be ghosted ) |
int ot::DA::getFaceNeighborValues | ( | unsigned int | eleID, |
const T * | in, | ||
T * | out, | ||
T * | coords, | ||
unsigned int * | neighID, | ||
unsigned int | face, | ||
NeighbourLevel & | level, | ||
unsigned int | dof | ||
) | 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. |
int ot::DA::getGlobalNodeIndices | ( | DendroIntL * | nodeIdx, |
unsigned int | ele | ||
) | const |
Returns the global node indices of a given element,.
[out] | nodeIdx | global node indices, |
[in] | ele | element id |
|
inline |
return the level of current octant
[in] | ele | elementID of the current octant |
|
inline |
computes the child number of the given octant
[in] | eleID | element ID returns the Morton child number of the eleID in the octant |
int ot::DA::getNodeIndices | ( | DendroIntL * | nodeIdx, |
unsigned int | ele, | ||
bool | isGhosted | ||
) | const |
Returns the local indices to the nodes of the current element.
nodes | Indices into the nodes of the given element. Should be allocated by the user prior to calling. |
|
inline |
return the TreeNode of the current octnat.
[in] | ele | elementID of the current octant |
void ot::DA::getOctreeBoundaryNodeIndices | ( | std::vector< unsigned int > & | bdyIndex, |
std::vector< double > & | coords, | ||
bool | isGhosted = false |
||
) |
Computes the octree writable boundary nodes.
T* ot::DA::getVecPointerToDof | ( | T * | in, |
unsigned int | dofInex, | ||
bool | isElemental = false , |
||
bool | isGhosted = false |
||
) | const |
returns a pointer to a dof index,
[in] | in | input vector pointer |
[in] | dofInex | dof index which is the pointer is needed, should be less than dof, value the vector created. |
[in] | isElemental | true if this is an elemental vector/ false otherwise |
[in] | isGhosted | true if this is a ghosted vector |
void ot::DA::ghostedNodalToNodalVec | ( | const T * | gVec, |
T *& | local, | ||
bool | isAllocated = false , |
||
unsigned int | dof = 1 |
||
) | const |
convert ghosted nodal vector to local vector (without ghosting)
[in] | gVec | ghosted vector |
[out] | local | local vector (assume an allocated vector) |
[in] | isAllocated | true if the out is allocated, false otherwise. |
[in] | dof | degrees of freedoms |
void ot::DA::intergridTransfer | ( | const T * | varIn, |
T *& | varOut, | ||
const ot::DA * | newDA, | ||
bool | isElemental = false , |
||
bool | isGhosted = false , |
||
unsigned int | dof = 1 |
||
) |
performs grid transfer operations after the remesh.
[in] | varIn | variable defined by oldDA |
[out] | varOut | variable defined by newDA. interpolate varOut from varIn. (Note: varOut allocated inside the function, no need to allocate outside) |
[in] | isElemental | true if it is an elemental vector |
[in] | isGhosted | true if allocated ghost vector |
[in] | dof | degrees of freedoms. |
bool ot::DA::isBoundaryOctant | ( | unsigned int | eleID | ) | const |
: returns true if specified eleID is a boundary element false if the eleID is local and not a boundary element.
[in] | eleID | element ID |
void ot::DA::nodalVecToGhostedNodal | ( | const T * | in, |
T *& | out, | ||
bool | isAllocated = false , |
||
unsigned int | dof = 1 |
||
) | const |
convert nodal local vector with ghosted buffer regions.
[in] | in | input vector (should be nodal and non ghosted) |
[out] | out | coverted nodal vector with ghost regions. |
[in] | isAllocated | true if the out is allocated, false otherwise. |
[in] | dof | degrees of freedoms |
ot::DA * ot::DA::remesh | ( | const DA_FLAGS::Refine * | flags, |
unsigned int | sz, | ||
unsigned int | grainSz = 100 , |
||
double | ld_bal = 0.3 , |
||
unsigned int | sfK = 2 , |
||
unsigned int(*)(const ot::TreeNode *) | getWeight = NULL |
||
) | const |
: Performs remesh based on the DA_FLAGS::Refine, which specifies no change, refine or coarsen.
[in] | flags | refinement flags. |
[in] | sz | size of the array flags (needs to be size of the local elements) |
[in] | grainSz | rougly the number of octants per core you need when you create the new da. |
[in] | ld_tol | load imbalance tolerance. |
[in] | sfK | splitter fix factor. better to be power of two. increase the value to 128 when running on > 64,000 cores |
void ot::DA::setVectorByFunction | ( | T * | local, |
std::function< void(T, T, T, T *)> | func, | ||
bool | isElemental = false , |
||
bool | isGhosted = false , |
||
unsigned int | dof = 1 |
||
) | const |
initialize a variable vector to a function depends on spatial coords.
void ot::DA::setVectorByScalar | ( | T * | local, |
const T * | value, | ||
bool | isElemental = false , |
||
bool | isGhosted = false , |
||
unsigned int | dof = 1 |
||
) | const |
initialize a variable vector to a function depends on spatial coords.
void ot::DA::vecTopvtu | ( | T * | local, |
const char * | fPrefix, | ||
char ** | nodalVarNames = NULL , |
||
bool | isElemental = false , |
||
bool | isGhosted = false , |
||
unsigned int | dof = 1 |
||
) |
write the vec to pvtu file
[in] | local | variable vector |
[in] | fPrefix | file name prefix |
[in] | isElemental | True if creating a elemental vector (cell data vector) false for a nodal vector |
[in] | isGhosted | True will allocate ghost nodal values as well, false will only allocate memory for local nodes. |
[in] | dof | degrees of freedoms |
void ot::DA::writeToGhostsEnd | ( | T * | vec, |
DA_FLAGS::WriteMode | mode, | ||
unsigned int | dof = 1 |
||
) |
Sync accumilation across ghost elements.
[in] | vec | vector pointer |
[in] | mode | mode of the write to ghost |
[in] | dof | degrees of freedoms. |