t8  UNKNOWN
t8code is a C library to manage a forest of adaptive space-trees of general element classes in parallel.
Typedefs | Enumerations | Functions
t8_forest_general.h File Reference

We define the forest of trees in this file. More...

#include <t8_cmesh.h>
#include <t8_element.h>
#include <t8_data/t8_containers.h>

Go to the source code of this file.

Typedefs

typedef struct t8_forestt8_forest_t
 Opaque pointer to a forest implementation.
 
typedef struct t8_treet8_tree_t
 
typedef void(* t8_generic_function_pointer) (void)
 This typedef is needed as a helper construct to properly be able to define a function that returns a pointer to a void fun(void) function. More...
 
typedef void(* t8_forest_replace_t) (t8_forest_t forest_old, t8_forest_t forest_new, t8_locidx_t which_tree, t8_eclass_scheme_c *ts, const int refine, const int num_outgoing, const t8_locidx_t first_outgoing, const int num_incoming, const t8_locidx_t first_incoming)
 Callback function prototype to replace one set of elements with another. More...
 
typedef int(* t8_forest_adapt_t) (t8_forest_t forest, t8_forest_t forest_from, t8_locidx_t which_tree, t8_locidx_t lelement_id, t8_eclass_scheme_c *ts, const int is_family, const int num_elements, t8_element_t *elements[])
 Callback function prototype to decide for refining and coarsening. More...
 

Enumerations

enum  t8_ghost_type_t {
  T8_GHOST_NONE = 0 ,
  T8_GHOST_FACES ,
  T8_GHOST_EDGES ,
  T8_GHOST_VERTICES
}
 This type controls, which neighbors count as ghost elements. More...
 

Functions

void t8_forest_init (t8_forest_t *pforest)
 Create a new forest with reference count one. More...
 
int t8_forest_is_initialized (t8_forest_t forest)
 Check whether a forest is not NULL, initialized and not committed. More...
 
int t8_forest_is_committed (t8_forest_t forest)
 Check whether a forest is not NULL, initialized and committed. More...
 
int t8_forest_no_overlap (t8_forest_t forest)
 Check whether the forest has local overlapping elements. More...
 
int t8_forest_is_equal (t8_forest_t forest_a, t8_forest_t forest_b)
 Check whether two committed forests have the same local elements. More...
 
void t8_forest_set_cmesh (t8_forest_t forest, t8_cmesh_t cmesh, sc_MPI_Comm comm)
 Set the cmesh associated to a forest. More...
 
void t8_forest_set_scheme (t8_forest_t forest, t8_scheme_cxx_t *scheme)
 Set the element scheme associated to a forest. More...
 
void t8_forest_set_level (t8_forest_t forest, int level)
 Set the initial refinement level to be used when forest is committed. More...
 
void t8_forest_set_copy (t8_forest_t forest, const t8_forest_t from)
 Set a forest as source for copying on committing. More...
 
void t8_forest_set_adapt (t8_forest_t forest, const t8_forest_t set_from, t8_forest_adapt_t adapt_fn, int recursive)
 Set a source forest with an adapt function to be adapted on committing. More...
 
void t8_forest_set_user_data (t8_forest_t forest, void *data)
 Set the user data of a forest. More...
 
void * t8_forest_get_user_data (const t8_forest_t forest)
 Return the user data pointer associated with a forest. More...
 
void t8_forest_set_user_function (t8_forest_t forest, t8_generic_function_pointer function)
 Set the user function pointer of a forest. More...
 
t8_generic_function_pointer t8_forest_get_user_function (const t8_forest_t forest)
 Return the user function pointer associated with a forest. More...
 
void t8_forest_set_partition (t8_forest_t forest, const t8_forest_t set_from, int set_for_coarsening)
 Set a source forest to be partitioned during commit. More...
 
void t8_forest_set_balance (t8_forest_t forest, const t8_forest_t set_from, int no_repartition)
 Set a source forest to be balanced during commit. More...
 
void t8_forest_set_ghost (t8_forest_t forest, int do_ghost, t8_ghost_type_t ghost_type)
 Enable or disable the creation of a layer of ghost elements. More...
 
void t8_forest_set_ghost_ext (t8_forest_t forest, int do_ghost, t8_ghost_type_t ghost_type, int ghost_version)
 Like t8_forest_set_ghost but with the additional options to change the ghost algorithm. More...
 
void t8_forest_set_load (t8_forest_t forest, const char *filename)
 
void t8_forest_comm_global_num_elements (t8_forest_t forest)
 Compute the global number of elements in a forest as the sum of the local element counts. More...
 
void t8_forest_commit (t8_forest_t forest)
 After allocating and adding properties to a forest, commit the changes. More...
 
int t8_forest_get_maxlevel (const t8_forest_t forest)
 Return the maximum allowed refinement level for any element in a forest. More...
 
t8_locidx_t t8_forest_get_local_num_elements (const t8_forest_t forest)
 Return the number of process local elements in the forest. More...
 
t8_gloidx_t t8_forest_get_global_num_elements (const t8_forest_t forest)
 Return the number of global elements in the forest. More...
 
t8_locidx_t t8_forest_get_num_ghosts (const t8_forest_t forest)
 Return the number of ghost elements of a forest. More...
 
t8_eclass_t t8_forest_get_eclass (const t8_forest_t forest, const t8_locidx_t ltreeid)
 Return the element class of a forest local tree. More...
 
t8_locidx_t t8_forest_get_local_id (const t8_forest_t forest, const t8_gloidx_t gtreeid)
 Given a global tree id compute the forest local id of this tree. More...
 
t8_locidx_t t8_forest_get_local_or_ghost_id (const t8_forest_t forest, const t8_gloidx_t gtreeid)
 Given a global tree id compute the forest local id of this tree. More...
 
t8_locidx_t t8_forest_ltreeid_to_cmesh_ltreeid (t8_forest_t forest, t8_locidx_t ltreeid)
 Given the local id of a tree in a forest, compute the tree's local id in the associated cmesh. More...
 
