Program Listing for File t8_cmesh.h

Return to documentation for file (src/t8_cmesh/t8_cmesh.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) 2015 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_CMESH_H
#define T8_CMESH_H

#include <t8.h>
#include <t8_data/t8_shmem.h>
#include <t8_cmesh/t8_cmesh_io/t8_cmesh_save.h>
#include <t8_element.h>
#include <t8_schemes/t8_scheme.h>

typedef struct t8_cmesh *t8_cmesh_t;

#include <t8_geometry/t8_geometry.h>

/* TODO: If including eclass were just for the cmesh_new routines, we should
 *       move them into a different file.
 *       However, when specifying the parent-child order in cmesh_reorder,
 *       we might keep the eclass interface for virtual functions.
 *       Actually, we need eclass in the type definition in cmesh.c.
 *       So we might as well use tree-related virtual functions there too.
 */
#include <t8_eclass.h>

/* TODO: make it legal to call cmesh_set functions multiple times,
 *       just overwrite the previous setting if no inconsistency can occur.
 *       edit: This should be achieved now.
 */

typedef struct t8_ctree *t8_ctree_t;
typedef struct t8_cghost *t8_cghost_t;

T8_EXTERN_C_BEGIN ();

void
t8_cmesh_init (t8_cmesh_t *pcmesh);

t8_cmesh_t
t8_cmesh_new ();

int
t8_cmesh_is_initialized (t8_cmesh_t cmesh);

int
t8_cmesh_is_committed (const t8_cmesh_t cmesh);

void
t8_cmesh_disable_negative_volume_check (t8_cmesh_t cmesh);

#if T8_ENABLE_DEBUG
int

t8_cmesh_validate_geometry (const t8_cmesh_t cmesh, const int check_for_negative_volume);
#endif

/* TODO: Currently it is not possible to destroy set_from before
 *       cmesh is destroyed. */
void
t8_cmesh_set_derive (t8_cmesh_t cmesh, t8_cmesh_t set_from);

t8_shmem_array_t
t8_cmesh_alloc_offsets (int mpisize, sc_MPI_Comm comm);

void
t8_cmesh_set_partition_range (t8_cmesh_t cmesh, int set_face_knowledge, t8_gloidx_t first_local_tree,
                              t8_gloidx_t last_local_tree);

void
t8_cmesh_set_partition_offsets (t8_cmesh_t cmesh, t8_shmem_array_t tree_offsets);

void
t8_cmesh_set_partition_uniform (t8_cmesh_t cmesh, const int element_level, const t8_scheme_c *scheme);

/* If level = 0  then no refinement is performed */
void
t8_cmesh_set_refine (t8_cmesh_t cmesh, const int level, const t8_scheme_c *scheme);

void
t8_cmesh_set_dimension (t8_cmesh_t cmesh, int dim);

void
t8_cmesh_set_tree_class (t8_cmesh_t cmesh, t8_gloidx_t gtree_id, t8_eclass_t tree_class);

void
t8_cmesh_set_attribute (t8_cmesh_t cmesh, t8_gloidx_t gtree_id, int package_id, int key, void *data, size_t data_size,
                        int data_persists);

void
t8_cmesh_set_attribute_string (t8_cmesh_t cmesh, t8_gloidx_t gtree_id, int package_id, int key, const char *string);

void
t8_cmesh_set_attribute_gloidx_array (t8_cmesh_t cmesh, t8_gloidx_t gtree_id, int package_id, int key,
                                     const t8_gloidx_t *data, const size_t data_count, int data_persists);

void
t8_cmesh_set_join (t8_cmesh_t cmesh, t8_gloidx_t gtree1, t8_gloidx_t gtree2, int face1, int face2, int orientation);

void
t8_cmesh_set_profiling (t8_cmesh_t cmesh, int set_profiling);

/* returns true if cmesh_a equals cmesh_b */
/* TODO: document
 * collective or serial */
int
t8_cmesh_is_equal (t8_cmesh_t cmesh_a, t8_cmesh_t cmesh_b);

int
t8_cmesh_is_empty (t8_cmesh_t cmesh);

t8_cmesh_t
t8_cmesh_bcast (t8_cmesh_t cmesh_in, int root, sc_MPI_Comm comm);

#if T8_ENABLE_METIS
/* TODO: document this. */
/* TODO: think about making this a pre-commit set_reorder function. */
void
t8_cmesh_reorder (t8_cmesh_t cmesh, sc_MPI_Comm comm);

/* TODO: think about a sensible interface for a parmetis reordering. */
#endif

void
t8_cmesh_register_geometry (t8_cmesh_t cmesh, t8_geometry_c *geometry);

void
t8_cmesh_set_tree_geometry (t8_cmesh_t cmesh, t8_gloidx_t gtreeid, const t8_geometry_c *geom);

void
t8_cmesh_commit (t8_cmesh_t cmesh, sc_MPI_Comm comm);

int
t8_cmesh_save (t8_cmesh_t cmesh, const char *fileprefix);

t8_cmesh_t
t8_cmesh_load (const char *filename, sc_MPI_Comm comm);

t8_cmesh_t
t8_cmesh_load_and_distribute (const char *fileprefix, int num_files, sc_MPI_Comm comm, t8_load_mode_t mode,
                              int procs_per_node);

int
t8_cmesh_comm_is_valid (t8_cmesh_t cmesh, sc_MPI_Comm comm);

