Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Static Public Attributes | List of all members
bitpit::ParaTree Class Reference

Para Tree is the user interface class. More...

#include <ParaTree.hpp>

Inheritance diagram for bitpit::ParaTree:
Inheritance graph
[legend]

Classes

struct  LoadBalanceRanges
 

Public Types

typedef std::unordered_map< int, std::array< uint32_t, 2 > > ExchangeRanges
 
enum  Operation {
  OP_NONE , OP_INIT , OP_PRE_ADAPT , OP_ADAPT_MAPPED ,
  OP_ADAPT_UNMAPPED , OP_LOADBALANCE_FIRST , OP_LOADBALANCE
}
 

Public Member Functions

 ParaTree (const ParaTree &other)
 
 ParaTree (const std::string &logfile=DEFAULT_LOG_FILE, MPI_Comm comm=MPI_COMM_WORLD)
 
 ParaTree (std::istream &stream, const std::string &logfile=DEFAULT_LOG_FILE, MPI_Comm comm=MPI_COMM_WORLD)
 
 ParaTree (uint8_t dim, const std::string &logfile=DEFAULT_LOG_FILE, MPI_Comm comm=MPI_COMM_WORLD)
 
virtual ~ParaTree ()
 
bool adapt (bool mapper_flag=false)
 
bool adaptGlobalCoarse (bool mapper_flag=false)
 
bool adaptGlobalRefine (bool mapper_flag=false)
 
bool check21Balance ()
 
bool checkToAdapt ()
 
void clearConnectivity (bool release=true)
 
template<class Impl >
void communicate (DataCommInterface< Impl > &userData)
 
void computeConnectivity ()
 
void computeIntersections ()
 
uint64_t computeNodePersistentKey (const Octant *oct, uint8_t node) const
 
uint64_t computeNodePersistentKey (uint32_t idx, uint8_t node) const
 
virtual void dump (std::ostream &stream, bool full=true)
 
LoadBalanceRanges evalLoadBalanceRanges (dvector *weights)
 
LoadBalanceRanges evalLoadBalanceRanges (uint8_t level, dvector *weights)
 
void expectedOctantAdapt (const Octant *octant, int8_t marker, octvector *result) const
 
void findAllCodimensionNeighbours (const Octant *oct, u32vector &neighbours, bvector &isghost)
 
void findAllCodimensionNeighbours (uint32_t idx, u32vector &neighbours, bvector &isghost)
 
void findAllNodeNeighbours (const Octant *oct, uint32_t node, u32vector &neighbours, bvector &isghost) const
 
void findAllNodeNeighbours (uint32_t idx, uint32_t node, u32vector &neighbours, bvector &isghost)
 
void findGhostAllCodimensionNeighbours (Octant *oct, u32vector &neighbours, bvector &isghost)
 
void findGhostAllCodimensionNeighbours (uint32_t idx, u32vector &neighbours, bvector &isghost)
 
void findGhostNeighbours (const Octant *oct, uint8_t face, uint8_t codim, u32vector &neighbours, bvector &isghost) const
 
void findGhostNeighbours (uint32_t idx, uint8_t face, uint8_t codim, u32vector &neighbours) const
 
void findGhostNeighbours (uint32_t idx, uint8_t face, uint8_t codim, u32vector &neighbours, bvector &isghost) const
 
void findNeighbours (const Octant *oct, uint8_t face, uint8_t codim, u32vector &neighbours, bvector &isghost) const
 
void findNeighbours (uint32_t idx, uint8_t face, uint8_t codim, u32vector &neighbours) const
 
void findNeighbours (uint32_t idx, uint8_t face, uint8_t codim, u32vector &neighbours, bvector &isghost) const
 
int findOwner (uint64_t morton) const
 
void freeComm ()
 
double getArea (const Intersection *inter) const
 
double getArea (const Octant *oct) const
 
double getArea (uint32_t idx) const
 
bool getBalance (const Octant *oct) const
 
bool getBalance (uint32_t idx) const
 
uint8_t getBalanceCodimension () const
 
const std::map< int, std::vector< uint32_t > > & getBordersPerProc () const
 
bool getBound (const Intersection *inter) const
 
bool getBound (const Octant *oct) const
 
bool getBound (const Octant *oct, uint8_t face) const
 
bool getBound (uint32_t idx) const
 
bool getBound (uint32_t idx, uint8_t face) const
 
darray3 getCenter (const Intersection *inter) const
 
darray3 getCenter (const Octant *oct) const
 
void getCenter (const Octant *oct, darray3 &centerCoords) const
 
darray3 getCenter (uint32_t idx) const
 
void getCenter (uint32_t idx, darray3 &centerCoords) const
 
MPI_Comm getComm () const
 
const u32vector2D & getConnectivity () const
 
const u32vector & getConnectivity (Octant *oct) const
 
const u32vector & getConnectivity (uint32_t idx) const
 
darray3 getCoordinates (const Octant *oct) const
 
darray3 getCoordinates (uint32_t idx) const
 
uint8_t getDim () const
 
virtual int getDumpVersion () const
 
const int8_t(* getEdgecoeffs () const)[3]
 
void getEdgecoeffs (int8_t edgecoeffs[12][3]) const
 
const uint8_t(* getEdgeface () const)[2]
 
void getEdgeface (uint8_t edgeface[12][2]) const
 
uint8_t getFace (const Intersection *inter) const
 
darray3 getFaceCenter (const Octant *oct, uint8_t face) const
 
void getFaceCenter (const Octant *oct, uint8_t face, darray3 &centerCoords) const
 
darray3 getFaceCenter (uint32_t idx, uint8_t face) const
 
void getFaceCenter (uint32_t idx, uint8_t face, darray3 &centerCoords) const
 
const uint8_t(* getFacenode () const)[4]
 
void getFacenode (uint8_t facenode[6][4]) const
 
uint8_t getFamilySplittingNode (const Octant *) const
 
bool getFiner (const Intersection *inter) const
 
uint64_t getFirstDescMorton () const
 
const u32vector2D & getGhostConnectivity () const
 
const u32vector & getGhostConnectivity (const Octant *oct) const
 
const u32vector & getGhostConnectivity (uint32_t idx) const
 
uint64_t getGhostGlobalIdx (uint32_t idx) const
 
int getGhostLayer (const Octant *oct) const
 
uint32_t getGhostLocalIdx (uint64_t gidx) const
 
OctantgetGhostOctant (uint32_t idx)
 
const OctantgetGhostOctant (uint32_t idx) const
 
uint64_t getGlobalIdx (const Octant *oct) const
 
uint64_t getGlobalIdx (uint32_t idx) const
 
uint64_t getGlobalNumOctants () const
 
uint32_t getIdx (const Octant *oct) const
 
uint32_t getIn (const Intersection *inter) const
 
octantIterator getInternalOctantsBegin ()
 
octantIterator getInternalOctantsEnd ()
 
IntersectiongetIntersection (uint32_t idx)
 
bool getIsGhost (const Intersection *inter) const
 
bool getIsGhost (const Octant *oct) const
 
bool getIsNewC (const Octant *oct) const
 
bool getIsNewC (uint32_t idx) const
 
bool getIsNewR (const Octant *oct) const
 
bool getIsNewR (uint32_t idx) const
 
double getL () const
 
uint64_t getLastDescMorton () const
 
uint64_t getLastDescMorton (const Octant *oct) const
 
uint64_t getLastDescMorton (uint32_t idx) const
 
Operation getLastOperation () const
 
uint8_t getLevel (const Intersection *inter) const
 
uint8_t getLevel (const Octant *oct) const
 
uint8_t getLevel (uint32_t idx) const
 
const LoadBalanceRangesgetLoadBalanceRanges () const
 
uint32_t getLocalIdx (uint64_t gidx) const
 
uint32_t getLocalIdx (uint64_t gidx, int rank) const
 
void getLocalIdx (uint64_t gidx, uint32_t &lidx, int &rank) const
 
uint8_t getLocalMaxDepth () const
 
double getLocalMaxSize () const
 
double getLocalMinSize () const
 
LoggergetLog ()
 
void getMapping (uint32_t idx, u32vector &mapper, bvector &isghost) const
 
void getMapping (uint32_t idx, u32vector &mapper, bvector &isghost, ivector &rank) const
 
int8_t getMarker (const Octant *oct) const
 
int8_t getMarker (uint32_t idx) const
 
int8_t getMaxDepth () const
 
uint32_t getMaxLength () const
 
int getMaxLevel () const
 
uint64_t getMorton (const Octant *oct) const
 
uint64_t getMorton (uint32_t idx) const
 
uint8_t getNchildren () const
 
uint8_t getNedges () const
 
uint8_t getNfaces () const
 
uint8_t getNnodes () const
 
uint8_t getNnodesperface () const
 
darray3 getNode (const Octant *oct, uint8_t node) const
 
void getNode (const Octant *oct, uint8_t node, darray3 &nodeCoords) const
 
darray3 getNode (uint32_t idx, uint8_t node) const
 
void getNode (uint32_t idx, uint8_t node, darray3 &nodeCoords) const
 
const int8_t(* getNodecoeffs () const)[3]
 
void getNodecoeffs (int8_t nodecoeffs[8][3]) const
 
darray3 getNodeCoordinates (uint32_t node) const
 
const uint8_t(* getNodeface () const)[3]
 
void getNodeface (uint8_t nodeface[8][3]) const
 
const u32array3 & getNodeLogicalCoordinates (uint32_t node) const
 
const u32arr3vector & getNodes () const
 
darr3vector getNodes (const Intersection *inter) const
 
darr3vector getNodes (const Octant *oct) const
 
void getNodes (const Octant *oct, darr3vector &nodesCoords) const
 
darr3vector getNodes (uint32_t idx) const
 
void getNodes (uint32_t idx, darr3vector &nodesCoords) const
 
std::size_t getNofGhostLayers () const
 
darray3 getNormal (const Intersection *inter) const
 
darray3 getNormal (const Octant *oct, uint8_t face) const
 
void getNormal (const Octant *oct, uint8_t face, darray3 &normal) const
 
darray3 getNormal (uint32_t idx, uint8_t face) const
 
void getNormal (uint32_t idx, uint8_t face, darray3 &normal) const
 
const int8_t(* getNormals () const)[3]
 
void getNormals (int8_t normals[6][3]) const
 
int getNproc () const
 
uint32_t getNumGhosts () const
 
uint32_t getNumIntersections () const
 
uint32_t getNumNodes () const
 
uint32_t getNumOctants () const
 
OctantgetOctant (uint32_t idx)
 
const OctantgetOctant (uint32_t idx) const
 
const uint8_t * getOppface () const
 
void getOppface (uint8_t oppface[6]) const
 
darray3 getOrigin () const
 
uint32_t getOut (const Intersection *inter) const
 
bool getOutIsGhost (const Intersection *inter) const
 
int getOwnerRank (uint64_t globalIdx) const
 
u32vector getOwners (const Intersection *inter) const
 
bool getParallel () const
 
const std::vector< uint64_t > & getPartitionFirstDesc () const
 
const std::vector< uint64_t > & getPartitionLastDesc () const
 
const std::vector< uint64_t > & getPartitionRangeGlobalIdx () const
 
bool getPbound (const Intersection *inter) const
 
bool getPbound (const Octant *oct) const
 
bool getPbound (const Octant *oct, uint8_t face) const
 
bool getPbound (uint32_t idx) const
 
bool getPbound (uint32_t idx, uint8_t face) const
 
octantIterator getPboundOctantsBegin ()
 
octantIterator getPboundOctantsEnd ()
 
bvector getPeriodic () const
 
bool getPeriodic (uint8_t i) const
 
octantID getPersistentIdx (const Octant *oct) const
 
octantID getPersistentIdx (uint32_t idx) const
 
OctantgetPointOwner (const darray3 &point)
 
OctantgetPointOwner (const darray3 &point, bool &isghost)
 
OctantgetPointOwner (const dvector &point)
 
OctantgetPointOwner (const dvector &point, bool &isghost)
 
uint32_t getPointOwnerIdx (const darray3 &point) const
 
uint32_t getPointOwnerIdx (const darray3 &point, bool &isghost) const
 
uint32_t getPointOwnerIdx (const double *point) const
 
uint32_t getPointOwnerIdx (const double *point, bool &isghost) const
 
uint32_t getPointOwnerIdx (const dvector &point) const
 