t8_locidx_t t8_forest_cmesh_ltreeid_to_ltreeid (t8_forest_t forest, t8_locidx_t lctreeid)
 Given the local id of a tree in the coarse mesh of a forest, compute the tree's local id in the forest. More...
 
t8_ctree_t t8_forest_get_coarse_tree (t8_forest_t forest, t8_locidx_t ltreeid)
 Given the local id of a tree in a forest, return the coarse tree of the cmesh that corresponds to this tree. More...
 
void t8_forest_leaf_face_neighbors (t8_forest_t forest, t8_locidx_t ltreeid, const t8_element_t *leaf, t8_element_t **pneighbor_leaves[], int face, int *dual_faces[], int *num_neighbors, t8_locidx_t **pelement_indices, t8_eclass_scheme_c **pneigh_scheme, int forest_is_balanced)
 Compute the leaf face neighbors of a forest. More...
 
void t8_forest_leaf_face_neighbors_ext (t8_forest_t forest, t8_locidx_t ltreeid, const t8_element_t *leaf, t8_element_t **pneighbor_leaves[], int face, int *dual_faces[], int *num_neighbors, t8_locidx_t **pelement_indices, t8_eclass_scheme_c **pneigh_scheme, int forest_is_balanced, t8_gloidx_t *gneigh_tree)
 
void t8_forest_ghost_exchange_data (t8_forest_t forest, sc_array_t *element_data)
 Exchange ghost information of user defined element data. More...
 
void t8_forest_ghost_print (t8_forest_t forest)
 Print the ghost structure of a forest. More...
 
void t8_forest_partition_cmesh (t8_forest_t forest, sc_MPI_Comm comm, int set_profiling)
 Change the cmesh associated to a forest to a partitioned cmesh that is partitioned according to the tree distribution in the forest. More...
 
sc_MPI_Comm t8_forest_get_mpicomm (const t8_forest_t forest)
 Return the mpi communicator associated to a forest. More...
 
t8_gloidx_t t8_forest_get_first_local_tree_id (const t8_forest_t forest)
 Return the global id of the first local tree of a forest. More...
 
t8_locidx_t t8_forest_get_num_local_trees (const t8_forest_t forest)
 Return the number of local trees of a given forest. More...
 
t8_locidx_t t8_forest_get_num_ghost_trees (const t8_forest_t forest)
 Return the number of ghost trees of a given forest. More...
 
t8_gloidx_t t8_forest_get_num_global_trees (const t8_forest_t forest)
 Return the number of global trees of a given forest. More...
 
t8_gloidx_t t8_forest_global_tree_id (const t8_forest_t forest, const t8_locidx_t ltreeid)
 Return the global id of a local tree or a ghost tree. More...
 
t8_tree_t t8_forest_get_tree (const t8_forest_t forest, const t8_locidx_t ltree_id)
 Return a pointer to a tree in a forest. More...
 
double * t8_forest_get_tree_vertices (t8_forest_t forest, t8_locidx_t ltreeid)
 Return a pointer to the vertex coordinates of a tree. More...
 
t8_element_array_tt8_forest_tree_get_leaves (const t8_forest_t forest, const t8_locidx_t ltree_id)
 Return the array of leaf elements of a local tree in a forest. More...
 
t8_cmesh_t t8_forest_get_cmesh (t8_forest_t forest)
 Return a cmesh associated to a forest. More...
 
t8_element_tt8_forest_get_element (t8_forest_t forest, t8_locidx_t lelement_id, t8_locidx_t *ltreeid)
 Return an element of the forest. More...
 
const t8_element_tt8_forest_get_element_in_tree (t8_forest_t forest, t8_locidx_t ltreeid, t8_locidx_t leid_in_tree)
 Return an element of a local tree in a forest. More...
 
t8_locidx_t t8_forest_get_tree_num_elements (t8_forest_t forest, t8_locidx_t ltreeid)
 Return the number of elements of a tree. More...
 
t8_locidx_t t8_forest_get_tree_element_offset (const t8_forest_t forest, const t8_locidx_t ltreeid)
 Return the element offset of a local tree, that is the number of elements in all trees with smaller local treeid. More...
 
t8_locidx_t t8_forest_get_tree_element_count (t8_tree_t tree)
 Return the number of elements of a tree. More...
 
t8_eclass_t t8_forest_get_tree_class (const t8_forest_t forest, const t8_locidx_t ltreeid)
 Return the eclass of a tree in a forest. More...
 
t8_gloidx_t t8_forest_get_first_local_element_id (t8_forest_t forest)
 Compute the global index of the first local element of a forest. More...
 
t8_scheme_cxx_tt8_forest_get_scheme (const t8_forest_t forest)
 Return the element scheme associated to a forest. More...
 
t8_eclass_scheme_ct8_forest_get_eclass_scheme (t8_forest_t forest, t8_eclass_t eclass)
 Return the eclass scheme of a given element class associated to a forest. More...
 
t8_eclass_t t8_forest_element_neighbor_eclass (t8_forest_t forest, t8_locidx_t ltreeid, const t8_element_t *elem, int face)
 Return the eclass of the tree in which a face neighbor of a given element lies. More...
 
t8_gloidx_t t8_forest_element_face_neighbor (t8_forest_t forest, t8_locidx_t ltreeid, const t8_element_t *elem, t8_element_t *neigh, t8_eclass_scheme_c *neigh_scheme, int face, int *neigh_face)
 Construct the face neighbor of an element, possibly across tree boundaries. More...
 
void t8_forest_iterate (t8_forest_t forest)
 
void t8_forest_element_points_inside (t8_forest_t forest, t8_locidx_t ltreeid, const t8_element_t *element, const double *points, int num_points, int *is_inside, const double tolerance)
 Query whether a batch of points lies inside an element. More...
 
t8_forest_t t8_forest_new_uniform (t8_cmesh_t cmesh, t8_scheme_cxx_t *scheme, const int level, const int do_face_ghost, sc_MPI_Comm comm)
 Build a uniformly refined forest on a coarse mesh. More...
 
