static void test_the_p4est (p4est_connectivity_t * conn, int N) { p4est_t *p4est; p4est_ghost_t *ghost; p4est_lnodes_t *lnodes; p4est = p4est_new (sc_MPI_COMM_WORLD, conn, 0, NULL, NULL); ghost = p4est_ghost_new (p4est, P4EST_CONNECT_FULL); lnodes = p4est_lnodes_new (p4est, ghost, N); p4est_lnodes_destroy (lnodes); p4est_ghost_destroy (ghost); p4est_destroy (p4est); }
p6est_profile_t * p6est_profile_new_local (p6est_t * p6est, p6est_ghost_t * ghost, p6est_profile_type_t ptype, p8est_connect_type_t btype, int degree) { p6est_profile_t *profile = P4EST_ALLOC (p6est_profile_t, 1); p4est_lnodes_t *lnodes; p4est_locidx_t nln, nle; p4est_topidx_t jt; p4est_t *columns = p6est->columns; p4est_tree_t *tree; sc_array_t *tquadrants; p4est_quadrant_t *col; p4est_qcoord_t diff = P4EST_ROOT_LEN - p6est->root_len; size_t first, last, count, zz, zy; p4est_locidx_t *en, (*lr)[2]; sc_array_t *lc; int i, j; p2est_quadrant_t *layer; sc_array_t *layers = p6est->layers; p4est_locidx_t nidx, enidx; p4est_connect_type_t hbtype; int8_t *c; sc_array_t *thisprof; sc_array_t *selfprof; sc_array_t *faceprof; sc_array_t *cornerprof; sc_array_t *work; sc_array_t oldprof; const int Nrp = degree + 1; P4EST_ASSERT (degree > 1); profile->ptype = ptype; profile->btype = btype; profile->lnode_changed[0] = NULL; profile->lnode_changed[1] = NULL; profile->enode_counts = NULL; profile->diff = diff; if (btype == P8EST_CONNECT_FACE) { hbtype = P4EST_CONNECT_FACE; } else { hbtype = P4EST_CONNECT_FULL; } if (ghost == NULL) { profile->cghost = p4est_ghost_new (p6est->columns, P4EST_CONNECT_FULL); profile->ghost_owned = 1; } else { P4EST_ASSERT (ghost->column_ghost->btype == P4EST_CONNECT_FULL); profile->cghost = ghost->column_ghost; profile->ghost_owned = 0; } if (ptype == P6EST_PROFILE_UNION) { P4EST_ASSERT (degree == 2); } profile->lnodes = lnodes = p4est_lnodes_new (p6est->columns, profile->cghost, degree); en = lnodes->element_nodes; nln = lnodes->num_local_nodes; nle = lnodes->num_local_elements; profile->lnode_ranges = P4EST_ALLOC_ZERO (p4est_locidx_t, 2 * nln); lr = (p4est_locidx_t (*)[2]) profile->lnode_ranges; profile->lnode_columns = lc = sc_array_new (sizeof (int8_t)); selfprof = sc_array_new (sizeof (int8_t)); work = sc_array_new (sizeof (int8_t)); faceprof = sc_array_new (sizeof (int8_t)); cornerprof = sc_array_new (sizeof (int8_t)); if (ptype == P6EST_PROFILE_UNION) { profile->lnode_changed[0] = P4EST_ALLOC (p4est_locidx_t, nln); profile->lnode_changed[1] = P4EST_ALLOC (p4est_locidx_t, nln); profile->enode_counts = P4EST_ALLOC (p4est_locidx_t, P4EST_INSUL * nle); profile->evenodd = 0; memset (profile->lnode_changed[0], -1, nln * sizeof (int)); } /* create the profiles for each node: layers are reduced to just their level * */ for (enidx = 0, jt = columns->first_local_tree; jt <= columns->last_local_tree; ++jt) { tree = p4est_tree_array_index (columns->trees, jt); tquadrants = &tree->quadrants; for (zz = 0; zz < tquadrants->elem_count; ++zz) { col = p4est_quadrant_array_index (tquadrants, zz); P6EST_COLUMN_GET_RANGE (col, &first, &last); count = last - first; sc_array_truncate (selfprof); c = (int8_t *) sc_array_push_count (selfprof, count); for (zy = first; zy < last; zy++) { layer = p2est_quadrant_array_index (layers, zy); *(c++) = layer->level; } if (ptype == P6EST_PROFILE_UNION) { p6est_profile_balance_self (selfprof, work); if (btype == P8EST_CONNECT_FACE) { p6est_profile_balance_face (selfprof, faceprof, work, diff); } else { p6est_profile_balance_full (selfprof, faceprof, work, diff); } if (btype == P8EST_CONNECT_EDGE) { p6est_profile_balance_face (selfprof, cornerprof, work, diff); } else if (btype == P8EST_CONNECT_FULL) { p6est_profile_balance_full (selfprof, cornerprof, work, diff); } } for (j = 0; j < Nrp; j++) { for (i = 0; i < Nrp; i++, enidx++) { nidx = en[enidx]; if (ptype == P6EST_PROFILE_UNION) { thisprof = NULL; if (!(i % degree) && !(j % degree)) { if (hbtype == P4EST_CONNECT_FACE) { /* skip corners if we don't need to balance them */ P4EST_ASSERT (!lr[nidx][0]); P4EST_ASSERT (!lr[nidx][1]); continue; } else { thisprof = cornerprof; } } else if ((i % degree) && (j % degree)) { thisprof = selfprof; } else { thisprof = faceprof; } count = thisprof->elem_count; profile->enode_counts[enidx] = count; if (!lr[nidx][1]) { /* if this node has not yet been initialized, initialize it */ lr[nidx][0] = lc->elem_count; lr[nidx][1] = count; c = (int8_t *) sc_array_push_count (lc, count); memcpy (c, thisprof->array, count * sizeof (int8_t)); } else { /* if this node has been initialized, combine the two profiles, * taking the finer layers from each */ sc_array_init_view (&oldprof, lc, lr[nidx][0], lr[nidx][1]); p6est_profile_union (thisprof, &oldprof, work); if (work->elem_count > oldprof.elem_count) { lr[nidx][0] = lc->elem_count; lr[nidx][1] = work->elem_count; c = (int8_t *) sc_array_push_count (lc, work->elem_count); memcpy (c, work->array, work->elem_count * work->elem_size); } } } else { count = selfprof->elem_count; if (!lr[nidx][1]) { /* if this node has not yet been initialized, initialize it */ lr[nidx][0] = lc->elem_count; lr[nidx][1] = count; c = (int8_t *) sc_array_push_count (lc, count); memcpy (c, selfprof->array, count * sizeof (int8_t)); } else { /* if this node has been initialized, combine the two profiles, * taking the coarser layers from each */ sc_array_init_view (&oldprof, lc, lr[nidx][0], lr[nidx][1]); p6est_profile_intersection (selfprof, &oldprof, work); P4EST_ASSERT (work->elem_count <= oldprof.elem_count); if (work->elem_count < oldprof.elem_count) { lr[nidx][1] = work->elem_count; memcpy (oldprof.array, work->array, work->elem_count * work->elem_size); } } } } } } } p6est_profile_compress (profile); sc_array_destroy (selfprof); sc_array_destroy (faceprof); sc_array_destroy (cornerprof); sc_array_destroy (work); return profile; }
int main (int argc, char *argv[]) { MPI_Comm comm = MPI_COMM_WORLD; p4est_t *p4est; p4est_connectivity_t *conn; p4est_ghost_t *ghost_layer; p4est_lnodes_t *lnodes; int rank; const int degree = 1; BFAM_MPI_CHECK(MPI_Init(&argc,&argv)); BFAM_MPI_CHECK(MPI_Comm_rank(comm, &rank)); bfam_log_init(rank, stdout, BFAM_LL_DEFAULT); bfam_signal_handler_set(); sc_init(comm, 0, 0, NULL, SC_LP_DEFAULT); p4est_init(NULL, SC_LP_DEFAULT); conn = p4est_connectivity_new_corner(); p4est = p4est_new_ext(comm, conn, 0, 0, 0, 0, NULL, NULL); refine_level = 1; p4est_refine(p4est, 1, refine_fn, NULL); p4est_balance(p4est, P4EST_CONNECT_FACE, NULL); p4est_partition(p4est, 1, NULL); p4est_vtk_write_file(p4est, NULL, "mesh"); ghost_layer = p4est_ghost_new(p4est, P4EST_CONNECT_FULL); lnodes = p4est_lnodes_new(p4est, ghost_layer, degree); /* * Output the mesh. It can be read using something like following command: * * mpirun -np 3 ./bfam_exam_p4est | grep MESH | sort -n -k 2 | sort -n -k 5 | gvim - */ fflush(stdout); BFAM_MPI_CHECK(MPI_Barrier(comm)); BFAM_ROOT_INFO("MESH 0 ------------ Mesh Begin ------------"); BFAM_ROOT_INFO("MESH 1 degree = %d", lnodes->degree); BFAM_ROOT_INFO("MESH 2 vnodes = %d", lnodes->vnodes); BFAM_INFO("MESH 3 num_local_elements = %jd", (intmax_t)lnodes->num_local_elements); BFAM_INFO("MESH 4 num_local_nodes = %jd", (intmax_t)lnodes->num_local_nodes); BFAM_INFO("MESH 5 owned_count = %jd", (intmax_t)lnodes->owned_count); BFAM_INFO("MESH 6 global_offset = %jd", (intmax_t)lnodes->global_offset); sc_array_t *global_nodes = sc_array_new(sizeof (p4est_gloidx_t)); sc_array_resize(global_nodes, lnodes->num_local_nodes); for(size_t zz = 0; zz < global_nodes->elem_count; ++zz) { *((p4est_gloidx_t *) sc_array_index(global_nodes, zz)) = p4est_lnodes_global_index(lnodes, zz); } p4est_lnodes_share_owned(global_nodes, lnodes); for(size_t zz = 0; zz < global_nodes->elem_count; ++zz) { const p4est_gloidx_t gn = *((p4est_gloidx_t *)sc_array_index(global_nodes, zz)); SC_CHECK_ABORT (gn == p4est_lnodes_global_index(lnodes, zz), "Lnodes: bad global index across procesors"); BFAM_INFO("MESH 7 global_nodes[%zu] = %jd", zz, (intmax_t)gn); } sc_array_destroy(global_nodes); p4est_topidx_t flt = p4est->first_local_tree; p4est_topidx_t llt = p4est->last_local_tree; p4est_locidx_t elid, elnid; p4est_topidx_t t; const double *v = conn->vertices; const p4est_topidx_t *tree_to_vertex = conn->tree_to_vertex; for(elid = 0, elnid = 0, t = flt; t <= llt; ++t) { p4est_tree_t *tree = p4est_tree_array_index(p4est->trees, t); const size_t count = tree->quadrants.elem_count; p4est_topidx_t vt[P4EST_CHILDREN]; for (int c = 0; c < P4EST_CHILDREN; ++c) { vt[c] = tree_to_vertex[t * P4EST_CHILDREN + c]; } for (size_t zz = 0; zz < count; ++zz, ++elid) { p4est_quadrant_t *q = p4est_quadrant_array_index(&tree->quadrants, zz); for(int jind = 0; jind < degree + 1; ++jind) { for(int iind = 0; iind < degree + 1; ++iind, ++elnid) { double xyz[3]; for (int j = 0; j < 3; ++j) { const p4est_qcoord_t len = P4EST_QUADRANT_LEN(q->level); const double rlen = (double) P4EST_ROOT_LEN; const double deg = (double) degree; const double qlen = ((double) len) / rlen; const double eta_x = ((double) q->x) / rlen + (((double) iind) / deg) * qlen; const double eta_y = ((double) q->y) / rlen + (((double) jind) / deg) * qlen; xyz[j] = ((1. - eta_y) * ((1. - eta_x) * v[3 * vt[0] + j] + eta_x * v[3 * vt[1] + j]) + eta_y * ((1. - eta_x) * v[3 * vt[2] + j] + eta_x * v[3 * vt[3] + j])); } const p4est_locidx_t nid = lnodes->element_nodes[elnid]; BFAM_INFO( "MESH 8 local_node[%03jd] = %03jd ( %25.16e %25.16e %25.16e )", (intmax_t)elnid, (intmax_t)nid, xyz[0], xyz[1], xyz[2]); } } } } BFAM_ROOT_INFO("MESH 9 ------------ Mesh End ------------"); p4est_lnodes_destroy(lnodes); p4est_ghost_destroy(ghost_layer); p4est_destroy(p4est); p4est_connectivity_destroy(conn); sc_finalize(); BFAM_MPI_CHECK(MPI_Finalize()); return EXIT_SUCCESS; }