#ifndef P4_TO_P8
#include <p4est_nodes.h>
#else
#include <p8est_nodes.h>
#endif
#include <sc_flops.h>
#include <sc_statistics.h>
#include <sc_options.h>
typedef enum
{
P4EST_CONFIG_NULL,
P4EST_CONFIG_UNIT,
P4EST_CONFIG_PERIODIC,
#ifndef P4_TO_P8
P4EST_CONFIG_THREE,
P4EST_CONFIG_MOEBIUS,
P4EST_CONFIG_STAR
#else
P4EST_CONFIG_ROTWRAP,
P4EST_CONFIG_TWOCUBES,
P4EST_CONFIG_ROTCUBES,
P4EST_CONFIG_SHELL
#endif
}
timings_config_t;
enum
{
TIMINGS_REFINE,
TIMINGS_BALANCE,
TIMINGS_BALANCE_A,
TIMINGS_BALANCE_COMM,
TIMINGS_BALANCE_B,
TIMINGS_BALANCE_A_COUNT_IN,
TIMINGS_BALANCE_A_COUNT_OUT,
TIMINGS_BALANCE_COMM_SENT,
TIMINGS_BALANCE_COMM_NZPEERS,
TIMINGS_BALANCE_B_COUNT_IN,
TIMINGS_BALANCE_B_COUNT_OUT,
TIMINGS_BALANCE_RANGES,
TIMINGS_BALANCE_NOTIFY,
TIMINGS_BALANCE_NOTIFY_ALLGATHER,
TIMINGS_BALANCE_A_ZERO_SENDS,
TIMINGS_BALANCE_A_ZERO_RECEIVES,
TIMINGS_BALANCE_B_ZERO_SENDS,
TIMINGS_BALANCE_B_ZERO_RECEIVES,
TIMINGS_REBALANCE,
TIMINGS_REBALANCE_A,
TIMINGS_REBALANCE_COMM,
TIMINGS_REBALANCE_B,
TIMINGS_REBALANCE_A_COUNT_IN,
TIMINGS_REBALANCE_A_COUNT_OUT,
TIMINGS_REBALANCE_COMM_SENT,
TIMINGS_REBALANCE_COMM_NZPEERS,
TIMINGS_REBALANCE_B_COUNT_IN,
TIMINGS_REBALANCE_B_COUNT_OUT,
TIMINGS_PARTITION,
TIMINGS_GHOSTS,
TIMINGS_NODES,
TIMINGS_TRILINEAR_OBSOLETE,
TIMINGS_REPARTITION,
TIMINGS_LNODES,
TIMINGS_LNODES3,
TIMINGS_LNODES7,
TIMINGS_NUM_STATS
};
typedef struct
{
timings_config_t config;
int mpisize;
int level;
unsigned checksum;
}
timings_regression_t;
typedef struct
{
sc_MPI_Comm mpicomm;
int mpisize;
int mpirank;
}
mpi_context_t;
static int refine_level = 0;
static int level_shift = 0;
static const timings_regression_t regression_oldschool[] =
{
#ifndef P4_TO_P8
{ P4EST_CONFIG_UNIT, 1, 10, 0x6e3e83c4U },
{ P4EST_CONFIG_UNIT, 1, 11, 0x334bc3deU },
{ P4EST_CONFIG_UNIT, 64, 14, 0xad908ce4U },
{ P4EST_CONFIG_UNIT, 256, 15, 0x9e7da646U },
{ P4EST_CONFIG_STAR, 1, 6, 0x14107b57U },
{ P4EST_CONFIG_STAR, 4, 6, 0x14107b57U },
{ P4EST_CONFIG_STAR, 52, 13, 0xc86c74d9U },
{ P4EST_CONFIG_STAR, 64, 13, 0xc86c74d9U },
#else
{ P4EST_CONFIG_UNIT, 1, 5, 0xe1ffa67bU },
{ P4EST_CONFIG_UNIT, 1, 6, 0x2cad814dU },
{ P4EST_CONFIG_UNIT, 3, 8, 0xeb252238U },
{ P4EST_CONFIG_PERIODIC, 1, 5, 0x99874fedU },
{ P4EST_CONFIG_PERIODIC, 2, 5, 0x575af6d5U },
{ P4EST_CONFIG_PERIODIC, 7, 6, 0xbc35524aU },
{ P4EST_CONFIG_ROTWRAP, 2, 6, 0x372f7402U },
{ P4EST_CONFIG_ROTWRAP, 7, 6, 0xa2f1ee48U },
{ P4EST_CONFIG_TWOCUBES, 5, 6, 0xa8b1f54eU },
{ P4EST_CONFIG_TWOCUBES, 8, 5, 0x98d3579dU },
{ P4EST_CONFIG_ROTCUBES, 1, 5, 0x404e4aa8U },
{ P4EST_CONFIG_ROTCUBES, 7, 6, 0x4c381706U },
{ P4EST_CONFIG_SHELL, 1, 4, 0x8c56f159U },
{ P4EST_CONFIG_SHELL, 3, 5, 0xafbc4f8cU },
{ P4EST_CONFIG_SHELL, 5, 6, 0xf6d9efb8U },
#endif
{ P4EST_CONFIG_NULL, 0, 0, 0 }
};
static const timings_regression_t regression_latest[] =
{
#ifndef P4_TO_P8
{ P4EST_CONFIG_UNIT, 1, 10, 0x6e3e83c4U },
{ P4EST_CONFIG_UNIT, 1, 11, 0x334bc3deU },
{ P4EST_CONFIG_UNIT, 64, 14, 0xad908ce4U },
{ P4EST_CONFIG_UNIT, 256, 15, 0x9e7da646U },
{ P4EST_CONFIG_STAR, 1, 6, 0x14107b57U },
{ P4EST_CONFIG_STAR, 4, 6, 0x14107b57U },
{ P4EST_CONFIG_STAR, 52, 13, 0xc86c74d9U },
{ P4EST_CONFIG_STAR, 64, 13, 0xc86c74d9U },
#else
{ P4EST_CONFIG_UNIT, 1, 5, 0xc2012e84U },
{ P4EST_CONFIG_UNIT, 1, 6, 0x2cad814dU },
{ P4EST_CONFIG_UNIT, 3, 8, 0xeb252238U },
{ P4EST_CONFIG_PERIODIC, 1, 5, 0x2776c9b7U },
{ P4EST_CONFIG_PERIODIC, 2, 5, 0x2776c9b7U },
{ P4EST_CONFIG_PERIODIC, 7, 6, 0x4f281079U },
{ P4EST_CONFIG_ROTWRAP, 2, 6, 0x372f7402U },
{ P4EST_CONFIG_ROTWRAP, 7, 6, 0x372f7402U },
{ P4EST_CONFIG_TWOCUBES, 5, 6, 0xa8b1f54eU },
{ P4EST_CONFIG_TWOCUBES, 8, 5, 0x0aad11d0U },
{ P4EST_CONFIG_ROTCUBES, 1, 5, 0x404e4aa8U },
{ P4EST_CONFIG_ROTCUBES, 7, 6, 0x4c381706U },
{ P4EST_CONFIG_SHELL, 1, 4, 0x8c56f159U },
{ P4EST_CONFIG_SHELL, 3, 5, 0xafbc4f8cU },
{ P4EST_CONFIG_SHELL, 5, 6, 0xf6d9efb8U },
#endif
{ P4EST_CONFIG_NULL, 0, 0, 0 }
};
static int
{
int qid;
if ((
int) q->
level >= refine_level) {
return 0;
}
if ((
int) q->
level < refine_level - level_shift) {
return 1;
}
return (qid == 0 || qid == 3
#ifdef P4_TO_P8
|| qid == 5 || qid == 6
#endif
);
}
int
main (int argc, char **argv)
{
int i;
int mpiret;
int wrongusage;
unsigned crc, gcrc;
const char *config_name;
const char *load_name;
const timings_regression_t *r, *regression;
timings_config_t config;
sc_statinfo_t stats[TIMINGS_NUM_STATS];
sc_flopinfo_t fi, snapshot;
mpi_context_t mpi_context, *mpi = &mpi_context;
sc_options_t *opt;
int overlap;
int subtree;
int borders;
int max_ranges;
int use_ranges, use_ranges_notify, use_balance_verify;
int oldschool, generate;
int first_argc;
int test_multiple_orders;
int skip_nodes, skip_lnodes;
int repartition_lnodes;
mpiret = sc_MPI_Init (&argc, &argv);
SC_CHECK_MPI (mpiret);
mpi->mpicomm = sc_MPI_COMM_WORLD;
mpiret = sc_MPI_Comm_size (mpi->mpicomm, &mpi->mpisize);
SC_CHECK_MPI (mpiret);
mpiret = sc_MPI_Comm_rank (mpi->mpicomm, &mpi->mpirank);
SC_CHECK_MPI (mpiret);
sc_init (mpi->mpicomm, 1, 1, NULL, SC_LP_DEFAULT);
#ifndef P4EST_ENABLE_DEBUG
sc_set_log_defaults (NULL, NULL, SC_LP_STATISTICS);
#endif
P4EST_GLOBAL_PRODUCTIONF (
"Size of %dtant: %lld bytes\n",
P4EST_DIM,
opt = sc_options_new (argv[0]);
sc_options_add_switch (opt, 'o', "new-balance-overlap", &overlap,
"use the new balance overlap algorithm");
sc_options_add_switch (opt, 's', "new-balance-subtree", &subtree,
"use the new balance subtree algorithm");
sc_options_add_switch (opt, 'b', "new-balance-borders", &borders,
"use borders in balance");
sc_options_add_int (opt, 'm', "max-ranges", &max_ranges, -1,
"override p4est_num_ranges");
sc_options_add_switch (opt, 'r', "ranges", &use_ranges,
"use ranges in balance");
sc_options_add_switch (opt, 't', "ranges-notify", &use_ranges_notify,
"use both ranges and notify");
sc_options_add_switch (opt, 'y', "balance-verify", &use_balance_verify,
"use verifications in balance");
sc_options_add_int (opt, 'l', "level", &refine_level, 0,
"initial refine level");
#ifndef P4_TO_P8
sc_options_add_string (opt, 'c', "configuration", &config_name, "unit",
"configuration: unit|periodic|three|moebius|star");
#else
sc_options_add_string (opt, 'c', "configuration", &config_name, "unit",
"configuration: unit|periodic|rotwrap|twocubes|rotcubes|shell");
#endif
sc_options_add_bool (opt, 0, "oldschool", &oldschool, 0,
"Use original p4est_new call");
sc_options_add_switch (opt, 0, "generate", &generate,
"Generate regression commands");
sc_options_add_string (opt, 'f', "load-forest", &load_name, NULL,
sc_options_add_switch (opt, 0, "multiple-lnodes",
&test_multiple_orders,
"Also time lnodes for orders 2, 4, and 8");
sc_options_add_switch (opt, 0, "skip-nodes", &skip_nodes, "Skip nodes");
sc_options_add_switch (opt, 0, "skip-lnodes", &skip_lnodes, "Skip lnodes");
sc_options_add_switch (opt, 0, "repartition-lnodes",
&repartition_lnodes,
"Repartition to load-balance lnodes");
opt, argc, argv);
if (first_argc < 0 || first_argc != argc) {
return 1;
}
if (max_ranges < -1 || max_ranges == 0) {
P4EST_GLOBAL_LERROR ("The -m / --max-ranges option must be positive\n");
return 1;
}
if (skip_lnodes) {
if (test_multiple_orders) {
SC_GLOBAL_PRODUCTION
("Warning: cannot test multiple lnode orders if --skip-lnodes is given.\n");
test_multiple_orders = 0;
}
}
wrongusage = 0;
config = P4EST_CONFIG_NULL;
if (!strcmp (config_name, "unit")) {
config = P4EST_CONFIG_UNIT;
}
else if (!strcmp (config_name, "periodic")) {
config = P4EST_CONFIG_PERIODIC;
}
#ifndef P4_TO_P8
else if (!strcmp (config_name, "three")) {
config = P4EST_CONFIG_THREE;
}
else if (!strcmp (config_name, "moebius")) {
config = P4EST_CONFIG_MOEBIUS;
}
else if (!strcmp (config_name, "star")) {
config = P4EST_CONFIG_STAR;
}
#else
else if (!strcmp (config_name, "rotwrap")) {
config = P4EST_CONFIG_ROTWRAP;
}
else if (!strcmp (config_name, "twocubes")) {
config = P4EST_CONFIG_TWOCUBES;
}
else if (!strcmp (config_name, "rotcubes")) {
config = P4EST_CONFIG_ROTCUBES;
}
else if (!strcmp (config_name, "shell")) {
config = P4EST_CONFIG_SHELL;
}
#endif
else if (load_name != NULL) {
config_name = load_name;
}
else {
wrongusage = 1;
}
if (wrongusage) {
P4EST_GLOBAL_LERRORF ("Wrong configuration name given: %s\n",
config_name);
sc_abort_collective ("Usage error");
}
level_shift = 4;
P4EST_GLOBAL_STATISTICSF
("Processors %d configuration %s level %d shift %d\n", mpi->mpisize,
config_name, refine_level, level_shift);
mpiret = sc_MPI_Barrier (mpi->mpicomm);
SC_CHECK_MPI (mpiret);
sc_flops_start (&fi);
regression = NULL;
if (load_name == NULL) {
#ifndef P4_TO_P8
if (config == P4EST_CONFIG_PERIODIC) {
}
else if (config == P4EST_CONFIG_THREE) {
}
else if (config == P4EST_CONFIG_MOEBIUS) {
}
else if (config == P4EST_CONFIG_STAR) {
}
else {
}
#else
if (config == P4EST_CONFIG_PERIODIC) {
}
else if (config == P4EST_CONFIG_ROTWRAP) {
}
else if (config == P4EST_CONFIG_TWOCUBES) {
}
else if (config == P4EST_CONFIG_ROTCUBES) {
}
else if (config == P4EST_CONFIG_SHELL) {
}
else {
}
#endif
if (oldschool) {
regression = regression_oldschool;
15, 0, 0, 0, NULL, NULL);
}
else {
regression = regression_latest;
0, refine_level - level_shift, 1, 0, NULL, NULL);
}
if (generate) {
const char *config_gname = NULL;
P4EST_GLOBAL_PRODUCTION ("Checksum regression tests available:\n");
for (r = regression; r->config != P4EST_CONFIG_NULL; ++r) {
switch (r->config) {
case P4EST_CONFIG_PERIODIC:
config_gname = "periodic";
break;
#ifndef P4_TO_P8
case P4EST_CONFIG_THREE:
config_gname = "three";
break;
case P4EST_CONFIG_MOEBIUS:
config_gname = "moebius";
break;
case P4EST_CONFIG_STAR:
config_gname = "star";
break;
#else
case P4EST_CONFIG_ROTWRAP:
config_gname = "rotwrap";
break;
case P4EST_CONFIG_TWOCUBES:
config_gname = "twocubes";
break;
case P4EST_CONFIG_ROTCUBES:
config_gname = "rotcubes";
break;
case P4EST_CONFIG_SHELL:
config_gname = "shell";
break;
#endif
default:
config_gname = "unit";
break;
}
P4EST_GLOBAL_PRODUCTIONF ("mpirun -np %3d %s%s -c %10s -l %2d\n",
oldschool ? " --oldschool" : "",
config_gname, r->level);
}
}
}
else {
}
P4EST_GLOBAL_STATISTICSF
("Balance: new overlap %d new subtree %d borders %d\n", overlap,
(overlap && subtree), (overlap && borders));
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_REFINE], snapshot.iwtime, "Refine");
#ifdef P4EST_TIMINGS_VTK
#endif
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_BALANCE], snapshot.iwtime, "Balance");
sc_stats_set1 (&stats[TIMINGS_BALANCE_A],
sc_stats_set1 (&stats[TIMINGS_BALANCE_COMM],
sc_stats_set1 (&stats[TIMINGS_BALANCE_B],
sc_stats_set1 (&stats[TIMINGS_BALANCE_A_COUNT_IN],
"Balance A count inlist");
sc_stats_set1 (&stats[TIMINGS_BALANCE_A_COUNT_OUT],
"Balance A count outlist");
sc_stats_set1 (&stats[TIMINGS_BALANCE_COMM_SENT],
"Balance sent second round");
sc_stats_set1 (&stats[TIMINGS_BALANCE_COMM_NZPEERS],
"Balance nonzero peers second round");
sc_stats_set1 (&stats[TIMINGS_BALANCE_B_COUNT_IN],
"Balance B count inlist");
sc_stats_set1 (&stats[TIMINGS_BALANCE_B_COUNT_OUT],
"Balance B count outlist");
sc_stats_set1 (&stats[TIMINGS_BALANCE_RANGES],
sc_stats_set1 (&stats[TIMINGS_BALANCE_NOTIFY],
sc_stats_set1 (&stats[TIMINGS_BALANCE_NOTIFY_ALLGATHER],
"Balance time for notify_allgather");
sc_stats_set1 (&stats[TIMINGS_BALANCE_A_ZERO_RECEIVES],
"Balance A zero receives");
sc_stats_set1 (&stats[TIMINGS_BALANCE_A_ZERO_SENDS],
"Balance A zero sends");
sc_stats_set1 (&stats[TIMINGS_BALANCE_B_ZERO_RECEIVES],
"Balance B zero receives");
sc_stats_set1 (&stats[TIMINGS_BALANCE_B_ZERO_SENDS],
"Balance B zero sends");
#ifdef P4EST_TIMINGS_VTK
#endif
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_REBALANCE], snapshot.iwtime, "Rebalance");
sc_stats_set1 (&stats[TIMINGS_REBALANCE_A],
sc_stats_set1 (&stats[TIMINGS_REBALANCE_COMM],
sc_stats_set1 (&stats[TIMINGS_REBALANCE_B],
sc_stats_set1 (&stats[TIMINGS_REBALANCE_A_COUNT_IN],
"Rebalance A count inlist");
sc_stats_set1 (&stats[TIMINGS_REBALANCE_A_COUNT_OUT],
"Rebalance A count outlist");
sc_stats_set1 (&stats[TIMINGS_REBALANCE_COMM_SENT],
"Rebalance sent second round");
sc_stats_set1 (&stats[TIMINGS_REBALANCE_COMM_NZPEERS],
"Rebalance nonzero peers second round");
sc_stats_set1 (&stats[TIMINGS_REBALANCE_B_COUNT_IN],
"Rebalance B count inlist");
sc_stats_set1 (&stats[TIMINGS_REBALANCE_B_COUNT_OUT],
"Rebalance B count outlist");
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_PARTITION], snapshot.iwtime, "Partition");
#ifdef P4EST_TIMINGS_VTK
#endif
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_GHOSTS], snapshot.iwtime, "Ghost layer");
if (!skip_nodes) {
sc_flops_snap (&fi, &snapshot);
nodes = p4est_nodes_new (
p4est, ghost);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_NODES], snapshot.iwtime, "Nodes");
}
else {
sc_stats_set1 (&stats[TIMINGS_NODES], 0., "Nodes");
}
sc_stats_set1 (&stats[TIMINGS_TRILINEAR_OBSOLETE], 0., "Unused");
if (!skip_nodes) {
p4est_nodes_destroy (nodes);
}
if (!skip_lnodes && repartition_lnodes) {
}
if (!skip_lnodes) {
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_LNODES], snapshot.iwtime, "L-Nodes");
}
else {
sc_stats_set1 (&stats[TIMINGS_LNODES], 0., "L-Nodes");
}
if (test_multiple_orders) {
if (repartition_lnodes) {
}
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_LNODES3], snapshot.iwtime, "L-Nodes 3");
if (repartition_lnodes) {
}
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_LNODES7], snapshot.iwtime, "L-Nodes 7");
}
else {
sc_stats_set1 (&stats[TIMINGS_LNODES3], 0., "L-Nodes 3");
sc_stats_set1 (&stats[TIMINGS_LNODES7], 0., "L-Nodes 7");
}
prev_quadrant = next_quadrant;
quadrant_counts[i] = (
p4est_locidx_t) (next_quadrant - prev_quadrant);
}
}
sc_flops_snap (&fi, &snapshot);
sc_flops_shot (&fi, &snapshot);
sc_stats_set1 (&stats[TIMINGS_REPARTITION], snapshot.iwtime, "Repartition");
P4EST_GLOBAL_PRODUCTIONF
(
"Done " P4EST_STRING "_partition_given shipped %lld quadrants %.3g%%\n",
(long long) global_shipped,
if (regression != NULL && mpi->mpirank == 0) {
for (r = regression; r->config != P4EST_CONFIG_NULL; ++r) {
if (r->config != config || r->mpisize != mpi->mpisize
|| r->level != refine_level)
continue;
SC_CHECK_ABORT (crc == r->checksum, "Checksum mismatch");
P4EST_GLOBAL_INFO ("Checksum regression OK\n");
break;
}
}
P4EST_GLOBAL_STATISTICSF ("Processors %d level %d shift %d"
" checksums 0x%08x 0x%08x\n",
mpi->mpisize, refine_level, level_shift,
crc, gcrc);
P4EST_GLOBAL_STATISTICSF ("Level %d refined to %lld balanced to %lld\n",
refine_level, (long long) count_refined,
(long long) count_balanced);
sc_stats_compute (mpi->mpicomm, TIMINGS_NUM_STATS, stats);
TIMINGS_NUM_STATS, stats, 1, 1);
sc_options_destroy (opt);
sc_finalize ();
mpiret = sc_MPI_Finalize ();
SC_CHECK_MPI (mpiret);
return 0;
}
void p4est_refine(p4est_t *p4est, int refine_recursive, p4est_refine_t refine_fn, p4est_init_t init_fn)
Refine a forest.
void p4est_destroy(p4est_t *p4est)
Destroy a p4est.
void p4est_partition(p4est_t *p4est, int allow_for_coarsening, p4est_weight_t weight_fn)
Equally partition the forest.
void p4est_balance(p4est_t *p4est, p4est_connect_type_t btype, p4est_init_t init_fn)
2:1 balance the size differences of neighboring elements in a forest.
p4est_t * p4est_load(const char *filename, sc_MPI_Comm mpicomm, size_t data_size, int load_data, void *user_pointer, p4est_connectivity_t **connectivity)
Load the complete connectivity/p4est structure from disk.
unsigned p4est_checksum(p4est_t *p4est)
Compute the checksum for a forest.
Routines for managing quadrants as elements of trees and subtrees.
p4est_gloidx_t p4est_partition_given(p4est_t *p4est, const p4est_locidx_t *num_quadrants_in_proc)
Partition p4est given the number of quadrants per proc.
#define P4EST_FREE(p)
free an allocated array
Definition: p4est_base.h:210
#define P4EST_ALLOC(t, n)
allocate a t-array with n elements
Definition: p4est_base.h:199
SC_DLL_PUBLIC int p4est_package_id
The package id for p4est within libsc.
int32_t p4est_topidx_t
Typedef for counting topological entities (trees, tree vertices).
Definition: p4est_base.h:93
int32_t p4est_locidx_t
Typedef for processor-local indexing of quadrants and nodes.
Definition: p4est_base.h:106
void p4est_init(sc_log_handler_t log_handler, int log_threshold)
Registers p4est with the SC Library and sets the logging behavior.
#define P4EST_ALLOC_ZERO(t, n)
allocate a t-array with n elements and zero
Definition: p4est_base.h:202
int64_t p4est_gloidx_t
Typedef for globally unique indexing of quadrants.
Definition: p4est_base.h:118
Routines for manipulating quadrants (neighbors, parents, children, etc.)
int p4est_quadrant_child_id(const p4est_quadrant_t *q)
Compute the position of this child within its siblings.
p4est_connectivity_t * p4est_connectivity_new_moebius(void)
Create a connectivity structure for a five-tree moebius band.
void p4est_connectivity_destroy(p4est_connectivity_t *connectivity)
Destroy a connectivity structure.
p4est_connectivity_t * p4est_connectivity_new_star(void)
Create a connectivity structure for a six-tree star.
#define P4EST_DIM
The spatial dimension.
Definition: p4est_connectivity.h:71
p4est_connectivity_t * p4est_connectivity_new_unitsquare(void)
Create a connectivity structure for the unit square.
#define P4EST_STRING
p4est identification string
Definition: p4est_connectivity.h:94
p4est_connectivity_t * p4est_connectivity_new_corner(void)
Create a connectivity structure for a three-tree mesh around a corner.
@ P4EST_CONNECT_FULL
= CORNER.
Definition: p4est_connectivity.h:119
p4est_connectivity_t * p4est_connectivity_new_periodic(void)
Create a connectivity structure for an all-periodic unit square.
Interface routines with extended capabilities.
p4est_t * p4est_new_ext(sc_MPI_Comm mpicomm, p4est_connectivity_t *connectivity, p4est_locidx_t min_quadrants, int min_level, int fill_uniform, size_t data_size, p4est_init_t init_fn, void *user_pointer)
Create a new forest.
Passing quadrants and data to neighboring processes.
p4est_ghost_t * p4est_ghost_new(p4est_t *p4est, p4est_connect_type_t btype)
Builds the ghost layer.
void p4est_ghost_destroy(p4est_ghost_t *ghost)
Frees all memory used for the ghost layer.
unsigned p4est_ghost_checksum(p4est_t *p4est, p4est_ghost_t *ghost)
Compute the parallel checksum of a ghost layer.
Generate Lobatto node numbers for any degree.
void p4est_partition_lnodes(p4est_t *p4est, p4est_ghost_t *ghost, int degree, int partition_for_coarsening)
Partition using weights based on the number of nodes assigned to each element in lnodes.
void p4est_lnodes_destroy(p4est_lnodes_t *lnodes)
Free all memory in a previously constructed lnodes structure.
p4est_lnodes_t * p4est_lnodes_new(p4est_t *p4est, p4est_ghost_t *ghost_layer, int degree)
Create a tensor-product Lobatto node structure for a given degree.
Routines for printing a forest and associated fields to VTK format.
void p4est_vtk_write_file(p4est_t *p4est, p4est_geometry_t *geom, const char *filename)
Write the p4est in VTK format.
Routines for managing quadrants as elements of trees and subtrees.
Routines for manipulating quadrants (neighbors, parents, children, etc.)
p8est_connectivity_t * p8est_connectivity_new_periodic(void)
Create a connectivity structure for an all-periodic unit cube.
p8est_connectivity_t * p8est_connectivity_new_unitcube(void)
Create a connectivity structure for the unit cube.
p8est_connectivity_t * p8est_connectivity_new_twocubes(void)
Create a connectivity structure that contains two cubes.
p8est_connectivity_t * p8est_connectivity_new_shell(void)
Create a connectivity structure that builds a spherical shell.
p8est_connectivity_t * p8est_connectivity_new_rotwrap(void)
Create a connectivity structure for a mostly periodic unit cube.
p8est_connectivity_t * p8est_connectivity_new_rotcubes(void)
Create a connectivity structure that contains a few cubes.
Interface routines with extended capabilities.
Passing quadrants and data to neighboring processes.
Generate Lobatto node numbers for any degree.
Routines for printing a forest and associated fields to VTK format.
This structure holds the 2D inter-tree connectivity information.
Definition: p4est_connectivity.h:190
Quadrants that neighbor the local domain.
Definition: p4est_ghost.h:46
Data pertaining to selecting, inspecting, and profiling algorithms.
Definition: p4est_extended.h:62
double balance_notify
time spent in sc_notify
Definition: p4est_extended.h:84
double balance_notify_allgather
time spent in sc_notify_allgather
Definition: p4est_extended.h:86
int balance_max_ranges
If positive and smaller than p4est_num ranges, overrides it.
Definition: p4est_extended.h:72
int use_balance_ranges_notify
If true, call both sc_ranges and sc_notify and verify consistency.
Definition: p4est_extended.h:68
double balance_ranges
time spent in sc_ranges
Definition: p4est_extended.h:83
int use_balance_ranges
Use sc_ranges to determine the asymmetric communication pattern.
Definition: p4est_extended.h:65
int use_balance_verify
Verify sc_ranges and/or sc_notify as applicable.
Definition: p4est_extended.h:70
Store a parallel numbering of Lobatto points of a given degree > 0.
Definition: p4est_lnodes.h:145
This structure holds complete parallel node information.
Definition: p4est_nodes.h:141
The 2D quadrant datatype.
Definition: p4est.h:76
int8_t level
level of refinement
Definition: p4est.h:80
The p4est forest datatype.
Definition: p4est.h:150
int mpisize
number of MPI processes
Definition: p4est.h:152
p4est_inspect_t * inspect
algorithmic switches
Definition: p4est.h:185
p4est_gloidx_t global_num_quadrants
number of quadrants on all trees on all processors
Definition: p4est.h:169