uint32_t getPointOwnerIdx (const dvector &point, bool &isghost) const
 
int getPointOwnerRank (darray3 point)
 
void getPreMapping (u32vector &idx, std::vector< int8_t > &markers, std::vector< bool > &isghost)
 
int8_t getPreMarker (Octant *oct)
 
int8_t getPreMarker (uint32_t idx)
 
int getRank () const
 
bool getSerial () const
 
double getSize (const Intersection *inter) const
 
double getSize (const Octant *oct) const
 
double getSize (uint32_t idx) const
 
uint64_t getStatus () const
 
double getTol () const
 
double getVolume (const Octant *oct) const
 
double getVolume (uint32_t idx) const
 
double getX (const Octant *oct) const
 
double getX (uint32_t idx) const
 
double getX0 () const
 
double getY (const Octant *oct) const
 
double getY (uint32_t idx) const
 
double getY0 () const
 
double getZ (const Octant *oct) const
 
double getZ (uint32_t idx) const
 
double getZ0 () const
 
bool isCommSet () const
 
bool isEdgeOnOctant (const Octant *edgeOctant, uint8_t edgeIndex, const Octant *octant) const
 
bool isFaceOnOctant (const Octant *faceOctant, uint8_t faceIndex, const Octant *octant) const
 
bool isNodeOnOctant (const Octant *nodeOctant, uint8_t nodeIndex, const Octant *octant) const
 
double levelToSize (uint8_t level) const
 
void loadBalance (const dvector *weight=NULL)
 
template<class Impl >
void loadBalance (DataLBInterface< Impl > &userData, dvector *weight=NULL)
 
template<class Impl >
void loadBalance (DataLBInterface< Impl > &userData, uint8_t &level, dvector *weight=NULL)
 
void loadBalance (uint8_t &level, const dvector *weight=NULL)
 
void preadapt ()
 
void printHeader ()
 
template<class Impl >
void privateLoadBalance (const uint32_t *partition, DataLBInterface< Impl > *userData=nullptr)
 
void replaceComm (MPI_Comm communicator, MPI_Comm *previousCommunicator=nullptr)
 
virtual void reset ()
 
virtual void restore (std::istream &stream)
 
void setBalance (Octant *oct, bool balance)
 
void setBalance (uint32_t idx, bool balance)
 
void setBalanceCodimension (uint8_t b21codim)
 
void setComm (MPI_Comm communicator)
 
void setMarker (Octant *oct, int8_t marker)
 
void setMarker (uint32_t idx, int8_t marker)
 
void setNofGhostLayers (std::size_t nofGhostLayers)
 
void setPeriodic (uint8_t i)
 
void settleMarkers ()
 
void setTol (double tol=1.0e-14)
 
void updateConnectivity ()
 
void write (const std::string &filename)
 
void writeTest (const std::string &filename, dvector data)
 

Static Public Attributes

static BITPIT_PUBLIC_API const std::string DEFAULT_LOG_FILE = "PABLO"
 

Detailed Description

Para Tree is the user interface class.

Date
17/dec/2015
Authors
Marco Cisternino
Edoardo Lombardi

The user should (read can...) work only with this class and its methods. The sizes are intended in reference physical domain with limits [0,1]. The transformation from the logical domain to the physical domain is defined by an internal mapping.

The partition of the octree is performed by following the Z-curve defined by the Morton index of the octants. By default it is a balanced partition over the number of octants for each process.

Class ParaTree has a dimensional parameter int dim and it accepts only two values: dim=2 and dim=3, for 2D and 3D respectively.

Definition at line 113 of file ParaTree.hpp.

Member Enumeration Documentation

◆ Operation

enum bitpit::ParaTree::Operation

Definition at line 123 of file ParaTree.hpp.

Constructor & Destructor Documentation

◆ ParaTree() [1/4]

bitpit::ParaTree::ParaTree ( const std::string & logfile = DEFAULT_LOG_FILE,
MPI_Comm comm = MPI_COMM_WORLD )

Default empty constructor of ParaTree.

Parameters
[in]logfileThe file name for the log of this object. PABLO.log is the default value.
[in]commThe MPI communicator used by the parallel octree. MPI_COMM_WORLD is the default value.

Definition at line 136 of file ParaTree.cpp.

◆ ParaTree() [2/4]

bitpit::ParaTree::ParaTree ( uint8_t dim,
const std::string & logfile = DEFAULT_LOG_FILE,
MPI_Comm comm = MPI_COMM_WORLD )

Default constructor of ParaTree. It builds one octant with node 0 in the Origin (0,0,0) and side of length 1.

Parameters
[in]dimThe space dimension of the m_octree.
[in]logfileThe file name for the log of this object. PABLO.log is the default value.
[in]commThe MPI communicator used by the parallel octree. MPI_COMM_WORLD is the default value.

Definition at line 161 of file ParaTree.cpp.

◆ ParaTree() [3/4]

bitpit::ParaTree::ParaTree ( std::istream & stream,
const std::string & logfile = DEFAULT_LOG_FILE,
MPI_Comm comm = MPI_COMM_WORLD )

Creates a new octree restoring the octree saved in the specified stream.

Parameters
[in]logfileThe file name for the log of this object. PABLO.log is the default value.
streamis the stream to read from
[in]commThe MPI communicator used by the parallel octree. MPI_COMM_WORLD is the default value.

Definition at line 191 of file ParaTree.cpp.

◆ ~ParaTree()

bitpit::ParaTree::~ParaTree ( )
virtual

Default Destructor of ParaTree.

Definition at line 211 of file ParaTree.cpp.

◆ ParaTree() [4/4]

bitpit::ParaTree::ParaTree ( const ParaTree & other)

Copy constructor of ParaTree.

Parameters
[in]otherclass of type ParaTree

Definition at line 223 of file ParaTree.cpp.

Member Function Documentation

◆ adapt()

bool bitpit::ParaTree::adapt ( bool mapper_flag = false)

Adapt the octree mesh with user setup for markers and 2:1 balancing conditions.

Parameters
[in]mapper_flagTrue/False if you want/don't want to track the changes in structure octant by a mapper.
NOTE: if mapper_flag = true the adapt method ends after a single level (refining/coarsening) adaptation.
The resulting markers will be increased/decreased by one.
Returns
Boolean true if adapt has done something.
NOTE: if pre-adapt method is called before an adapt call the adapt method do not perform pre-adapt process (no 2:1 balance check).
Examples
PABLO_example_00001.cpp, PABLO_example_00002.cpp, PABLO_example_00004.cpp, PABLO_example_00005.cpp, PABLO_example_00007.cpp, PABLO_example_00008.cpp, PABLO_example_00009.cpp, and PABLO_example_00010.cpp.

Definition at line 3698 of file ParaTree.cpp.

◆ adaptGlobalCoarse()

bool bitpit::ParaTree::adaptGlobalCoarse ( bool mapper_flag = false)

Adapt the octree mesh coarsening all the octants by one level. Optionally track the changes in structure octant by a mapper.

Parameters
[in]mapper_flagTrue/false for tracking/not tracking the changes in structure octant .
Examples
PABLO_example_00002.cpp, and PABLO_example_00005.cpp.

Definition at line 3794 of file ParaTree.cpp.

◆ adaptGlobalRefine()

bool bitpit::ParaTree::adaptGlobalRefine ( bool mapper_flag = false)

Adapt the octree mesh refining all the octants by one level. Optionally track the changes in structure octant by a mapper.

Parameters
[in]mapper_flagTrue/false for tracking/not tracking the changes in structure octant .
Examples
PABLO_example_00001.cpp, PABLO_example_00002.cpp, PABLO_example_00003.cpp, PABLO_example_00004.cpp, PABLO_example_00005.cpp, PABLO_example_00007.cpp, PABLO_example_00008.cpp, PABLO_example_00009.cpp, PABLO_example_00010.cpp, and PABLO_example_00011.cpp.

Definition at line 3713 of file ParaTree.cpp.

◆ check21Balance()

bool bitpit::ParaTree::check21Balance ( )

Check if the grid is 2:1 balanced across intersection of balanceCodim codimension

Returns
true if the grid is 2:1 balanced

Definition at line 4065 of file ParaTree.cpp.

◆ checkToAdapt()

bool bitpit::ParaTree::checkToAdapt ( )

Check to control if the tree has to be adapted.

Returns
Boolean true if at least one octant has marker not zero.

Definition at line 3668 of file ParaTree.cpp.

◆ clearConnectivity()

void bitpit::ParaTree::clearConnectivity ( bool release = true)

Clear the connectivity of octants.

Parameters
[in]releaseif it's true the memory hold by the connectivity will be released, otherwise the connectivity will be cleared but its memory will not be released

Definition at line 3965 of file ParaTree.cpp.

◆ communicate()

template<class Impl >
void bitpit::ParaTree::communicate ( DataCommInterface< Impl > & userData)
inline

Communicate data provided by the user between the processes.

Definition at line 600 of file ParaTree.hpp.

◆ computeConnectivity()

void bitpit::ParaTree::computeConnectivity ( )

Compute the connectivity of octants and store the coordinates of nodes.

Examples
PABLO_example_00001.cpp, PABLO_example_00002.cpp, PABLO_example_00005.cpp, and PABLO_example_00011.cpp.

Definition at line 3955 of file ParaTree.cpp.

◆ computeIntersections()

void bitpit::ParaTree::computeIntersections ( )

Compute the intersection between octants (local, ghost, boundary).

Definition at line 4496 of file ParaTree.cpp.

◆ computeNodePersistentKey() [1/2]

uint64_t bitpit::ParaTree::computeNodePersistentKey ( const Octant * oct,
uint8_t node ) const

Compute the persistent XYZ key of the specified node of an octant (without level).

Parameters
[in]octPointer to the target octant
[in]nodeIndex of the target node.
Returns
persistent XYZ key of the node.

Definition at line 1961 of file ParaTree.cpp.

◆ computeNodePersistentKey() [2/2]

uint64_t bitpit::ParaTree::computeNodePersistentKey ( uint32_t idx,
uint8_t node ) const

Compute the persistent XYZ key of the specified node of an octant (without level).

Parameters
[in]idxLocal index of the target octant.
[in]nodeIndex of the target node.
Returns
persistent XYZ key of the node.

Definition at line 1501 of file ParaTree.cpp.

◆ dump()

void bitpit::ParaTree::dump ( std::ostream & stream,
bool full = true )
virtual

Write the octree to the specified stream.

Parameters
streamis the stream to write to
fullis the flag for a complete dump with mapping structureof last operation of the tree

Reimplemented in bitpit::PabloUniform.

Definition at line 502 of file ParaTree.cpp.

◆ evalLoadBalanceRanges() [1/2]

ParaTree::LoadBalanceRanges bitpit::ParaTree::evalLoadBalanceRanges ( dvector * weights)

Evaluate the elements of the current partition that will be exchanged with other processes during the load balance.

Parameters
[in]weightsare the weights of the local octants (if a null pointer is given a uniform distribution is used)
Returns
The ranges of local ids that will be exchanged with other processes.

Definition at line 4229 of file ParaTree.cpp.

◆ evalLoadBalanceRanges() [2/2]

ParaTree::LoadBalanceRanges bitpit::ParaTree::evalLoadBalanceRanges ( uint8_t level,
dvector * weights )

Evaluate the elements of the current partition that will be exchanged with other processes during the load balance.

The families of octants of a desired level are retained compact on the same process.

Parameters
[in]levelis the level of the families that will be retained compact on the same process
[in]weightsare the weights of the local octants (if a null pointer is given a uniform distribution is used)
Returns
The ranges of local ids that will be exchanged with other processes.

Definition at line 4267 of file ParaTree.cpp.

◆ expectedOctantAdapt()

void bitpit::ParaTree::expectedOctantAdapt ( const Octant * oct,
int8_t marker,
octvector * result ) const

Given an input testing marker, it gets the expected resulting octants of an adaption of a target octant.

Parameters
[in]octPointer to target octant
[in]markerTesting adaptation marker
[out]resultpointer to vector of expcted octants resulting from adaptation [marker<0 : one-coarsening octant father; marker=0, copy of the current input octant itself; marker>0 one-refinement octant children]

Definition at line 3327 of file ParaTree.cpp.

◆ findAllCodimensionNeighbours() [1/2]