t8_forest_t t8_forest_new_adapt (t8_forest_t forest_from, t8_forest_adapt_t adapt_fn, int recursive, int do_face_ghost, void *user_data)
 Build a adapted forest from another forest. More...
 
void t8_forest_ref (t8_forest_t forest)
 Increase the reference counter of a forest. More...
 
void t8_forest_unref (t8_forest_t *pforest)
 Decrease the reference counter of a forest. More...
 

Detailed Description

We define the forest of trees in this file.

Typedef Documentation

◆ t8_forest_adapt_t

typedef int(* t8_forest_adapt_t) (t8_forest_t forest, t8_forest_t forest_from, t8_locidx_t which_tree, t8_locidx_t lelement_id, t8_eclass_scheme_c *ts, const int is_family, const int num_elements, t8_element_t *elements[])

Callback function prototype to decide for refining and coarsening.

If is_family equals 1, the first num_elements in elements form a family and we decide whether this family should be coarsened or only the first element should be refined. Otherwise is_family must equal zero and we consider the first entry of the element array for refinement. Entries of the element array beyond the first num_elements are undefined.

Parameters
[in]forestthe forest to which the new elements belong
[in]forest_fromthe forest that is adapted.
[in]which_treethe local tree containing elements
[in]lelement_idthe local element id in forest_old in the tree of the current element
[in]tsthe eclass scheme of the tree
[in]is_familyif 1, the first num_elements entries in elements form a family. If 0, they do not.
[in]num_elementsthe number of entries in elements that are defined
[in]elementsPointers to a family or, if is_family is zero, pointer to one element.
Returns
1 if the first entry in elements should be refined, -1 if the family elements shall be coarsened, -2 if the first entry in elements should be removed, 0 else.

◆ t8_forest_replace_t

typedef void(* t8_forest_replace_t) (t8_forest_t forest_old, t8_forest_t forest_new, t8_locidx_t which_tree, t8_eclass_scheme_c *ts, const int refine, const int num_outgoing, const t8_locidx_t first_outgoing, const int num_incoming, const t8_locidx_t first_incoming)

Callback function prototype to replace one set of elements with another.

This is used by the replace routine which can be called after adapt, when the elements of an existing, valid forest are changed. The callback allows the user to make changes to the elements of the new forest that are either refined, coarsened or the same as elements in the old forest.

Parameters
[in]forest_oldThe forest that is adapted
[in]forest_newThe forest that is newly constructed from forest_old
[in]which_treeThe local tree containing first_outgoing and first_incoming
[in]tsThe eclass scheme of the tree
[in]refine-1 if family in forest_old got coarsened, 0 if element has not been touched, 1 if element got refined and -2 if element got removed. See return of t8_forest_adapt_t.
[in]num_outgoingThe number of outgoing elements.
[in]first_outgoingThe tree local index of the first outgoing element. 0 <= first_outgoing < which_tree->num_elements
[in]num_incomingThe number of incoming elements.
[in]first_incomingThe tree local index of the first incoming element. 0 <= first_incom < new_which_tree->num_elements

If an element is being refined, refine and num_outgoing will be 1 and num_incoming will be the number of children. If a family is being coarsened, refine will be -1, num_outgoing will be the number of family members and num_incoming will be 1. If an element is being removed, refine and num_outgoing will be 1 and num_incoming will be 0. Else refine will be 0 and num_outgoing and num_incoming will both be 1.

See also
t8_forest_iterate_replace

◆ t8_generic_function_pointer

typedef void(* t8_generic_function_pointer) (void)

This typedef is needed as a helper construct to properly be able to define a function that returns a pointer to a void fun(void) function.

See also
t8_forest_get_user_function.

Enumeration Type Documentation

◆ t8_ghost_type_t

This type controls, which neighbors count as ghost elements.

Currently, we support face-neighbors. Vertex and edge neighbors will eventually be added.

Enumerator
T8_GHOST_NONE 

Do not create ghost layer.

T8_GHOST_FACES 

Consider all face (codimension 1) neighbors.

T8_GHOST_EDGES 

Consider all edge (codimension 2) and face neighbors.

T8_GHOST_VERTICES 

Consider all vertex (codimension 3) and edge and face neighbors.

Function Documentation

◆ t8_forest_cmesh_ltreeid_to_ltreeid()

t8_locidx_t t8_forest_cmesh_ltreeid_to_ltreeid ( t8_forest_t  forest,
t8_locidx_t  lctreeid 
)

Given the local id of a tree in the coarse mesh of a forest, compute the tree's local id in the forest.

Parameters
[in]forestThe forest.
[in]ltreeidThe local id of a tree in the coarse mesh of forest.
Returns
The local id of the tree in the forest. -1 if the tree is not forest local. forest must be committed before calling this function.
Note
For forest local trees, this is the inverse function of t8_forest_ltreeid_to_cmesh_ltreeid.
See also
https://github.com/DLR-AMR/t8code/wiki/Tree-indexing for more details about tree indexing.

◆ t8_forest_comm_global_num_elements()

void t8_forest_comm_global_num_elements ( t8_forest_t  forest)

Compute the global number of elements in a forest as the sum of the local element counts.

Parameters
[in]forestThe forest.

◆ t8_forest_commit()

void t8_forest_commit ( t8_forest_t  forest)

After allocating and adding properties to a forest, commit the changes.

This call sets up the internal state of the forest.

Parameters
[in,out]forestMust be created with t8_forest_init and specialized with t8_forest_set_* calls first.

◆ t8_forest_element_face_neighbor()

t8_gloidx_t t8_forest_element_face_neighbor ( t8_forest_t  forest,
t8_locidx_t  ltreeid,
const t8_element_t elem,
t8_element_t neigh,
t8_eclass_scheme_c neigh_scheme,
int  face,
int *  neigh_face 
)

Construct the face neighbor of an element, possibly across tree boundaries.

Returns the global tree-id of the tree in which the neighbor element lies in.

