/** Coarsen by the L2 error estimate of the initial condition. * * Given the maximum global error, we enforce that each quadrant's portion of * the error must not exceed is fraction of the total volume of the domain * (which is 1). * * \param [in] p4est the forest * \param [in] which_tree the tree in the forest containing \a children * \param [in] children a family of quadrants * * \return 1 if \a children should be coarsened, 0 otherwise. */ static int step3_coarsen_initial_condition (p4est_t * p4est, p4est_topidx_t which_tree, p4est_quadrant_t * children[]) { p4est_quadrant_t parent; step3_ctx_t *ctx = (step3_ctx_t *) p4est->user_pointer; double global_err = ctx->max_err; double global_err2 = global_err * global_err; double h; step3_data_t parentdata; double parentmidpoint[3]; double vol, err2; /* get the parent of the first child (the parent of all children) */ p4est_quadrant_parent (children[0], &parent); step3_get_midpoint (p4est, which_tree, &parent, parentmidpoint); parentdata.u = step3_initial_condition (parentmidpoint, parentdata.du, ctx); h = (double) P4EST_QUADRANT_LEN (parent.level) / (double) P4EST_ROOT_LEN; /* the quadrant's volume is also its volume fraction */ #ifdef P4_TO_P8 vol = h * h * h; #else vol = h * h; #endif parent.p.user_data = (void *) (&parentdata); err2 = step3_error_sqr_estimate (&parent); if (err2 < global_err2 * vol) { return 1; } else { return 0; } }
static void check_linear_id (const p4est_quadrant_t * q1, const p4est_quadrant_t * q2) { int l; int comp = p4est_quadrant_compare (q1, q2); int level = (int) SC_MIN (q1->level, q2->level); uint64_t id1 = p4est_quadrant_linear_id (q1, level); uint64_t id2 = p4est_quadrant_linear_id (q2, level); p4est_quadrant_t quad, par, anc; /* test linear id */ if (p4est_quadrant_is_ancestor (q1, q2)) { SC_CHECK_ABORT (id1 == id2 && comp < 0, "Ancestor 1"); } else if (p4est_quadrant_is_ancestor (q2, q1)) { SC_CHECK_ABORT (id1 == id2 && comp > 0, "Ancestor 2"); } else { SC_CHECK_ABORT ((comp == 0 && id1 == id2) || (comp < 0 && id1 < id2) || (comp > 0 && id1 > id2), "compare"); } /* test ancestor and parent functions */ par = quad = *q1; for (l = quad.level - 1; l >= 0; --l) { p4est_quadrant_parent (&par, &par); p4est_quadrant_ancestor (&quad, l, &anc); SC_CHECK_ABORT (p4est_quadrant_is_equal (&par, &anc), "Ancestor test A"); } par = quad = *q2; for (l = quad.level - 1; l >= 0; --l) { p4est_quadrant_parent (&par, &par); p4est_quadrant_ancestor (&quad, l, &anc); SC_CHECK_ABORT (p4est_quadrant_is_equal (&par, &anc), "Ancestor test B"); } }
static int coarsen_fn (p4est_t * p4est, p4est_topidx_t which_tree, p4est_quadrant_t * q[]) { int pid; p4est_quadrant_t p; SC_CHECK_ABORT (p4est_quadrant_is_familypv (q), "Coarsen invocation"); if (q[0]->level <= 2) return 0; p4est_quadrant_parent (q[0], &p); pid = p4est_quadrant_child_id (&p); return pid == 3; }
static void p4est_coarsen_old (p4est_t * p4est, int coarsen_recursive, p4est_coarsen_t coarsen_fn, p4est_init_t init_fn) { #ifdef P4EST_ENABLE_DEBUG size_t data_pool_size; #endif int i, maxlevel; int couldbegood; size_t zz; size_t incount, removed; size_t cidz, first, last, rest, before; p4est_locidx_t num_quadrants, prev_offset; p4est_topidx_t jt; p4est_tree_t *tree; p4est_quadrant_t *c[P4EST_CHILDREN]; p4est_quadrant_t *cfirst, *clast; sc_array_t *tquadrants; P4EST_GLOBAL_PRODUCTIONF ("Into " P4EST_STRING "_coarsen_old with %lld total quadrants\n", (long long) p4est->global_num_quadrants); p4est_log_indent_push (); P4EST_ASSERT (p4est_is_valid (p4est)); /* loop over all local trees */ prev_offset = 0; for (jt = p4est->first_local_tree; jt <= p4est->last_local_tree; ++jt) { tree = p4est_tree_array_index (p4est->trees, jt); tquadrants = &tree->quadrants; #ifdef P4EST_ENABLE_DEBUG data_pool_size = 0; if (p4est->user_data_pool != NULL) { data_pool_size = p4est->user_data_pool->elem_count; } #endif removed = 0; /* initial log message for this tree */ P4EST_VERBOSEF ("Into coarsen tree %lld with %llu\n", (long long) jt, (unsigned long long) tquadrants->elem_count); /* Initialize array indices. If children are coarsened, the array will have an empty window. first index of the first child to be considered last index of the last child before the hole in the array before number of children before the hole in the array rest index of the first child after the hole in the array */ first = last = 0; before = rest = 1; /* run through the array and coarsen recursively */ incount = tquadrants->elem_count; while (rest + P4EST_CHILDREN - 1 - before < incount) { couldbegood = 1; for (zz = 0; zz < P4EST_CHILDREN; ++zz) { if (zz < before) { c[zz] = p4est_quadrant_array_index (tquadrants, first + zz); if (zz != (size_t) p4est_quadrant_child_id (c[zz])) { couldbegood = 0; break; } } else { c[zz] = p4est_quadrant_array_index (tquadrants, rest + zz - before); } } if (couldbegood && p4est_quadrant_is_familypv (c) && coarsen_fn (p4est, jt, c)) { /* coarsen now */ for (zz = 0; zz < P4EST_CHILDREN; ++zz) { p4est_quadrant_free_data (p4est, c[zz]); } tree->quadrants_per_level[c[0]->level] -= P4EST_CHILDREN; cfirst = c[0]; p4est_quadrant_parent (c[0], cfirst); p4est_quadrant_init_data (p4est, jt, cfirst, init_fn); tree->quadrants_per_level[cfirst->level] += 1; p4est->local_num_quadrants -= P4EST_CHILDREN - 1; removed += P4EST_CHILDREN - 1; rest += P4EST_CHILDREN - before; if (coarsen_recursive) { last = first; cidz = (size_t) p4est_quadrant_child_id (cfirst); if (cidz > first) first = 0; else first -= cidz; } else { /* don't coarsen again, move the counters and the hole */ P4EST_ASSERT (first == last && before == 1); if (rest < incount) { ++first; cfirst = p4est_quadrant_array_index (tquadrants, first); clast = p4est_quadrant_array_index (tquadrants, rest); *cfirst = *clast; last = first; ++rest; } } } else { /* do nothing, just move the counters and the hole */ ++first; if (first > last) { if (first != rest) { cfirst = p4est_quadrant_array_index (tquadrants, first); clast = p4est_quadrant_array_index (tquadrants, rest); *cfirst = *clast; } last = first; ++rest; } } before = last - first + 1; } /* adjust final array size */ first = last; if (first + 1 < rest) { while (rest < incount) { ++first; cfirst = p4est_quadrant_array_index (tquadrants, first); clast = p4est_quadrant_array_index (tquadrants, rest); *cfirst = *clast; ++rest; } sc_array_resize (tquadrants, first + 1); } /* compute maximum level */ maxlevel = 0; num_quadrants = 0; for (i = 0; i <= P4EST_QMAXLEVEL; ++i) { P4EST_ASSERT (tree->quadrants_per_level[i] >= 0); num_quadrants += tree->quadrants_per_level[i]; /* same type */ if (tree->quadrants_per_level[i] > 0) { maxlevel = i; } } tree->maxlevel = (int8_t) maxlevel; tree->quadrants_offset = prev_offset; prev_offset += num_quadrants; /* do some sanity checks */ P4EST_ASSERT (num_quadrants == (p4est_locidx_t) tquadrants->elem_count); P4EST_ASSERT (tquadrants->elem_count == incount - removed); if (p4est->user_data_pool != NULL) { P4EST_ASSERT (data_pool_size - removed == p4est->user_data_pool->elem_count); } P4EST_ASSERT (p4est_tree_is_sorted (tree)); P4EST_ASSERT (p4est_tree_is_complete (tree)); /* final log message for this tree */ P4EST_VERBOSEF ("Done coarsen tree %lld now %llu\n", (long long) jt, (unsigned long long) tquadrants->elem_count); } if (p4est->last_local_tree >= 0) { for (; jt < p4est->connectivity->num_trees; ++jt) { tree = p4est_tree_array_index (p4est->trees, jt); tree->quadrants_offset = p4est->local_num_quadrants; } } /* compute global number of quadrants */ p4est_comm_count_quadrants (p4est); P4EST_ASSERT (p4est_is_valid (p4est)); p4est_log_indent_pop (); P4EST_GLOBAL_PRODUCTIONF ("Done " P4EST_STRING "_coarsen_old with %lld total quadrants\n", (long long) p4est->global_num_quadrants); }
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) { 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; }
int check_balance_seeds (p4est_quadrant_t * q, p4est_quadrant_t * p, p4est_connect_type_t b, sc_array_t * seeds) { int ib; int level = q->level; p4est_quadrant_t *s, *t; sc_array_t *thislevel = sc_array_new (sizeof (p4est_quadrant_t)); sc_array_t *nextlevel = sc_array_new (sizeof (p4est_quadrant_t)); sc_array_t *temparray; p4est_quadrant_t temp1, temp2; int f, c; #ifdef P4_TO_P8 int e; #endif int stop = 0; sc_array_resize (seeds, 0); s = (p4est_quadrant_t *) sc_array_push (thislevel); p4est_quadrant_sibling (q, s, 0); #ifndef P4_TO_P8 if (b == P4EST_CONNECT_FACE) { ib = 0; } else { ib = 1; } #else if (b == P8EST_CONNECT_FACE) { ib = 0; } else if (b == P8EST_CONNECT_EDGE) { ib = 1; } else { ib = 2; } #endif while (level > p->level + 1) { size_t nlast = thislevel->elem_count; size_t zz; stop = 0; for (zz = 0; zz < nlast; zz++) { s = p4est_quadrant_array_index (thislevel, zz); P4EST_ASSERT (p4est_quadrant_child_id (s) == 0); p4est_quadrant_parent (s, &temp1); for (f = 0; f < P4EST_FACES; f++) { p4est_quadrant_face_neighbor (&temp1, f, &temp2); if (is_farther (&temp1, p, &temp2)) { continue; } if (p4est_quadrant_is_ancestor (p, &temp2)) { stop = 1; sc_array_resize (seeds, seeds->elem_count + 1); t = p4est_quadrant_array_index (seeds, seeds->elem_count - 1); p4est_quadrant_sibling (&temp2, t, 0); } else if (p4est_quadrant_is_inside_root (&temp2)) { t = (p4est_quadrant_t *) sc_array_push (nextlevel); p4est_quadrant_sibling (&temp2, t, 0); } } if (ib == 0) { continue; } #ifdef P4_TO_P8 for (e = 0; e < P8EST_EDGES; e++) { p8est_quadrant_edge_neighbor (&temp1, e, &temp2); if (is_farther (&temp1, p, &temp2)) { continue; } if (p4est_quadrant_is_ancestor (p, &temp2)) { stop = 1; sc_array_resize (seeds, seeds->elem_count + 1); t = p4est_quadrant_array_index (seeds, seeds->elem_count - 1); p4est_quadrant_sibling (&temp2, t, 0); } else if (p4est_quadrant_is_inside_root (&temp2)) { t = (p4est_quadrant_t *) sc_array_push (nextlevel); p4est_quadrant_sibling (&temp2, t, 0); } } if (ib == 1) { continue; } #endif for (c = 0; c < P4EST_CHILDREN; c++) { p4est_quadrant_corner_neighbor (&temp1, c, &temp2); if (is_farther (&temp1, p, &temp2)) { continue; } if (p4est_quadrant_is_ancestor (p, &temp2)) { stop = 1; sc_array_resize (seeds, seeds->elem_count + 1); t = p4est_quadrant_array_index (seeds, seeds->elem_count - 1); p4est_quadrant_sibling (&temp2, t, 0); } else if (p4est_quadrant_is_inside_root (&temp2)) { t = (p4est_quadrant_t *) sc_array_push (nextlevel); p4est_quadrant_sibling (&temp2, t, 0); } } } if (stop) { sc_array_sort (seeds, p4est_quadrant_compare); sc_array_uniq (seeds, p4est_quadrant_compare); #ifdef P4_TO_P8 if (!ib && seeds->elem_count == 1) { sc_array_sort (nextlevel, p4est_quadrant_compare); sc_array_uniq (nextlevel, p4est_quadrant_compare); temparray = thislevel; thislevel = nextlevel; nextlevel = temparray; sc_array_reset (nextlevel); level--; nlast = thislevel->elem_count; for (zz = 0; zz < nlast; zz++) { s = p4est_quadrant_array_index (thislevel, zz); P4EST_ASSERT (p4est_quadrant_child_id (s) == 0); p4est_quadrant_parent (s, &temp1); for (f = 0; f < P4EST_FACES; f++) { p4est_quadrant_face_neighbor (&temp1, f, &temp2); if (p4est_quadrant_is_ancestor (p, &temp2)) { int f2; p4est_quadrant_t a; p4est_quadrant_t u; t = p4est_quadrant_array_index (seeds, 0); p8est_quadrant_parent (t, &a); for (f2 = 0; f2 < P8EST_FACES; f2++) { if (f2 / 2 == f / 2) { continue; } p8est_quadrant_face_neighbor (&a, f2, &u); if (p8est_quadrant_is_equal (&temp2, &u) || p8est_quadrant_is_sibling (&temp2, &u)) { break; } } if (f2 == P8EST_FACES) { sc_array_resize (seeds, seeds->elem_count + 1); t = p4est_quadrant_array_index (seeds, seeds->elem_count - 1); p4est_quadrant_sibling (&temp2, t, 0); } } } } } #endif sc_array_sort (seeds, p4est_quadrant_compare); sc_array_uniq (seeds, p4est_quadrant_compare); break; } sc_array_sort (nextlevel, p4est_quadrant_compare); sc_array_uniq (nextlevel, p4est_quadrant_compare); temparray = thislevel; thislevel = nextlevel; nextlevel = temparray; sc_array_reset (nextlevel); level--; } sc_array_destroy (thislevel); sc_array_destroy (nextlevel); return stop; }