void bitpit::ParaTree::findAllCodimensionNeighbours ( const Octant * oct,
u32vector & neighbours,
bvector & isghost )

Finds all the neighbours of an internal octant through all its boundaries of any codimension. Returns a vector with the index of neighbours in their structure (octants or ghosts) and sets isghost[i] = true if the i-th neighbour is ghost in the local tree. Neighbours are not sorted by Morton.

Parameters
[in]octpointer to the current octant
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2886 of file ParaTree.cpp.

◆ findAllCodimensionNeighbours() [2/2]

void bitpit::ParaTree::findAllCodimensionNeighbours ( uint32_t idx,
u32vector & neighbours,
bvector & isghost )

Finds all the neighbours of an internal octant through all its boundaries of any codimension. Returns a vector with the index of neighbours in their structure (octants or ghosts) and sets isghost[i] = true if the i-th neighbour is ghost in the local tree.

Parameters
[in]idxIndex of current octant
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2872 of file ParaTree.cpp.

◆ findAllNodeNeighbours() [1/2]

void bitpit::ParaTree::findAllNodeNeighbours ( const Octant * oct,
uint32_t node,
u32vector & neighbours,
bvector & isghost ) const

Finds all the neighbours of a node

Parameters
[in]octPointer to current octant
[in]nodeIndex of node passed through for neighbours finding
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2773 of file ParaTree.cpp.

◆ findAllNodeNeighbours() [2/2]

void bitpit::ParaTree::findAllNodeNeighbours ( uint32_t idx,
uint32_t node,
u32vector & neighbours,
bvector & isghost )

Finds all the neighbours of a node

Parameters
[in]idxIndex of current octant
[in]nodeIndex of node passed through for neighbours finding
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2857 of file ParaTree.cpp.

◆ findGhostAllCodimensionNeighbours() [1/2]

void bitpit::ParaTree::findGhostAllCodimensionNeighbours ( Octant * oct,
u32vector & neighbours,
bvector & isghost )

Finds all the neighbours of a ghost octant through all its boundaries of any codimension. Returns a vector with the index of neighbours in their structure (octants or ghosts) and sets isghost[i] = true if the i-th neighbour is ghost in the local tree. Neighbours are not sorted by Morton.

Parameters
[in]octpointer to the current octant
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2932 of file ParaTree.cpp.

◆ findGhostAllCodimensionNeighbours() [2/2]

void bitpit::ParaTree::findGhostAllCodimensionNeighbours ( uint32_t idx,
u32vector & neighbours,
bvector & isghost )

Finds all the neighbours of a ghost octant through all its boundaries of any codimension. Returns a vector with the index of neighbours in their structure (octants or ghosts) and sets isghost[i] = true if the i-th neighbour is ghost in the local tree.

Parameters
[in]idxIndex of current octant
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2918 of file ParaTree.cpp.

◆ findGhostNeighbours() [1/3]

void bitpit::ParaTree::findGhostNeighbours ( const Octant * oct,
uint8_t entityIdx,
uint8_t entityCodim,
u32vector & neighbours,
bvector & isghost ) const

Finds all the neighbours of ghost octant through the specified entity (face/edge/node). Returns a vector with the index of the neighbours in their container and sets isghost[i] = true if the i-th neighbour is ghost in the local tree.

Parameters
[in]octPointer to current ghost octant
[in]entityIdxIndex of face/edge/node passed through for neighbours finding
[in]entityCodimCodimension of the entity (1=face, 2=edge and 3=vertex for 3D trees, 1=face, 2=vertex for 2D trees)
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2760 of file ParaTree.cpp.

◆ findGhostNeighbours() [2/3]

void bitpit::ParaTree::findGhostNeighbours ( uint32_t idx,
uint8_t entityIdx,
uint8_t entityCodim,
u32vector & neighbours ) const

Finds the internal neighbours of the octant through the specified entity (face/edge/node). Returns a vector with the index of the neighbours in their container.

Parameters
[in]idxIndex of current octant
[in]entityIdxIndex of face/edge/node passed through for neighbours finding
[in]entityCodimCodimension of the entity (1=face, 2=edge and 3=vertex for 3D trees, 1=face, 2=vertex for 2D trees)
[out]neighboursVector with the index of the neighbours in their container

Definition at line 2725 of file ParaTree.cpp.

◆ findGhostNeighbours() [3/3]

void bitpit::ParaTree::findGhostNeighbours ( uint32_t idx,
uint8_t entityIdx,
uint8_t entityCodim,
u32vector & neighbours,
bvector & isghost ) const

Finds the ghost neighbours of the octant through the specified entity (face/edge/node). Returns a vector with the index of the neighbours in their container and sets isghost[i] = true if the i-th neighbour is ghost in the local tree.

Parameters
[in]idxIndex of current octant
[in]entityIdxIndex of face/edge/node passed through for neighbours finding
[in]entityCodimCodimension of the entity (1=face, 2=edge and 3=vertex for 3D trees, 1=face, 2=vertex for 2D trees)
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2743 of file ParaTree.cpp.

◆ findNeighbours() [1/3]

void bitpit::ParaTree::findNeighbours ( const Octant * oct,
uint8_t entityIdx,
uint8_t entityCodim,
u32vector & neighbours,
bvector & isghost ) const

Finds all the neighbours of an octant through the specified entity (face/edge/node). Returns a vector with the index of the neighbours in their container and sets isghost[i] = true if the i-th neighbour is ghost in the local tree.

Parameters
[in]octPointer to current octant
[in]entityIdxIndex of face/edge/node passed through for neighbours finding
[in]entityCodimCodimension of the entity (1=face, 2=edge and 3=vertex for 3D trees, 1=face, 2=vertex for 2D trees)
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.

Definition at line 2711 of file ParaTree.cpp.

◆ findNeighbours() [2/3]

void bitpit::ParaTree::findNeighbours ( uint32_t idx,
uint8_t entityIdx,
uint8_t entityCodim,
u32vector & neighbours ) const

Finds all the internal neighbours of a local octant through the specified entity (face/edge/node). Returns a vector with the index of the neighbours in their container.

Parameters
[in]idxIndex of current octant
[in]entityIdxIndex of face/edge/node passed through for neighbours finding
[in]entityCodimCodimension of the entity (1=face, 2=edge and 3=vertex for 3D trees, 1=face, 2=vertex for 2D trees)
[out]neighboursVector of internal neighbours indices in octants/ghosts structure

Definition at line 2694 of file ParaTree.cpp.

◆ findNeighbours() [3/3]

void bitpit::ParaTree::findNeighbours ( uint32_t idx,
uint8_t entityIdx,
uint8_t entityCodim,
u32vector & neighbours,
bvector & isghost ) const

Finds the neighbours (both local and ghost ones) of the octant through the specified entity (face/edge/node). Returns a vector with the index of the neighbours in their container and sets isghost[i] = true if the i-th neighbour is ghost in the local tree.

Parameters
[in]idxIndex of current octant
[in]entityIdxIndex of face/edge/node passed through for neighbours finding
[in]entityCodimCodimension of the entity (1=face, 2=edge and 3=vertex for 3D trees, 1=face, 2=vertex for 2D trees)
[out]neighboursVector with the index of the neighbours in their container
[out]isghostVector with boolean flag; true if the respective octant in neighbours is a ghost octant. Can be ignored in serial runs.
Examples
PABLO_example_00003.cpp, and PABLO_example_00011.cpp.

Definition at line 2679 of file ParaTree.cpp.

◆ findOwner()

int bitpit::ParaTree::findOwner ( uint64_t morton) const

It finds the process owning the element definded by the Morton number passed as argument The Morton number can be obtained using the method getMorton() of Octant.

Parameters
[in]mortonMorton number of the element you want find the owner of
Returns
Rank of the process owning the element

Definition at line 3893 of file ParaTree.cpp.

◆ freeComm()

void bitpit::ParaTree::freeComm ( )

Free the MPI communicator

Definition at line 916 of file ParaTree.cpp.

◆ getArea() [1/3]

double bitpit::ParaTree::getArea ( const Intersection * inter) const

Get the area of an intersection (for 2D case the same value of getSize).

Parameters
[in]interPointer to target intersection.
Returns
Area of intersection.

Definition at line 2384 of file ParaTree.cpp.

◆ getArea() [2/3]

double bitpit::ParaTree::getArea ( const Octant * oct) const

Get the area of an octant (for 2D case the same value of getSize).

Parameters
[in]octPointer to the target octant
Returns
Area of octant.

Definition at line 1771 of file ParaTree.cpp.

◆ getArea() [3/3]

double bitpit::ParaTree::getArea ( uint32_t idx) const

Get the area of an octant (for 2D case the same value of getSize).

Parameters
[in]idxLocal index of target octant.
Returns
Area of octant.

Definition at line 1320 of file ParaTree.cpp.

◆ getBalance() [1/2]

bool bitpit::ParaTree::getBalance ( const Octant * oct) const

Get the balancing condition of an octant.

Parameters
[in]octPointer to the target octant
Returns
Has octant to be balanced?

Definition at line 1970 of file ParaTree.cpp.

◆ getBalance() [2/2]

bool bitpit::ParaTree::getBalance ( uint32_t idx) const

Get the balancing condition of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Has octant to be balanced?

Definition at line 1510 of file ParaTree.cpp.

◆ getBalanceCodimension()

uint8_t bitpit::ParaTree::getBalanceCodimension ( ) const

Get the codimension for 2:1 balancing

Returns
Maximum codimension of the entity through which the 2:1 balance is performed.

Definition at line 2174 of file ParaTree.cpp.

◆ getBordersPerProc()

const std::map< int, u32vector > & bitpit::ParaTree::getBordersPerProc ( ) const

Get a map of border octants per process

Returns
A map of border octants per process

Definition at line 2560 of file ParaTree.cpp.

◆ getBound() [1/5]

bool bitpit::ParaTree::getBound ( const Intersection * inter) const

Get if an intersection is a boundary domain intersection.

Parameters
[in]interPointer to target intersection.
Returns
Boundary or not boundary?.

Definition at line 2295 of file ParaTree.cpp.

◆ getBound() [2/5]

bool bitpit::ParaTree::getBound ( const Octant * oct) const

Get the bound condition of the octant

Parameters
[in]octPointer to the target octant
Returns
Is the octant a boundary octant?

Definition at line 1989 of file ParaTree.cpp.

◆ getBound() [3/5]

bool bitpit::ParaTree::getBound ( const Octant * oct,
uint8_t face ) const

Get the bound condition of the face of the octant

Parameters
[in]octPointer to the target octant
[in]faceIndex of the face
Returns
Is the face a boundary face?

Definition at line 1980 of file ParaTree.cpp.

◆ getBound() [4/5]

bool bitpit::ParaTree::getBound ( uint32_t idx) const

Get the bound condition of the face of the octant

Parameters
[in]idxLocal index of the target octant
Returns
Is the octant a boundary octant?

Definition at line 1529 of file ParaTree.cpp.

◆ getBound() [5/5]

bool bitpit::ParaTree::getBound ( uint32_t idx,
uint8_t face ) const

Get the bound condition of the face of the octant

Parameters
[in]idxLocal index of the target octant
[in]faceIndex of the face
Returns
Is the face a boundary face?

Definition at line 1520 of file ParaTree.cpp.

◆ getCenter() [1/5]

darray3 bitpit::ParaTree::getCenter ( const Intersection * inter) const

Get the coordinates of the center of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
Coordinates of the center of intersection.

Definition at line 2398 of file ParaTree.cpp.

◆ getCenter() [2/5]

darray3 bitpit::ParaTree::getCenter ( const Octant * oct) const

Get the coordinates of the center of an octant.

Parameters
[in]octPointer to the target octant
Returns
Coordinates of the center of octant.

Definition at line 1799 of file ParaTree.cpp.

◆ getCenter() [3/5]

void bitpit::ParaTree::getCenter ( const Octant * oct,
darray3 & centerCoords ) const

Get the coordinates of the center of an octant.

Parameters
[in]octPointer to the target octant
[out]centerCoordsCoordinates of the center of octant.

Definition at line 1789 of file ParaTree.cpp.

◆ getCenter() [4/5]

darray3 bitpit::ParaTree::getCenter ( uint32_t idx) const

Get the coordinates of the center of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Coordinates of the center of octant.

Definition at line 1348 of file ParaTree.cpp.

◆ getCenter() [5/5]