Parameters
[in]elemThe element to be considered.
[in,out]neighOn input an allocated element of the scheme of the face_neighbors eclass. On output, this element's data is filled with the data of the face neighbor. If the neighbor does not exist the data could be modified arbitrarily.
[in]neigh_schemeThe eclass scheme of neigh.
[in]faceThe number of the face along which the neighbor should be constructed.
[out]neigh_faceThe number of the face viewed from perspective of neigh.
Returns
The global tree-id of the tree in which neigh is in. -1 if there exists no neighbor across that face.

◆ t8_forest_element_neighbor_eclass()

t8_eclass_t t8_forest_element_neighbor_eclass ( t8_forest_t  forest,
t8_locidx_t  ltreeid,
const t8_element_t elem,
int  face 
)

Return the eclass of the tree in which a face neighbor of a given element lies.

Parameters
[in]forest.A committed forest.
[in]ltreeid.The local tree in which the element lies.
[in]elem.An element in the tree ltreeid.
[in]face.A face number of elem.
Returns
The local tree id of the tree in which the face neighbor of elem across face lies.

◆ t8_forest_element_points_inside()

void t8_forest_element_points_inside ( t8_forest_t  forest,
t8_locidx_t  ltreeid,
const t8_element_t element,
const double *  points,
int  num_points,
int *  is_inside,
const double  tolerance 
)

Query whether a batch of points lies inside an element.

For bilinearly interpolated elements.

Note
For 2D quadrilateral elements this function is only an approximation. It is correct if the four vertices lie in the same plane, but it may produce only approximate results if the vertices do not lie in the same plane.
Parameters
[in]forestThe forest.
[in]ltree_idThe forest local id of the tree in which the element is.
[in]elementThe element.
[in]points3-dimensional coordinates of the points to check
[in]num_pointsThe number of points to check
[in,out]is_insideAn array of length num_points, filled with 0/1 on output. True (non-zero) if a point lies within an element, false otherwise. The return value is also true if the point lies on the element boundary. Thus, this function may return true for different leaf elements, if they are neighbors and the point lies on the common boundary.
[in]toleranceTolerance that we allow the point to not exactly match the element. If this value is larger we detect more points. If it is zero we probably do not detect points even if they are inside due to rounding errors.

◆ t8_forest_get_cmesh()

t8_cmesh_t t8_forest_get_cmesh ( t8_forest_t  forest)

Return a cmesh associated to a forest.

Parameters
[in]forestThe forest.
Returns
The cmesh associated to the forest.

◆ t8_forest_get_coarse_tree()

t8_ctree_t t8_forest_get_coarse_tree ( t8_forest_t  forest,
t8_locidx_t  ltreeid 
)

Given the local id of a tree in a forest, return the coarse tree of the cmesh that corresponds to this tree.

Parameters
[in]forestThe forest.
[in]ltreeidThe local id of a tree in the forest.
Returns
The coarse tree that matches the forest tree with local id ltreeid.

◆ t8_forest_get_eclass()

t8_eclass_t t8_forest_get_eclass ( const t8_forest_t  forest,
const t8_locidx_t  ltreeid 
)

Return the element class of a forest local tree.

Parameters
[in]forestThe forest.
[in]ltreeidThe local id of a tree in forest.
Returns
The element class of the tree ltreeid. forest must be committed before calling this function.

◆ t8_forest_get_eclass_scheme()

t8_eclass_scheme_c* t8_forest_get_eclass_scheme ( t8_forest_t  forest,
t8_eclass_t  eclass 
)

Return the eclass scheme of a given element class associated to a forest.

Parameters
[in]forest.A committed forest.
[in]eclass.An element class.
Returns
The eclass scheme of eclass associated to forest.
See also
t8_forest_set_scheme
Note
The forest is not required to have trees of class eclass.

◆ t8_forest_get_element()

t8_element_t* t8_forest_get_element ( t8_forest_t  forest,
t8_locidx_t  lelement_id,
t8_locidx_t ltreeid 
)

Return an element of the forest.

Parameters
[in]forestThe forest.
[in]lelement_idThe local id of an element in forest.
[out]ltreeidIf not NULL, on output the local tree id of the tree in which the element lies in.
Returns
A pointer to the element. NULL if this element does not exist.
Note
This function performs a binary search. For constant access, use t8_forest_get_element_in_tree forest must be committed before calling this function.

◆ t8_forest_get_element_in_tree()

const t8_element_t* t8_forest_get_element_in_tree ( t8_forest_t  forest,
t8_locidx_t  ltreeid,
t8_locidx_t  leid_in_tree 
)

Return an element of a local tree in a forest.

Parameters
[in]forestThe forest.
[in]ltreeidAn id of a local tree in the forest.
[in]leid_in_treeThe index of an element in the tree.
Returns
A pointer to the element.
Note
If the tree id is know, this function should be preferred over t8_forest_get_element. forest must be committed before calling this function.

◆ t8_forest_get_first_local_element_id()

t8_gloidx_t t8_forest_get_first_local_element_id ( t8_forest_t  forest)

Compute the global index of the first local element of a forest.

This function is collective.

Parameters
[in]forestA committed forest, whose first element's index is computed.
Returns
The global index of forest's first local element. Forest must be committed when calling this function. This function is collective and must be called on each process.

◆ t8_forest_get_first_local_tree_id()

t8_gloidx_t t8_forest_get_first_local_tree_id ( const t8_forest_t  forest)

Return the global id of the first local tree of a forest.

Parameters
[in]forestThe forest.
Returns
The global id of the first local tree in forest.

◆ t8_forest_get_global_num_elements()

t8_gloidx_t t8_forest_get_global_num_elements ( const t8_forest_t  forest)

Return the number of global elements in the forest.

Parameters
[in]forestA forest.
Returns
The number of elements (summed over all processes) in forest. forest must be committed before calling this function.

◆ t8_forest_get_local_id()

t8_locidx_t t8_forest_get_local_id ( const t8_forest_t  forest,
const t8_gloidx_t  gtreeid 
)

Given a global tree id compute the forest local id of this tree.

