Program Listing for File t8_scheme.h

Return to documentation for file (src/t8_schemes/t8_scheme.h)

/*
  This file is part of t8code.
  t8code is a C library to manage a collection (a forest) of multiple
  connected adaptive space-trees of general element classes in parallel.

  Copyright (C) 2024 the developers

  t8code is free software; you can redistribute it and/or modify
  it under the terms of the GNU General Public License as published by
  the Free Software Foundation; either version 2 of the License, or
  (at your option) any later version.

  t8code is distributed in the hope that it will be useful,
  but WITHOUT ANY WARRANTY; without even the implied warranty of
  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
  GNU General Public License for more details.

  You should have received a copy of the GNU General Public License
  along with t8code; if not, write to the Free Software Foundation, Inc.,
  51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/

#ifndef T8_SCHEME_H
#define T8_SCHEME_H

#include <t8.h>
#include <t8_element/t8_element.h>

typedef struct t8_scheme t8_scheme_c;

T8_EXTERN_C_BEGIN ();

void
t8_scheme_ref (t8_scheme_c *scheme);

void
t8_scheme_unref (t8_scheme_c **pscheme);

size_t
t8_element_get_element_size (const t8_scheme_c *scheme, const t8_eclass_t tree_class);

int
t8_element_refines_irregular (const t8_scheme_c *scheme, const t8_eclass_t tree_class);

int
t8_element_get_maxlevel (const t8_scheme_c *scheme, const t8_eclass_t tree_class);

int
t8_element_get_level (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

void
t8_element_copy (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *source,
                 t8_element_t *dest);

int
t8_element_compare (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *elem1,
                    const t8_element_t *elem2);

int
t8_element_is_equal (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *elem1,
                     const t8_element_t *elem2);

int
element_is_refinable (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

void
t8_element_get_parent (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                       t8_element_t *parent);

int
t8_element_get_num_siblings (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

void
t8_element_get_sibling (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *elem,
                        const int sibid, t8_element_t *sibling);

int
t8_element_get_num_corners (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

int
t8_element_get_num_faces (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

int
t8_element_get_max_num_faces (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

int
t8_element_get_num_children (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

int
t8_get_max_num_children (const t8_scheme_c *scheme, const t8_eclass_t tree_class);

int
t8_element_get_num_face_children (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                                  const int face);

int
t8_element_get_face_corner (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                            const int face, const int corner);

int
t8_element_get_corner_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                            const int corner, const int face);

void
t8_element_get_child (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                      const int childid, t8_element_t *child);

void
t8_element_get_children (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                         const int length, t8_element_t *c[]);

int
t8_element_get_child_id (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

int
t8_element_get_ancestor_id (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                            const int level);

int
t8_elements_are_family (const t8_scheme_c *scheme, const t8_eclass_t tree_class, t8_element_t *const *fam);

void
t8_element_get_nca (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *elem1,
                    const t8_element_t *elem2, t8_element_t *nca);

t8_element_shape_t
t8_element_get_face_shape (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                           const int face);

void
t8_element_get_children_at_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                                 const int face, t8_element_t *children[], const int num_children, int *child_indices);

int
t8_element_face_get_child_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                                const int face, const int face_child);

int
t8_element_face_get_parent_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                                 const int face);

int
t8_element_get_tree_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                          const int face);

void
t8_element_transform_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *elem1,
                           t8_element_t *elem2, const int orientation, const int sign, const int is_smaller_face);

int
t8_element_extrude_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *face,
                         t8_element_t *element, int root_face);

void
t8_element_get_boundary_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                              const int face, t8_element_t *boundary);

void
t8_element_get_first_descendant_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class,
                                      const t8_element_t *element, const int face, t8_element_t *first_desc,
                                      const int level);

void
t8_element_get_last_descendant_face (const t8_scheme_c *scheme, const t8_eclass_t tree_class,
                                     const t8_element_t *element, const int face, t8_element_t *last_desc,
                                     const int level);

int
t8_element_is_root_boundary (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                             const int face);

int
t8_element_get_face_neighbor_inside (const t8_scheme_c *scheme, const t8_eclass_t tree_class,
                                     const t8_element_t *element, t8_element_t *neigh, const int face, int *neigh_face);

t8_element_shape_t
t8_element_get_shape (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

void
t8_element_set_linear_id (const t8_scheme_c *scheme, const t8_eclass_t tree_class, t8_element_t *element,
                          const int level, const t8_linearidx_t id);

t8_linearidx_t
t8_element_get_linear_id (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                          const int level);

void
t8_element_get_first_descendant (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                                 t8_element_t *desc, const int level);

void
t8_element_get_last_descendant (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                                t8_element_t *desc, const int level);

void
t8_element_get_successor (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *elem1,
                          t8_element_t *elem2);

void
t8_element_get_vertex_reference_coords (const t8_scheme_c *scheme, const t8_eclass_t tree_class,
                                        const t8_element_t *element, const int vertex, double coords[]);

void
t8_element_get_reference_coords (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                                 const double *ref_coords, const size_t num_coords, double out_coords[]);

t8_gloidx_t
t8_element_count_leaves (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                         const int level);

t8_gloidx_t
t8_element_count_leaves_from_root (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const int level);

#if T8_ENABLE_DEBUG
int
t8_element_is_valid (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

void
t8_element_debug_print (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element);

#endif
void
t8_element_to_string (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const t8_element_t *element,
                      char *debug_string, const int string_size);

void
t8_element_new (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const int length, t8_element_t **elems);

void
t8_element_init (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const int length, t8_element_t *elem);

void
t8_element_deinit (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const int length, t8_element_t *elems);

void
t8_element_destroy (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const int length, t8_element_t **elems);

void
t8_element_set_to_root (const t8_scheme_c *scheme, const t8_eclass_t tree_class, t8_element_t *element);

void
t8_element_MPI_Pack (const t8_scheme_c *scheme, const t8_eclass_t tree_class, t8_element_t **const elements,
                     const unsigned int count, void *send_buffer, const int buffer_size, int *position,
                     sc_MPI_Comm comm);

void
t8_element_MPI_Pack_size (const t8_scheme_c *scheme, const t8_eclass_t tree_class, const unsigned int count,
                          sc_MPI_Comm comm, int *pack_size);

void
t8_element_MPI_Unpack (const t8_scheme_c *scheme, const t8_eclass_t tree_class, void *recvbuf, const int buffer_size,
                       int *position, t8_element_t **elements, const unsigned int count, sc_MPI_Comm comm);

T8_EXTERN_C_END ();

#endif /* !T8_SCHEME_H */