int
t8_cmesh_is_partitioned (t8_cmesh_t cmesh);

int
t8_cmesh_get_dimension (const t8_cmesh_t cmesh);

t8_gloidx_t
t8_cmesh_get_num_trees (t8_cmesh_t cmesh);

t8_locidx_t
t8_cmesh_get_num_local_trees (t8_cmesh_t cmesh);

t8_locidx_t
t8_cmesh_get_num_ghosts (t8_cmesh_t cmesh);

t8_gloidx_t
t8_cmesh_get_first_treeid (t8_cmesh_t cmesh);

const t8_geometry_c *
t8_cmesh_get_tree_geometry (t8_cmesh_t cmesh, t8_gloidx_t gtreeid);

int
t8_cmesh_treeid_is_local_tree (const t8_cmesh_t cmesh, const t8_locidx_t ltreeid);

int
t8_cmesh_treeid_is_ghost (const t8_cmesh_t cmesh, const t8_locidx_t ltreeid);

t8_locidx_t
t8_cmesh_ltreeid_to_ghostid (const t8_cmesh_t cmesh, const t8_locidx_t ltreeid);

/* TODO: Replace this iterator with a new one that does not need the
 *        treeid to be part of the ctree struct */
/* TODO: should this and the next function be part of the interface? */
t8_ctree_t
t8_cmesh_get_first_tree (t8_cmesh_t cmesh);

/* TODO: should this function behave like first_tree if tree argument is NULL? */
t8_ctree_t
t8_cmesh_get_next_tree (t8_cmesh_t cmesh, t8_ctree_t tree);

t8_ctree_t
t8_cmesh_get_tree (t8_cmesh_t cmesh, t8_locidx_t ltree_id);

t8_eclass_t
t8_cmesh_get_tree_class (t8_cmesh_t cmesh, t8_locidx_t ltree_id);

int
t8_cmesh_tree_face_is_boundary (t8_cmesh_t cmesh, t8_locidx_t ltree_id, int face);

t8_eclass_t
t8_cmesh_get_ghost_class (t8_cmesh_t cmesh, t8_locidx_t lghost_id);

t8_gloidx_t
t8_cmesh_get_global_id (t8_cmesh_t cmesh, t8_locidx_t local_id);

t8_locidx_t
t8_cmesh_get_local_id (t8_cmesh_t cmesh, t8_gloidx_t global_id);

t8_locidx_t
t8_cmesh_get_face_neighbor (const t8_cmesh_t cmesh, const t8_locidx_t ltreeid, const int face, int *dual_face,
                            int *orientation);

t8_eclass_t
t8_cmesh_get_tree_face_neighbor_eclass (const t8_cmesh_t cmesh, const t8_locidx_t ltreeid, const int face);

void
t8_cmesh_print_profile (t8_cmesh_t cmesh);

double *
t8_cmesh_get_tree_vertices (t8_cmesh_t cmesh, t8_locidx_t ltreeid);

void *
t8_cmesh_get_attribute (const t8_cmesh_t cmesh, const int package_id, const int key, const t8_locidx_t ltree_id);

t8_gloidx_t *
t8_cmesh_get_attribute_gloidx_array (const t8_cmesh_t cmesh, const int package_id, const int key,
                                     const t8_locidx_t ltree_id, const size_t data_count);

t8_shmem_array_t
t8_cmesh_get_partition_table (t8_cmesh_t cmesh);

/* TODO: remove get_ when there is no risk of confusion? Convention?
 *       Update: use get throughout for access functions that do not change the object.
 * */

void
t8_cmesh_uniform_bounds_equal_element_count (t8_cmesh_t cmesh, const int level, const t8_scheme_c *tree_scheme,
                                             t8_gloidx_t *first_local_tree, t8_gloidx_t *child_in_tree_begin,
                                             t8_gloidx_t *last_local_tree, t8_gloidx_t *child_in_tree_end,
                                             int8_t *first_tree_shared);

void
t8_cmesh_uniform_bounds_for_irregular_refinement (const t8_cmesh_t cmesh, const int level, const t8_scheme_c *scheme,
                                                  t8_gloidx_t *first_local_tree, t8_gloidx_t *child_in_tree_begin,
                                                  t8_gloidx_t *last_local_tree, t8_gloidx_t *child_in_tree_end,
                                                  int8_t *first_tree_shared, sc_MPI_Comm comm);

void
t8_cmesh_ref (t8_cmesh_t cmesh);

void
t8_cmesh_unref (t8_cmesh_t *pcmesh);

void
t8_cmesh_destroy (t8_cmesh_t *pcmesh);

void
t8_cmesh_coords_axb (const double *coords_in, double *coords_out, int num_vertices, double alpha, const double b[3]);

void
t8_cmesh_translate_coordinates (const double *coords_in, double *coords_out, const int num_vertices,
                                const double translate[3]);

void
t8_cmesh_new_translate_vertices_to_attributes (const t8_locidx_t *tvertices, const double *vertices,
                                               double *attr_vertices, const int num_vertices);

void
t8_cmesh_debug_print_trees (const t8_cmesh_t cmesh, sc_MPI_Comm comm);

int
t8_cmesh_get_local_bounding_box (const t8_cmesh_t cmesh, double bounds[6]);
T8_EXTERN_C_END ();

#endif /* !T8_CMESH_H */