If the tree is a local tree, then the local id is between 0 and the number of local trees. If the tree is not a local tree, a negative number is returned.

Parameters
[in]forestThe forest.
[in]gtreeidThe global id of a tree.
Returns
The tree's local id in forest, if it is a local tree. A negative number if not.
See also
https://github.com/DLR-AMR/t8code/wiki/Tree-indexing for more details about tree indexing.

◆ t8_forest_get_local_num_elements()

t8_locidx_t t8_forest_get_local_num_elements ( const t8_forest_t  forest)

Return the number of process local elements in the forest.

Parameters
[in]forestA forest.
Returns
The number of elements on this process in forest. forest must be committed before calling this function.

◆ t8_forest_get_local_or_ghost_id()

t8_locidx_t t8_forest_get_local_or_ghost_id ( const t8_forest_t  forest,
const t8_gloidx_t  gtreeid 
)

Given a global tree id compute the forest local id of this tree.

If the tree is a local tree, then the local id is between 0 and the number of local trees. If the tree is a ghost, then the local id is between num_local_trees and num_local_trees + num_ghost_trees. If the tree is neither a local tree nor a ghost tree, a negative number is returned.

Parameters
[in]forestThe forest.
[in]gtreeidThe global id of a tree.
Returns
The tree's local id in forest, if it is a local tree. num_local_trees + the ghosts id, if it is a ghost tree. A negative number if not.
See also
https://github.com/DLR-AMR/t8code/wiki/Tree-indexing for more details about tree indexing.

◆ t8_forest_get_maxlevel()

int t8_forest_get_maxlevel ( const t8_forest_t  forest)

Return the maximum allowed refinement level for any element in a forest.

Parameters
[in]forestA forest.
Returns
The maximum level of refinement that is allowed for an element in this forest. It is guaranteed that any tree in forest can be refined this many times and it is not allowed to refine further. forest must be committed before calling this function. For forest with a single element class (non-hybrid) maxlevel is the maximum refinement level of this element class, whilst for hybrid forests the maxlevel is the minimum of all maxlevels of the element classes in this forest.

◆ t8_forest_get_mpicomm()

sc_MPI_Comm t8_forest_get_mpicomm ( const t8_forest_t  forest)

Return the mpi communicator associated to a forest.

Parameters
[in]forestThe forest.
Returns
The mpi communicator of forest. forest must be committed before calling this function.

◆ t8_forest_get_num_ghost_trees()

t8_locidx_t t8_forest_get_num_ghost_trees ( const t8_forest_t  forest)

Return the number of ghost trees of a given forest.

Parameters
[in]forestThe forest.
Returns
The number of ghost trees of that forest.

◆ t8_forest_get_num_ghosts()

t8_locidx_t t8_forest_get_num_ghosts ( const t8_forest_t  forest)

Return the number of ghost elements of a forest.

Parameters
[in]forestThe forest.
Returns
The number of ghost elements stored in the ghost structure of forest. 0 if no ghosts were constructed.
See also
t8_forest_set_ghost forest must be committed before calling this function.

◆ t8_forest_get_num_global_trees()

t8_gloidx_t t8_forest_get_num_global_trees ( const t8_forest_t  forest)

Return the number of global trees of a given forest.

Parameters
[in]forestThe forest.
Returns
The number of global trees of that forest.

◆ t8_forest_get_num_local_trees()

t8_locidx_t t8_forest_get_num_local_trees ( const t8_forest_t  forest)

Return the number of local trees of a given forest.

Parameters
[in]forestThe forest.
Returns
The number of local trees of that forest.

◆ t8_forest_get_scheme()

t8_scheme_cxx_t* t8_forest_get_scheme ( const t8_forest_t  forest)

Return the element scheme associated to a forest.

Parameters
[in]forest.A committed forest.
Returns
The element scheme of the forest.
See also
t8_forest_set_scheme

◆ t8_forest_get_tree()

t8_tree_t t8_forest_get_tree ( const t8_forest_t  forest,
const t8_locidx_t  ltree_id 
)

Return a pointer to a tree in a forest.

Parameters
[in]forestThe forest.
[in]ltree_idThe local id of the tree.
Returns
A pointer to the tree with local id ltree_id. forest must be committed before calling this function.

◆ t8_forest_get_tree_class()

t8_eclass_t t8_forest_get_tree_class ( const t8_forest_t  forest,
const t8_locidx_t  ltreeid 
)

Return the eclass of a tree in a forest.

Parameters
[in]forestThe forest.
[in]ltreeidThe local id of a tree (local or ghost) in forest.
Returns
The element class of the tree with local id ltreeid.

◆ t8_forest_get_tree_element_count()

t8_locidx_t t8_forest_get_tree_element_count ( t8_tree_t  tree)

Return the number of elements of a tree.

Parameters
[in]treeA tree in a forest.
Returns
The number of elements of that tree.

◆ t8_forest_get_tree_element_offset()

t8_locidx_t t8_forest_get_tree_element_offset ( const t8_forest_t  forest,
const t8_locidx_t  ltreeid 
)

Return the element offset of a local tree, that is the number of elements in all trees with smaller local treeid.

Parameters
[in]forestThe forest.
[in]ltreeidA local id of a tree.
Returns
The number of leaf elements on all local tree with id < ltreeid.
Note
forest must be committed before calling this function.

◆ t8_forest_get_tree_num_elements()

t8_locidx_t t8_forest_get_tree_num_elements ( t8_forest_t  forest,
t8_locidx_t  ltreeid 
)

Return the number of elements of a tree.

Parameters
[in]forestThe forest.
[in]ltreeidA local id of a tree.
Returns
The number of elements in the local tree ltreeid.

◆ t8_forest_get_tree_vertices()

double* t8_forest_get_tree_vertices ( t8_forest_t  forest,
t8_locidx_t  ltreeid 
)

Return a pointer to the vertex coordinates of a tree.

Parameters
[in]forestThe forest.
[in]ltreeidThe id of a local tree.
Returns
If stored, a pointer to the vertex coordinates of tree. If no coordinates for this tree are found, NULL.