void bitpit::ParaTree::getCenter ( uint32_t idx,
darray3 & centerCoords ) const

Get the coordinates of the center of an octant.

Parameters
[in]idxLocal index of target octant.
[out]centerCoordsCoordinates of the center of octant.

Definition at line 1338 of file ParaTree.cpp.

◆ getComm()

MPI_Comm bitpit::ParaTree::getComm ( ) const

Get thecommunicator used by octree between processes.

Returns
MPI Communicator.

Definition at line 944 of file ParaTree.cpp.

◆ getConnectivity() [1/3]

const u32vector2D & bitpit::ParaTree::getConnectivity ( ) const

Get the connectivity of the octants

Returns
Constant reference to the connectivity matrix of noctants*nnodes with the connectivity of each octant (4/8 indices of nodes for 2D/3D case).

Definition at line 3981 of file ParaTree.cpp.

◆ getConnectivity() [2/3]

const u32vector & bitpit::ParaTree::getConnectivity ( Octant * oct) const

Get the local connectivity of an octant

Parameters
[in]octPointer to an octant
Returns
Constant reference to the connectivity of the octant (4/8 indices of nodes for 2D/3D case).

Definition at line 4000 of file ParaTree.cpp.

◆ getConnectivity() [3/3]

const u32vector & bitpit::ParaTree::getConnectivity ( uint32_t idx) const

Get the local connectivity of an octant

Parameters
[in]idxLocal index of octant
Returns
Constant reference to the connectivity of the octant (4/8 indices of nodes for 2D/3D case).

Definition at line 3991 of file ParaTree.cpp.

◆ getCoordinates() [1/2]

darray3 bitpit::ParaTree::getCoordinates ( const Octant * oct) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]octPointer to the target octant
Returns
Coordinates of node 0.

Definition at line 1726 of file ParaTree.cpp.

◆ getCoordinates() [2/2]

darray3 bitpit::ParaTree::getCoordinates ( uint32_t idx) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]idxLocal index of target octant.
Returns
Coordinates X,Y,Z of node 0.

Definition at line 1275 of file ParaTree.cpp.

◆ getDim()

uint8_t bitpit::ParaTree::getDim ( ) const

Get the dimension of the octree.

Returns
Dimension of the octree (2D/3D).

Definition at line 780 of file ParaTree.cpp.

◆ getDumpVersion()

int bitpit::ParaTree::getDumpVersion ( ) const
virtual

Get the version associated to the binary dumps.

Returns
The version associated to the binary dumps.

Reimplemented in bitpit::PabloUniform.

Definition at line 488 of file ParaTree.cpp.

◆ getEdgecoeffs() [1/2]

const int8_t(* bitpit::ParaTree::getEdgecoeffs ( ) )[3]

Get the normals per edge (in 2D case not to be considered at all).

Returns
Pointer to components (x,y,z) of the "normals" per edge.

Definition at line 1219 of file ParaTree.cpp.

◆ getEdgecoeffs() [2/2]

void bitpit::ParaTree::getEdgecoeffs ( int8_t edgecoeffs[12][3]) const

Get the normals per edge (in 2D case not to be considered at all).

Parameters
[out]edgecoeffsComponents (x,y,z) of the "normals" per edge.

Definition at line 1155 of file ParaTree.cpp.

◆ getEdgeface() [1/2]

const uint8_t(* bitpit::ParaTree::getEdgeface ( ) )[2]

Get the edge-face connectivity for 12 edge (in 2D case not to be considered at all).

Returns
Pointer to connectivity edge-face. edgeface[i][0:1] = local indices of faces sharing the i-th edge of an octant.

Definition at line 1203 of file ParaTree.cpp.

◆ getEdgeface() [2/2]

void bitpit::ParaTree::getEdgeface ( uint8_t edgeface[12][2]) const

Get the edge-face connectivity for 12 edge (in 2D case not to be considered at all).

Parameters
[out]edgefaceConnectivity edge-face. edgeface[i][0:1] = local indices of faces sharing the i-th edge of an octant.

Definition at line 1130 of file ParaTree.cpp.

◆ getFace()

uint8_t bitpit::ParaTree::getFace ( const Intersection * inter) const

Get the face index of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
Face index of the finer octant of intersection (owners[getFiner(inter)]).

Definition at line 2322 of file ParaTree.cpp.

◆ getFaceCenter() [1/4]

darray3 bitpit::ParaTree::getFaceCenter ( const Octant * oct,
uint8_t face ) const

Get the coordinates of the center of a face of an octant.

Parameters
[in]octPointer to the target octant
[in]faceIndex of the target face.
Returns
Coordinates of the center of the face-th face af octant.

Definition at line 1812 of file ParaTree.cpp.

◆ getFaceCenter() [2/4]

void bitpit::ParaTree::getFaceCenter ( const Octant * oct,
uint8_t face,
darray3 & centerCoords ) const

Get the coordinates of the center of a face of an octant.

Parameters
[in]octPointer to the target octant
[in]faceIndex of the target face.
[out]centerCoordsCoordinates of the center of the face-th face af octant.

Definition at line 1825 of file ParaTree.cpp.

◆ getFaceCenter() [3/4]

darray3 bitpit::ParaTree::getFaceCenter ( uint32_t idx,
uint8_t face ) const

Get the coordinates of the center of a face of an octant.

Parameters
[in]idxLocal index of target octant.
[in]faceIndex of the target face.
Returns
Coordinates of the center of the face-th face af octant.

Definition at line 1361 of file ParaTree.cpp.

◆ getFaceCenter() [4/4]

void bitpit::ParaTree::getFaceCenter ( uint32_t idx,
uint8_t face,
darray3 & centerCoords ) const

Get the coordinates of the center of a face of an octant.

Parameters
[in]idxLocal index of target octant.
[in]faceIndex of the target face.
[out]centerCoordsCoordinates of the center of the face-th face af octant.

Definition at line 1374 of file ParaTree.cpp.

◆ getFacenode() [1/2]

const uint8_t(* bitpit::ParaTree::getFacenode ( ) )[4]

Get the face-node connectivity for 6 faces (in 2D case consider only the first 4 terms).

Returns
Pointer to connectivity face-node. facenode[i][0:1] = local indices of nodes of the i-th face of an octant.

Definition at line 1185 of file ParaTree.cpp.

◆ getFacenode() [2/2]

void bitpit::ParaTree::getFacenode ( uint8_t facenode[6][4]) const

Get the face-node connectivity for 6 faces (in 2D case consider only the first 4 terms).

Parameters
[out]facenodeConnectivity face-node. facenode[i][0:1] = local indices of nodes of the i-th face of an octant.

Definition at line 1104 of file ParaTree.cpp.

◆ getFamilySplittingNode()

uint8_t bitpit::ParaTree::getFamilySplittingNode ( const Octant * oct) const

Get the local index of the node of a target octant, corresponding to the splitting node of its family; i.e. the index of the local node coincident with the center point of its father.

Parameters
[in]octPointer to target octant
Returns
Local index of octant node corresponding to the splitting family node.

Definition at line 3317 of file ParaTree.cpp.

◆ getFiner()

bool bitpit::ParaTree::getFiner ( const Intersection * inter) const

Get the finer owner octant of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
The finer octant of the owners of intersection (false/true = 0/1).

Definition at line 2286 of file ParaTree.cpp.

◆ getFirstDescMorton()

uint64_t bitpit::ParaTree::getFirstDescMorton ( ) const

Get the first possible descendant with maximum refinement level of the local tree.

Returns
Constant reference to the first finest descendant of the local tree.

Definition at line 2181 of file ParaTree.cpp.

◆ getGhostConnectivity() [1/3]

const u32vector2D & bitpit::ParaTree::getGhostConnectivity ( ) const

Get the connectivity of the ghost octants

Returns
Constant reference to connectivity matrix [nghostoctants*nnodes] with the connectivity of each octant (4/8 indices of nodes for 2D/3D case).

Definition at line 4035 of file ParaTree.cpp.

◆ getGhostConnectivity() [2/3]

const u32vector & bitpit::ParaTree::getGhostConnectivity ( const Octant * oct) const

Get the local connectivity of a ghost octant

Parameters
[in]octPointer to a ghost octant
Returns
Constant reference to the connectivity of the ghost octant (4/8 indices of nodes for 2D/3D case).

Definition at line 4055 of file ParaTree.cpp.

◆ getGhostConnectivity() [3/3]

const u32vector & bitpit::ParaTree::getGhostConnectivity ( uint32_t idx) const

Get the local connectivity of a ghost octant

Parameters
[in]idxLocal index of ghost octant
Returns
Constant reference to the connectivity of the ghost octant (4/8 indices of nodes for 2D/3D case).

Definition at line 4045 of file ParaTree.cpp.

◆ getGhostGlobalIdx()

uint64_t bitpit::ParaTree::getGhostGlobalIdx ( uint32_t idx) const

Get the global index of a ghost octant.

Parameters
[in]idxLocal index of target ghost octant.
Returns
Global index of ghost octant.

Definition at line 1650 of file ParaTree.cpp.

◆ getGhostLayer()

int bitpit::ParaTree::getGhostLayer ( const Octant * oct) const

Get the layer number of the ghost halo an octant belong to.

Parameters
[in]octPointer to target octant.
Returns
the layer in the ghost halo. For internal (non-ghost) octants the function will return -1.

Definition at line 2518 of file ParaTree.cpp.

◆ getGhostLocalIdx()

uint32_t bitpit::ParaTree::getGhostLocalIdx ( uint64_t gidx) const

Get the local index of a ghost octant.

Parameters
[in]gidxGlobal index of target octant.
Returns
Local index of the ghost octant.

Definition at line 1630 of file ParaTree.cpp.

◆ getGhostOctant() [1/2]

Octant * bitpit::ParaTree::getGhostOctant ( uint32_t idx)

Get a ghost octant as pointer to the target octant. NOTE: no checks will be performed on the ghost octant index.

Parameters
[in]idxLocal index (in ghosts structure) of target ghost octant.
Returns
Pointer to target ghost octant.

Definition at line 2489 of file ParaTree.cpp.

◆ getGhostOctant() [2/2]

const Octant * bitpit::ParaTree::getGhostOctant ( uint32_t idx) const

Get a ghost octant as constant pointer to the target octant. NOTE: no checks will be performed on the ghost octant index.

Parameters
[in]idxLocal index (in ghosts structure) of target ghost octant.
Returns
Constant pointer to target ghost octant.

Definition at line 2499 of file ParaTree.cpp.

◆ getGlobalIdx() [1/2]

uint64_t bitpit::ParaTree::getGlobalIdx ( const Octant * oct) const

Get the global index of an octant.

Parameters
[in]octPointer to target octant.
Returns
Global index of octant.

Definition at line 2049 of file ParaTree.cpp.

◆ getGlobalIdx() [2/2]

uint64_t bitpit::ParaTree::getGlobalIdx ( uint32_t idx) const

Get the global index of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Global index of octant.

Definition at line 1575 of file ParaTree.cpp.

◆ getGlobalNumOctants()

uint64_t bitpit::ParaTree::getGlobalNumOctants ( ) const

Get the global number of octants.

Returns
Global number of octants.

Definition at line 788 of file ParaTree.cpp.

◆ getIdx()

uint32_t bitpit::ParaTree::getIdx ( const Octant * oct) const

Get the local index of an octant.

Parameters
[in]octPointer to target octant.
Returns
Local index of octant.

Definition at line 2035 of file ParaTree.cpp.

◆ getIn()

uint32_t bitpit::ParaTree::getIn ( const Intersection * inter) const

Get the owner octant of an intersection with inner normal.

Parameters
[in]interPointer to target intersection.
Returns
Index of the octant owner with inner normal.

Definition at line 2343 of file ParaTree.cpp.

◆ getInternalOctantsBegin()

octantIterator bitpit::ParaTree::getInternalOctantsBegin ( )

Get the begin position for the iterator of the local internal octants.

Returns
Iterator begin of the local internal octants (dereferencing results in a pointer to an octant).

Definition at line 2205 of file ParaTree.cpp.

◆ getInternalOctantsEnd()

octantIterator bitpit::ParaTree::getInternalOctantsEnd ( )

Get the end position for the iterator of the local internal octants.

Returns
Iterator end of the local internal octants (dereferencing results in a pointer to an octant).

