p4est  2.8.643-dbc7-dirty
p4est is a software library for parallel adaptive mesh refinement.
p4est_bits.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 
32 #ifndef P4EST_BITS_H
33 #define P4EST_BITS_H
34 
35 #include <p4est.h>
36 #include <sc_random.h>
37 
38 SC_EXTERN_C_BEGIN;
39 
44 
50 void p4est_quadrant_print (int log_priority,
51  const p4est_quadrant_t * q);
52 
57  const p4est_quadrant_t * q2);
58 
65  p4est_quadrant_t * copy);
66 
72  const p4est_quadrant_t * q2);
73 
79  q1,
80  const p4est_quadrant_t *
81  q2);
82 
89 int p4est_quadrant_compare (const void *v1, const void *v2);
90 
101  const p4est_qcoord_t v2[]);
102 
109 int p4est_quadrant_disjoint (const void *v1, const void *v2);
110 
117 int p4est_quadrant_compare_piggy (const void *v1,
118  const void *v2);
119 
126  const void *v2);
127 
131 int p4est_quadrant_equal_fn (const void *v1, const void *v2,
132  const void *u);
133 
136 unsigned p4est_quadrant_hash_fn (const void *v, const void *u);
137 
143 int p4est_node_equal_piggy_fn (const void *v1,
144  const void *v2, const void *u);
145 
150 unsigned p4est_node_hash_piggy_fn (const void *v, const void *u);
151 
157  p4est_quadrant_t * r);
158 
163 
170  int level, p4est_quadrant_t * q);
171 
177  const p4est_quadrant_t * n);
178 
184  int level);
185 
190 
196  coord[]);
197 
203  q);
204 
210 
216  q);
217 
223  q);
224 
232  int inside);
233 
240  int level);
241 
247 
253 
260  const p4est_quadrant_t * q2);
261 
267  const p4est_quadrant_t * q2);
268 
272  const p4est_quadrant_t * q1,
273  const p4est_quadrant_t * q2,
274  const p4est_quadrant_t * q3);
275 
280 
285 
292  const p4est_quadrant_t * r);
293 
299  const p4est_quadrant_t * r);
300 
307  const p4est_quadrant_t * r);
308 
315  const p4est_quadrant_t * r);
316 
324  const p4est_quadrant_t * r);
325 
331  const p4est_quadrant_t * r);
332 
336  const p4est_quadrant_t * q);
337 
341  const p4est_quadrant_t *
342  q);
343 
354  const p4est_quadrant_t * l,
355  const p4est_quadrant_t * a);
356 
366  p4est_quadrant_t * q);
367 
377  p4est_quadrant_t * q);
378 
388  int level, p4est_quadrant_t * r);
389 
399  p4est_quadrant_t * r);
400 
408  p4est_quadrant_t * r,
409  int sibling_id);
410 
418  p4est_quadrant_t * r, int child_id);
419 
427  int face,
428  p4est_quadrant_t * r);
429 
447  * q, p4est_topidx_t t,
448  int face,
449  p4est_quadrant_t * r,
450  int *nface,
452  conn);
453 
468  * q, int face,
469  p4est_quadrant_t n[],
471  nur[]);
472 
491  * q, int face,
492  p4est_quadrant_t n[]);
493 
501  q, int corner,
502  p4est_quadrant_t * r);
503 
520  q, p4est_locidx_t t,
521  int corner,
522  sc_array_t * quads,
523  sc_array_t *
524  treeids,
525  sc_array_t *
526  ncorners,
528  * conn);
529 
537  p4est_quadrant_t * q,
538  int corner,
540  r);
541 
549  int corner,
550  p4est_quadrant_t * r);
551 
559  p4est_quadrant_t * c0,
560  p4est_quadrant_t * c1,
561  p4est_quadrant_t * c2,
562  p4est_quadrant_t * c3);
563 
571  p4est_quadrant_t c[]);
572 
580  p4est_quadrant_t * c[]);
581 
588  q, p4est_quadrant_t * fd,
589  int level);
590 
597  q, p4est_quadrant_t * ld,
598  int level);
599 
608  q, p4est_quadrant_t * r,
609  int c, int level);
610 
620  q1,
621  const p4est_quadrant_t *
622  q2, p4est_quadrant_t * r);
623 
629  q1,
630  const p4est_quadrant_t *
631  q2,
632  p4est_quadrant_t * r);
633 
645  p4est_quadrant_t * r,
646  const int ftransform[]);
647 
658  coords_in[],
660  coords_out[],
661  const int ftransform[]);
662 
666  int corner, int inside);
667 
674  int icorner, int inside);
675 
682  p4est_quadrant_t * r,
683  int corner);
684 
701  quadrant, int level);
702 
711  int level, uint64_t id);
712 
720  quadrant,
721  p4est_quadrant_t * result);
722 
730  quadrant,
731  p4est_quadrant_t * result);
732 
739  sc_rand_state_t * rstate);
740 
752  (const p4est_neighbor_transform_t * nt,
753  const p4est_quadrant_t * self_quad, p4est_quadrant_t * neigh_quad);
754 
766  (const p4est_neighbor_transform_t * nt,
767  const p4est_quadrant_t * neigh_quad, p4est_quadrant_t * self_quad);
768 SC_EXTERN_C_END;
769 
770 #endif /* !P4EST_BITS_H */
The top-level 2D p4est interface.
int32_t p4est_qcoord_t
Typedef for quadrant coordinates.
Definition: p4est_base.h:81
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
void p4est_quadrant_half_corner_neighbor(const p4est_quadrant_t *q, int corner, p4est_quadrant_t *r)
Compute the half size corner neighbor of a quadrant.
void p4est_nearest_common_ancestor_D(const p4est_quadrant_t *q1, const p4est_quadrant_t *q2, p4est_quadrant_t *r)
Computes the nearest common ancestor of two quadrants in the same tree.
int p4est_quadrant_is_next(const p4est_quadrant_t *q, const p4est_quadrant_t *r)
Test if two quadrants follow each other in the tree with no holes.
void p4est_quadrant_half_face_neighbors(const p4est_quadrant_t *q, int face, p4est_quadrant_t n[], p4est_quadrant_t nur[])
Get the smaller face neighbors of q.
int p4est_coordinates_compare(const p4est_qcoord_t v1[], const p4est_qcoord_t v2[])
Compare two sets of coordinates in their Morton ordering.
void p4est_quadrant_last_descendant(const p4est_quadrant_t *q, p4est_quadrant_t *ld, int level)
Compute the last descendant of a quadrant on a given level.
int p4est_node_equal_piggy_fn(const void *v1, const void *v2, const void *u)
Test if two nodes are in the same tree and have equal Morton indices.
int p4est_quadrant_is_node(const p4est_quadrant_t *q, int inside)
Test if a quadrant is used to represent a mesh node.
void p4est_quadrant_set_morton(p4est_quadrant_t *quadrant, int level, uint64_t id)
Set quadrant Morton indices based on linear position in uniform grid.
int p4est_quadrant_overlaps_tree(p4est_tree_t *tree, const p4est_quadrant_t *q)
Test if a quadrant has at least partial overlap with a tree.
unsigned p4est_node_hash_piggy_fn(const void *v, const void *u)
Compute hash value of a node based on its tree and Morton index.
int p4est_quadrant_is_familyv(const p4est_quadrant_t q[])
Test if 4 quadrants are siblings in Morton ordering, array version.
int p4est_quadrant_is_ancestor(const p4est_quadrant_t *q, const p4est_quadrant_t *r)
Test if a quadrant is an ancestor of another quadrant.
int p4est_quadrant_child_id(const p4est_quadrant_t *q)
Compute the position of this child within its siblings.
unsigned p4est_quadrant_hash_fn(const void *v, const void *u)
Computes a hash value for a quadrant by the lookup3 method.
void p4est_quadrant_face_neighbor(const p4est_quadrant_t *q, int face, p4est_quadrant_t *r)
Compute the face neighbor of a quadrant.
void p4est_quadrant_pad(p4est_quadrant_t *q)
Write -1 into the pad8 and pad16 members of a quadrant.
int p4est_quadrant_compare_piggy(const void *v1, const void *v2)
Compare two quadrants in their Morton ordering and the which_tree member.
int p4est_quadrant_disjoint(const void *v1, const void *v2)
Compare two quadrants in their Morton ordering, with equivalence if the two quadrants overlap.
void p4est_nearest_common_ancestor(const p4est_quadrant_t *q1, const p4est_quadrant_t *q2, p4est_quadrant_t *r)
Computes the nearest common ancestor of two quadrants in the same tree.
int p4est_quadrant_is_inside_tree(p4est_tree_t *tree, const p4est_quadrant_t *q)
Test if a quadrant is completely contained within a tree.
void p4est_quadrant_corner_neighbor(const p4est_quadrant_t *q, int corner, p4est_quadrant_t *r)
Compute the corner neighbor of a quadrant.
int p4est_quadrant_is_outside_face(const p4est_quadrant_t *q)
Test if a quadrant is outside a tree face boundary (no corner).
int p4est_quadrant_is_outside_corner(const p4est_quadrant_t *q)
Test if a quadrant is outside a tree corner boundary.
void p4est_neighbor_transform_quadrant(const p4est_neighbor_transform_t *nt, const p4est_quadrant_t *self_quad, p4est_quadrant_t *neigh_quad)
Transform a quadrant from self's coordinate system to neighbor's coordinate system.
int p4est_quadrant_compare_local_num(const void *v1, const void *v2)
Compare two quadrants with respect to their local_num in the piggy3 member.
void p4est_quadrant_first_descendant(const p4est_quadrant_t *q, p4est_quadrant_t *fd, int level)
Compute the first descendant of a quadrant on a given level.
int p4est_coordinates_is_inside_root(const p4est_qcoord_t coord[])
Test if Morton indices are inside the unit tree.
void p4est_node_unclamp(p4est_quadrant_t *n)
Move a clamped node out on the border.
void p4est_node_to_quadrant(const p4est_quadrant_t *n, int level, p4est_quadrant_t *q)
Find the enclosing quadrant of a given node at a given level.
int p4est_quadrant_is_ancestor_D(const p4est_quadrant_t *q, const p4est_quadrant_t *r)
Test if a quadrant is an ancestor of another quadrant.
int p4est_quadrant_is_equal(const p4est_quadrant_t *q1, const p4est_quadrant_t *q2)
Test if two quadrants have equal Morton indices.
void p4est_quadrant_children(const p4est_quadrant_t *q, p4est_quadrant_t *c0, p4est_quadrant_t *c1, p4est_quadrant_t *c2, p4est_quadrant_t *c3)
Compute the 4 children of a quadrant.
void p4est_quadrant_shift_corner(const p4est_quadrant_t *q, p4est_quadrant_t *r, int corner)
Shifts a quadrant until it touches the specified corner from the inside.
void p4est_quadrant_ancestor(const p4est_quadrant_t *q, int level, p4est_quadrant_t *r)
Compute the ancestor of a quadrant at a given level.
void p4est_quadrant_enlarge_first(const p4est_quadrant_t *a, p4est_quadrant_t *q)
Enlarge a quadrant as long as its first corner stays the same.
void p4est_quadrant_successor(const p4est_quadrant_t *quadrant, p4est_quadrant_t *result)
Compute the successor according to the Morton index in a uniform mesh.
int p4est_coordinates_is_valid(const p4est_qcoord_t coord[], int level)
Test if Morton indices are valid and are inside the unit tree.
int p4est_quadrant_is_valid(const p4est_quadrant_t *q)
Test if a quadrant has valid Morton indices and is inside the unit tree.
int p4est_quadrant_is_family(const p4est_quadrant_t *q0, const p4est_quadrant_t *q1, const p4est_quadrant_t *q2, const p4est_quadrant_t *q3)
Test if 4 quadrants are siblings in Morton ordering.
int p4est_quadrant_contains_node(const p4est_quadrant_t *q, const p4est_quadrant_t *n)
Decide if a node is completely contained within a quadrant.
void p4est_quadrant_enlarge_last(const p4est_quadrant_t *a, p4est_quadrant_t *q)
Enlarge a quadrant as long as its last corner stays the same.
int p4est_quadrant_is_next_D(const p4est_quadrant_t *q, const p4est_quadrant_t *r)
Test if two quadrants follow each other in the tree with no holes.
void p4est_quadrant_copy(const p4est_quadrant_t *q, p4est_quadrant_t *copy)
Copy the Morton indices of the quadrant q.
uint64_t p4est_quadrant_linear_id(const p4est_quadrant_t *quadrant, int level)
Computes the linear position of a quadrant in a uniform grid.
void p4est_neighbor_transform_quadrant_reverse(const p4est_neighbor_transform_t *nt, const p4est_quadrant_t *neigh_quad, p4est_quadrant_t *self_quad)
Transform a quadrant from a neighbors's coordinate system to self's coordinate system.
int p4est_quadrant_is_parent(const p4est_quadrant_t *q, const p4est_quadrant_t *r)
Test if a quadrant is the parent of another quadrant.
int p4est_quadrant_equal_fn(const void *v1, const void *v2, const void *u)
Test if two quadrants have equal Morton indices, callback version.
p4est_topidx_t p4est_quadrant_face_neighbor_extra(const p4est_quadrant_t *q, p4est_topidx_t t, int face, p4est_quadrant_t *r, int *nface, p4est_connectivity_t *conn)
Compute the face neighbor of a quadrant, transforming across tree boundaries if necessary.
int p4est_quadrant_is_inside_root(const p4est_quadrant_t *q)
Test if a quadrant is inside the unit tree.
void p4est_quadrant_all_face_neighbors(const p4est_quadrant_t *q, int face, p4est_quadrant_t n[])
Create all possible face neighbors of q.
void p4est_quadrant_parent(const p4est_quadrant_t *q, p4est_quadrant_t *r)
Compute the parent of a quadrant.
int p4est_quadrant_overlaps(const p4est_quadrant_t *q1, const p4est_quadrant_t *q2)
Test if two quadrants overlap.
void p4est_quadrant_corner_descendant(const p4est_quadrant_t *q, p4est_quadrant_t *r, int c, int level)
Compute the descendant of a quadrant touching a given corner.
void p4est_quadrant_transform_corner(p4est_quadrant_t *q, int icorner, int inside)
Move a quadrant inside or diagonally outside a corner position.
void p4est_quadrant_sibling(const p4est_quadrant_t *q, p4est_quadrant_t *r, int sibling_id)
Compute a specific sibling of a quadrant.
void p4est_quadrant_child(const p4est_quadrant_t *q, p4est_quadrant_t *r, int child_id)
Compute a specific child of a quadrant.
void p4est_quadrant_corner_neighbor_extra(const p4est_quadrant_t *q, p4est_locidx_t t, int corner, sc_array_t *quads, sc_array_t *treeids, sc_array_t *ncorners, p4est_connectivity_t *conn)
Compute the corner neighbors of a quadrant, transforming across tree boundaries if necessary.
void p4est_quadrant_srand(const p4est_quadrant_t *q, sc_rand_state_t *rstate)
Initialize a random number generator by quadrant coordinates.
void p4est_quadrant_childrenv(const p4est_quadrant_t *q, p4est_quadrant_t c[])
Compute the 4 children of a quadrant, array version.
int p4est_quadrant_is_parent_D(const p4est_quadrant_t *q, const p4est_quadrant_t *r)
Test if a quadrant is the parent of another quadrant.
void p4est_node_clamp_inside(const p4est_quadrant_t *n, p4est_quadrant_t *r)
Clamp a node inside the unit tree if it sits on a high border.
int p4est_quadrant_compare(const void *v1, const void *v2)
Compare two quadrants in their Morton ordering.
int p4est_quadrant_is_sibling(const p4est_quadrant_t *q1, const p4est_quadrant_t *q2)
Test if two quadrants are siblings.
void p4est_quadrant_corner_node(const p4est_quadrant_t *q, int corner, p4est_quadrant_t *r)
Compute the corner node of a quadrant.
void p4est_quadrant_print(int log_priority, const p4est_quadrant_t *q)
Prints one line with quadrant's x, y and level.
int p4est_quadrant_is_familypv(p4est_quadrant_t *q[])
Test if 4 quadrants are siblings in Morton ordering, array version.
int p4est_quadrant_is_inside_3x3(const p4est_quadrant_t *q)
Test if a quadrant is inside the 3x3 box around the root tree.
void p4est_quadrant_predecessor(const p4est_quadrant_t *quadrant, p4est_quadrant_t *result)
Compute the predecessor according to the Morton index in a uniform mesh.
int p4est_quadrant_is_sibling_D(const p4est_quadrant_t *q1, const p4est_quadrant_t *q2)
Test if two quadrants are siblings.
void p4est_quadrant_transform_face(const p4est_quadrant_t *q, p4est_quadrant_t *r, const int ftransform[])
Transforms a quadrant/node across a face between trees.
void p4est_coordinates_transform_face(const p4est_qcoord_t coords_in[], p4est_qcoord_t coords_out[], const int ftransform[])
Transforms coordinates across a face between trees.
int p4est_quadrant_is_extended(const p4est_quadrant_t *q)
Test if a quadrant has valid Morton indices in the 3x3 box around root.
int p4est_quadrant_is_equal_piggy(const p4est_quadrant_t *q1, const p4est_quadrant_t *q2)
Test if two quadrants have equal Morton indices and the same tree id.
int p4est_quadrant_ancestor_id(const p4est_quadrant_t *q, int level)
Compute the position of the ancestor of this child at level level within its siblings.
int p4est_quadrant_is_first_last(const p4est_quadrant_t *f, const p4est_quadrant_t *l, const p4est_quadrant_t *a)
Whether two descendants of a quadrant are first and last, up to size.
int p4est_quadrant_touches_corner(const p4est_quadrant_t *q, int corner, int inside)
Checks if a quadrant touches a corner (diagonally inside or outside).
void p4est_quadrant_childrenpv(const p4est_quadrant_t *q, p4est_quadrant_t *c[])
Compute the 4 children of a quadrant, array version.
uint64_t sc_rand_state_t
This structure holds the 2D inter-tree connectivity information.
Definition: p4est_connectivity.h:153
Generic interface for transformations beteen a tree and any of its neighbors.
Definition: p4est_connectivity.h:207
The 2D quadrant datatype.
Definition: p4est.h:72
The p4est tree datatype.
Definition: p4est.h:115