◆ t8_forest_get_user_data()

void* t8_forest_get_user_data ( const t8_forest_t  forest)

Return the user data pointer associated with a forest.

Parameters
[in]forestThe forest.
Returns
The user data pointer of forest. The forest does not need be committed before calling this function.
See also
t8_forest_set_user_data

◆ t8_forest_get_user_function()

t8_generic_function_pointer t8_forest_get_user_function ( const t8_forest_t  forest)

Return the user function pointer associated with a forest.

Parameters
[in]forestThe forest.
Returns
The user function pointer of forest. The forest does not need be committed before calling this function.
See also
t8_forest_set_user_function

◆ t8_forest_ghost_exchange_data()

void t8_forest_ghost_exchange_data ( t8_forest_t  forest,
sc_array_t *  element_data 
)

Exchange ghost information of user defined element data.

Parameters
[in]forestThe forest. Must be committed.
[in]element_dataAn array of length num_local_elements + num_ghosts storing one value for each local element and ghost in forest. After calling this function the entries for the ghost elements are update with the entries in the element_data array of the corresponding owning process.
Note
This function is collective and hence must be called by all processes in the forest's MPI Communicator.

◆ t8_forest_ghost_print()

void t8_forest_ghost_print ( t8_forest_t  forest)

Print the ghost structure of a forest.

Only used for debugging.

◆ t8_forest_global_tree_id()

t8_gloidx_t t8_forest_global_tree_id ( const t8_forest_t  forest,
const t8_locidx_t  ltreeid 
)

Return the global id of a local tree or a ghost tree.

Parameters
[in]forestThe forest.
[in]ltreeidAn id 0 <= ltreeid < num_local_trees + num_ghosts specifying a local tree or ghost tree.
Returns
The global id corresponding to the tree with local id ltreeid. forest must be committed before calling this function.
See also
https://github.com/DLR-AMR/t8code/wiki/Tree-indexing for more details about tree indexing.

◆ t8_forest_init()

void t8_forest_init ( t8_forest_t pforest)

Create a new forest with reference count one.

This forest needs to be specialized with the t8_forest_set_* calls. Currently it is manatory to either call the functions t8_forest_set_mpicomm, t8_forest_set_cmesh, and t8_forest_set_scheme, or to call one of t8_forest_set_copy, t8_forest_set_adapt, or t8_forest_set_partition. It is illegal to mix these calls, or to call more than one of the three latter functions Then it needs to be set up with t8_forest_commit.

Parameters
[in,out]pforestOn input, this pointer must be non-NULL. On return, this pointer set to the new forest.

◆ t8_forest_is_committed()

int t8_forest_is_committed ( t8_forest_t  forest)

Check whether a forest is not NULL, initialized and committed.

In addition, it asserts that the forest is consistent as much as possible.

Parameters
[in]forestThis forest is examined. May be NULL.
Returns
True if forest is not NULL and t8_forest_init has been called on it as well as t8_forest_commit. False otherwise.

◆ t8_forest_is_equal()

int t8_forest_is_equal ( t8_forest_t  forest_a,
t8_forest_t  forest_b 
)

Check whether two committed forests have the same local elements.

Parameters
[in]forest_aThe first forest.
[in]forest_bThe second forest.
Returns
True if forest_a and forest_b do have the same number of local trees and each local tree has the same elements, that is t8_element_equal returns true for each pair of elements of forest_a and forest_b.
Note
This function is not collective. It only returns the state on the current rank.

◆ t8_forest_is_initialized()

int t8_forest_is_initialized ( t8_forest_t  forest)

Check whether a forest is not NULL, initialized and not committed.

In addition, it asserts that the forest is consistent as much as possible.

Parameters
[in]forestThis forest is examined. May be NULL.
Returns
True if forest is not NULL, t8_forest_init has been called on it, but not t8_forest_commit. False otherwise.

◆ t8_forest_leaf_face_neighbors()

void t8_forest_leaf_face_neighbors ( t8_forest_t  forest,
t8_locidx_t  ltreeid,
const t8_element_t leaf,
t8_element_t **  pneighbor_leaves[],
int  face,
int *  dual_faces[],
int *  num_neighbors,
t8_locidx_t **  pelement_indices,
t8_eclass_scheme_c **  pneigh_scheme,
int  forest_is_balanced 
)

Compute the leaf face neighbors of a forest.

Parameters
[in]forestThe forest. Must have a valid ghost layer.
[in]ltreeidA local tree id.
[in]leafA leaf in tree ltreeid of forest.
[out]neighbor_leavesUnallocated on input. On output the neighbor leaves are stored here.
[in]faceThe index of the face across which the face neighbors are searched.
[out]dual_faceOn output the face id's of the neighboring elements' faces.
[out]num_neighborsOn output the number of neighbor leaves.
[out]pelement_indicesUnallocated on input. On output the element indices of the neighbor leaves are stored here. 0, 1, ... num_local_el - 1 for local leaves and num_local_el , ... , num_local_el + num_ghosts - 1 for ghosts.
[out]pneigh_schemeOn output the eclass scheme of the neighbor elements.
[in]forest_is_balancedTrue if we know that forest is balanced, false otherwise.
Note
If there are no face neighbors, then *neighbor_leaves = NULL, num_neighbors = 0, and *pelement_indices = NULL on output.
Currently forest must be balanced.
forest must be committed before calling this function.

◆ t8_forest_ltreeid_to_cmesh_ltreeid()

t8_locidx_t t8_forest_ltreeid_to_cmesh_ltreeid ( t8_forest_t  forest,
t8_locidx_t  ltreeid 
)

Given the local id of a tree in a forest, compute the tree's local id in the associated cmesh.

Parameters
[in]forestThe forest.
[in]ltreeidThe local id of a tree or ghost in the forest.
Returns
The local id of the tree in the cmesh associated with the forest. forest must be committed before calling this function.
Note
For forest local trees, this is the inverse function of t8_forest_cmesh_ltreeid_to_ltreeid.
See also
https://github.com/DLR-AMR/t8code/wiki/Tree-indexing for more details about tree indexing.

