p4est 2.8.6
p4est is a software library for parallel adaptive mesh refinement.
p8est_extended.h
Go to the documentation of this file.
1/*
2 This file is part of p4est.
3 p4est is a C library to manage a collection (a forest) of multiple
4 connected adaptive quadtrees or octrees in parallel.
5
6 Copyright (C) 2010 The University of Texas System
7 Additional copyright (C) 2011 individual authors
8 Written by Carsten Burstedde, Lucas C. Wilcox, and Tobin Isaac
9
10 p4est is free software; you can redistribute it and/or modify
11 it under the terms of the GNU General Public License as published by
12 the Free Software Foundation; either version 2 of the License, or
13 (at your option) any later version.
14
15 p4est is distributed in the hope that it will be useful,
16 but WITHOUT ANY WARRANTY; without even the implied warranty of
17 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
18 GNU General Public License for more details.
19
20 You should have received a copy of the GNU General Public License
21 along with p4est; if not, write to the Free Software Foundation, Inc.,
22 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
23*/
24
25/********************************************************************
26 * IMPORTANT NOTE *
27 * *
28 * These interfaces are intended for those who like finer control. *
29 * The API offers extended versions of some basic p4est functions. *
30 * The API may change without notice. *
31 ********************************************************************/
32
40#ifndef P8EST_EXTENDED_H
41#define P8EST_EXTENDED_H
42
43#include <sc_uint128.h>
44#include <p8est_mesh.h>
45#include <p8est_iterate.h>
46#include <p8est_lnodes.h>
47#include <p8est_io.h>
48
49SC_EXTERN_C_BEGIN;
50
56typedef sc_uint128_t p8est_lid_t;
57
58/* Data pertaining to selecting, inspecting, and profiling algorithms.
59 * A pointer to this structure is hooked into the p8est main structure.
60 *
61 * TODO: Describe the purpose of various switches, counters, and timings.
62 *
63 * The balance_ranges and balance_notify* times are collected
64 * whenever an inspect structure is present in p8est.
65 */
67{
78 size_t balance_A_count_in;
79 size_t balance_A_count_out;
80 size_t balance_comm_sent;
81 size_t balance_comm_nzpeers;
82 size_t balance_B_count_in;
83 size_t balance_B_count_out;
84 size_t balance_zero_sends[2], balance_zero_receives[2];
85 double balance_A;
86 double balance_comm;
87 double balance_B;
92 int use_B;
93};
94
114typedef void (*p8est_replace_t) (p8est_t * p8est,
115 p4est_topidx_t which_tree,
116 int num_outgoing,
117 p8est_quadrant_t * outgoing[],
118 int num_incoming,
119 p8est_quadrant_t * incoming[]);
120
129 const p8est_lid_t * b);
130
138 const p8est_lid_t * b);
139
146void p8est_lid_init (p8est_lid_t * input, uint64_t high,
147 uint64_t low);
148
153
158
162void p8est_lid_set_uint64 (p8est_lid_t * input, uint64_t u);
163
173 int bit_number);
174
182void p8est_lid_set_bit (p8est_lid_t * input, int bit_number);
183
191void p8est_lid_copy (const p8est_lid_t * input,
192 p8est_lid_t * output);
193
203 const p8est_lid_t * b,
204 p8est_lid_t * result);
205
216 const p8est_lid_t * b,
217 p8est_lid_t * result);
218
227 p8est_lid_t * result);
228
239 const p8est_lid_t * b,
240 p8est_lid_t * result);
241
252 const p8est_lid_t * b,
253 p8est_lid_t * result);
254
266 unsigned shift_count,
267 p8est_lid_t * result);
268
280 unsigned shift_count,
281 p8est_lid_t * result);
282
290 const p8est_lid_t * b);
291
300 const p8est_lid_t * b);
301
309 const p8est_lid_t * b);
310
318 const p8est_lid_t * b);
337 quadrant, int level,
338 p8est_lid_t * id);
339
348 quadrant, int level,
349 const p8est_lid_t * id);
350
381p8est_t *p8est_new_ext (sc_MPI_Comm mpicomm,
382 p8est_connectivity_t * connectivity,
383 p4est_locidx_t min_quadrants,
384 int min_level, int fill_uniform,
385 size_t data_size, p8est_init_t init_fn,
386 void *user_pointer);
387
403 p8est_ghost_t * ghost,
404 int compute_tree_index,
405 int compute_level_lists,
407
422p8est_t *p8est_copy_ext (p8est_t * input, int copy_data,
423 int duplicate_mpicomm);
424
447 int refine_recursive, int maxlevel,
448 p8est_refine_t refine_fn,
449 p8est_init_t init_fn,
450 p8est_replace_t replace_fn);
451
473void p8est_coarsen_ext (p8est_t * p8est, int coarsen_recursive,
474 int callback_orphans,
475 p8est_coarsen_t coarsen_fn,
476 p8est_init_t init_fn,
477 p8est_replace_t replace_fn);
478
492 p8est_init_t init_fn,
493 p8est_replace_t replace_fn);
494
495void p8est_balance_subtree_ext (p8est_t * p8est,
497 p4est_topidx_t which_tree,
498 p8est_init_t init_fn,
499 p8est_replace_t replace_fn);
500
518 int partition_for_coarsening,
519 p8est_weight_t weight_fn);
520
529 num_quadrants_in_proc);
530
537 p8est_ghost_t * ghost_layer,
538 void *user_data,
539 p8est_iter_volume_t iter_volume,
540 p8est_iter_face_t iter_face,
541 p8est_iter_edge_t iter_edge,
542 p8est_iter_corner_t iter_corner,
543 int remote);
544
562void p8est_save_ext (const char *filename, p8est_t * p8est,
563 int save_data, int save_partition);
564
586p8est_t *p8est_load_ext (const char *filename, sc_MPI_Comm mpicomm,
587 size_t data_size, int load_data,
588 int autopartition, int broadcasthead,
589 void *user_pointer,
590 p8est_connectivity_t ** connectivity);
591
595p8est_t *p8est_source_ext (sc_io_source_t * src,
596 sc_MPI_Comm mpicomm, size_t data_size,
597 int load_data, int autopartition,
598 int broadcasthead, void *user_pointer,
599 p8est_connectivity_t ** connectivity);
600
601#ifdef P4EST_ENABLE_FILE_DEPRECATED
602
612p8est_file_context_t *p8est_file_open_read_ext (sc_MPI_Comm mpicomm,
613 const char *filename,
614 char *user_string,
616 global_num_quadrants,
617 int *errcode);
618
627p8est_file_context_t *p8est_file_read_field_ext (p8est_file_context_t * fc,
628 p4est_gloidx_t * gfq,
629 size_t quadrant_size,
630 sc_array_t * quadrant_data,
631 char *user_string,
632 int *errcode);
633
634#endif /* P4EST_ENABLE_FILE_DEPRECATED */
635
679 p8est_ghost_t ** ghost,
680 p8est_lnodes_t ** lnodes,
682 int overlap,
684 first_local_quad,
685 sc_array_t * out_points_per_dim,
686 sc_array_t * out_cone_sizes,
687 sc_array_t * out_cones,
688 sc_array_t *
689 out_cone_orientations,
690 sc_array_t * out_vertex_coords,
691 sc_array_t * out_children,
692 sc_array_t * out_parents,
693 sc_array_t * out_childids,
694 sc_array_t * out_leaves,
695 sc_array_t * out_remotes,
696 int custom_numbering);
697
698SC_EXTERN_C_END;
699
700#endif /* !P8EST_EXTENDED_H */
int32_t p4est_topidx_t
Typedef for counting topological entities (trees, tree vertices).
Definition: p4est_base.h:93
int32_t p4est_locidx_t
Typedef for processor-local indexing of quadrants and nodes.
Definition: p4est_base.h:106
int64_t p4est_gloidx_t
Typedef for globally unique indexing of quadrants.
Definition: p4est_base.h:118
int(* p8est_weight_t)(p8est_t *p8est, p4est_topidx_t which_tree, p8est_quadrant_t *quadrant)
Callback function prototype to calculate weights for partitioning.
Definition: p8est.h:229
int(* p8est_coarsen_t)(p8est_t *p8est, p4est_topidx_t which_tree, p8est_quadrant_t *quadrants[])
Callback function prototype to decide for coarsening.
Definition: p8est.h:219
int(* p8est_refine_t)(p8est_t *p8est, p4est_topidx_t which_tree, p8est_quadrant_t *quadrant)
Callback function prototype to decide for refinement.
Definition: p8est.h:209
void(* p8est_init_t)(p8est_t *p8est, p4est_topidx_t which_tree, p8est_quadrant_t *quadrant)
Callback function prototype to initialize the quadrant's user data.
Definition: p8est.h:199
p8est_connect_type_t
Characterize a type of adjacency.
Definition: p8est_connectivity.h:119
void p8est_lid_bitwise_and_inplace(p8est_lid_t *a, const p8est_lid_t *b)
Calculates the bitwise and of the uint128_t a and the uint128_t b.
int p8est_lid_chk_bit(const p8est_lid_t *input, int bit_number)
Returns the bit_number-th bit of input.
void p8est_coarsen_ext(p8est_t *p8est, int coarsen_recursive, int callback_orphans, p8est_coarsen_t coarsen_fn, p8est_init_t init_fn, p8est_replace_t replace_fn)
Coarsen a forest.
void p8est_quadrant_linear_id_ext128(const p8est_quadrant_t *quadrant, int level, p8est_lid_t *id)
Computes the linear position as p8est_lid_t of a quadrant in a uniform grid.
void p8est_refine_ext(p8est_t *p8est, int refine_recursive, int maxlevel, p8est_refine_t refine_fn, p8est_init_t init_fn, p8est_replace_t replace_fn)
Refine a forest with a bounded refinement level and a replace option.
void p8est_get_plex_data_ext(p8est_t *p8est, p8est_ghost_t **ghost, p8est_lnodes_t **lnodes, p8est_connect_type_t ctype, int overlap, p4est_locidx_t *first_local_quad, sc_array_t *out_points_per_dim, sc_array_t *out_cone_sizes, sc_array_t *out_cones, sc_array_t *out_cone_orientations, sc_array_t *out_vertex_coords, sc_array_t *out_children, sc_array_t *out_parents, sc_array_t *out_childids, sc_array_t *out_leaves, sc_array_t *out_remotes, int custom_numbering)
Create the data necessary to create a PETsc DMPLEX representation of a forest, as well as the accompa...
void p8est_lid_bitwise_or_inplace(p8est_lid_t *a, const p8est_lid_t *b)
Calculates the bitwise or of the uint128_t a and the uint128_t b.
p4est_gloidx_t p8est_partition_ext(p8est_t *p8est, int partition_for_coarsening, p8est_weight_t weight_fn)
Repartition the forest.
void p8est_quadrant_set_morton_ext128(p8est_quadrant_t *quadrant, int level, const p8est_lid_t *id)
Set quadrant Morton indices based on linear position given as p8est_lid_t in uniform grid.
void p8est_lid_set_bit(p8est_lid_t *input, int bit_number)
Sets the exponent-th bit of input to one.
void p8est_lid_set_uint64(p8est_lid_t *input, uint64_t u)
Initializes a linear index to an unsigned 64 bit integer.
sc_uint128_t p8est_lid_t
A datatype to handle the linear id in 3D.
Definition: p8est_extended.h:56
int p8est_lid_is_equal(const p8est_lid_t *a, const p8est_lid_t *b)
Checks if the p8est_lid_t a and the p8est_lid_t b are equal.
void p8est_lid_init(p8est_lid_t *input, uint64_t high, uint64_t low)
Initializes a linear index to a given value.
void p8est_lid_bitwise_and(const p8est_lid_t *a, const p8est_lid_t *b, p8est_lid_t *result)
Calculates the bitwise and of the uint128_t a and the uint128_t b.
void p8est_lid_sub(const p8est_lid_t *a, const p8est_lid_t *b, p8est_lid_t *result)
Subtracts the p8est_lid_t b from the p8est_lid_t a.
p4est_gloidx_t p8est_partition_for_coarsening(p8est_t *p8est, p4est_locidx_t *num_quadrants_in_proc)
Correct partition to allow one level of coarsening.
void p8est_lid_shift_right(const p8est_lid_t *input, unsigned shift_count, p8est_lid_t *result)
Calculates the bit right shift of uint128_t input by shift_count bits.
void p8est_save_ext(const char *filename, p8est_t *p8est, int save_data, int save_partition)
Save the complete connectivity/p8est data to disk.
p8est_mesh_t * p8est_mesh_new_ext(p8est_t *p8est, p8est_ghost_t *ghost, int compute_tree_index, int compute_level_lists, p8est_connect_type_t btype)
Create a new mesh.
void p8est_lid_shift_left(const p8est_lid_t *input, unsigned shift_count, p8est_lid_t *result)
Calculates the bit left shift of uint128_t input by shift_count bits.
p8est_t * p8est_new_ext(sc_MPI_Comm mpicomm, p8est_connectivity_t *connectivity, p4est_locidx_t min_quadrants, int min_level, int fill_uniform, size_t data_size, p8est_init_t init_fn, void *user_pointer)
Create a new forest.
void p8est_iterate_ext(p8est_t *p8est, p8est_ghost_t *ghost_layer, void *user_data, p8est_iter_volume_t iter_volume, p8est_iter_face_t iter_face, p8est_iter_edge_t iter_edge, p8est_iter_corner_t iter_corner, int remote)
p8est_iterate_ext adds the option remote: if this is false, then it is the same as p8est_iterate; if ...
void p8est_lid_add_inplace(p8est_lid_t *a, const p8est_lid_t *b)
Adds the p8est_lid_t b to the p8est_lid_t a.
void p8est_lid_set_zero(p8est_lid_t *input)
Initializes a linear index to zero.
void p8est_lid_set_one(p8est_lid_t *input)
Initializes a linear index to one.
p8est_t * p8est_source_ext(sc_io_source_t *src, sc_MPI_Comm mpicomm, size_t data_size, int load_data, int autopartition, int broadcasthead, void *user_pointer, p8est_connectivity_t **connectivity)
The same as p8est_load_ext, but reading the connectivity/p8est from an open sc_io_source_t stream.
void p8est_lid_copy(const p8est_lid_t *input, p8est_lid_t *output)
Copies an initialized p8est_lid_t to a p8est_lid_t.
void(* p8est_replace_t)(p8est_t *p8est, p4est_topidx_t which_tree, int num_outgoing, p8est_quadrant_t *outgoing[], int num_incoming, p8est_quadrant_t *incoming[])
Callback function prototype to replace one set of quadrants with another.
Definition: p8est_extended.h:114
void p8est_lid_sub_inplace(p8est_lid_t *a, const p8est_lid_t *b)
Subtracts the uint128_t b from the uint128_t a.
void p8est_lid_add(const p8est_lid_t *a, const p8est_lid_t *b, p8est_lid_t *result)
Adds the uint128_t b to the uint128_t a.
void p8est_lid_bitwise_neg(const p8est_lid_t *a, p8est_lid_t *result)
Calculates the bitwise negation of the uint128_t a.
void p8est_balance_ext(p8est_t *p8est, p8est_connect_type_t btype, p8est_init_t init_fn, p8est_replace_t replace_fn)
2:1 balance the size differences of neighboring elements in a forest.
int p8est_lid_compare(const p8est_lid_t *a, const p8est_lid_t *b)
Compare the p8est_lid_t a and the p8est_lid_t b.
p8est_t * p8est_load_ext(const char *filename, sc_MPI_Comm mpicomm, size_t data_size, int load_data, int autopartition, int broadcasthead, void *user_pointer, p8est_connectivity_t **connectivity)
Load the complete connectivity/p4est structure from disk.
p8est_t * p8est_copy_ext(p8est_t *input, int copy_data, int duplicate_mpicomm)
Make a deep copy of a p8est.
void p8est_lid_bitwise_or(const p8est_lid_t *a, const p8est_lid_t *b, p8est_lid_t *result)
Calculates the bitwise or of the uint128_t a and b.
Provide functions to serialize/deserialize a forest.
Iteration over mesh topology via callbacks.
void(* p8est_iter_edge_t)(p8est_iter_edge_info_t *info, void *user_data)
The prototype for a function that p8est_iterate will execute wherever the edge is an edge of all quad...
Definition: p8est_iterate.h:217
void(* p8est_iter_corner_t)(p8est_iter_corner_info_t *info, void *user_data)
The prototype for a function that p8est_iterate will execute wherever the corner is a corner for all ...
Definition: p8est_iterate.h:273
void(* p8est_iter_volume_t)(p8est_iter_volume_info_t *info, void *user_data)
The prototype for a function that p8est_iterate() will execute at every quadrant local to the current...
Definition: p8est_iterate.h:62
void(* p8est_iter_face_t)(p8est_iter_face_info_t *info, void *user_data)
The prototype for a function that p8est_iterate() will execute wherever two quadrants share a face: t...
Definition: p8est_iterate.h:136
Forest topology in a conventional mesh format.
This structure holds the 3D inter-tree connectivity information.
Definition: p8est_connectivity.h:215
quadrants that neighbor the local domain
Definition: p8est_ghost.h:41
Definition: p8est_extended.h:67
int use_balance_verify
Verify sc_ranges and/or sc_notify as applicable.
Definition: p8est_extended.h:75
double balance_ranges
time spent in sc_ranges
Definition: p8est_extended.h:88
double balance_notify_allgather
time spent in sc_notify_allgather
Definition: p8est_extended.h:91
int use_balance_ranges_notify
If true, call both sc_ranges and sc_notify and verify consistency.
Definition: p8est_extended.h:73
int use_balance_ranges
Use sc_ranges to determine the asymmetric communication pattern.
Definition: p8est_extended.h:70
double balance_notify
time spent in sc_notify
Definition: p8est_extended.h:89
int balance_max_ranges
If positive and smaller than p8est_num ranges, overrides it.
Definition: p8est_extended.h:77
Store a parallel numbering of Lobatto points of a given degree > 0.
Definition: p8est_lnodes.h:112
This structure contains complete mesh information on a 2:1 balanced forest.
Definition: p8est_mesh.h:170
The 3D quadrant (i.e., octant) datatype.
Definition: p8est.h:68
The p8est forest datatype.
Definition: p8est.h:132