sc_keyvalue_t * sc_keyvalue_new () { sc_keyvalue_t *kv; kv = SC_ALLOC (sc_keyvalue_t, 1); kv->hash = sc_hash_new (sc_keyvalue_entry_hash, sc_keyvalue_entry_equal, NULL, NULL); kv->value_allocator = sc_mempool_new (sizeof (sc_keyvalue_entry_t)); return kv; }
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; }
trilinear_mesh_t * p8est_trilinear_mesh_new_from_nodes (p4est_t * p4est, p4est_nodes_t * nodes) { const int num_procs = p4est->mpisize; const int rank = p4est->mpirank; int mpiret; int k, owner; #ifdef P4EST_DEBUG int prev_owner = 0; int64_t prev_fvnid = -1; #endif int *sharers; size_t current, zz, num_sharers; int32_t e, n, local_owned_end; int64_t global_borrowed, global_shared; int64_t local_counts[5], global_counts[5]; int32link_t *lynk, **tail; p4est_topidx_t which_tree; p4est_locidx_t *local_node, *shared_offsets; p4est_tree_t *tree; p4est_quadrant_t *q; p4est_indep_t *in; p8est_hang4_t *fh; p8est_hang2_t *eh; trilinear_elem_t *elem; trilinear_anode_t *anode; trilinear_dnode_t *dnode; trilinear_mesh_t *mesh; trilinear_mesh_pid_t *elem_pids; trilinear_mesh_pid_t *node_pids; sc_recycle_array_t *rarr; P4EST_GLOBAL_PRODUCTIONF ("Into trilinear_mesh_extract with %lld total elements\n", (long long) p4est->global_num_quadrants); /* Allocate output data structure. */ mesh = P4EST_ALLOC_ZERO (trilinear_mesh_t, 1); memset (mesh, -1, sizeof (*mesh)); shared_offsets = nodes->shared_offsets; /* Assign local counts. */ P4EST_ASSERT (nodes->num_local_quadrants == p4est->local_num_quadrants); mesh->local_elem_num = p4est->local_num_quadrants; mesh->local_anode_num = nodes->indep_nodes.elem_count; mesh->local_dnode_num = nodes->face_hangings.elem_count + nodes->edge_hangings.elem_count; mesh->local_onode_num = nodes->num_owned_indeps; mesh->local_owned_offset = nodes->offset_owned_indeps; mesh->local_node_num = mesh->local_anode_num + mesh->local_dnode_num; local_owned_end = mesh->local_owned_offset + mesh->local_onode_num; /* Communicate global counts. */ local_counts[0] = mesh->local_elem_num; local_counts[1] = mesh->local_anode_num; local_counts[2] = mesh->local_onode_num; local_counts[3] = mesh->local_dnode_num; local_counts[4] = nodes->num_owned_shared; mpiret = MPI_Allreduce (local_counts, global_counts, 5, MPI_LONG_LONG_INT, MPI_SUM, p4est->mpicomm); SC_CHECK_MPI (mpiret); P4EST_ASSERT (global_counts[0] == p4est->global_num_quadrants); mesh->total_elem_num = global_counts[0]; global_borrowed = global_counts[1] - global_counts[2]; mesh->total_anode_num = global_counts[2]; mesh->total_dnode_num = global_counts[3]; global_shared = global_counts[4]; mesh->total_node_num = mesh->total_anode_num + mesh->total_dnode_num; /* Allocate the mesh memory. */ mesh->elem_table = P4EST_ALLOC (trilinear_elem_t, mesh->local_elem_num); mesh->node_table = P4EST_ALLOC (trilinear_node_t, mesh->local_node_num); mesh->fvnid_count_table = P4EST_ALLOC (int64_t, num_procs + 1); mesh->fvnid_interval_table = P4EST_ALLOC (int64_t, num_procs + 1); mesh->all_fvnid_start = mesh->fvnid_interval_table; mesh->sharer_pool = sc_mempool_new (sizeof (int32link_t)); mesh->elem_pids = P4EST_ALLOC (trilinear_mesh_pid_t, mesh->local_node_num); mesh->node_pids = P4EST_ALLOC (trilinear_mesh_pid_t, mesh->local_node_num); /* Assign global free variable information. */ mesh->fvnid_interval_table[0] = 0; for (k = 0; k < num_procs; ++k) { mesh->fvnid_interval_table[k + 1] = mesh->fvnid_interval_table[k] + (mesh->fvnid_count_table[k] = nodes->global_owned_indeps[k]); } mesh->fvnid_count_table[num_procs] = -1; mesh->global_fvnid_num = mesh->fvnid_interval_table[num_procs]; mesh->global_fvnid_start = 0; mesh->global_fvnid_end = mesh->global_fvnid_num - 1; P4EST_ASSERT (mesh->global_fvnid_num == mesh->total_anode_num); /* Assign element information. */ local_node = nodes->local_nodes; which_tree = p4est->first_local_tree; elem_pids = mesh->elem_pids; if (which_tree >= 0) { tree = p4est_tree_array_index (p4est->trees, which_tree); current = 0; for (e = 0; e < mesh->local_elem_num; ++e) { if (current == tree->quadrants.elem_count) { ++which_tree; tree = p4est_tree_array_index (p4est->trees, which_tree); current = 0; } q = p4est_quadrant_array_index (&tree->quadrants, current); elem = mesh->elem_table + e; for (k = 0; k < P4EST_CHILDREN; ++k) { elem->local_node_id[k] = *local_node++; } elem->lx = (tick_t) q->x; elem->ly = (tick_t) q->y; elem->lz = (tick_t) q->z; elem->size = P4EST_QUADRANT_LEN (q->level); elem->data = q->p.user_data; elem_pids[e] = (trilinear_mesh_pid_t) which_tree; ++current; } P4EST_ASSERT (which_tree == p4est->last_local_tree); P4EST_ASSERT (current == tree->quadrants.elem_count); } /* Assign anchored node information. */ mesh->anode_table = mesh->node_table; mesh->onode_table = mesh->node_table + mesh->local_owned_offset; mesh->dnode_table = mesh->node_table + mesh->local_anode_num; node_pids = mesh->node_pids; for (n = 0; n < mesh->local_anode_num; ++n) { anode = &mesh->node_table[n].anchored; in = (p4est_indep_t *) sc_array_index (&nodes->indep_nodes, (size_t) n); anode->point.x = in->x; anode->point.y = in->y; anode->point.z = in->z; node_pids[n] = (trilinear_mesh_pid_t) in->p.piggy3.which_tree; if (n < mesh->local_owned_offset) { owner = nodes->nonlocal_ranks[n]; P4EST_ASSERT (owner < rank && owner >= prev_owner); } else if (n >= local_owned_end) { owner = nodes->nonlocal_ranks[n - mesh->local_onode_num]; P4EST_ASSERT (owner > rank && owner >= prev_owner); } else { owner = rank; } anode->fvnid = mesh->all_fvnid_start[owner] + in->p.piggy3.local_num; P4EST_ASSERT (anode->fvnid > prev_fvnid); if (in->pad8 == 0) { P4EST_ASSERT (in->pad16 == -1); P4EST_ASSERT (shared_offsets == NULL || shared_offsets[n] == -1); anode->share = NULL; } else { P4EST_ASSERT (in->pad8 > 0); num_sharers = (size_t) in->pad8; rarr = (sc_recycle_array_t *) sc_array_index (&nodes->shared_indeps, num_sharers - 1); if (nodes->shared_offsets == NULL) { P4EST_ASSERT (in->pad16 >= 0); zz = (size_t) in->pad16; } else { P4EST_ASSERT (in->pad16 == -1); zz = (size_t) shared_offsets[n]; } sharers = (int *) sc_array_index (&rarr->a, zz); tail = &anode->share; for (zz = 0; zz < num_sharers; ++zz) { *tail = lynk = (int32link_t *) sc_mempool_alloc (mesh->sharer_pool); lynk->id = (int32_t) sharers[zz]; tail = &lynk->next; } *tail = NULL; } #ifdef P4EST_DEBUG prev_owner = owner; prev_fvnid = anode->fvnid; #endif } /* Assign face hanging node information. */ for (zz = 0; zz < nodes->face_hangings.elem_count; ++n, ++zz) { dnode = &mesh->node_table[n].dangling; fh = (p8est_hang4_t *) sc_array_index (&nodes->face_hangings, zz); dnode->point.x = fh->x; dnode->point.y = fh->y; dnode->point.z = fh->z; dnode->type = 0; /* Not used in Rhea. */ dnode->local_anode_id[0] = fh->p.piggy.depends[0]; dnode->local_anode_id[1] = fh->p.piggy.depends[1]; dnode->local_anode_id[2] = fh->p.piggy.depends[2]; dnode->local_anode_id[3] = fh->p.piggy.depends[3]; node_pids[n] = (trilinear_mesh_pid_t) fh->p.piggy.which_tree; } /* Assign edge hanging node information. */ for (zz = 0; zz < nodes->edge_hangings.elem_count; ++n, ++zz) { dnode = &mesh->node_table[n].dangling; eh = (p8est_hang2_t *) sc_array_index (&nodes->edge_hangings, zz); dnode->point.x = eh->x; dnode->point.y = eh->y; dnode->point.z = eh->z; dnode->type = 0; /* Not used in Rhea. */ dnode->local_anode_id[0] = eh->p.piggy.depends[0]; dnode->local_anode_id[1] = eh->p.piggy.depends[1]; dnode->local_anode_id[2] = dnode->local_anode_id[3] = -1; node_pids[n] = (trilinear_mesh_pid_t) eh->p.piggy.which_tree; } P4EST_ASSERT (n == mesh->local_node_num); /* Assign the remaining variables. */ mesh->mpicomm = p4est->mpicomm; mesh->mpisize = (int32_t) num_procs; mesh->mpirank = (int32_t) rank; mesh->recsize = (int32_t) p4est->data_size; mesh->destructor = p8est_trilinear_mesh_destroy; /* These members are incomplete and need to be filled later. */ memset (mesh->bounds, 0, 6 * sizeof (int)); memset (mesh->sizes, 0, 3 * sizeof (int)); mesh->minsize = mesh->maxsize = 0; mesh->ticksize = 0.; mesh->extra_info = NULL; mesh->gid = -1; /* We are done */ P4EST_GLOBAL_PRODUCTIONF ("Done trilinear_mesh_extract" " with %lld anodes %lld %lld\n", (long long) mesh->total_anode_num, (long long) global_borrowed, (long long) global_shared); return mesh; }