Definition at line 2213 of file ParaTree.cpp.

◆ getIntersection()

Intersection * bitpit::ParaTree::getIntersection ( uint32_t idx)

Get a pointer to target intersection.

Parameters
[in]idxLocal index of intersection.
Returns
Pointer to target intersection.

Definition at line 2262 of file ParaTree.cpp.

◆ getIsGhost() [1/2]

bool bitpit::ParaTree::getIsGhost ( const Intersection * inter) const

Get if an intersection is an intersection between an internal and a ghost element.

Parameters
[in]interPointer to target intersection.
Returns
Ghost or not ghost?.

Definition at line 2304 of file ParaTree.cpp.

◆ getIsGhost() [2/2]

bool bitpit::ParaTree::getIsGhost ( const Octant * oct) const

Get the nature of an octant.

Parameters
[in]octPointer to target octant.
Returns
Is octant ghost?

Definition at line 2508 of file ParaTree.cpp.

◆ getIsNewC() [1/2]

bool bitpit::ParaTree::getIsNewC ( const Octant * oct) const

Get if the octant is new after coarsening.

Parameters
[in]octPointer to the target octant
Returns
Is octant new?

Definition at line 2026 of file ParaTree.cpp.

◆ getIsNewC() [2/2]

bool bitpit::ParaTree::getIsNewC ( uint32_t idx) const

Get if the octant is new after coarsening.

Parameters
[in]idxLocal index of target octant.
Returns
Is octant new?
Examples
PABLO_example_00004.cpp, PABLO_example_00007.cpp, PABLO_example_00008.cpp, and PABLO_example_00009.cpp.

Definition at line 1566 of file ParaTree.cpp.

◆ getIsNewR() [1/2]

bool bitpit::ParaTree::getIsNewR ( const Octant * oct) const

Get if the octant is new after refinement.

Parameters
[in]octPointer to the target octant
Returns
Is octant new?

Definition at line 2017 of file ParaTree.cpp.

◆ getIsNewR() [2/2]

bool bitpit::ParaTree::getIsNewR ( uint32_t idx) const

Get if the octant is new after refinement.

Parameters
[in]idxLocal index of target octant.
Returns
Is octant new?
Examples
PABLO_example_00004.cpp, PABLO_example_00007.cpp, and PABLO_example_00009.cpp.

Definition at line 1557 of file ParaTree.cpp.

◆ getL()

double bitpit::ParaTree::getL ( ) const

Get the length of the domain.

Returns
Length of the octree.

Definition at line 1015 of file ParaTree.cpp.

◆ getLastDescMorton() [1/3]

uint64_t bitpit::ParaTree::getLastDescMorton ( ) const

Get the last possible descendant with maximum refinement level of the local tree.

Returns
Constant reference to the last finest descendant of the local tree.

Definition at line 2188 of file ParaTree.cpp.

◆ getLastDescMorton() [2/3]

uint64_t bitpit::ParaTree::getLastDescMorton ( const Octant * oct) const

Get the morton index of the last possible descendant with maximum refinement level of a target octant.

Parameters
[in]octPointer to the target octant
Returns
Morton index of the last finest descendant of the target octant.

Definition at line 1950 of file ParaTree.cpp.

◆ getLastDescMorton() [3/3]

uint64_t bitpit::ParaTree::getLastDescMorton ( uint32_t idx) const

Get the morton index of the last possible descendant with maximum refinement level of a target octant.

Parameters
[in]idxLocal index of the target octant.
Returns
Morton index of the last finest descendant of the target octant.

Definition at line 2197 of file ParaTree.cpp.

◆ getLastOperation()

ParaTree::Operation bitpit::ParaTree::getLastOperation ( ) const