◆ t8_forest_new_adapt()

t8_forest_t t8_forest_new_adapt ( t8_forest_t  forest_from,
t8_forest_adapt_t  adapt_fn,
int  recursive,
int  do_face_ghost,
void *  user_data 
)

Build a adapted forest from another forest.

Parameters
[in]forest_fromThe forest to refine
[in]adapt_fnAdapt function to use
[in]replace_fnReplace function to use
[in]recursiveIf true adptation is recursive
[in]do_face_ghostIf true, a layer of ghost elements is created for the forest.
[in]user_dataIf not NULL, the user data pointer of the forest is set to this value.
Returns
A new forest that is adapted from forest_from.
Note
This is equivalent to calling t8_forest_init, t8_forest_set_adapt, t8_forest_set_ghost, and t8_forest_commit

◆ t8_forest_new_uniform()

t8_forest_t t8_forest_new_uniform ( t8_cmesh_t  cmesh,
t8_scheme_cxx_t scheme,
const int  level,
const int  do_face_ghost,
sc_MPI_Comm  comm 
)

Build a uniformly refined forest on a coarse mesh.

Parameters
[in]cmeshA coarse mesh.
[in]schemeAn eclass scheme.
[in]levelAn initial uniform refinement level.
[in]do_face_ghostIf true, a layer of ghost elements is created for the forest.
[in]commMPI communicator to use.
Returns
A uniform forest with coarse mesh cmesh, eclass_scheme scheme and refinement level level.
Note
This is equivalent to calling t8_forest_init, t8_forest_set_cmesh, t8_forest_set_scheme, t8_forest_set_level, and t8_forest_commit.

◆ t8_forest_no_overlap()

int t8_forest_no_overlap ( t8_forest_t  forest)

Check whether the forest has local overlapping elements.

Parameters
[in]forestThe forest to consider.
Returns
True if forest has no elements which are inside each other.
Note
This function is collective, but only checks local overlapping on each process.
See also
t8_forest_partition_test_boundary_element if you also want to test for global overlap across the process boundaries.

◆ t8_forest_partition_cmesh()

void t8_forest_partition_cmesh ( t8_forest_t  forest,
sc_MPI_Comm  comm,
int  set_profiling 
)

Change the cmesh associated to a forest to a partitioned cmesh that is partitioned according to the tree distribution in the forest.

Parameters
[in,out]forestThe forest.
[in]commThe MPI communicator that is used to partition and commit the cmesh.
[in]set_profilingIf true, profiling for the new cmesh will be enabled.
See also
t8_cmesh_set_profiling,
t8_cmesh_print_profile
t8_cmesh.h

◆ t8_forest_ref()

void t8_forest_ref ( t8_forest_t  forest)

Increase the reference counter of a forest.

Parameters
[in,out]forestOn input, this forest must exist with positive reference count. It may be in any state.

◆ t8_forest_set_adapt()

void t8_forest_set_adapt ( t8_forest_t  forest,
const t8_forest_t  set_from,
t8_forest_adapt_t  adapt_fn,
int  recursive 
)

Set a source forest with an adapt function to be adapted on committing.

By default, the forest takes ownership of the source set_from such that it will be destroyed on calling t8_forest_commit. To keep ownership of set_from, call t8_forest_ref before passing it into this function. This means that it is ILLEGAL to continue using set_from or dereferencing it UNLESS it is referenced directly before passing it into this function.

Parameters
[in,out]forestThe forest
[in]set_fromThe source forest from which forest will be adapted. We take ownership. This can be prevented by referencing set_from. If NULL, a previously (or later) set forest will be taken (t8_forest_set_partition, t8_forest_set_balance).
[in]adapt_fnThe adapt function used on committing.
[in]recursiveA flag specifying whether adaptation is to be done recursively or not. If the value is zero, adaptation is not recursive and it is recursive otherwise.
Note
This setting can be combined with t8_forest_set_partition and t8_forest_set_balance. The order in which these operations are executed is always 1) Adapt 2) Balance 3) Partition
This setting may not be combined with t8_forest_set_copy and overwrites this setting.

◆ t8_forest_set_balance()

void t8_forest_set_balance ( t8_forest_t  forest,
const t8_forest_t  set_from,
int  no_repartition 
)

Set a source forest to be balanced during commit.

A forest is said to be balanced if each element has face neighbors of level at most +1 or -1 of the element's level.

Parameters
[in,out]forestThe forest.
[in]set_fromA second forest that should be balanced. We take ownership. This can be prevented by referencing set_from. If NULL, a previously (or later) set forest will be taken (t8_forest_set_adapt, t8_forest_set_partition)
[in]no_repartitionBalance constructs several intermediate forest that are refined from each other. In order to maintain a balanced load these forest are repartitioned in each round and the resulting forest is load-balanced per default. If this behaviour is not desired, no_repartition should be set to true. If no_repartition is false, an additional call of t8_forest_set_partition is not necessary.
Note
This setting can be combined with t8_forest_set_adapt and t8_forest_set_balance. The order in which these operations are executed is always 1) Adapt 2) Balance 3) Partition.
This setting may not be combined with t8_forest_set_copy and overwrites this setting.

◆ t8_forest_set_cmesh()

void t8_forest_set_cmesh ( t8_forest_t  forest,
t8_cmesh_t  cmesh,
sc_MPI_Comm  comm 
)

Set the cmesh associated to a forest.

By default, the forest takes ownership of the cmesh such that it will be destroyed when the forest is destroyed. To keep ownership of the cmesh, call t8_cmesh_ref before passing it to t8_forest_set_cmesh. This means that it is ILLEGAL to continue using cmesh or dereferencing it UNLESS it is referenced directly before passing it into this function.

Parameters
[in,out]forestThe forest whose cmesh variable will be set.
[in]cmeshThe cmesh to be set. We take ownership. This can be prevented by referencing cmesh.

