p4est_wrap_t * p4est_wrap_new_copy (p4est_wrap_t * source, size_t data_size, p4est_replace_t replace_fn, void *user_pointer) { p4est_wrap_t *pp; P4EST_ASSERT (source != NULL); pp = P4EST_ALLOC_ZERO (p4est_wrap_t, 1); pp->hollow = 1; sc_refcount_init_invalid (&pp->conn_rc); pp->conn_owner = (source->conn_owner != NULL ? source->conn_owner : source); pp->conn = pp->conn_owner->conn; sc_refcount_ref (&pp->conn_owner->conn_rc); pp->p4est_dim = P4EST_DIM; pp->p4est_half = P4EST_HALF; pp->p4est_faces = P4EST_FACES; pp->p4est_children = P4EST_CHILDREN; pp->btype = source->btype; pp->replace_fn = replace_fn; pp->p4est = p4est_copy (source->p4est, 0); if (data_size > 0) { p4est_reset_data (pp->p4est, data_size, NULL, NULL); } pp->weight_exponent = 0; /* keep this even though using ALLOC_ZERO */ pp->p4est->user_pointer = pp; pp->user_pointer = user_pointer; return pp; }
static void p4est_coarsen_both (p4est_t * p4est, int coarsen_recursive, p4est_coarsen_t coarsen_fn, p4est_init_t init_fn) { int success; p4est_t *copy; copy = p4est_copy (p4est, 1); p4est_coarsen_old (copy, coarsen_recursive, coarsen_fn, init_fn); coarsen_callback_count = 0; p4est_coarsen_ext (p4est, coarsen_recursive, 1, coarsen_fn, init_fn, NULL); SC_CHECK_ABORT (coarsen_recursive || coarsen_callback_count == (int) p4est->local_num_quadrants, "Coarsen count"); success = p4est_is_equal (p4est, copy, 1); SC_CHECK_ABORT (success, "Coarsen mismatch"); p4est_destroy (copy); }
int main (int argc, char **argv) { sc_MPI_Comm mpicomm; int mpiret; int size, rank; unsigned crcF, crcC; p4est_connectivity_t *connectivity; p4est_t *p4est; p4est_t *p4estF, *p4estC; #ifdef P4_TO_P8 unsigned crcE; p4est_t *p4estE; #endif /* initialize */ mpiret = sc_MPI_Init (&argc, &argv); SC_CHECK_MPI (mpiret); mpicomm = sc_MPI_COMM_WORLD; mpiret = sc_MPI_Comm_size (mpicomm, &size); SC_CHECK_MPI (mpiret); mpiret = sc_MPI_Comm_rank (mpicomm, &rank); SC_CHECK_MPI (mpiret); sc_init (mpicomm, 1, 1, NULL, SC_LP_DEFAULT); p4est_init (NULL, SC_LP_DEFAULT); /* create forest and refine */ #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, 0, NULL, NULL); p4est_refine (p4est, 1, refine_fn, NULL); /* test face balance */ p4estF = p4est_copy (p4est, 0); #ifndef P4_TO_P8 p4est_balance (p4estF, P4EST_CONNECT_FACE, NULL); #else p4est_balance (p4estF, P8EST_CONNECT_FACE, NULL); #endif crcF = p4est_checksum (p4estF); P4EST_GLOBAL_INFOF ("Face balance with %lld quadrants and crc 0x%08x\n", (long long) p4estF->global_num_quadrants, crcF); #ifdef P4_TO_P8 /* test edge balance */ p4estE = p4est_copy (p4est, 1); p4est_balance (p4estF, P8EST_CONNECT_EDGE, NULL); p4est_balance (p4estE, P8EST_CONNECT_EDGE, NULL); crcE = p4est_checksum (p4estE); SC_CHECK_ABORT (crcE == p4est_checksum (p4estF), "mismatch A"); P4EST_GLOBAL_INFOF ("Edge balance with %lld quadrants and crc 0x%08x\n", (long long) p4estE->global_num_quadrants, crcE); #endif /* test corner balance */ p4estC = p4est_copy (p4est, 1); #ifndef P4_TO_P8 p4est_balance (p4estF, P4EST_CONNECT_CORNER, NULL); p4est_balance (p4estC, P4EST_CONNECT_CORNER, NULL); #else p4est_balance (p4estF, P8EST_CONNECT_CORNER, NULL); p4est_balance (p4estC, P8EST_CONNECT_CORNER, NULL); #endif crcC = p4est_checksum (p4estC); SC_CHECK_ABORT (crcC == p4est_checksum (p4estF), "mismatch B"); P4EST_GLOBAL_INFOF ("Corner balance with %lld quadrants and crc 0x%08x\n", (long long) p4estC->global_num_quadrants, crcC); /* destroy forests and connectivity */ p4est_destroy (p4est); p4est_destroy (p4estF); #ifdef P4_TO_P8 p4est_destroy (p4estE); #endif p4est_destroy (p4estC); p4est_connectivity_destroy (connectivity); /* clean up and exit */ sc_finalize (); mpiret = sc_MPI_Finalize (); SC_CHECK_MPI (mpiret); return 0; }
static void test_build_local (sc_MPI_Comm mpicomm) { sc_array_t *points; p4est_connectivity_t *conn; p4est_t *p4est, *built, *copy; test_build_t stb, *tb = &stb; /* 0. prepare data that we will reuse */ tb->maxlevel = 7 - P4EST_DIM; tb->counter = -1; tb->wrapper = 3; tb->init_default = -1; tb->init_add = -1; tb->count_add = -1; tb->last_tree = -1; tb->build = NULL; #ifndef P4_TO_P8 conn = p4est_connectivity_new_moebius (); #else conn = p8est_connectivity_new_rotcubes (); #endif /* P4_TO_P8 */ p4est = p4est_new_ext (mpicomm, conn, 0, 0, 2, 0, NULL, tb); p4est_refine (p4est, 1, test_build_refine, NULL); p4est_partition (p4est, 0, NULL); /* TODO: enrich tests with quadrant data */ /* 1. Create a p4est that shall be identical to the old one. */ tb->build = p4est_build_new (p4est, 0, NULL, NULL); p4est_search_local (p4est, 0, test_search_local_1, NULL, NULL); built = p4est_build_complete (tb->build); SC_CHECK_ABORT (p4est_is_equal (p4est, built, 0), "Mismatch build_local 1"); p4est_destroy (built); /* 2. Create a p4est that is as coarse as possible. * Coarsen recursively, compare. */ tb->build = p4est_build_new (p4est, 4, NULL, NULL); p4est_search_local (p4est, 0, test_search_local_2, NULL, NULL); built = p4est_build_complete (tb->build); copy = p4est_copy (p4est, 0); p4est_coarsen (copy, 1, test_build_coarsen, NULL); SC_CHECK_ABORT (p4est_is_equal (copy, built, 0), "Mismatch build_local 2"); p4est_destroy (copy); p4est_destroy (built); /* 3. Create a p4est with some random pattern for demonstration */ tb->init_default = 0; tb->init_add = 0; tb->count_add = 0; tb->build = p4est_build_new (p4est, 0, test_search_init_3, tb); p4est_build_init_add (tb->build, test_search_init_add_3); p4est_search_local (p4est, 1, test_search_local_3, NULL, NULL); built = p4est_build_complete (tb->build); p4est_build_verify_3 (built); SC_CHECK_ABORT (p4est_is_valid (built), "Invalid build_local 3"); p4est_destroy (built); /* 4. Create a p4est from a search with one quadrant per tree */ tb->init_default = 0; tb->init_add = 0; tb->count_add = 0; tb->last_tree = -1; tb->build = p4est_build_new (p4est, sizeof (long), test_search_init_4, tb); p4est_build_init_add (tb->build, test_search_init_add_4); p4est_search_local (p4est, 0, test_search_local_4, NULL, NULL); built = p4est_build_complete (tb->build); p4est_build_verify_4 (built); SC_CHECK_ABORT (p4est_is_valid (built), "Invalid build_local 4"); p4est_destroy (built); /* 5. Create a p4est from a multiple-item search */ points = sc_array_new_size (sizeof (int8_t), 2); *(int8_t *) sc_array_index (points, 0) = 0; *(int8_t *) sc_array_index (points, 1) = 1; tb->wrapper = 5; tb->init_default = 0; tb->init_add = 0; tb->build = p4est_build_new (p4est, 0, NULL, tb); p4est_search_local (p4est, 0, NULL, test_search_point_5, points); built = p4est_build_complete (tb->build); #if 0 p4est_build_verify_5 (built); #endif SC_CHECK_ABORT (p4est_is_valid (built), "Invalid build_local 5"); p4est_destroy (built); sc_array_destroy (points); /* clean up */ p4est_destroy (p4est); p4est_connectivity_destroy (conn); }
int main (int argc, char **argv) { int rank; int num_procs; int mpiret; sc_MPI_Comm mpicomm; p4est_t *p4est, *copy; p4est_connectivity_t *connectivity; int i; p4est_topidx_t t; size_t qz; p4est_locidx_t num_quadrants_on_last; p4est_locidx_t *num_quadrants_in_proc; p4est_gloidx_t *pertree1, *pertree2; p4est_quadrant_t *quad; p4est_tree_t *tree; user_data_t *user_data; int64_t sum; unsigned crc; mpiret = sc_MPI_Init (&argc, &argv); SC_CHECK_MPI (mpiret); mpicomm = sc_MPI_COMM_WORLD; mpiret = sc_MPI_Comm_rank (mpicomm, &rank); SC_CHECK_MPI (mpiret); sc_init (mpicomm, 1, 1, NULL, SC_LP_DEFAULT); /* create connectivity and forest structures */ #ifdef P4_TO_P8 connectivity = p8est_connectivity_new_twocubes (); #else connectivity = p4est_connectivity_new_corner (); #endif p4est = p4est_new_ext (mpicomm, connectivity, 15, 0, 0, sizeof (user_data_t), init_fn, NULL); pertree1 = P4EST_ALLOC (p4est_gloidx_t, p4est->connectivity->num_trees + 1); pertree2 = P4EST_ALLOC (p4est_gloidx_t, p4est->connectivity->num_trees + 1); num_procs = p4est->mpisize; num_quadrants_in_proc = P4EST_ALLOC (p4est_locidx_t, num_procs); /* refine and balance to make the number of elements interesting */ test_pertree (p4est, NULL, pertree1); p4est_refine (p4est, 1, refine_fn, init_fn); test_pertree (p4est, NULL, pertree1); /* Set an arbitrary partition. * * Since this is just a test we assume the global number of * quadrants will fit in an int32_t */ num_quadrants_on_last = (p4est_locidx_t) p4est->global_num_quadrants; for (i = 0; i < num_procs - 1; ++i) { num_quadrants_in_proc[i] = (p4est_locidx_t) i + 1; /* type ok */ num_quadrants_on_last -= (p4est_locidx_t) i + 1; /* type ok */ } num_quadrants_in_proc[num_procs - 1] = num_quadrants_on_last; SC_CHECK_ABORT (num_quadrants_on_last > 0, "Negative number of quadrants on the last processor"); /* Save a checksum of the original forest */ crc = p4est_checksum (p4est); /* partition the forest */ (void) p4est_partition_given (p4est, num_quadrants_in_proc); test_pertree (p4est, pertree1, pertree2); /* Double check that we didn't loose any quads */ SC_CHECK_ABORT (crc == p4est_checksum (p4est), "bad checksum, missing a quad"); /* count the actual number of quadrants per proc */ SC_CHECK_ABORT (num_quadrants_in_proc[rank] == p4est->local_num_quadrants, "partition failed, wrong number of quadrants"); /* check user data content */ for (t = p4est->first_local_tree; t <= p4est->last_local_tree; ++t) { tree = p4est_tree_array_index (p4est->trees, t); for (qz = 0; qz < tree->quadrants.elem_count; ++qz) { quad = p4est_quadrant_array_index (&tree->quadrants, qz); user_data = (user_data_t *) quad->p.user_data; sum = quad->x + quad->y + quad->level; SC_CHECK_ABORT (user_data->a == t, "bad user_data, a"); SC_CHECK_ABORT (user_data->sum == sum, "bad user_data, sum"); } } /* do a weighted partition with uniform weights */ p4est_partition (p4est, 0, weight_one); test_pertree (p4est, pertree1, pertree2); SC_CHECK_ABORT (crc == p4est_checksum (p4est), "bad checksum after uniformly weighted partition"); /* copy the p4est */ copy = p4est_copy (p4est, 1); SC_CHECK_ABORT (crc == p4est_checksum (copy), "bad checksum after copy"); /* do a weighted partition with many zero weights */ weight_counter = 0; weight_index = (rank == 1) ? 1342 : 0; p4est_partition (copy, 0, weight_once); test_pertree (copy, pertree1, pertree2); SC_CHECK_ABORT (crc == p4est_checksum (copy), "bad checksum after unevenly weighted partition 1"); /* do a weighted partition with many zero weights */ weight_counter = 0; weight_index = 0; p4est_partition (copy, 0, weight_once); test_pertree (copy, pertree1, pertree2); SC_CHECK_ABORT (crc == p4est_checksum (copy), "bad checksum after unevenly weighted partition 2"); /* do a weighted partition with many zero weights * * Since this is just a test we assume the local number of * quadrants will fit in an int */ weight_counter = 0; weight_index = (rank == num_procs - 1) ? ((int) copy->local_num_quadrants - 1) : 0; p4est_partition (copy, 0, weight_once); test_pertree (copy, pertree1, pertree2); SC_CHECK_ABORT (crc == p4est_checksum (copy), "bad checksum after unevenly weighted partition 3"); /* check user data content */ for (t = copy->first_local_tree; t <= copy->last_local_tree; ++t) { tree = p4est_tree_array_index (copy->trees, t); for (qz = 0; qz < tree->quadrants.elem_count; ++qz) { quad = p4est_quadrant_array_index (&tree->quadrants, qz); user_data = (user_data_t *) quad->p.user_data; sum = quad->x + quad->y + quad->level; SC_CHECK_ABORT (user_data->a == t, "bad user_data, a"); SC_CHECK_ABORT (user_data->sum == sum, "bad user_data, sum"); } } /* Add another test. Overwrites pertree1, pertree2 */ test_partition_circle (mpicomm, connectivity, pertree1, pertree2); /* clean up and exit */ P4EST_FREE (pertree1); P4EST_FREE (pertree2); P4EST_FREE (num_quadrants_in_proc); p4est_destroy (p4est); p4est_destroy (copy); p4est_connectivity_destroy (connectivity); sc_finalize (); mpiret = sc_MPI_Finalize (); SC_CHECK_MPI (mpiret); return 0; }
static void test_partition_circle (sc_MPI_Comm mpicomm, p4est_connectivity_t * connectivity, p4est_gloidx_t * pertree1, p4est_gloidx_t * pertree2) { int i, j; int num_procs; int empty_proc1, empty_proc2; unsigned crc1, crc2; p4est_gloidx_t global_num; p4est_locidx_t *new_counts; p4est_t *p4est, *copy; /* Create a forest and make a copy */ circle_count = 0; p4est = p4est_new_ext (mpicomm, connectivity, 0, 3, 1, sizeof (int), circle_init, NULL); num_procs = p4est->mpisize; test_pertree (p4est, NULL, pertree1); global_num = p4est->global_num_quadrants; crc1 = p4est_checksum (p4est); copy = p4est_copy (p4est, 1); P4EST_ASSERT (p4est_checksum (copy) == crc1); new_counts = P4EST_ALLOC (p4est_locidx_t, num_procs); /* Partition with one empty processor */ if (num_procs > 1) { P4EST_GLOBAL_INFO ("First circle partition\n"); empty_proc1 = num_procs / 3; j = 0; for (i = 0; i < num_procs; ++i) { if (i == empty_proc1) { new_counts[i] = 0; } else { new_counts[i] = p4est_partition_cut_gloidx (global_num, j + 1, num_procs - 1) - p4est_partition_cut_gloidx (global_num, j, num_procs - 1); P4EST_ASSERT (new_counts[i] >= 0); ++j; } } P4EST_ASSERT (j == num_procs - 1); p4est_partition_given (p4est, new_counts); test_pertree (p4est, pertree1, pertree2); crc2 = p4est_checksum (p4est); SC_CHECK_ABORT (crc1 == crc2, "First checksum mismatch"); } /* Partition with two empty processors */ if (num_procs > 2) { P4EST_GLOBAL_INFO ("Second circle partition\n"); empty_proc1 = (2 * num_procs) / 3 - 2; empty_proc2 = (2 * num_procs) / 3; j = 0; for (i = 0; i < num_procs; ++i) { if (i == empty_proc1 || i == empty_proc2) { new_counts[i] = 0; } else { new_counts[i] = p4est_partition_cut_gloidx (global_num, j + 1, num_procs - 2) - p4est_partition_cut_gloidx (global_num, j, num_procs - 2); P4EST_ASSERT (new_counts[i] >= 0); ++j; } } P4EST_ASSERT (j == num_procs - 2); p4est_partition_given (p4est, new_counts); test_pertree (p4est, pertree1, pertree2); crc2 = p4est_checksum (p4est); SC_CHECK_ABORT (crc1 == crc2, "Second checksum mismatch"); } /* Uniform partition */ P4EST_GLOBAL_INFO ("Third circle partition\n"); p4est_partition (p4est, 0, NULL); test_pertree (p4est, pertree1, pertree2); crc2 = p4est_checksum (p4est); SC_CHECK_ABORT (crc1 == crc2, "Third checksum mismatch"); SC_CHECK_ABORT (p4est_is_equal (p4est, copy, 1), "Forest mismatch"); P4EST_FREE (new_counts); p4est_destroy (copy); p4est_destroy (p4est); }