Get the last operation perforfmed by the octree (initialization, adapt (mapped or unmapped), loadbalance (first or not).

Returns
Last operation performed by the octree.

Definition at line 836 of file ParaTree.cpp.

◆ getLevel() [1/3]

uint8_t bitpit::ParaTree::getLevel ( const Intersection * inter) const

Get the level of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
Level of intersection.

Definition at line 2274 of file ParaTree.cpp.

◆ getLevel() [2/3]

uint8_t bitpit::ParaTree::getLevel ( const Octant * oct) const

Get the level of an octant.

Parameters
[in]octPointer to the target octant
Returns
Level of octant.

Definition at line 1932 of file ParaTree.cpp.

◆ getLevel() [3/3]

uint8_t bitpit::ParaTree::getLevel ( uint32_t idx) const

Get the level of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Level of octant.
Examples
PABLO_example_00001.cpp, and PABLO_example_00010.cpp.

Definition at line 1481 of file ParaTree.cpp.

◆ getLoadBalanceRanges()

const ParaTree::LoadBalanceRanges & bitpit::ParaTree::getLoadBalanceRanges ( ) const

Get a map of elements sent to the other processes during load balance

Returns
an unordered map associating rank to sent elements given by index extremes of a chunck

Definition at line 2526 of file ParaTree.cpp.

◆ getLocalIdx() [1/3]

uint32_t bitpit::ParaTree::getLocalIdx ( uint64_t gidx) const

Get the local index of an octant.

Parameters
[in]gidxGlobal index of target octant.
Returns
Local index of octant.

Definition at line 1616 of file ParaTree.cpp.

◆ getLocalIdx() [2/3]

uint32_t bitpit::ParaTree::getLocalIdx ( uint64_t gidx,
int rank ) const

Get the local index of an octant.

Parameters
[in]gidxGlobal index of target octant.
[in]rankRank of the process owning the octant as a local one. getOwnerRank can be used to get the rank of the process owning the octant locally.
Returns
Local index of octant.

Definition at line 1590 of file ParaTree.cpp.

◆ getLocalIdx() [3/3]

void bitpit::ParaTree::getLocalIdx ( uint64_t gidx,
uint32_t & lidx,
int & rank ) const

Get the local index of an octant and the rank owning the octant.

Parameters
[in]gidxGlobal index of target octant.
[out]lidxIndex of the octant.
[out]rankRank of the process owning the octant as a local one.

Definition at line 1605 of file ParaTree.cpp.

◆ getLocalMaxDepth()

uint8_t bitpit::ParaTree::getLocalMaxDepth ( ) const

Get the local depth of the octree.

Returns
Local depth of the octree.

Definition at line 2142 of file ParaTree.cpp.

◆ getLocalMaxSize()

double bitpit::ParaTree::getLocalMaxSize ( ) const

Get the local current maximum size of the octree.

Returns
Local current maximum size of the local partition of the octree.

Definition at line 2159 of file ParaTree.cpp.

◆ getLocalMinSize()

double bitpit::ParaTree::getLocalMinSize ( ) const

Get the local current minimum size reached by the octree.

Returns
Local current minimum size of the local partition of the octree.

Definition at line 2150 of file ParaTree.cpp.

◆ getLog()

Logger & bitpit::ParaTree::getLog ( )

Get the logger.

Returns
Reference to logger object.

Definition at line 828 of file ParaTree.cpp.

◆ getMapping() [1/2]

void bitpit::ParaTree::getMapping ( uint32_t idx,
u32vector & mapper,
bvector & isghost ) const

Get mapping info of an octant after an adapting with tracking changes.

Parameters
[in]idxIndex of new octant.
[out]mapperMapper from new octants to old octants. I.e. mapper[i] = j -> the i-th octant after adapt was in the j-th position before adapt; if the i-th octant is new after refinement the j-th old octant was the father of the new octant; if the i-th octant is new after coarsening the j-th old octant was a child of the new octant (mapper size = 4).
[out]isghostInfo on ghostness of old octants. I.e. isghost[i] = true/false -> the mapper[i] = j-th old octant was a local/ghost octant.
Examples
PABLO_example_00004.cpp, PABLO_example_00007.cpp, PABLO_example_00008.cpp, and PABLO_example_00009.cpp.

Definition at line 3355 of file ParaTree.cpp.

◆ getMapping() [2/2]

void bitpit::ParaTree::getMapping ( uint32_t idx,
u32vector & mapper,
bvector & isghost,
ivector & rank ) const

Get mapping info of an octant after an adapting or loadbalance with tracking changes.

Parameters
[in]idxIndex of new octant.
[out]mapperMapper from new octants to old octants. I.e. mapper[i] = j -> the i-th octant after adapt was in the j-th position before adapt; if the i-th octant is new after refinement or loadbalance the j-th old octant was the father of the new octant or the same octant respectively; if the i-th octant is new after coarsening the j-th old octant was a child of the new octant (mapper size = 4).
[out]isghostInfo on ghostness of old octants.
[out]rankProcess where the octant was located before the adapt/loadbalance. I.e. isghost[i] = true/false -> the mapper[i] = j-th old octant was a local/ghost octant.

Definition at line 3405 of file ParaTree.cpp.

◆ getMarker() [1/2]

int8_t bitpit::ParaTree::getMarker ( const Octant * oct) const

Get the refinement marker of an octant.

Parameters
[in]octPointer to the target octant
Returns
Marker of octant.

Definition at line 1909 of file ParaTree.cpp.

◆ getMarker() [2/2]

int8_t bitpit::ParaTree::getMarker ( uint32_t idx) const

Get the refinement marker of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Marker of octant.

Definition at line 1458 of file ParaTree.cpp.

◆ getMaxDepth()

int8_t bitpit::ParaTree::getMaxDepth ( ) const

Get the current maximum size of the octree. If the tree is empty a negative number is returned.

Returns
Current maximum size of the octree.

Definition at line 3883 of file ParaTree.cpp.

◆ getMaxLength()

uint32_t bitpit::ParaTree::getMaxLength ( ) const

Get the length of the domain in logical domain.

Returns
Length of the octree in logical domain.

Definition at line 1031 of file ParaTree.cpp.

◆ getMaxLevel()

int bitpit::ParaTree::getMaxLevel ( ) const

Get the maximum level of refinement allowed for this octree.

Returns
Maximum refinement level for the octree.

Definition at line 1023 of file ParaTree.cpp.

◆ getMorton() [1/2]

uint64_t bitpit::ParaTree::getMorton ( const Octant * oct) const

Compute the Morton index of an octant (without level).

Parameters
[in]octPointer to the target octant
Returns
morton Morton index of the octant.

Definition at line 1941 of file ParaTree.cpp.

◆ getMorton() [2/2]

uint64_t bitpit::ParaTree::getMorton ( uint32_t idx) const

Compute the Morton index of an octant (without level).

Parameters
[in]idxLocal index of target octant.
Returns
morton Morton index of the octant.

Definition at line 1490 of file ParaTree.cpp.

◆ getNchildren()

uint8_t bitpit::ParaTree::getNchildren ( ) const

Get the number of possible children for each octant (4 for 2D case, 8 for 3D case)

Returns
Number of children for octant.

Definition at line 1063 of file ParaTree.cpp.

◆ getNedges()

uint8_t bitpit::ParaTree::getNedges ( ) const

Get the number of edges for each octant (0 for 2D case, 12 for 3D case)

Returns
Number of edges for octant.
Examples
PABLO_example_00011.cpp.

Definition at line 1055 of file ParaTree.cpp.

◆ getNfaces()

uint8_t bitpit::ParaTree::getNfaces ( ) const

Get the number of faces for each octant (4 for 2D case, 6 for 3D case)

Returns
Number of faces for octant.
Examples
PABLO_example_00011.cpp.

Definition at line 1047 of file ParaTree.cpp.

◆ getNnodes()

uint8_t bitpit::ParaTree::getNnodes ( ) const

Get the number of nodes for each octant (4 for 2D case, 8 for 3D case)

Returns
Number of nodes for octant.
Examples
PABLO_example_00011.cpp.

Definition at line 1039 of file ParaTree.cpp.

◆ getNnodesperface()

uint8_t bitpit::ParaTree::getNnodesperface ( ) const

Get the number of nodes for each face of an octant (2 for 2D case, 4 for 3D case)

Returns
Number of nodes for face for an octant.

Definition at line 1071 of file ParaTree.cpp.

◆ getNode() [1/4]

darray3 bitpit::ParaTree::getNode ( const Octant * oct,
uint8_t node ) const

Get the coordinates of single node of an octant.

Parameters
[in]octPointer to the target octant
[in]nodeIndex of the target node.
Returns
Coordinates of the node-th node of octant.

Definition at line 1836 of file ParaTree.cpp.

◆ getNode() [2/4]

void bitpit::ParaTree::getNode ( const Octant * oct,
uint8_t node,
darray3 & nodeCoords ) const

Get the coordinates of the center of a face of an octant.

Parameters
[in]octPointer to the target octant
[in]nodeIndex of the target node.
[out]nodeCoordsCoordinates of the node-th node of octant.

Definition at line 1849 of file ParaTree.cpp.

◆ getNode() [3/4]

darray3 bitpit::ParaTree::getNode ( uint32_t idx,
uint8_t node ) const

Get the coordinates of single node of an octant.

Parameters
[in]idxLocal index of target octant.
[in]nodeIndex of the target node.
Returns
Coordinates of the node-th node of octant.

Definition at line 1385 of file ParaTree.cpp.

◆ getNode() [4/4]

void bitpit::ParaTree::getNode ( uint32_t idx,
uint8_t node,
darray3 & nodeCoords ) const

Get the coordinates of the center of a face of an octant.

Parameters
[in]idxLocal index of target octant.
[in]nodeIndex of the target node.
[out]nodeCoordsCoordinates of the node-th node of octant.

Definition at line 1398 of file ParaTree.cpp.

◆ getNodecoeffs() [1/2]

const int8_t(* bitpit::ParaTree::getNodecoeffs ( ) )[3]

Get the normals of the nodes (in 2D case consider only the first 4).

Returns
Pointer to components (x,y,z) of the "normals" of the nodes.

Definition at line 1211 of file ParaTree.cpp.

◆ getNodecoeffs() [2/2]

void bitpit::ParaTree::getNodecoeffs ( int8_t nodecoeffs[8][3]) const

Get the normals of the nodes (in 2D case consider only the first 4).

Parameters
[out]nodecoeffsComponents (x,y,z) of the "normals" of the nodes.

Definition at line 1142 of file ParaTree.cpp.

◆ getNodeCoordinates()

darray3 bitpit::ParaTree::getNodeCoordinates ( uint32_t node) const

Get the physical coordinates of a node

Parameters
[in]nodeLocal index of node
Returns
Vector with the coordinates of the node.

Definition at line 4026 of file ParaTree.cpp.

◆ getNodeface() [1/2]

const uint8_t(* bitpit::ParaTree::getNodeface ( ) )[3]

Get the node-face connectivity for 8 nodes (in 2D case consider only the first 4 terms).

Returns
Pointer to connectivity node-face. nodeface[i][0:1] = local indices of faces sharing the i-th node of an octant.

Definition at line 1194 of file ParaTree.cpp.

◆ getNodeface() [2/2]

void bitpit::ParaTree::getNodeface ( uint8_t nodeface[8][3]) const

Get the node-face connectivity for 8 nodes (in 2D case consider only the first 4 terms).

Parameters
[out]nodefaceConnectivity node-face. nodeface[i][0:1] = local indices of faces sharing the i-th node of an octant.

Definition at line 1117 of file ParaTree.cpp.

◆ getNodeLogicalCoordinates()

const u32array3 & bitpit::ParaTree::getNodeLogicalCoordinates ( uint32_t node) const

Get the logical coordinates of a node

Parameters
[in]nodeLocal index of node
Returns
Constant reference to a vector containing the coordinates of the node.

Definition at line 4017 of file ParaTree.cpp.

◆ getNodes() [1/6]

const u32arr3vector & bitpit::ParaTree::getNodes ( ) const

Get the logical coordinates of the nodes

Returns
Constant reference to the nodes matrix [nnodes*3] with the coordinates of the nodes.

Definition at line 4008 of file ParaTree.cpp.

◆ getNodes() [2/6]

darr3vector bitpit::ParaTree::getNodes ( const Intersection * inter) const

Get the coordinates of the nodes of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
Coordinates of the nodes of intersection.

Definition at line 2419 of file ParaTree.cpp.

◆ getNodes() [3/6]

darr3vector bitpit::ParaTree::getNodes ( const Octant * oct) const

Get the coordinates of the nodes of an octant.

Parameters
[in]octPointer to the target octant
Returns
Coordinates of the nodes of octant.

Definition at line 1870 of file ParaTree.cpp.

◆ getNodes() [4/6]

void bitpit::ParaTree::getNodes ( const Octant * oct,
darr3vector & nodes ) const

Get the coordinates of the nodes of an octant.

Parameters
[in]octPointer to the target octant
[out]nodesCoordinates of the nodes of octant.

Definition at line 1859 of file ParaTree.cpp.

◆ getNodes() [5/6]

darr3vector bitpit::ParaTree::getNodes ( uint32_t idx) const

Get the coordinates of the nodes of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Coordinates of the nodes of octant.

Definition at line 1419 of file ParaTree.cpp.

◆ getNodes() [6/6]

void bitpit::ParaTree::getNodes ( uint32_t idx,
darr3vector & nodesCoords ) const

Get the coordinates of the nodes of an octant.

Parameters
[in]idxLocal index of target octant.
[out]nodesCoordinates of the nodes of octant.

Definition at line 1408 of file ParaTree.cpp.

◆ getNofGhostLayers()

std::size_t bitpit::ParaTree::getNofGhostLayers ( ) const

Get the number of ghosts layers

Returns
the number of ghosts layers

Definition at line 2534 of file ParaTree.cpp.

◆ getNormal() [1/5]

darray3 bitpit::ParaTree::getNormal ( const Intersection * inter) const

Get the normal of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
Coordinates of the normal of intersection.

Definition at line 2444 of file ParaTree.cpp.

◆ getNormal() [2/5]

darray3 bitpit::ParaTree::getNormal ( const Octant * oct,
uint8_t face ) const

Get the normal of a face of an octant.

Parameters
[in]octPointer to the target octant
[in]faceIndex of the face for normal computing.
Returns
Normal of the face.

Definition at line 1896 of file ParaTree.cpp.

◆ getNormal() [3/5]

void bitpit::ParaTree::getNormal ( const Octant * oct,
uint8_t face,
darray3 & normal ) const

Get the normal of a face of an octant.

Parameters
[in]octPointer to the target octant
[in]faceIndex of the face for normal computing.
[out]normalCoordinates of the normal of face.

Definition at line 1884 of file ParaTree.cpp.

◆ getNormal() [4/5]

darray3 bitpit::ParaTree::getNormal ( uint32_t idx,
uint8_t face ) const

Get the normal of a face of an octant.

Parameters
[in]idxLocal index of target octant.
[in]faceIndex of the face for normal computing.
Returns
Normal of the face.

Definition at line 1445 of file ParaTree.cpp.

◆ getNormal() [5/5]

void bitpit::ParaTree::getNormal ( uint32_t idx,
uint8_t face,
darray3 & normal ) const

Get the normal of a face of an octant.

Parameters
[in]idxLocal index of target octant.
[in]faceIndex of the face for normal computing.
[out]normalCoordinates of the normal of face.

Definition at line 1433 of file ParaTree.cpp.

◆ getNormals() [1/2]

const int8_t(* bitpit::ParaTree::getNormals ( ) )[3]

Get the components (in logical domain) of the 6 normals to the faces of an octant (for the 2D case consider only the first 4)

Returns
Pointer to normals array[6][3] to the faces of an octant.

Definition at line 1166 of file ParaTree.cpp.

◆ getNormals() [2/2]

void bitpit::ParaTree::getNormals ( int8_t normals[6][3]) const

Get the components (in logical domain) of the 6 normals to the faces of an octant (for the 2D case consider only the first 4)

Parameters
[out]normalsNormals array[6][3] to the faces of an octant.

Definition at line 1079 of file ParaTree.cpp.

◆ getNproc()

int bitpit::ParaTree::getNproc ( ) const

Get the total number of processes.

Returns
Number of processes.

Definition at line 820 of file ParaTree.cpp.

◆ getNumGhosts()

uint32_t bitpit::ParaTree::getNumGhosts ( ) const

Get the local number of ghost octants.

Returns
Local number of ghost octants.
Examples
PABLO_example_00007.cpp, PABLO_example_00008.cpp, and PABLO_example_00009.cpp.

Definition at line 2126 of file ParaTree.cpp.

◆ getNumIntersections()

uint32_t bitpit::ParaTree::getNumIntersections ( ) const

Get the local number of intersections.

Returns
Local number of intersections.

Definition at line 2253 of file ParaTree.cpp.

◆ getNumNodes()

uint32_t bitpit::ParaTree::getNumNodes ( ) const

Get the local number of nodes.

Returns
Local total number of nodes.

Definition at line 2134 of file ParaTree.cpp.

◆ getNumOctants()

uint32_t bitpit::ParaTree::getNumOctants ( ) const

◆ getOctant() [1/2]

Octant * bitpit::ParaTree::getOctant ( uint32_t idx)

Get an octant as pointer to the target octant. NOTE: no checks will be performed on the octant index.

Parameters
[in]idxLocal index of target octant.
Returns
Pointer to target octant.
Examples
PABLO_example_00001.cpp.

Definition at line 2469 of file ParaTree.cpp.

◆ getOctant() [2/2]

const Octant * bitpit::ParaTree::getOctant ( uint32_t idx) const

Get an octant as constant pointer to the target octant. NOTE: no checks will be performed on the octant index.

Parameters
[in]idxLocal index of target octant.
Returns
Constant pointer to target octant.

Definition at line 2479 of file ParaTree.cpp.

◆ getOppface() [1/2]

const uint8_t * bitpit::ParaTree::getOppface ( ) const

Get the indices of the faces of virtual octants opposed to the 6 faces of an octant (for the 2D case consider only the first 4 terms).

Returns
Pointer to opposed faces array[6] to the faces of an octant (oppface[i] = Index of the face of an octant neighbour through the i-th face of the current octant).

Definition at line 1177 of file ParaTree.cpp.

◆ getOppface() [2/2]

void bitpit::ParaTree::getOppface ( uint8_t oppface[6]) const

Get the indices of the faces of virtual octants opposed to the 6 faces of an octant (for the 2D case consider only the first 4 terms).

Parameters
[out]oppfaceOpposed faces array[6] to the faces of an octant (oppface[i] = Index of the face of an octant neighbour through the i-th face of the current octant).

Definition at line 1093 of file ParaTree.cpp.

◆ getOrigin()

darray3 bitpit::ParaTree::getOrigin ( ) const

Get the coordinates of the origin of the octree.

Returns
Coordinates of the origin.

Definition at line 983 of file ParaTree.cpp.

◆ getOut()

uint32_t bitpit::ParaTree::getOut ( const Intersection * inter) const

Get the owner octant of an intersection with outer normal.

Parameters
[in]interPointer to target intersection.
Returns
Index of the octant owner with outer normal.

Definition at line 2352 of file ParaTree.cpp.

◆ getOutIsGhost()

bool bitpit::ParaTree::getOutIsGhost ( const Intersection * inter) const

Get if the owner octant with outer normal is a ghost octant.

Parameters
[in]interPointer to target intersection.
Returns
Is the octant owner with outer normal a ghost octant?.

Definition at line 2361 of file ParaTree.cpp.

◆ getOwnerRank()

int bitpit::ParaTree::getOwnerRank ( uint64_t globalIndex) const

It finds the process owning the element definded by the global index passed as argument The global index can be computed using the methods getGlobalIdx or getGhostGlobalIdx.

Parameters
[in]globalIndexindex of the element you want find the owner of
Returns
Rank of the process owning the element

Definition at line 3937 of file ParaTree.cpp.

◆ getOwners()

u32vector bitpit::ParaTree::getOwners ( const Intersection * inter) const

Get the owner octants of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
A couple of octants owners of intersection.

Definition at line 2331 of file ParaTree.cpp.

◆ getParallel()

bool bitpit::ParaTree::getParallel ( ) const

Get if the octree is parallel.

Returns
Is the octree distributed?.

Definition at line 804 of file ParaTree.cpp.

◆ getPartitionFirstDesc()

const std::vector< uint64_t > & bitpit::ParaTree::getPartitionFirstDesc ( ) const

Get a constant reference to the vector containing the Morton number of the first octant on each process.

Returns
Constant reference to the vector containing the Morton number of the first octant on each process.

Definition at line 965 of file ParaTree.cpp.

◆ getPartitionLastDesc()

const std::vector< uint64_t > & bitpit::ParaTree::getPartitionLastDesc ( ) const

Get a constant reference to the vector containing the Morton number of the last possible octant on each process.

Returns
Constant reference to the vector containing the Morton number of the last possible octant on each process.

Definition at line 975 of file ParaTree.cpp.

◆ getPartitionRangeGlobalIdx()

const std::vector< uint64_t > & bitpit::ParaTree::getPartitionRangeGlobalIdx ( ) const

Get a constant reference to the vector containing the global index of the last existing octant in each process.

Returns
A constant reference to the vector containing the global index of the last existing octant in each process.

Definition at line 955 of file ParaTree.cpp.

◆ getPbound() [1/5]

bool bitpit::ParaTree::getPbound ( const Intersection * inter) const

Get if an intersection is a boundary intersection for a process.

Parameters
[in]interPointer to target intersection.
Returns
Process boundary or not boundary?.

Definition at line 2313 of file ParaTree.cpp.

◆ getPbound() [2/5]

bool bitpit::ParaTree::getPbound ( const Octant * oct) const

Get the partition bound condition of the face of the octant

Parameters
[in]octPointer to the target octant
Returns
Is the octant a partition boundary octant?

Definition at line 2008 of file ParaTree.cpp.

◆ getPbound() [3/5]

bool bitpit::ParaTree::getPbound ( const Octant * oct,
uint8_t face ) const

Get the partition bound condition of the face of the octant

Parameters
[in]octPointer to the target octant
[in]faceIndex of the face
Returns
Is the face a partition boundary face?

Definition at line 1999 of file ParaTree.cpp.

◆ getPbound() [4/5]

bool bitpit::ParaTree::getPbound ( uint32_t idx) const

Get the partition bound condition of the face of the octant

Parameters
[in]idxLocal index of the target octant
Returns
Is the octant a partition boundary octant?

Definition at line 1548 of file ParaTree.cpp.

◆ getPbound() [5/5]

bool bitpit::ParaTree::getPbound ( uint32_t idx,
uint8_t face ) const

Get the partition bound condition of the face of the octant

Parameters
[in]idxLocal index of the target octant
[in]faceIndex of the face
Returns
Is the face a partition boundary face?

Definition at line 1539 of file ParaTree.cpp.

◆ getPboundOctantsBegin()

octantIterator bitpit::ParaTree::getPboundOctantsBegin ( )

Get the begin position for the iterator of the local border of process octants.

Returns
Iterator begin of the local border of process octants (dereferencing results in a pointer to an octant).

Definition at line 2221 of file ParaTree.cpp.

◆ getPboundOctantsEnd()

octantIterator bitpit::ParaTree::getPboundOctantsEnd ( )

Get the end position for the iterator of the local border of process octants.

Returns
Iterator end of the local border of process octants (dereferencing results in a pointer to an octant).

Definition at line 2229 of file ParaTree.cpp.

◆ getPeriodic() [1/2]

bvector bitpit::ParaTree::getPeriodic ( ) const

Get the periodic condition of the boundaries.

Returns
Vector with the periodic conditions (true/false) of each boundary.

Definition at line 1228 of file ParaTree.cpp.

◆ getPeriodic() [2/2]

bool bitpit::ParaTree::getPeriodic ( uint8_t i) const

Get the periodic condition of a target boundary.

Parameters
[in]iIndex of the target boundary face.
Returns
Boolean with the periodic conditions (true/false) of the target boundary.

Definition at line 1237 of file ParaTree.cpp.

◆ getPersistentIdx() [1/2]

bitset< 72 > bitpit::ParaTree::getPersistentIdx ( const Octant * oct) const

Get the persistent index of an octant.

Parameters
[in]octPointer to the target octant
Returns
Persistent index of octant, i.e. a bitset composed by Morton index and level of octant.

Definition at line 2068 of file ParaTree.cpp.

◆ getPersistentIdx() [2/2]

bitset< 72 > bitpit::ParaTree::getPersistentIdx ( uint32_t idx) const

Get the persistent index of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Persistent index of octant, i.e. a bitset composed by Morton index and level of octant.

Definition at line 1683 of file ParaTree.cpp.

◆ getPointOwner() [1/4]

Octant * bitpit::ParaTree::getPointOwner ( const darray3 & point)

Get the internal octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
Returns
Pointer to octant owner of target point (=NULL if point is outside of the domain).

Definition at line 2972 of file ParaTree.cpp.

◆ getPointOwner() [2/4]

Octant * bitpit::ParaTree::getPointOwner ( const darray3 & point,
bool & isghost )

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
[out]isghostBoolean flag, true if the octant found is ghost
Returns
Pointer to octant owner of target point (=NULL if point is outside of the ghostd domain).

Definition at line 2986 of file ParaTree.cpp.

◆ getPointOwner() [3/4]

Octant * bitpit::ParaTree::getPointOwner ( const dvector & point)

Get the internal octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
Returns
Pointer to octant owner of target point (=NULL if point is outside of the domain).

Definition at line 2941 of file ParaTree.cpp.

◆ getPointOwner() [4/4]

Octant * bitpit::ParaTree::getPointOwner ( const dvector & point,
bool & isghost )

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
[out]isghostBoolean flag, true if the octant found is ghost
Returns
Pointer to octant owner of target point (=NULL if point is outside of the ghosted domain).

Definition at line 2955 of file ParaTree.cpp.

◆ getPointOwnerIdx() [1/6]

uint32_t bitpit::ParaTree::getPointOwnerIdx ( const darray3 & point) const

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
Returns
Index of octant owner of target point (max uint32_t representable if point outside of the domain).

Definition at line 3002 of file ParaTree.cpp.

◆ getPointOwnerIdx() [2/6]

uint32_t bitpit::ParaTree::getPointOwnerIdx ( const darray3 & point,
bool & isghost ) const

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
[out]isghostBoolean flag, true if the octant found is ghost
Returns
Index of octant owner of target point (max uint32_t representable if point outside of the ghosted domain).

Definition at line 3012 of file ParaTree.cpp.

◆ getPointOwnerIdx() [3/6]

uint32_t bitpit::ParaTree::getPointOwnerIdx ( const double * point) const

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
Returns
Index of octant owner of target point (max uint32_t representable if point outside of the domain).

Definition at line 3043 of file ParaTree.cpp.

◆ getPointOwnerIdx() [4/6]

uint32_t bitpit::ParaTree::getPointOwnerIdx ( const double * point,
bool & isghost ) const

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
[out]isghostBoolean flag, true if the octant found is ghost
Returns
Index of octant owner of target point (max uint32_t representable if point outside of the ghosted domain).

end ghosts search

Definition at line 3122 of file ParaTree.cpp.

◆ getPointOwnerIdx() [5/6]

uint32_t bitpit::ParaTree::getPointOwnerIdx ( const dvector & point) const

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
Returns
Index of octant owner of target point (max uint32_t representable if point outside of the domain).

Definition at line 3021 of file ParaTree.cpp.

◆ getPointOwnerIdx() [6/6]

uint32_t bitpit::ParaTree::getPointOwnerIdx ( const dvector & point,
bool & isghost ) const

Get the octant owner of an input point.

Parameters
[in]pointCoordinates of target point.
[out]isghostBoolean flag, true if the octant found is ghost
Returns
Index of octant owner of target point (max uint32_t representable if point outside of the ghosted domain).

Definition at line 3032 of file ParaTree.cpp.

◆ getPointOwnerRank()

int bitpit::ParaTree::getPointOwnerRank ( darray3 point)

Get the octant owner rank of an input point.

Parameters
[in]pointCoordinates of target point.
Returns
Owner rank of target point (negative if out of global domain).

Definition at line 3266 of file ParaTree.cpp.

◆ getPreMapping()

void bitpit::ParaTree::getPreMapping ( u32vector & idx,
std::vector< int8_t > & markers,
std::vector< bool > & isghost )

Get octants with marker different from zero and the related markers and ghostness info. The methods has to be called after a apredapt, otherwise it calls preadapt method.

Parameters
[out]idxVector of local indices of octants with marker different from zero.
[out]markersVector with markers related to octants in the idx list.
[out]isghostInfo on ghostness of octants.

Definition at line 3441 of file ParaTree.cpp.

◆ getPreMarker() [1/2]

int8_t bitpit::ParaTree::getPreMarker ( Octant * oct)

Get the refinement marker of an octant after a preadapt.

Parameters
[in]octPointer to the target octant
Returns
Marker of octant. NOTE: if a last operation is not preadapt, it calls preadapt method.

Definition at line 1919 of file ParaTree.cpp.

◆ getPreMarker() [2/2]

int8_t bitpit::ParaTree::getPreMarker ( uint32_t idx)

Get the refinement marker of an octant after a preadapt.

Parameters
[in]idxLocal index of target octant.
Returns
Marker of octant. NOTE: if a last operation is not preadapt, it calls preadapt method.

Definition at line 1468 of file ParaTree.cpp.

◆ getRank()

int bitpit::ParaTree::getRank ( ) const

Get the rank of local process.

Returns
Rank of local process.

Definition at line 812 of file ParaTree.cpp.

◆ getSerial()

bool bitpit::ParaTree::getSerial ( ) const

Get if the octree is serial.

Returns
Is the octree serial?.

Definition at line 796 of file ParaTree.cpp.

◆ getSize() [1/3]

double bitpit::ParaTree::getSize ( const Intersection * inter) const

Get the size of an intersection.

Parameters
[in]interPointer to target intersection.
Returns
Size of intersection.

Definition at line 2370 of file ParaTree.cpp.

◆ getSize() [2/3]

double bitpit::ParaTree::getSize ( const Octant * oct) const

Get the size of an octant, i.e. the side length.

Parameters
[in]octPointer to the target octant
Returns
Size of octant.

Definition at line 1762 of file ParaTree.cpp.

◆ getSize() [3/3]

double bitpit::ParaTree::getSize ( uint32_t idx) const

Get the size of an octant, i.e. the side length.

Parameters
[in]idxLocal index of target octant.
Returns
Size of octant.

Definition at line 1311 of file ParaTree.cpp.

◆ getStatus()

uint64_t bitpit::ParaTree::getStatus ( ) const

Get the status label of the octree.

Returns
Status label of the octree.

Definition at line 2110 of file ParaTree.cpp.

◆ getTol()

double bitpit::ParaTree::getTol ( ) const

Get the tolerance used in geometric operations.

Definition at line 1244 of file ParaTree.cpp.

◆ getVolume() [1/2]

double bitpit::ParaTree::getVolume ( const Octant * oct) const

Get the volume of an octant.

Parameters
[in]octPointer to the target octant
Returns
Volume of octant.

Definition at line 1780 of file ParaTree.cpp.

◆ getVolume() [2/2]

double bitpit::ParaTree::getVolume ( uint32_t idx) const

Get the volume of an octant.

Parameters
[in]idxLocal index of target octant.
Returns
Volume of octant.

Definition at line 1329 of file ParaTree.cpp.

◆ getX() [1/2]

double bitpit::ParaTree::getX ( const Octant * oct) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]octPointer to the target octant
Returns
Coordinate X of node 0.

