p4est_t * p4est_new_points (sc_MPI_Comm mpicomm, p4est_connectivity_t * connectivity, int maxlevel, p4est_quadrant_t * points, p4est_locidx_t num_points, p4est_locidx_t max_points, size_t data_size, p4est_init_t init_fn, void *user_pointer) { int mpiret; int num_procs, rank; int i, isizet; size_t lcount; size_t *nmemb; #ifdef P4EST_ENABLE_DEBUG size_t zz; #endif p4est_topidx_t jt, num_trees; p4est_topidx_t first_tree, last_tree, next_tree; p4est_quadrant_t *first_quad, *next_quad, *quad; p4est_quadrant_t a, b, c, f, l, n; p4est_tree_t *tree; p4est_t *p4est; p4est_points_state_t ppstate; P4EST_GLOBAL_PRODUCTIONF ("Into " P4EST_STRING "_new_points with max level %d max points %lld\n", maxlevel, (long long) max_points); p4est_log_indent_push (); P4EST_ASSERT (p4est_connectivity_is_valid (connectivity)); P4EST_ASSERT (max_points >= -1); /* retrieve MPI information */ mpiret = sc_MPI_Comm_size (mpicomm, &num_procs); SC_CHECK_MPI (mpiret); mpiret = sc_MPI_Comm_rank (mpicomm, &rank); SC_CHECK_MPI (mpiret); /* This implementation runs in O(P/p * maxlevel) * with P the total number of points, p the number of processors. * Two optimizations are possible: * 1. At startup remove points that lead to duplicate quadrants. * 2. Use complete_region between successive points instead of * the call to refine. This should give O(N/p) * maxlevel * with N the total number of quadrants. */ /* parallel sort the incoming points */ lcount = (size_t) num_points; nmemb = P4EST_ALLOC_ZERO (size_t, num_procs); isizet = (int) sizeof (size_t); mpiret = sc_MPI_Allgather (&lcount, isizet, sc_MPI_BYTE, nmemb, isizet, sc_MPI_BYTE, mpicomm); SC_CHECK_MPI (mpiret); sc_psort (mpicomm, points, nmemb, sizeof (p4est_quadrant_t), p4est_quadrant_compare_piggy); P4EST_FREE (nmemb); #ifdef P4EST_ENABLE_DEBUG first_quad = points; for (zz = 1; zz < lcount; ++zz) { next_quad = points + zz; P4EST_ASSERT (p4est_quadrant_compare_piggy (first_quad, next_quad) <= 0); first_quad = next_quad; } #endif /* create the p4est */ p4est = P4EST_ALLOC_ZERO (p4est_t, 1); ppstate.points = points; ppstate.num_points = num_points; ppstate.max_points = max_points; ppstate.current = 0; ppstate.maxlevel = maxlevel; /* assign some data members */ p4est->data_size = 2 * sizeof (p4est_locidx_t); /* temporary */ p4est->user_pointer = &ppstate; p4est->connectivity = connectivity; num_trees = connectivity->num_trees; /* create parallel environment */ p4est_comm_parallel_env_create (p4est, mpicomm); /* allocate memory pools */ p4est->user_data_pool = sc_mempool_new (p4est->data_size); p4est->quadrant_pool = sc_mempool_new (sizeof (p4est_quadrant_t)); P4EST_GLOBAL_PRODUCTIONF ("New " P4EST_STRING " with %lld trees on %d processors\n", (long long) num_trees, num_procs); /* allocate trees */ p4est->trees = sc_array_new (sizeof (p4est_tree_t)); sc_array_resize (p4est->trees, num_trees); for (jt = 0; jt < num_trees; ++jt) { tree = p4est_tree_array_index (p4est->trees, jt); sc_array_init (&tree->quadrants, sizeof (p4est_quadrant_t)); P4EST_QUADRANT_INIT (&tree->first_desc); P4EST_QUADRANT_INIT (&tree->last_desc); tree->quadrants_offset = 0; for (i = 0; i <= P4EST_QMAXLEVEL; ++i) { tree->quadrants_per_level[i] = 0; } for (; i <= P4EST_MAXLEVEL; ++i) { tree->quadrants_per_level[i] = -1; } tree->maxlevel = 0; } p4est->local_num_quadrants = 0; p4est->global_num_quadrants = 0; /* create point based partition */ P4EST_QUADRANT_INIT (&f); p4est->global_first_position = P4EST_ALLOC_ZERO (p4est_quadrant_t, num_procs + 1); if (num_points == 0) { P4EST_VERBOSE ("Empty processor"); first_tree = p4est->first_local_tree = -1; first_quad = NULL; } else { /* we are probably not empty */ if (rank == 0) { first_tree = p4est->first_local_tree = 0; p4est_quadrant_set_morton (&f, maxlevel, 0); } else { first_tree = p4est->first_local_tree = points->p.which_tree; p4est_node_to_quadrant (points, maxlevel, &f); } first_quad = &f; } last_tree = p4est->last_local_tree = -2; p4est_comm_global_partition (p4est, first_quad); first_quad = p4est->global_first_position + rank; next_quad = p4est->global_first_position + (rank + 1); next_tree = next_quad->p.which_tree; if (first_tree >= 0 && p4est_quadrant_is_equal (first_quad, next_quad) && first_quad->p.which_tree == next_quad->p.which_tree) { /* if all our points are consumed by the next processor we are empty */ first_tree = p4est->first_local_tree = -1; } if (first_tree >= 0) { /* we are definitely not empty */ if (next_quad->x == 0 && next_quad->y == 0 #ifdef P4_TO_P8 && next_quad->z == 0 #endif ) { last_tree = p4est->last_local_tree = next_tree - 1; } else { last_tree = p4est->last_local_tree = next_tree; } P4EST_ASSERT (first_tree <= last_tree); } /* fill the local trees */ P4EST_QUADRANT_INIT (&a); P4EST_QUADRANT_INIT (&b); P4EST_QUADRANT_INIT (&c); P4EST_QUADRANT_INIT (&l); n = *next_quad; n.level = (int8_t) maxlevel; for (jt = first_tree; jt <= last_tree; ++jt) { int onlyone = 0; int includeb = 0; tree = p4est_tree_array_index (p4est->trees, jt); /* determine first local quadrant of this tree */ if (jt == first_tree) { a = *first_quad; a.level = (int8_t) maxlevel; first_quad = next_quad = NULL; /* free to use further down */ P4EST_ASSERT (p4est_quadrant_is_valid (&a)); } else { p4est_quadrant_set_morton (&a, maxlevel, 0); P4EST_ASSERT (jt < next_tree || p4est_quadrant_compare (&a, &n) < 0); } /* enlarge first local quadrant if possible */ if (jt < next_tree) { while (p4est_quadrant_child_id (&a) == 0 && a.level > 0) { p4est_quadrant_parent (&a, &a); } P4EST_ASSERT (jt == first_tree || a.level == 0); } else { for (c = a; p4est_quadrant_child_id (&c) == 0; a = c) { p4est_quadrant_parent (&c, &c); p4est_quadrant_last_descendant (&c, &l, maxlevel); if (p4est_quadrant_compare (&l, &n) >= 0) { break; } } P4EST_ASSERT (a.level > 0); P4EST_ASSERT ((p4est_quadrant_last_descendant (&a, &l, maxlevel), p4est_quadrant_compare (&l, &n) < 0)); } p4est_quadrant_first_descendant (&a, &tree->first_desc, P4EST_QMAXLEVEL); /* determine largest possible last quadrant of this tree */ if (jt < next_tree) { p4est_quadrant_last_descendant (&a, &l, maxlevel); p4est_quadrant_set_morton (&b, 0, 0); p4est_quadrant_last_descendant (&b, &b, maxlevel); if (p4est_quadrant_is_equal (&l, &b)) { onlyone = 1; } else { includeb = 1; for (c = b; p4est_quadrant_child_id (&c) == P4EST_CHILDREN - 1; b = c) { p4est_quadrant_parent (&c, &c); p4est_quadrant_first_descendant (&c, &f, maxlevel); if (p4est_quadrant_compare (&l, &f) >= 0) { break; } } } } else { b = n; } /* create a complete tree */ if (onlyone) { quad = p4est_quadrant_array_push (&tree->quadrants); *quad = a; p4est_quadrant_init_data (p4est, jt, quad, p4est_points_init); tree->maxlevel = a.level; ++tree->quadrants_per_level[a.level]; } else { p4est_complete_region (p4est, &a, 1, &b, includeb, tree, jt, p4est_points_init); quad = p4est_quadrant_array_index (&tree->quadrants, tree->quadrants.elem_count - 1); } tree->quadrants_offset = p4est->local_num_quadrants; p4est->local_num_quadrants += tree->quadrants.elem_count; p4est_quadrant_last_descendant (quad, &tree->last_desc, P4EST_QMAXLEVEL); /* verification */ #ifdef P4EST_ENABLE_DEBUG first_quad = p4est_quadrant_array_index (&tree->quadrants, 0); for (zz = 1; zz < tree->quadrants.elem_count; ++zz) { next_quad = p4est_quadrant_array_index (&tree->quadrants, zz); P4EST_ASSERT (((p4est_locidx_t *) first_quad->p.user_data)[1] == ((p4est_locidx_t *) next_quad->p.user_data)[0]); first_quad = next_quad; } #endif } if (last_tree >= 0) { for (; jt < num_trees; ++jt) { tree = p4est_tree_array_index (p4est->trees, jt); tree->quadrants_offset = p4est->local_num_quadrants; } } /* compute some member variables */ p4est->global_first_quadrant = P4EST_ALLOC (p4est_gloidx_t, num_procs + 1); p4est_comm_count_quadrants (p4est); /* print more statistics */ P4EST_VERBOSEF ("total local quadrants %lld\n", (long long) p4est->local_num_quadrants); P4EST_ASSERT (p4est_is_valid (p4est)); p4est_log_indent_pop (); P4EST_GLOBAL_PRODUCTIONF ("Done " P4EST_STRING "_new_points with %lld total quadrants\n", (long long) p4est->global_num_quadrants); /* refine to have one point per quadrant */ if (max_points >= 0) { p4est_refine_ext (p4est, 1, maxlevel, p4est_points_refine, p4est_points_init, NULL); #ifdef P4EST_ENABLE_DEBUG for (jt = first_tree; jt <= last_tree; ++jt) { tree = p4est_tree_array_index (p4est->trees, jt); first_quad = p4est_quadrant_array_index (&tree->quadrants, 0); for (zz = 1; zz < tree->quadrants.elem_count; ++zz) { next_quad = p4est_quadrant_array_index (&tree->quadrants, zz); P4EST_ASSERT (((p4est_locidx_t *) first_quad->p.user_data)[1] == ((p4est_locidx_t *) next_quad->p.user_data)[0]); first_quad = next_quad; } } #endif } /* initialize user pointer and data size */ p4est_reset_data (p4est, data_size, init_fn, user_pointer); return p4est; }
int main (int argc, char **argv) { sc_MPI_Comm mpicomm; int mpiret; int mpisize, mpirank; unsigned crc; #ifndef P4_TO_P8 size_t kz; int8_t l; p4est_quadrant_t *q; p4est_tree_t stree, *tree = &stree; #endif p4est_t *p4est; p4est_connectivity_t *connectivity; /* initialize MPI */ mpiret = sc_MPI_Init (&argc, &argv); SC_CHECK_MPI (mpiret); mpicomm = sc_MPI_COMM_WORLD; mpiret = sc_MPI_Comm_size (mpicomm, &mpisize); SC_CHECK_MPI (mpiret); mpiret = sc_MPI_Comm_rank (mpicomm, &mpirank); SC_CHECK_MPI (mpiret); sc_init (mpicomm, 1, 1, NULL, SC_LP_DEFAULT); p4est_init (NULL, SC_LP_DEFAULT); #ifndef P4_TO_P8 connectivity = p4est_connectivity_new_star (); #else connectivity = p8est_connectivity_new_rotcubes (); #endif p4est = p4est_new_ext (mpicomm, connectivity, 0, 0, 0, 4, NULL, NULL); #ifndef P4_TO_P8 /* build empty tree */ sc_array_init (&tree->quadrants, sizeof (p4est_quadrant_t)); for (l = 0; l <= P4EST_MAXLEVEL; ++l) { tree->quadrants_per_level[l] = 0; } tree->maxlevel = 0; /* insert two quadrants */ sc_array_resize (&tree->quadrants, 4); q = p4est_quadrant_array_index (&tree->quadrants, 0); p4est_quadrant_set_morton (q, 3, 13); q = p4est_quadrant_array_index (&tree->quadrants, 1); p4est_quadrant_set_morton (q, 1, 1); q = p4est_quadrant_array_index (&tree->quadrants, 2); p4est_quadrant_set_morton (q, 1, 2); q = p4est_quadrant_array_index (&tree->quadrants, 3); p4est_quadrant_set_morton (q, 1, 3); for (kz = 0; kz < tree->quadrants.elem_count; ++kz) { q = p4est_quadrant_array_index (&tree->quadrants, kz); q->p.user_data = sc_mempool_alloc (p4est->user_data_pool); ++tree->quadrants_per_level[q->level]; tree->maxlevel = (int8_t) SC_MAX (tree->maxlevel, q->level); } /* balance the tree, print and destroy */ #if 0 p4est_balance_subtree (p4est, P4EST_CONNECT_FULL, 0, NULL); p4est_tree_print (SC_LP_INFO, tree); #endif for (kz = 0; kz < tree->quadrants.elem_count; ++kz) { q = p4est_quadrant_array_index (&tree->quadrants, kz); sc_mempool_free (p4est->user_data_pool, q->p.user_data); } sc_array_reset (&tree->quadrants); #endif /* !P4_TO_P8 */ /* check reset data function */ p4est_reset_data (p4est, 0, init_fn, NULL); p4est_reset_data (p4est, 0, NULL, NULL); /* refine and balance the forest */ SC_CHECK_ABORT (p4est_is_balanced (p4est, P4EST_CONNECT_FULL), "Balance 1"); p4est_refine (p4est, 1, refine_fn, NULL); SC_CHECK_ABORT (!p4est_is_balanced (p4est, P4EST_CONNECT_FULL), "Balance 2"); p4est_balance (p4est, P4EST_CONNECT_FULL, NULL); SC_CHECK_ABORT (p4est_is_balanced (p4est, P4EST_CONNECT_FULL), "Balance 3"); /* check reset data function */ p4est_reset_data (p4est, 17, NULL, NULL); p4est_reset_data (p4est, 8, init_fn, NULL); /* checksum and partition */ crc = p4est_checksum (p4est); p4est_partition (p4est, 0, NULL); SC_CHECK_ABORT (p4est_checksum (p4est) == crc, "Partition"); SC_CHECK_ABORT (p4est_is_balanced (p4est, P4EST_CONNECT_FULL), "Balance 4"); /* check reset data function */ p4est_reset_data (p4est, 3, NULL, NULL); p4est_reset_data (p4est, 3, NULL, NULL); /* checksum and rebalance */ crc = p4est_checksum (p4est); p4est_balance (p4est, P4EST_CONNECT_FULL, NULL); SC_CHECK_ABORT (p4est_checksum (p4est) == crc, "Rebalance"); /* clean up and exit */ P4EST_ASSERT (p4est->user_data_pool->elem_count == (size_t) p4est->local_num_quadrants); p4est_destroy (p4est); p4est_connectivity_destroy (connectivity); sc_finalize (); mpiret = sc_MPI_Finalize (); SC_CHECK_MPI (mpiret); return 0; }
static void test_periodic (p8est_connectivity_t * conn) { int i; int iface, nface; int iedge, icorner; int ft[9]; size_t zz; p4est_topidx_t itree, ntree; p8est_edge_info_t ei; p8est_edge_transform_t *et; p8est_corner_info_t ci; p8est_corner_transform_t *ct; sc_array_t *eta, *cta; itree = 0; for (iface = 0; iface < 6; ++iface) { ntree = p8est_find_face_transform (conn, itree, iface, ft); nface = p8est_face_dual[iface]; SC_CHECK_ABORT (ntree == itree, "PF tree"); for (i = 0; i < 3; ++i) { SC_CHECK_ABORTF (ft[i] == ft[i + 3], "PF axis %d", i); } SC_CHECK_ABORT (ft[6] == 0 && ft[7] == 0, "PF reverse"); SC_CHECK_ABORT (ft[8] == 2 * (iface % 2) + nface % 2, "PF code"); } eta = &ei.edge_transforms; sc_array_init (eta, sizeof (p8est_edge_transform_t)); for (iedge = 0; iedge < 12; ++iedge) { p8est_find_edge_transform (conn, itree, iedge, &ei); SC_CHECK_ABORT ((int) ei.iedge == iedge, "PE ei"); SC_CHECK_ABORT (eta->elem_count == 1, "PE count"); for (zz = 0; zz < eta->elem_count; ++zz) { et = p8est_edge_array_index (eta, zz); SC_CHECK_ABORT (et->ntree == itree, "PE tree"); SC_CHECK_ABORT ((int) et->nedge + iedge == 8 * (iedge / 4) + 3, "PE edge"); SC_CHECK_ABORT (et->nflip == 0, "PE flip"); SC_CHECK_ABORT (et->corners == et->nedge % 4, "PE corners"); SC_CHECK_ABORT ((int) et->naxis[0] == iedge / 4 && et->naxis[1] < et->naxis[2] && et->naxis[0] + et->naxis[1] + et->naxis[2] == 3, "PE axis"); } } sc_array_reset (eta); cta = &ci.corner_transforms; sc_array_init (cta, sizeof (p8est_corner_transform_t)); for (icorner = 0; icorner < 8; ++icorner) { p8est_find_corner_transform (conn, itree, icorner, &ci); SC_CHECK_ABORT ((int) ci.icorner == icorner, "PC ci"); SC_CHECK_ABORT (cta->elem_count == 1, "PC count"); for (zz = 0; zz < cta->elem_count; ++zz) { ct = p8est_corner_array_index (cta, zz); SC_CHECK_ABORT (ct->ntree == itree, "PC tree"); SC_CHECK_ABORT (ct->ncorner + icorner == 7, "PC corner"); } } sc_array_reset (cta); }
int main (int argc, char **argv) { const p4est_qcoord_t qone = 1; int mpiret; int k; int level, mid, cid; int id0, id1, id2, id3; int64_t index1, index2; size_t iz, jz, incount; p4est_qcoord_t mh = P4EST_QUADRANT_LEN (P4EST_QMAXLEVEL); p4est_connectivity_t *connectivity; p4est_t *p4est1; p4est_t *p4est2; p4est_tree_t *t1, *t2, tree; p4est_quadrant_t *p, *q1, *q2; p4est_quadrant_t r, s; p4est_quadrant_t c0, c1, c2, c3; p4est_quadrant_t cv[P4EST_CHILDREN], *cp[P4EST_CHILDREN]; p4est_quadrant_t A, B, C, D, E, F, G, H, I, P, Q; p4est_quadrant_t a, f, g, h; uint64_t Aid, Fid; /* initialize MPI */ mpiret = sc_MPI_Init (&argc, &argv); SC_CHECK_MPI (mpiret); /* create connectivity and forest structures */ connectivity = p4est_connectivity_new_unitsquare (); p4est1 = p4est_new_ext (sc_MPI_COMM_SELF, connectivity, 15, 0, 0, 0, NULL, NULL); p4est2 = p4est_new_ext (sc_MPI_COMM_SELF, connectivity, 15, 0, 0, 8, NULL, NULL); /* refine the second tree to a uniform level */ p4est_refine (p4est1, 1, refine_none, NULL); p4est_refine (p4est2, 1, refine_some, NULL); t1 = p4est_tree_array_index (p4est1->trees, 0); t2 = p4est_tree_array_index (p4est2->trees, 0); SC_CHECK_ABORT (p4est_tree_is_sorted (t1), "is_sorted"); SC_CHECK_ABORT (p4est_tree_is_sorted (t2), "is_sorted"); /* run a bunch of cross-tests */ p = NULL; for (iz = 0; iz < t1->quadrants.elem_count; ++iz) { q1 = p4est_quadrant_array_index (&t1->quadrants, iz); /* test the index conversion */ index1 = p4est_quadrant_linear_id (q1, (int) q1->level); p4est_quadrant_set_morton (&r, (int) q1->level, index1); index2 = p4est_quadrant_linear_id (&r, (int) r.level); SC_CHECK_ABORT (index1 == index2, "index conversion"); level = (int) q1->level - 1; if (level >= 0) { index1 = p4est_quadrant_linear_id (q1, level); p4est_quadrant_set_morton (&r, level, index1); index2 = p4est_quadrant_linear_id (&r, level); SC_CHECK_ABORT (index1 == index2, "index conversion"); } /* test the is_next function */ if (p != NULL) { SC_CHECK_ABORT (p4est_quadrant_is_next (p, q1), "is_next"); } p = q1; /* test the is_family function */ p4est_quadrant_children (q1, &c0, &c1, &c2, &c3); SC_CHECK_ABORT (p4est_quadrant_is_family (&c0, &c1, &c2, &c3), "is_family"); SC_CHECK_ABORT (!p4est_quadrant_is_family (&c1, &c0, &c2, &c3), "is_family"); SC_CHECK_ABORT (!p4est_quadrant_is_family (&c0, &c0, &c1, &c2), "is_family"); p4est_quadrant_childrenv (q1, cv); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c0, &cv[0]), "is_family"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c1, &cv[1]), "is_family"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c2, &cv[2]), "is_family"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c3, &cv[3]), "is_family"); SC_CHECK_ABORT (p4est_quadrant_is_family (&cv[0], &cv[1], &cv[2], &cv[3]), "is_family"); cp[0] = &cv[0]; cp[1] = &cv[1]; cp[2] = &cv[2]; cp[3] = &cv[3]; SC_CHECK_ABORT (p4est_quadrant_is_familypv (cp), "is_family"); cv[1] = cv[0]; SC_CHECK_ABORT (!p4est_quadrant_is_familyv (cv), "is_family"); cp[1] = &c1; SC_CHECK_ABORT (p4est_quadrant_is_familypv (cp), "is_family"); cp[2] = &c3; SC_CHECK_ABORT (!p4est_quadrant_is_familypv (cp), "is_family"); /* test the sibling function */ mid = p4est_quadrant_child_id (q1); for (cid = 0; cid < 4; ++cid) { p4est_quadrant_sibling (q1, &r, cid); if (cid != mid) { SC_CHECK_ABORT (p4est_quadrant_is_sibling (q1, &r), "sibling"); } else { SC_CHECK_ABORT (p4est_quadrant_is_equal (q1, &r), "sibling"); } } /* test t1 against itself */ for (jz = 0; jz < t1->quadrants.elem_count; ++jz) { q2 = p4est_quadrant_array_index (&t1->quadrants, jz); /* test the comparison function */ SC_CHECK_ABORT (p4est_quadrant_compare (q1, q2) == -p4est_quadrant_compare (q2, q1), "compare"); SC_CHECK_ABORT ((p4est_quadrant_compare (q1, q2) == 0) == p4est_quadrant_is_equal (q1, q2), "is_equal"); /* test the descriptive versions of functions */ SC_CHECK_ABORT (p4est_quadrant_is_sibling_D (q1, q2) == p4est_quadrant_is_sibling (q1, q2), "is_sibling"); SC_CHECK_ABORT (p4est_quadrant_is_parent_D (q1, q2) == p4est_quadrant_is_parent (q1, q2), "is_parent"); SC_CHECK_ABORT (p4est_quadrant_is_parent_D (q2, q1) == p4est_quadrant_is_parent (q2, q1), "is_parent"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor_D (q1, q2) == p4est_quadrant_is_ancestor (q1, q2), "is_ancestor"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor_D (q2, q1) == p4est_quadrant_is_ancestor (q2, q1), "is_ancestor"); SC_CHECK_ABORT (p4est_quadrant_is_next_D (q1, q2) == p4est_quadrant_is_next (q1, q2), "is_next"); SC_CHECK_ABORT (p4est_quadrant_is_next_D (q2, q1) == p4est_quadrant_is_next (q2, q1), "is_next"); p4est_nearest_common_ancestor_D (q1, q2, &r); p4est_nearest_common_ancestor (q1, q2, &s); SC_CHECK_ABORT (p4est_quadrant_is_equal (&r, &s), "common_ancestor"); p4est_nearest_common_ancestor_D (q2, q1, &r); p4est_nearest_common_ancestor (q2, q1, &s); SC_CHECK_ABORT (p4est_quadrant_is_equal (&r, &s), "common_ancestor"); } /* test t1 against t2 */ for (jz = 0; jz < t2->quadrants.elem_count; ++jz) { q2 = p4est_quadrant_array_index (&t2->quadrants, jz); /* test the comparison function */ SC_CHECK_ABORT (p4est_quadrant_compare (q1, q2) == -p4est_quadrant_compare (q2, q1), "compare"); SC_CHECK_ABORT ((p4est_quadrant_compare (q1, q2) == 0) == p4est_quadrant_is_equal (q1, q2), "is_equal"); /* test the descriptive versions of functions */ SC_CHECK_ABORT (p4est_quadrant_is_sibling_D (q1, q2) == p4est_quadrant_is_sibling (q1, q2), "is_sibling"); SC_CHECK_ABORT (p4est_quadrant_is_parent_D (q1, q2) == p4est_quadrant_is_parent (q1, q2), "is_parent"); SC_CHECK_ABORT (p4est_quadrant_is_parent_D (q2, q1) == p4est_quadrant_is_parent (q2, q1), "is_parent"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor_D (q1, q2) == p4est_quadrant_is_ancestor (q1, q2), "is_ancestor"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor_D (q2, q1) == p4est_quadrant_is_ancestor (q2, q1), "is_ancestor"); SC_CHECK_ABORT (p4est_quadrant_is_next_D (q1, q2) == p4est_quadrant_is_next (q1, q2), "is_next"); SC_CHECK_ABORT (p4est_quadrant_is_next_D (q2, q1) == p4est_quadrant_is_next (q2, q1), "is_next"); p4est_nearest_common_ancestor_D (q1, q2, &r); p4est_nearest_common_ancestor (q1, q2, &s); SC_CHECK_ABORT (p4est_quadrant_is_equal (&r, &s), "common_ancestor"); p4est_nearest_common_ancestor_D (q2, q1, &r); p4est_nearest_common_ancestor (q2, q1, &s); SC_CHECK_ABORT (p4est_quadrant_is_equal (&r, &s), "common_ancestor"); } } p = NULL; for (iz = 0; iz < t2->quadrants.elem_count; ++iz) { q1 = p4est_quadrant_array_index (&t2->quadrants, iz); /* test the is_next function */ if (p != NULL) { SC_CHECK_ABORT (p4est_quadrant_is_next (p, q1), "is_next"); } p = q1; } /* test the coarsen function */ p4est_coarsen (p4est1, 1, coarsen_none, NULL); p4est_coarsen (p4est1, 1, coarsen_all, NULL); p4est_coarsen (p4est2, 1, coarsen_some, NULL); /* test the linearize algorithm */ incount = t2->quadrants.elem_count; (void) p4est_linearize_tree (p4est2, t2); SC_CHECK_ABORT (incount == t2->quadrants.elem_count, "linearize"); /* this is user_data neutral only when p4est1->data_size == 0 */ sc_array_init (&tree.quadrants, sizeof (p4est_quadrant_t)); sc_array_resize (&tree.quadrants, 18); q1 = p4est_quadrant_array_index (&tree.quadrants, 0); q2 = p4est_quadrant_array_index (&t2->quadrants, 0); *q1 = *q2; q2 = p4est_quadrant_array_index (&t2->quadrants, 1); for (k = 0; k < 3; ++k) { q1 = p4est_quadrant_array_index (&tree.quadrants, (size_t) (k + 1)); *q1 = *q2; q1->level = (int8_t) (q1->level + k); } for (k = 0; k < 10; ++k) { q1 = p4est_quadrant_array_index (&tree.quadrants, (size_t) (k + 4)); q2 = p4est_quadrant_array_index (&t2->quadrants, (size_t) (k + 3)); *q1 = *q2; q1->level = (int8_t) (q1->level + k); } for (k = 0; k < 4; ++k) { q1 = p4est_quadrant_array_index (&tree.quadrants, (size_t) (k + 14)); q2 = p4est_quadrant_array_index (&t2->quadrants, (size_t) (k + 12)); *q1 = *q2; q1->level = (int8_t) (q1->level + 10 + k); } tree.maxlevel = 0; for (k = 0; k <= P4EST_QMAXLEVEL; ++k) { tree.quadrants_per_level[k] = 0; } for (; k <= P4EST_MAXLEVEL; ++k) { tree.quadrants_per_level[k] = -1; } incount = tree.quadrants.elem_count; for (iz = 0; iz < incount; ++iz) { q1 = p4est_quadrant_array_index (&tree.quadrants, iz); ++tree.quadrants_per_level[q1->level]; tree.maxlevel = (int8_t) SC_MAX (tree.maxlevel, q1->level); } SC_CHECK_ABORT (!p4est_tree_is_linear (&tree), "is_linear"); (void) p4est_linearize_tree (p4est1, &tree); SC_CHECK_ABORT (incount - 3 == tree.quadrants.elem_count, "linearize"); sc_array_reset (&tree.quadrants); /* create a partial tree and check overlap */ sc_array_resize (&tree.quadrants, 3); q1 = p4est_quadrant_array_index (&tree.quadrants, 0); p4est_quadrant_set_morton (q1, 1, 1); q1 = p4est_quadrant_array_index (&tree.quadrants, 1); p4est_quadrant_set_morton (q1, 2, 8); q1 = p4est_quadrant_array_index (&tree.quadrants, 2); p4est_quadrant_set_morton (q1, 2, 9); for (k = 0; k <= P4EST_QMAXLEVEL; ++k) { tree.quadrants_per_level[k] = 0; } for (; k <= P4EST_MAXLEVEL; ++k) { tree.quadrants_per_level[k] = -1; } tree.quadrants_per_level[1] = 1; tree.quadrants_per_level[2] = 2; tree.maxlevel = 2; p4est_quadrant_first_descendant (p4est_quadrant_array_index (&tree.quadrants, 0), &tree.first_desc, P4EST_QMAXLEVEL); p4est_quadrant_last_descendant (p4est_quadrant_array_index (&tree.quadrants, tree.quadrants.elem_count - 1), &tree.last_desc, P4EST_QMAXLEVEL); SC_CHECK_ABORT (p4est_tree_is_complete (&tree), "is_complete"); p4est_quadrant_set_morton (&D, 0, 0); SC_CHECK_ABORT (p4est_quadrant_overlaps_tree (&tree, &D), "overlaps 0"); p4est_quadrant_set_morton (&A, 1, 0); SC_CHECK_ABORT (!p4est_quadrant_overlaps_tree (&tree, &A), "overlaps 1"); p4est_quadrant_set_morton (&A, 1, 1); SC_CHECK_ABORT (p4est_quadrant_overlaps_tree (&tree, &A), "overlaps 2"); p4est_quadrant_set_morton (&A, 1, 2); SC_CHECK_ABORT (p4est_quadrant_overlaps_tree (&tree, &A), "overlaps 3"); p4est_quadrant_set_morton (&A, 1, 3); SC_CHECK_ABORT (!p4est_quadrant_overlaps_tree (&tree, &A), "overlaps 4"); p4est_quadrant_set_morton (&B, 3, 13); SC_CHECK_ABORT (!p4est_quadrant_overlaps_tree (&tree, &B), "overlaps 5"); p4est_quadrant_set_morton (&B, 3, 25); SC_CHECK_ABORT (p4est_quadrant_overlaps_tree (&tree, &B), "overlaps 6"); p4est_quadrant_set_morton (&B, 3, 39); SC_CHECK_ABORT (p4est_quadrant_overlaps_tree (&tree, &B), "overlaps 7"); p4est_quadrant_set_morton (&B, 3, 40); SC_CHECK_ABORT (!p4est_quadrant_overlaps_tree (&tree, &B), "overlaps 8"); p4est_quadrant_set_morton (&C, 4, 219); SC_CHECK_ABORT (!p4est_quadrant_overlaps_tree (&tree, &C), "overlaps 9"); sc_array_reset (&tree.quadrants); /* destroy the p4est and its connectivity structure */ p4est_destroy (p4est1); p4est_destroy (p4est2); p4est_connectivity_destroy (connectivity); /* This will test the ability to address negative quadrants */ P4EST_QUADRANT_INIT (&A); P4EST_QUADRANT_INIT (&B); P4EST_QUADRANT_INIT (&C); P4EST_QUADRANT_INIT (&D); P4EST_QUADRANT_INIT (&E); P4EST_QUADRANT_INIT (&F); P4EST_QUADRANT_INIT (&G); P4EST_QUADRANT_INIT (&H); P4EST_QUADRANT_INIT (&I); P4EST_QUADRANT_INIT (&P); P4EST_QUADRANT_INIT (&Q); A.x = -qone << P4EST_MAXLEVEL; A.y = -qone << P4EST_MAXLEVEL; A.level = 0; B.x = qone << P4EST_MAXLEVEL; B.y = -qone << P4EST_MAXLEVEL; B.level = 0; C.x = -qone << P4EST_MAXLEVEL; C.y = qone << P4EST_MAXLEVEL; C.level = 0; D.x = qone << P4EST_MAXLEVEL; D.y = qone << P4EST_MAXLEVEL; D.level = 0; /* this one is outside the 3x3 box */ E.x = -qone << (P4EST_MAXLEVEL + 1); E.y = -qone; E.level = 0; F.x = P4EST_ROOT_LEN + (P4EST_ROOT_LEN - mh); F.y = P4EST_ROOT_LEN + (P4EST_ROOT_LEN - mh); F.level = P4EST_QMAXLEVEL; G.x = -mh; G.y = -mh; G.level = P4EST_QMAXLEVEL; H.x = -qone << (P4EST_MAXLEVEL - 1); H.y = -qone << (P4EST_MAXLEVEL - 1); H.level = 1; I.x = -qone << P4EST_MAXLEVEL; I.y = -qone << (P4EST_MAXLEVEL - 1); I.level = 1; check_linear_id (&A, &A); check_linear_id (&A, &B); check_linear_id (&A, &C); check_linear_id (&A, &D); /* check_linear_id (&A, &E); */ check_linear_id (&A, &F); check_linear_id (&A, &G); check_linear_id (&A, &H); check_linear_id (&A, &I); check_linear_id (&B, &A); check_linear_id (&B, &B); check_linear_id (&B, &C); check_linear_id (&B, &D); /* check_linear_id (&B, &E); */ check_linear_id (&B, &F); check_linear_id (&B, &G); check_linear_id (&B, &H); check_linear_id (&B, &I); check_linear_id (&D, &A); check_linear_id (&D, &B); check_linear_id (&D, &C); check_linear_id (&D, &D); /* check_linear_id (&D, &E); */ check_linear_id (&D, &F); check_linear_id (&D, &G); check_linear_id (&D, &H); check_linear_id (&D, &I); check_linear_id (&G, &A); check_linear_id (&G, &B); check_linear_id (&G, &C); check_linear_id (&G, &D); /* check_linear_id (&G, &E); */ check_linear_id (&G, &F); check_linear_id (&G, &G); check_linear_id (&G, &H); check_linear_id (&G, &I); check_linear_id (&I, &A); check_linear_id (&I, &B); check_linear_id (&I, &C); check_linear_id (&I, &D); /* check_linear_id (&I, &E); */ check_linear_id (&I, &F); check_linear_id (&I, &G); check_linear_id (&I, &H); check_linear_id (&I, &I); SC_CHECK_ABORT (p4est_quadrant_is_extended (&A) == 1, "is_extended A"); SC_CHECK_ABORT (p4est_quadrant_is_extended (&B) == 1, "is_extended B"); SC_CHECK_ABORT (p4est_quadrant_is_extended (&C) == 1, "is_extended C"); SC_CHECK_ABORT (p4est_quadrant_is_extended (&D) == 1, "is_extended D"); SC_CHECK_ABORT (!p4est_quadrant_is_extended (&E) == 1, "!is_extended E"); SC_CHECK_ABORT (p4est_quadrant_is_extended (&F) == 1, "is_extended F"); SC_CHECK_ABORT (p4est_quadrant_is_extended (&G) == 1, "is_extended G"); SC_CHECK_ABORT (p4est_quadrant_compare (&A, &A) == 0, "compare"); SC_CHECK_ABORT (p4est_quadrant_compare (&A, &B) > 0, "compare"); SC_CHECK_ABORT (p4est_quadrant_compare (&B, &A) < 0, "compare"); SC_CHECK_ABORT (p4est_quadrant_compare (&F, &F) == 0, "compare"); SC_CHECK_ABORT (p4est_quadrant_compare (&G, &F) > 0, "compare"); SC_CHECK_ABORT (p4est_quadrant_compare (&F, &G) < 0, "compare"); A.p.which_tree = 0; B.p.piggy1.which_tree = 0; SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&A, &A) == 0, "compare_piggy"); SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&A, &B) > 0, "compare_piggy"); SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&B, &A) < 0, "compare_piggy"); F.p.piggy2.which_tree = 0; G.p.which_tree = 0; SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&F, &F) == 0, "compare_piggy"); SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&G, &F) > 0, "compare_piggy"); SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&F, &G) < 0, "compare_piggy"); F.p.piggy1.which_tree = (p4est_topidx_t) P4EST_TOPIDX_MAX - 3; G.p.piggy2.which_tree = (p4est_topidx_t) P4EST_TOPIDX_MAX / 2; SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&F, &F) == 0, "compare_piggy"); SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&G, &F) < 0, "compare_piggy"); SC_CHECK_ABORT (p4est_quadrant_compare_piggy (&F, &G) > 0, "compare_piggy"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&A, &A) == 1, "is_equal"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&F, &F) == 1, "is_equal"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&G, &G) == 1, "is_equal"); /* Not sure if these make sense because D, O and A are all level 0 */ #if 0 SC_CHECK_ABORT (p4est_quadrant_is_sibling (&D, &O) == 1, "is_sibling"); SC_CHECK_ABORT (p4est_quadrant_is_sibling (&D, &A) == 0, "is_sibling"); SC_CHECK_ABORT (p4est_quadrant_is_sibling_D (&D, &O) == 1, "is_sibling_D"); SC_CHECK_ABORT (p4est_quadrant_is_sibling_D (&D, &A) == 0, "is_sibling_D"); #endif SC_CHECK_ABORT (p4est_quadrant_is_sibling (&I, &H) == 1, "is_sibling"); SC_CHECK_ABORT (p4est_quadrant_is_sibling (&I, &G) == 0, "is_sibling"); SC_CHECK_ABORT (p4est_quadrant_is_sibling_D (&I, &H) == 1, "is_sibling_D"); SC_CHECK_ABORT (p4est_quadrant_is_sibling_D (&I, &G) == 0, "is_sibling_D"); SC_CHECK_ABORT (p4est_quadrant_is_parent (&A, &H) == 1, "is_parent"); SC_CHECK_ABORT (p4est_quadrant_is_parent (&H, &A) == 0, "is_parent"); SC_CHECK_ABORT (p4est_quadrant_is_parent (&A, &D) == 0, "is_parent"); SC_CHECK_ABORT (p4est_quadrant_is_parent_D (&A, &H) == 1, "is_parent_D"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor (&A, &G) == 1, "is_ancestor"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor (&G, &A) == 0, "is_ancestor"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor_D (&A, &G) == 1, "is_ancestor_D"); SC_CHECK_ABORT (p4est_quadrant_is_ancestor_D (&G, &A) == 0, "is_ancestor_D"); /* SC_CHECK_ABORT (p4est_quadrant_is_next (&F, &E) == 1, "is_next"); */ SC_CHECK_ABORT (p4est_quadrant_is_next (&A, &H) == 0, "is_next"); /* SC_CHECK_ABORT (p4est_quadrant_is_next_D (&F, &E) == 1, "is_next_D"); */ SC_CHECK_ABORT (p4est_quadrant_is_next_D (&A, &H) == 0, "is_next_D"); p4est_quadrant_parent (&H, &a); SC_CHECK_ABORT (p4est_quadrant_is_equal (&A, &a) == 1, "parent"); p4est_quadrant_sibling (&I, &h, 3); SC_CHECK_ABORT (p4est_quadrant_is_equal (&H, &h) == 1, "sibling"); p4est_quadrant_children (&A, &c0, &c1, &c2, &c3); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c2, &I) == 1, "children"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c3, &H) == 1, "children"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c3, &G) == 0, "children"); SC_CHECK_ABORT (p4est_quadrant_is_family (&c0, &c1, &c2, &c3) == 1, "is_family"); id0 = p4est_quadrant_child_id (&c0); id1 = p4est_quadrant_child_id (&c1); id2 = p4est_quadrant_child_id (&c2); id3 = p4est_quadrant_child_id (&c3); SC_CHECK_ABORT (id0 == 0 && id1 == 1 && id2 == 2 && id3 == 3, "child_id"); SC_CHECK_ABORT (p4est_quadrant_child_id (&G) == 3, "child_id"); p4est_quadrant_first_descendant (&A, &c1, 1); SC_CHECK_ABORT (p4est_quadrant_is_equal (&c0, &c1) == 1, "first_descendant"); p4est_quadrant_last_descendant (&A, &g, P4EST_QMAXLEVEL); SC_CHECK_ABORT (p4est_quadrant_is_equal (&G, &g) == 1, "last_descendant"); Fid = p4est_quadrant_linear_id (&F, P4EST_QMAXLEVEL); p4est_quadrant_set_morton (&f, P4EST_QMAXLEVEL, Fid); SC_CHECK_ABORT (p4est_quadrant_is_equal (&F, &f) == 1, "set_morton/linear_id"); Aid = p4est_quadrant_linear_id (&A, 0); p4est_quadrant_set_morton (&a, 0, Aid); SC_CHECK_ABORT (Aid == 15, "linear_id"); SC_CHECK_ABORT (p4est_quadrant_is_equal (&A, &a) == 1, "set_morton/linear_id"); p4est_nearest_common_ancestor (&I, &H, &a); SC_CHECK_ABORT (p4est_quadrant_is_equal (&A, &a) == 1, "ancestor"); p4est_nearest_common_ancestor_D (&I, &H, &a); SC_CHECK_ABORT (p4est_quadrant_is_equal (&A, &a) == 1, "ancestor_D"); for (k = 0; k < 16; ++k) { if (k != 4 && k != 6 && k != 8 && k != 9 && k != 12 && k != 13 && k != 14) { p4est_quadrant_set_morton (&E, 0, (uint64_t) k); } } p4est_quadrant_set_morton (&P, 0, 10); p4est_quadrant_set_morton (&Q, 0, 11); SC_CHECK_ABORT (p4est_quadrant_is_next (&P, &Q), "is_next"); SC_CHECK_ABORT (!p4est_quadrant_is_next (&A, &Q), "is_next"); sc_finalize (); mpiret = sc_MPI_Finalize (); SC_CHECK_MPI (mpiret); return 0; }
static void test_weird (void) { const p4est_topidx_t num_edges = 1, num_ett = 2; const p4est_topidx_t num_corners = 1, num_ctt = 4; int i; size_t zz; p8est_edge_info_t ei; p8est_edge_transform_t *et; p8est_corner_info_t ci; p8est_corner_transform_t *ct; sc_array_t *eta, *cta; p8est_connectivity_t *conn; conn = p8est_connectivity_new (0, 1, num_edges, num_ett, num_corners, num_ctt); for (i = 0; i < 6; ++i) { conn->tree_to_tree[i] = 0; conn->tree_to_face[i] = (int8_t) i; } conn->tree_to_face[4] = 5; conn->tree_to_face[5] = 4; for (i = 0; i < 12; ++i) { conn->tree_to_edge[i] = -1; } conn->tree_to_edge[5] = 0; conn->tree_to_edge[7] = 0; conn->edge_to_tree[0] = 0; conn->edge_to_tree[1] = 0; conn->edge_to_edge[0] = 5; conn->edge_to_edge[1] = 19; conn->ett_offset[0] = 0; for (i = 0; i < 8; ++i) { conn->tree_to_corner[i] = -1; } conn->tree_to_corner[0] = 0; conn->tree_to_corner[1] = 0; conn->tree_to_corner[4] = 0; conn->tree_to_corner[5] = 0; conn->corner_to_tree[0] = 0; conn->corner_to_tree[1] = 0; conn->corner_to_tree[2] = 0; conn->corner_to_tree[3] = 0; conn->corner_to_corner[0] = 0; conn->corner_to_corner[1] = 1; conn->corner_to_corner[2] = 4; conn->corner_to_corner[3] = 5; conn->ctt_offset[0] = 0; P4EST_ASSERT (p8est_connectivity_is_valid (conn)); eta = &ei.edge_transforms; sc_array_init (eta, sizeof (p8est_edge_transform_t)); for (i = 0; i < 2; ++i) { p8est_find_edge_transform (conn, 0, weird_edges[i][0], &ei); SC_CHECK_ABORT ((int) ei.iedge == weird_edges[i][0], "WE ei"); SC_CHECK_ABORT (eta->elem_count == 1, "WE count A"); for (zz = 0; zz < eta->elem_count; ++zz) { et = p8est_edge_array_index (eta, zz); SC_CHECK_ABORT (et->ntree == 0, "WE tree"); SC_CHECK_ABORT ((int) et->nedge == weird_edges[i][1], "WE edge"); SC_CHECK_ABORT (et->nflip == 1, "WE flip"); SC_CHECK_ABORT (et->corners == et->nedge % 4, "WE corners"); SC_CHECK_ABORT (et->naxis[0] == 1 && et->naxis[1] == 0 && et->naxis[2] == 2, "WE axis"); } } sc_array_reset (eta); cta = &ci.corner_transforms; sc_array_init (cta, sizeof (p8est_corner_transform_t)); for (i = 0; i < 8; ++i) { p8est_find_corner_transform (conn, 0, i, &ci); SC_CHECK_ABORT ((int) ci.icorner == i, "WC ci"); SC_CHECK_ABORT ((int) cta->elem_count == 2 - (i & 0x02), "WC count"); for (zz = 0; zz < cta->elem_count; ++zz) { ct = p8est_corner_array_index (cta, zz); SC_CHECK_ABORT (ct->ntree == 0, "WC tree"); SC_CHECK_ABORT ((size_t) ct->ncorner == 4 * zz + !(i % 2), "WC corner"); } } sc_array_reset (cta); p8est_connectivity_destroy (conn); }
static void test_rotwrap (p8est_connectivity_t * conn) { int i; int iface, nface; int iedge, icorner; int ft[9]; size_t zz; p4est_topidx_t itree, ntree; p8est_edge_info_t ei; p8est_edge_transform_t *et; p8est_corner_info_t ci; p8est_corner_transform_t *ct; sc_array_t *eta, *cta; itree = 0; for (iface = 0; iface < 6; ++iface) { ntree = p8est_find_face_transform (conn, itree, iface, ft); if (iface == 2 || iface == 3) { SC_CHECK_ABORT (ntree == -1, "RF tree"); continue; } nface = p8est_face_dual[iface]; SC_CHECK_ABORT (ntree == itree, "RF tree"); if (iface == 0 || iface == 1) { for (i = 0; i < 3; ++i) { SC_CHECK_ABORTF (ft[i] == (i + 1) % 3, "RFA axis A%d", i); SC_CHECK_ABORTF (ft[i] == ft[i + 3], "RFA axis B%d", i); } SC_CHECK_ABORT (ft[6] == 0 && ft[7] == 0, "RFA reverse"); } else { for (i = 0; i < 3; ++i) { SC_CHECK_ABORTF (ft[i] == i, "RFB axis A%d", i); } SC_CHECK_ABORT (ft[0] == ft[4], "RFB axis B0"); SC_CHECK_ABORT (ft[1] == ft[3], "RFB axis B1"); SC_CHECK_ABORT (ft[2] == ft[5], "RFB axis B2"); SC_CHECK_ABORT (ft[6] == (iface != 4) && ft[7] == (iface != 5), "RFB reverse"); } SC_CHECK_ABORT (ft[8] == 2 * (iface % 2) + nface % 2, "RF code"); } eta = &ei.edge_transforms; sc_array_init (eta, sizeof (p8est_edge_transform_t)); for (iedge = 0; iedge < 12; ++iedge) { p8est_find_edge_transform (conn, itree, iedge, &ei); SC_CHECK_ABORT ((int) ei.iedge == iedge, "RE ei"); SC_CHECK_ABORT ((int) eta->elem_count == 2 - (iedge / 4), "RE count AB"); for (zz = 0; zz < eta->elem_count; ++zz) { et = p8est_edge_array_index (eta, zz); SC_CHECK_ABORT (et->ntree == itree, "RE tree"); SC_CHECK_ABORT ((int) et->nedge == rotwrap_edges[iedge][zz], "RE edge"); SC_CHECK_ABORT ((int) et->nflip == rotwrap_flip[iedge][zz], "RE flip"); SC_CHECK_ABORT (et->corners == et->nedge % 4, "RE corners"); SC_CHECK_ABORT ((int) et->naxis[0] == rotwrap_axes[iedge][zz] && et->naxis[1] < et->naxis[2] && et->naxis[0] + et->naxis[1] + et->naxis[2] == 3, "RE axis"); } } sc_array_reset (eta); cta = &ci.corner_transforms; sc_array_init (cta, sizeof (p8est_corner_transform_t)); for (icorner = 0; icorner < 8; ++icorner) { p8est_find_corner_transform (conn, itree, icorner, &ci); SC_CHECK_ABORT ((int) ci.icorner == icorner, "RC ci"); SC_CHECK_ABORT (cta->elem_count == 2, "RC count"); for (zz = 0; zz < cta->elem_count; ++zz) { ct = p8est_corner_array_index (cta, zz); SC_CHECK_ABORT (ct->ntree == itree, "RC tree"); SC_CHECK_ABORT ((int) ct->ncorner == rotwrap_corners[icorner][zz], "RC corner"); } } sc_array_reset (cta); }
p6est_lnodes_t * p6est_lnodes_new (p6est_t * p6est, p6est_ghost_t * ghost, int degree) { p6est_lnodes_t *lnodes; p6est_profile_t *profile; p4est_lnodes_t *clnodes; int nperelem = (degree + 1) * (degree + 1) * (degree + 1); /* int nperface = (degree - 1) * (degree - 1); */ /* int nperedge = (degree - 1); */ p4est_locidx_t ncid, cid, enid, *en; p4est_locidx_t nnodecols; p4est_locidx_t nelemcols; p4est_locidx_t nll; p4est_locidx_t nlayers; p4est_locidx_t *layernodecount; p4est_locidx_t *layernodeoffsets; p4est_locidx_t (*lr)[2]; p4est_locidx_t ncolnodes; p4est_locidx_t *global_owned_count; p4est_locidx_t num_owned, num_local; p4est_gloidx_t gnum_owned, offset; p4est_gloidx_t *owned_offsets; int i, j, k; int mpisize = p6est->mpisize; int mpiret; sc_array_t lnoview; size_t zz, nsharers; int Nrp = degree + 1; if (degree == 1) { p4est_locidx_t eid, nid, enid2, nid2; p4est_locidx_t *newnum, newlocal, newowned; P4EST_GLOBAL_PRODUCTION ("Into adapt p6est_lnodes_new for degree = 1\n"); p4est_log_indent_push (); /* adapt 2 to 1 */ lnodes = p6est_lnodes_new (p6est, ghost, 2); nll = p6est->layers->elem_count; num_local = lnodes->num_local_nodes; num_owned = lnodes->owned_count; en = lnodes->element_nodes; newnum = P4EST_ALLOC (p4est_locidx_t, P8EST_INSUL * nll); memset (newnum, -1, P8EST_INSUL * nll * sizeof (p4est_locidx_t)); for (enid = 0, eid = 0; eid < nll; eid++) { for (k = 0; k < 3; k++) { for (j = 0; j < 3; j++) { for (i = 0; i < 3; i++, enid++) { if (k != 1 && j != 1 && i != 1) { newnum[en[enid]] = 0; } } } } } newlocal = 0; newowned = 0; for (nid = 0; nid < num_local; nid++) { if (newnum[nid] >= 0) { newnum[nid] = newlocal++; if (nid < num_owned) { newowned++; } } } /* compress en */ enid2 = 0; for (enid = 0, eid = 0; eid < nll; eid++) { for (k = 0; k < 3; k++) { for (j = 0; j < 3; j++) { for (i = 0; i < 3; i++, enid++) { if (k != 1 && j != 1 && i != 1) { en[enid2++] = newnum[en[enid]]; } } } } } P4EST_ASSERT (enid2 == P8EST_CHILDREN * nll); lnodes->element_nodes = P4EST_REALLOC (en, p4est_locidx_t, P8EST_CHILDREN * nll); owned_offsets = P4EST_ALLOC (p4est_gloidx_t, mpisize + 1); mpiret = sc_MPI_Allgather (&newowned, 1, P4EST_MPI_LOCIDX, lnodes->global_owned_count, 1, P4EST_MPI_LOCIDX, p6est->mpicomm); owned_offsets[0] = 0; for (i = 0; i < mpisize; i++) { owned_offsets[i + 1] = owned_offsets[i] + lnodes->global_owned_count[i]; } lnodes->global_offset = owned_offsets[p6est->mpirank]; lnodes->num_local_nodes = newlocal; lnodes->owned_count = newowned; lnodes->degree = 1; lnodes->vnodes = P8EST_CHILDREN; lnodes->nonlocal_nodes = P4EST_REALLOC (lnodes->nonlocal_nodes, p4est_gloidx_t, newlocal - newowned); nsharers = lnodes->sharers->elem_count; for (zz = 0; zz < nsharers; zz++) { size_t nshared, zy, zw; p6est_lnodes_rank_t *rank = p6est_lnodes_rank_array_index (lnodes->sharers, zz); if (rank->owned_count) { if (rank->rank != p6est->mpirank) { p4est_locidx_t newrankowned = 0; p4est_locidx_t newrankoffset = -1; for (nid = rank->owned_offset; nid < rank->owned_offset + rank->owned_count; nid++) { if (newnum[nid] >= 0) { lnodes->nonlocal_nodes[newnum[nid] - newowned] = owned_offsets[rank->rank]; newrankowned++; if (newrankoffset < 0) { newrankoffset = newnum[nid]; } } } rank->owned_offset = newrankoffset; rank->owned_count = newrankowned; } else { rank->owned_offset = 0; rank->owned_count = newowned; } } rank->shared_mine_count = 0; rank->shared_mine_offset = -1; zw = 0; nshared = rank->shared_nodes.elem_count; for (zy = 0; zy < nshared; zy++) { nid = *((p4est_locidx_t *) sc_array_index (&rank->shared_nodes, zy)); if (newnum[nid] >= 0) { p4est_locidx_t *lp; lp = (p4est_locidx_t *) sc_array_index (&rank->shared_nodes, zw++); *lp = newnum[nid]; if (newnum[nid] < newowned) { rank->shared_mine_count++; if (rank->shared_mine_offset == -1) { rank->shared_mine_offset = zw - 1; } } } } sc_array_resize (&rank->shared_nodes, zw); } /* send local numbers to others */ { sc_array_t view; sc_array_init_data (&view, newnum, sizeof (p4est_locidx_t), newlocal); p6est_lnodes_share_owned (&view, lnodes); } nid2 = 0; for (nid = num_owned; nid < num_local; nid++) { if (newnum[nid] >= 0) { lnodes->nonlocal_nodes[nid2++] += (p4est_gloidx_t) newnum[nid]; } } P4EST_ASSERT (nid2 == newlocal - newowned); P4EST_FREE (owned_offsets); P4EST_FREE (newnum); p4est_log_indent_pop (); P4EST_GLOBAL_PRODUCTION ("Done adapt p6est_lnodes_new for degree = 1\n"); return lnodes; } P4EST_GLOBAL_PRODUCTION ("Into p6est_lnodes_new\n"); p4est_log_indent_push (); P4EST_ASSERT (degree >= 1); lnodes = P4EST_ALLOC (p6est_lnodes_t, 1); /* first get the profile */ profile = p6est_profile_new_local (p6est, ghost, P6EST_PROFILE_INTERSECTION, P8EST_CONNECT_FULL, degree); p6est_profile_sync (profile); lr = (p4est_locidx_t (*)[2]) profile->lnode_ranges; clnodes = profile->lnodes; nnodecols = clnodes->num_local_nodes; nelemcols = clnodes->num_local_elements; en = clnodes->element_nodes; layernodecount = P4EST_ALLOC_ZERO (p4est_locidx_t, nnodecols); layernodeoffsets = P4EST_ALLOC_ZERO (p4est_locidx_t, nnodecols + 1); for (cid = 0, enid = 0; cid < nelemcols; cid++) { for (j = 0; j < Nrp; j++) { for (i = 0; i < Nrp; i++, enid++) { ncid = en[enid]; nlayers = lr[ncid][1]; P4EST_ASSERT (nlayers); ncolnodes = nlayers * degree + 1; layernodecount[ncid] = ncolnodes; } } } num_owned = 0; num_local = 0; for (ncid = 0; ncid < nnodecols; ncid++) { num_local += layernodecount[ncid]; if (ncid < clnodes->owned_count) { num_owned += layernodecount[ncid]; } } P4EST_VERBOSEF ("p6est_lnodes: %d owned %d local\n", num_owned, num_local); if (nnodecols) { layernodeoffsets[0] = 0; for (ncid = 0; ncid < nnodecols; ncid++) { layernodeoffsets[ncid + 1] = layernodeoffsets[ncid] + layernodecount[ncid]; } } gnum_owned = num_owned; owned_offsets = P4EST_ALLOC (p4est_gloidx_t, mpisize + 1); global_owned_count = P4EST_ALLOC (p4est_locidx_t, mpisize); mpiret = sc_MPI_Allgather (&gnum_owned, 1, P4EST_MPI_GLOIDX, owned_offsets, 1, P4EST_MPI_GLOIDX, p6est->mpicomm); SC_CHECK_MPI (mpiret); offset = 0; for (i = 0; i < mpisize; i++) { global_owned_count[i] = (p4est_locidx_t) owned_offsets[i]; gnum_owned = owned_offsets[i]; owned_offsets[i] = offset; offset += gnum_owned; } owned_offsets[mpisize] = offset; nll = p6est->layers->elem_count; nsharers = clnodes->sharers->elem_count; lnodes->mpicomm = p6est->mpicomm; lnodes->num_local_nodes = num_local; lnodes->owned_count = num_owned; lnodes->global_offset = owned_offsets[p6est->mpirank]; lnodes->nonlocal_nodes = P4EST_ALLOC (p4est_gloidx_t, num_local - num_owned); lnodes->sharers = sc_array_new_size (sizeof (p6est_lnodes_rank_t), nsharers); lnodes->global_owned_count = global_owned_count; lnodes->degree = degree; lnodes->vnodes = nperelem; lnodes->num_local_elements = nll; lnodes->face_code = P4EST_ALLOC (p6est_lnodes_code_t, nll); lnodes->element_nodes = P4EST_ALLOC (p4est_locidx_t, nperelem * nll); p6est_profile_element_to_node (p6est, profile, layernodeoffsets, lnodes->element_nodes, lnodes->face_code); for (zz = 0; zz < nsharers; zz++) { p4est_lnodes_rank_t *crank = p4est_lnodes_rank_array_index (clnodes->sharers, zz); p6est_lnodes_rank_t *rank = p6est_lnodes_rank_array_index (lnodes->sharers, zz); size_t zy; size_t nshared; rank->rank = crank->rank; sc_array_init (&rank->shared_nodes, sizeof (p4est_locidx_t)); nshared = crank->shared_nodes.elem_count; rank->owned_offset = -1; rank->owned_count = 0; rank->shared_mine_count = 0; rank->shared_mine_offset = -1; for (zy = 0; zy < nshared; zy++) { p4est_locidx_t cnid = *((p4est_locidx_t *) sc_array_index (&crank->shared_nodes, zy)); p4est_locidx_t *lp; p4est_locidx_t nthis, il; p4est_locidx_t old_count = rank->shared_nodes.elem_count; nthis = layernodecount[cnid]; lp = (p4est_locidx_t *) sc_array_push_count (&rank->shared_nodes, nthis); for (il = 0; il < nthis; il++) { lp[il] = layernodeoffsets[cnid] + il; if (zy >= (size_t) crank->shared_mine_offset && (p4est_locidx_t) zy - crank->shared_mine_offset < crank->shared_mine_count) { rank->shared_mine_count++; if (rank->shared_mine_offset == -1) { rank->shared_mine_offset = old_count + il; } } if (cnid >= crank->owned_offset && cnid - crank->owned_offset < crank->owned_count) { rank->owned_count++; if (rank->owned_offset == -1) { rank->owned_offset = lp[il]; } } } } if (rank->rank == p6est->mpirank) { rank->owned_offset = 0; rank->owned_count = num_owned; } } memcpy (layernodecount, layernodeoffsets, nnodecols * sizeof (p4est_locidx_t)); sc_array_init_data (&lnoview, layernodecount, sizeof (p4est_locidx_t), (size_t) nnodecols); p4est_lnodes_share_owned (&lnoview, clnodes); for (zz = 0; zz < nsharers; zz++) { p4est_lnodes_rank_t *crank = p4est_lnodes_rank_array_index (clnodes->sharers, zz); if (crank->rank == p6est->mpirank) { continue; } for (ncid = crank->owned_offset; ncid < crank->owned_offset + crank->owned_count; ncid++) { p4est_gloidx_t owners_offset; p4est_locidx_t nid; P4EST_ASSERT (ncid >= clnodes->owned_count); owners_offset = owned_offsets[crank->rank] + layernodecount[ncid]; for (nid = layernodeoffsets[ncid]; nid < layernodeoffsets[ncid + 1]; nid++) { P4EST_ASSERT (nid >= num_owned); P4EST_ASSERT (nid < num_local); lnodes->nonlocal_nodes[nid - num_owned] = owners_offset++; } } } p6est_profile_destroy (profile); P4EST_FREE (owned_offsets); P4EST_FREE (layernodecount); P4EST_FREE (layernodeoffsets); p4est_log_indent_pop (); P4EST_GLOBAL_PRODUCTION ("Done p6est_lnodes_new\n"); return lnodes; }
void p8est_quadrant_edge_neighbor_extra (const p4est_quadrant_t * q, p4est_topidx_t t, int edge, sc_array_t * quads, sc_array_t * treeids, p4est_connectivity_t * conn) { p4est_quadrant_t temp; p4est_quadrant_t *qp; p4est_topidx_t *tp; int face; p8est_edge_info_t ei; p8est_edge_transform_t *et; sc_array_t *eta; size_t etree; eta = &ei.edge_transforms; P4EST_ASSERT (SC_ARRAY_IS_OWNER (quads)); P4EST_ASSERT (quads->elem_count == 0); P4EST_ASSERT (quads->elem_size == sizeof (p4est_quadrant_t)); P4EST_ASSERT (SC_ARRAY_IS_OWNER (treeids)); P4EST_ASSERT (treeids->elem_count == 0); P4EST_ASSERT (treeids->elem_size == sizeof (p4est_topidx_t)); p8est_quadrant_edge_neighbor (q, edge, &temp); if (p4est_quadrant_is_inside_root (&temp)) { qp = p4est_quadrant_array_push (quads); *qp = temp; tp = (p4est_topidx_t *) sc_array_push (treeids); *tp = t; return; } if (!p8est_quadrant_is_outside_edge (&temp)) { qp = p4est_quadrant_array_push (quads); tp = (p4est_topidx_t *) sc_array_push (treeids); face = p8est_edge_faces[edge][0]; p4est_quadrant_face_neighbor (q, face, &temp); if (p4est_quadrant_is_inside_root (&temp)) { face = p8est_edge_faces[edge][1]; *tp = p8est_quadrant_face_neighbor_extra (&temp, t, face, qp, conn); if (*tp == -1) { qp = (p4est_quadrant_t *) sc_array_pop (quads); tp = (p4est_topidx_t *) sc_array_pop (treeids); } return; } face = p8est_edge_faces[edge][1]; p4est_quadrant_face_neighbor (q, face, &temp); P4EST_ASSERT (p4est_quadrant_is_inside_root (&temp)); face = p8est_edge_faces[edge][0]; *tp = p8est_quadrant_face_neighbor_extra (&temp, t, face, qp, conn); if (*tp == -1) { qp = (p4est_quadrant_t *) sc_array_pop (quads); tp = (p4est_topidx_t *) sc_array_pop (treeids); } return; } sc_array_init (eta, sizeof (p8est_edge_transform_t)); p8est_find_edge_transform (conn, t, edge, &ei); sc_array_resize (quads, eta->elem_count); sc_array_resize (treeids, eta->elem_count); for (etree = 0; etree < eta->elem_count; etree++) { qp = p4est_quadrant_array_index (quads, etree); tp = (p4est_topidx_t *) sc_array_index (treeids, etree); et = p8est_edge_array_index (eta, etree); p8est_quadrant_transform_edge (&temp, qp, &ei, et, 1); *tp = et->ntree; } sc_array_reset (eta); }
p4est_mesh_t * p4est_mesh_new_ext (p4est_t * p4est, p4est_ghost_t * ghost, int compute_tree_index, int compute_level_lists, p4est_connect_type_t btype) { int do_corner = 0; int do_volume = 0; int rank; p4est_locidx_t lq, ng; p4est_locidx_t jl; p4est_mesh_t *mesh; P4EST_ASSERT (p4est_is_balanced (p4est, btype)); mesh = P4EST_ALLOC_ZERO (p4est_mesh_t, 1); lq = mesh->local_num_quadrants = p4est->local_num_quadrants; ng = mesh->ghost_num_quadrants = (p4est_locidx_t) ghost->ghosts.elem_count; if (btype == P4EST_CONNECT_FULL) { do_corner = 1; } do_volume = (compute_tree_index || compute_level_lists ? 1 : 0); /* Optional map of tree index for each quadrant */ if (compute_tree_index) { mesh->quad_to_tree = P4EST_ALLOC (p4est_topidx_t, lq); } mesh->ghost_to_proc = P4EST_ALLOC (int, ng); mesh->quad_to_quad = P4EST_ALLOC (p4est_locidx_t, P4EST_FACES * lq); mesh->quad_to_face = P4EST_ALLOC (int8_t, P4EST_FACES * lq); mesh->quad_to_half = sc_array_new (P4EST_HALF * sizeof (p4est_locidx_t)); /* Optional per-level lists of quadrants */ if (compute_level_lists) { mesh->quad_level = P4EST_ALLOC (sc_array_t, P4EST_QMAXLEVEL + 1); for (jl = 0; jl <= P4EST_QMAXLEVEL; ++jl) { sc_array_init (mesh->quad_level + jl, sizeof (p4est_locidx_t)); } } /* Populate ghost information */ rank = 0; for (jl = 0; jl < ng; ++jl) { while (ghost->proc_offsets[rank + 1] <= jl) { ++rank; P4EST_ASSERT (rank < p4est->mpisize); } mesh->ghost_to_proc[jl] = rank; } /* Fill face arrays with default values */ memset (mesh->quad_to_quad, -1, P4EST_FACES * lq * sizeof (p4est_locidx_t)); memset (mesh->quad_to_face, -25, P4EST_FACES * lq * sizeof (int8_t)); if (do_corner) { /* Initialize corner information to a consistent state */ mesh->quad_to_corner = P4EST_ALLOC (p4est_locidx_t, P4EST_CHILDREN * lq); memset (mesh->quad_to_corner, -1, P4EST_CHILDREN * lq * sizeof (p4est_locidx_t)); mesh->corner_offset = sc_array_new (sizeof (p4est_locidx_t)); *(p4est_locidx_t *) sc_array_push (mesh->corner_offset) = 0; mesh->corner_quad = sc_array_new (sizeof (p4est_locidx_t)); mesh->corner_corner = sc_array_new (sizeof (int8_t)); } /* Call the forest iterator to collect face connectivity */ p4est_iterate (p4est, ghost, mesh, (do_volume ? mesh_iter_volume : NULL), mesh_iter_face, #ifdef P4_TO_P8 NULL, #endif do_corner ? mesh_iter_corner : NULL); return mesh; }