◆ t8_forest_set_copy()

void t8_forest_set_copy ( t8_forest_t  forest,
const t8_forest_t  from 
)

Set a forest as source for copying on committing.

By default, the forest takes ownership of the source from such that it will be destroyed on calling t8_forest_commit. To keep ownership of from, call t8_forest_ref before passing it into this function. This means that it is ILLEGAL to continue using from or dereferencing it UNLESS it is referenced directly before passing it into this function.

Parameters
[in,out]forestThe forest.
[in]fromA second forest from which forest will be copied in t8_forest_commit.
Note
This setting cannot be combined with t8_forest_set_adapt, t8_forest_set_partition, or t8_forest_set_balance and overwrites these settings.

◆ t8_forest_set_ghost()

void t8_forest_set_ghost ( t8_forest_t  forest,
int  do_ghost,
t8_ghost_type_t  ghost_type 
)

Enable or disable the creation of a layer of ghost elements.

On default no ghosts are created.

Parameters
[in]forestThe forest.
[in]do_ghostIf non-zero a ghost layer will be created.
[in]ghost_typeControls which neighbors count as ghost elements, currently only T8_GHOST_FACES is supported. This value is ignored if do_ghost = 0.

◆ t8_forest_set_ghost_ext()

void t8_forest_set_ghost_ext ( t8_forest_t  forest,
int  do_ghost,
t8_ghost_type_t  ghost_type,
int  ghost_version 
)

Like t8_forest_set_ghost but with the additional options to change the ghost algorithm.

This is used for debugging and timing the algorithm. An application should almost always use t8_forest_set_ghost.

Parameters
[in]ghost_versionIf 1, the iterative ghost algorithm for balanced forests is used. If 2, the iterative algorithm for unbalanced forests. If 3, the top-down search algorithm for unbalanced forests.
See also
t8_forest_set_ghost

◆ t8_forest_set_level()

void t8_forest_set_level ( t8_forest_t  forest,
int  level 
)

Set the initial refinement level to be used when forest is committed.

Parameters
[in,out]forestThe forest whose level will be set.
[in]levelThe initial refinement level of forest, when it is committed.
Note
This setting cannot be combined with any of the derived forest methods (t8_forest_set_copy, t8_forest_set_adapt, t8_forest_set_partition, and t8_forest_set_balance) and overwrites any of these settings. If this function is used, then the forest is created from scratch as a uniform refinement of the specified cmesh (t8_forest_set_cmesh, t8_forest_set_scheme).

◆ t8_forest_set_partition()

void t8_forest_set_partition ( t8_forest_t  forest,
const t8_forest_t  set_from,
int  set_for_coarsening 
)

Set a source forest to be partitioned during commit.

The partitioning is done according to the SFC and each rank is assigned the same (maybe +1) number of elements.

Parameters
[in,out]forestThe forest.
[in]set_fromA second forest that should be partitioned. We take ownership. This can be prevented by referencing set_from. If NULL, a previously (or later) set forest will be taken (t8_forest_set_adapt, t8_forest_set_balance).
[in]set_for_coarseningCURRENTLY DISABLED. If true, then the partitions are choose such that coarsening an element once is a process local operation.
Note
This setting can be combined with t8_forest_set_adapt and t8_forest_set_balance. The order in which these operations are executed is always 1) Adapt 2) Balance 3) Partition If t8_forest_set_balance is called with the no_repartition parameter set as false, it is not necessary to call t8_forest_set_partition additionally.
This setting may not be combined with t8_forest_set_copy and overwrites this setting.

◆ t8_forest_set_scheme()

void t8_forest_set_scheme ( t8_forest_t  forest,
t8_scheme_cxx_t scheme 
)

Set the element scheme associated to a forest.

By default, the forest takes ownership of the scheme such that it will be destroyed when the forest is destroyed. To keep ownership of the scheme, call t8_scheme_ref before passing it to t8_forest_set_scheme. This means that it is ILLEGAL to continue using scheme or dereferencing it UNLESS it is referenced directly before passing it into this function.

Parameters
[in,out]forestThe forest whose scheme variable will be set.
[in]schemeThe scheme to be set. We take ownership. This can be prevented by referencing scheme.

◆ t8_forest_set_user_data()

void t8_forest_set_user_data ( t8_forest_t  forest,
void *  data 
)

Set the user data of a forest.

This can i.e. be used to pass user defined arguments to the adapt routine.

Parameters
[in,out]forestThe forest
[in]dataA pointer to user data. t8code will never touch the data. The forest does not need be committed before calling this function.
See also
t8_forest_get_user_data

◆ t8_forest_set_user_function()

void t8_forest_set_user_function ( t8_forest_t  forest,
t8_generic_function_pointer  function 
)

Set the user function pointer of a forest.

This can i.e. be used to pass user defined functions to the adapt routine.

Parameters
[in,out]forestThe forest
[in]functionA pointer to a user defined function. t8code will never touch the function. The forest does not need be committed before calling this function.
Note
function can be an arbitrary function with return value and parameters of your choice. When accessing it with t8_forest_get_user_function you should cast it into the proper type.
See also
t8_forest_get_user_function

◆ t8_forest_tree_get_leaves()

t8_element_array_t* t8_forest_tree_get_leaves ( const t8_forest_t  forest,
const t8_locidx_t  ltree_id 
)

Return the array of leaf elements of a local tree in a forest.

Parameters
[in]forestThe forest.
[in]ltree_idThe local id of a local tree of forest.
Returns
An array of t8_element_t * storing all leaf elements of this tree.

◆ t8_forest_unref()

void t8_forest_unref ( t8_forest_t pforest)

Decrease the reference counter of a forest.

If the counter reaches zero, this forest is destroyed. In this case, the forest dereferences its cmesh and scheme members.

Parameters
[in,out]pforestOn input, the forest pointed to must exist with positive reference count. It may be in any state. If the reference count reaches zero, the forest is destroyed and this pointer set to NULL. Otherwise, the pointer is not changed and the forest is not modified in other ways.