Definition at line 1735 of file ParaTree.cpp.

◆ getX() [2/2]

double bitpit::ParaTree::getX ( uint32_t idx) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]idxLocal index of target octant.
Returns
Coordinate X of node 0.

Definition at line 1284 of file ParaTree.cpp.

◆ getX0()

double bitpit::ParaTree::getX0 ( ) const

Get the coordinate X of the origin of the octree.

Returns
Coordinate X of the origin.

Definition at line 991 of file ParaTree.cpp.

◆ getY() [1/2]

double bitpit::ParaTree::getY ( const Octant * oct) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]octPointer to the target octant
Returns
Coordinate Y of node 0.

Definition at line 1744 of file ParaTree.cpp.

◆ getY() [2/2]

double bitpit::ParaTree::getY ( uint32_t idx) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]idxLocal index of target octant.
Returns
Coordinate Y of node 0.

Definition at line 1293 of file ParaTree.cpp.

◆ getY0()

double bitpit::ParaTree::getY0 ( ) const

Get the coordinate Y of the origin of the octree.

Returns
Coordinate Y of the origin.

Definition at line 999 of file ParaTree.cpp.

◆ getZ() [1/2]

double bitpit::ParaTree::getZ ( const Octant * oct) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]octPointer to the target octant
Returns
Coordinate Z of node 0.

Definition at line 1753 of file ParaTree.cpp.

◆ getZ() [2/2]

double bitpit::ParaTree::getZ ( uint32_t idx) const

Get the coordinates of an octant, i.e. the coordinates of its node 0.

Parameters
[in]idxLocal index of target octant.
Returns
Coordinate Z of node 0.

Definition at line 1302 of file ParaTree.cpp.

◆ getZ0()

double bitpit::ParaTree::getZ0 ( ) const

Get the coordinate Z of the origin of the octree.

Returns
Coordinate Z of the origin.

Definition at line 1007 of file ParaTree.cpp.

◆ isCommSet()

bool bitpit::ParaTree::isCommSet ( ) const

Check if the communicator to be used for parallel communications has already been set.

Returns
Returns true if the communicator has been set, false otherwise.

Definition at line 936 of file ParaTree.cpp.

◆ isEdgeOnOctant()

bool bitpit::ParaTree::isEdgeOnOctant ( const Octant * edgeOctant,
uint8_t edgeIndex,
const Octant * octant ) const

Check if an edge lies on the specified octant.

Parameters
[in]edgeOctantPointer to the octant owning the edge
[in]edgeIndexLocal index of the edge
[in]octantPointer to the octant for which the check has to be berformed

Definition at line 3538 of file ParaTree.cpp.

◆ isFaceOnOctant()

bool bitpit::ParaTree::isFaceOnOctant ( const Octant * faceOctant,
uint8_t faceIndex,
const Octant * octant ) const

Check if a face lies on the specified octant.

Parameters
[in]faceOctantPointer to the octant owning the face
[in]faceIndexLocal index of the face
[in]octantPointer to the octant for which the check has to be berformed

Definition at line 3584 of file ParaTree.cpp.

◆ isNodeOnOctant()

bool bitpit::ParaTree::isNodeOnOctant ( const Octant * nodeOctant,
uint8_t nodeIndex,
const Octant * octant ) const

Check if a node lies on the specified octant.

Parameters
[in]nodeOctantPointer to the octant owning the node
[in]nodeIndexLocal index of the node
[in]octantPointer to the octant for which the check has to be berformed

Definition at line 3497 of file ParaTree.cpp.

◆ levelToSize()

double bitpit::ParaTree::levelToSize ( uint8_t level) const

Get the size of an octant corresponding to a target level.

Parameters
[in]levelInput level.
Returns
Size of an octant of input level.

Definition at line 4484 of file ParaTree.cpp.

◆ loadBalance() [1/4]

void bitpit::ParaTree::loadBalance ( const dvector * weight = NULL)

Distribute Load-Balancing the octants (with user defined weights) of the whole tree over the processes of the job following the Morton order. Until loadBalance is not called for the first time the mesh is serial.

Parameters
[in]weightPointer to a vector of weights of the local octants (weight=NULL is uniform distribution).
Examples
PABLO_example_00005.cpp, PABLO_example_00007.cpp, PABLO_example_00008.cpp, PABLO_example_00009.cpp, and PABLO_example_00010.cpp.

Definition at line 4137 of file ParaTree.cpp.

◆ loadBalance() [2/4]

template<class Impl >
void bitpit::ParaTree::loadBalance ( DataLBInterface< Impl > & userData,
dvector * weight = NULL )
inline

Distribute Load-Balancing the octants (with user defined weights) of the whole tree and data provided by the user over the processes of the job following the Morton order. Until loadBalance is not called for the first time the mesh is serial. Even distribute data provided by the user between the processes.

Parameters
[in]userDataUser interface to distribute the data during loadBalance.
[in]weightPointer to a vector of weights of the local octants (weight=NULL is uniform distribution).

Definition at line 664 of file ParaTree.hpp.

◆ loadBalance() [3/4]

template<class Impl >
void bitpit::ParaTree::loadBalance ( DataLBInterface< Impl > & userData,
uint8_t & level,
dvector * weight = NULL )
inline

Distribute Load-Balanced the octants (with user defined weights) of the whole tree and data provided by the user over the processes of the job. Until loadBalance is not called for the first time the mesh is serial. The families of octants of a desired level are retained compact on the same process. Even distribute data provided by the user between the processes.

Parameters
[in]userDataUser interface to distribute the data during loadBalance.
[in]levelNumber of level over the max depth reached in the tree at which families of octants are fixed compact on the same process (level=0 is classic LoadBalance).
[in]weightPointer to a vector of weights of the local octants (weight=NULL is uniform distribution).

Definition at line 715 of file ParaTree.hpp.

◆ loadBalance() [4/4]

void bitpit::ParaTree::loadBalance ( uint8_t & level,
const dvector * weight = NULL )

Distribute Load-Balanced the octants (with user defined weights) of the whole tree over the processes of the job. Until loadBalance is not called for the first time the mesh is serial. The families of octants of a desired level are retained compact on the same process.

Parameters
[in]levelNumber of level over the max depth reached in the tree at which families of octants are fixed compact on the same process (level=0 is classic LoadBalance).
[in]weightPointer to a vector of weights of the local octants (weight=NULL is uniform distribution).

Definition at line 4184 of file ParaTree.cpp.

◆ preadapt()

void bitpit::ParaTree::preadapt ( )

Pre-adapt the octree mesh with user setup for markers and 2:1 balancing conditions.

The user can call pre-adapt and then adapt or only adapt, however after the pre-adapt function has been called the adapt function is mandatory.

Examples
PABLO_example_00010.cpp.

Definition at line 3650 of file ParaTree.cpp.

◆ printHeader()

void bitpit::ParaTree::printHeader ( )

Print the initial PABLO header.

Definition at line 755 of file ParaTree.cpp.

◆ privateLoadBalance()

template<class Impl >
void bitpit::ParaTree::privateLoadBalance ( const uint32_t * partition,
DataLBInterface< Impl > * userData = nullptr )
inline

Distribute Load-Balancing octants and user data of the whole tree over the processes of the job following a given partition distribution. Until loadBalance is not called for the first time the mesh is serial.

Parameters
[in]partitionTarget distribution of octants over processes.
[in,out]userDataUser data that will be distributed among the processes.

Definition at line 763 of file ParaTree.hpp.

◆ replaceComm()

void bitpit::ParaTree::replaceComm ( MPI_Comm communicator,
MPI_Comm * previousCommunicator = nullptr )

Set the MPI communicator to be used for parallel communications. If the communicator is already set, it will be replaced with the new one only if the two communicators are equivalent, i.e., the rank of the processes have to be the same in both communicators. The tree will not use the specified communicator directly, instead a duplicates is ceated. Previous communicator will be freed or returned depending on the received arguments.

Parameters
communicatoris the communicator to be used for parallel communications.
[in,out]previousCommunicatorif a valid pointer is received, on output this parameter will contain the previous communicator. If this patameter is a null pointer, the previous communicator will be freed.

Definition at line 878 of file ParaTree.cpp.

◆ reset()

void bitpit::ParaTree::reset ( )
virtual

Reset the octree

Reimplemented in bitpit::PabloUniform.

Definition at line 448 of file ParaTree.cpp.

◆ restore()

void bitpit::ParaTree::restore ( std::istream & stream)
virtual

Restore the octree from the specified stream.

Parameters
streamStream to read from.

Reimplemented in bitpit::PabloUniform.

Definition at line 588 of file ParaTree.cpp.

◆ setBalance() [1/2]

void bitpit::ParaTree::setBalance ( Octant * oct,
bool balance )

Set the balancing condition of an octant.

Parameters
[in]octPointer to the target octant
[in]balanceHas octant to be 2:1 balanced in adapting procedure?

Definition at line 2094 of file ParaTree.cpp.

◆ setBalance() [2/2]

void bitpit::ParaTree::setBalance ( uint32_t idx,
bool balance )

Set the balancing condition of an octant.

Parameters
[in]idxLocal index of target octant.
[in]balanceHas octant to be 2:1 balanced in adapting procedure?
Examples
PABLO_example_00001.cpp, PABLO_example_00002.cpp, PABLO_example_00004.cpp, PABLO_example_00005.cpp, PABLO_example_00007.cpp, PABLO_example_00008.cpp, PABLO_example_00009.cpp, and PABLO_example_00010.cpp.

Definition at line 1709 of file ParaTree.cpp.

◆ setBalanceCodimension()

void bitpit::ParaTree::setBalanceCodimension ( uint8_t b21codim)

Set the codimension for 2:1 balancing

Parameters
[in]b21codimMaximum codimension of the entity through which the 2:1 balance is performed (1 = 2:1 balance through edges (default); 2 = 2:1 balance through nodes and edges).
Examples
PABLO_example_00001.cpp.

Definition at line 2237 of file ParaTree.cpp.

◆ setComm()

void bitpit::ParaTree::setComm ( MPI_Comm communicator)

Set the MPI communicator to be used for parallel communications. The tree will not use the specified communicator directly, instead a duplicates is ceated.

Parameters
communicatoris the communicator to be used for parallel communications.

Definition at line 848 of file ParaTree.cpp.

◆ setMarker() [1/2]

void bitpit::ParaTree::setMarker ( Octant * oct,
int8_t marker )

Set the refinement marker of an octant.

Parameters
[in]octPointer to the target octant
[in]markerRefinement marker of octant (n=n refinement in adapt, -n=n coarsening in adapt, default=0).

Definition at line 2081 of file ParaTree.cpp.

◆ setMarker() [2/2]

void bitpit::ParaTree::setMarker ( uint32_t idx,
int8_t marker )

Set the refinement marker of an octant.

Parameters
[in]idxLocal index of target octant.
[in]markerRefinement marker of octant (n=n refinement in adapt, -n=n coarsening in adapt, default=0).
Examples
PABLO_example_00001.cpp, PABLO_example_00002.cpp, PABLO_example_00004.cpp, PABLO_example_00005.cpp, PABLO_example_00007.cpp, PABLO_example_00008.cpp, PABLO_example_00009.cpp, and PABLO_example_00010.cpp.

Definition at line 1696 of file ParaTree.cpp.

◆ setNofGhostLayers()

void bitpit::ParaTree::setNofGhostLayers ( std::size_t nofGhostLayers)

Set the number of ghosts layers

Parameters
[in]nofGhostLayersThe number of ghost layers in the ghost halo

Definition at line 2542 of file ParaTree.cpp.

◆ setPeriodic()

void bitpit::ParaTree::setPeriodic ( uint8_t i)

Set the periodic condition of a target boundary (implicitly set the periodic face).

Parameters
[in]iIndex of the target boundary face (even the opp face will be set).
Examples
PABLO_example_00011.cpp.

Definition at line 1252 of file ParaTree.cpp.

◆ settleMarkers()

void bitpit::ParaTree::settleMarkers ( )

Rearrange the octree markers with user setup for markers and 2:1 balancing conditions.

This function is analogous to pre-adapt method, but no function is mandatory after a settleMarkers call.

Note: the last operation tag is not changed after a settleMarkers call.

Definition at line 3631 of file ParaTree.cpp.

◆ setTol()

void bitpit::ParaTree::setTol ( double tol = 1.0e-14)

Set the tolerance used in geometric operations.

Parameters
[in]tolDesired tolerance.

Definition at line 1262 of file ParaTree.cpp.

◆ updateConnectivity()

void bitpit::ParaTree::updateConnectivity ( )

◆ write()

void bitpit::ParaTree::write ( const std::string & filename)

Write the physical octree mesh in .vtu format in a user-defined file. If the connectivity is not stored, the method temporary computes it. If the connectivity of ghost octants is already computed, the method writes the ghosts on file.

Parameters
[in]filenameName of output file (PABLO will add the total number of processes p000# and the current rank s000#).

Definition at line 5856 of file ParaTree.cpp.

◆ writeTest()

void bitpit::ParaTree::writeTest ( const std::string & filename,
dvector data )

Write the physical octree mesh in .vtu format with data for test in a user-defined file. If the connectivity is not stored, the method temporary computes it. The method doesn't write the ghosts on file.

Parameters
[in]filenameName of output file (PABLO will add the total number of processes p000# and the current rank s000#).
[in]dataVector of double with user data.

Definition at line 6011 of file ParaTree.cpp.

Member Data Documentation

◆ DEFAULT_LOG_FILE

const std::string bitpit::ParaTree::DEFAULT_LOG_FILE = "PABLO"
static

Default name of logger file.

Definition at line 119 of file ParaTree.hpp.


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