12 #include "../config/config.hpp"
18 #include "../general/binaryio.hpp"
26 using namespace bin_io;
50 int &curved,
int &is_nc)
51 :
NCMesh(input, version, curved, is_nc)
57 MPI_Comm_rank(
MyComm, &my_rank);
66 "Parallel mesh file doesn't seem to match current MPI setup. "
67 "Loading a parallel NC mesh with a non-matching communicator "
68 "size is not supported.");
71 MPI_Allreduce(&iso, &
Iso, 1, MPI_C_BOOL, MPI_LAND,
MyComm);
79 , MyComm(other.MyComm)
80 , NRanks(other.NRanks)
102 for (
int i = 0; i < 3; i++)
123 int f_index =
faces[face].index;
126 owner = std::min(owner, el.
rank);
137 el_loc = (el.
index << 4) | local;
190 int e_index =
nodes[enode].edge_index;
193 owner = std::min(owner, el.
rank);
204 el_loc = (el.
index << 4) | local;
251 int v_index =
nodes[vnode].vert_index;
254 owner = std::min(owner, el.
rank);
265 el_loc = (el.
index << 4) | local;
304 for (
int i = 0; i < num; i++)
316 for (
int i = 0; i < list.
masters.Size(); i++)
320 char master_old_flag = master_flag;
324 int si = list.
slaves[j].index;
328 master_flag |= slave_flag;
329 slave_flag |= master_old_flag;
343 for (
int i = 0; i < list.
conforming.Size(); i++)
350 for (
int i = 0; i < list.
masters.Size(); i++)
357 for (
int i = 0; i < list.
slaves.Size(); i++)
359 int si = list.
slaves[i].index;
369 if (lhs.size() == rhs.size())
371 for (
unsigned i = 0; i < lhs.size(); i++)
373 if (lhs[i] < rhs[i]) {
return true; }
377 return lhs.size() < rhs.size();
383 for (
unsigned i = 1; i < group.size(); i++)
385 if (group[i] <= group[i-1]) {
return false; }
393 if (group.size() == 1 && group[0] ==
MyRank)
397 MFEM_ASSERT(group_sorted(group),
"invalid group");
409 MFEM_ASSERT(rank != INT_MAX,
"invalid rank");
410 static std::vector<int> group;
420 for (
unsigned i = 0; i < group.size(); i++)
422 if (group[i] == rank) {
return true; }
433 entity_group.
SetSize(nentities);
439 int begin = 0, end = 0;
440 while (begin < index_rank.
Size())
442 int index = index_rank[begin].from;
443 if (index >= nentities)
447 while (end < index_rank.
Size() && index_rank[end].from ==
index)
451 group.resize(end - begin);
452 for (
int i = begin; i < end; i++)
454 group[i - begin] = index_rank[i].to;
463 for (
int i = 0; i < ranks.
Size(); i++)
476 int v[4], e[4], eo[4];
497 for (
int j = 0; j < 2; j++)
521 for (
int j = 0; j < nfv; j++)
536 for (
int i = 0; i < 3; i++)
549 for (
int i = 0; i < 2; i++)
558 if (local) { local[i] = lf; }
562 for (
int j = 0; j < 4; j++)
564 ids[i][j] = e[i]->
node[fv[j]];
574 if (
Dim < 3) {
return; }
583 for (
auto face =
faces.begin(); face !=
faces.end(); ++face)
585 if (face->elem[0] >= 0 && face->elem[1] >= 0 && face->index <
NFaces)
590 if (e1->
rank == e2->
rank) {
continue; }
591 if (e1->
rank > e2->
rank) { std::swap(e1, e2); }
606 for (i = j = 0; i < bdr_vertices.
Size(); i++)
608 if (bdr_vertices[i] <
NVertices) { bdr_vertices[j++] = bdr_vertices[i]; }
613 for (i = j = 0; i < bdr_edges.
Size(); i++)
615 if (bdr_edges[i] <
NEdges) { bdr_edges[j++] = bdr_edges[i]; }
630 for (
int i = 0; i < nleaves; i++)
645 for (
int i = 0; i < nleaves; i++)
653 else if (boundary_set[i] && etype)
670 for (
int i = 0; i < 8 && el.
child[i] >= 0; i++)
700 static void set_to_array(
const std::set<T> &set,
Array<T> &array)
704 for (std::set<int>::iterator it = set.begin(); it != set.end(); ++it)
721 set_to_array(ranks, neighbors);
733 group_shared.
MakeI(ngroups-1);
737 for (
int i = 0; i < conf_group.
Size(); i++)
741 if (entity_geom && (*entity_geom)[i] !=
geom) {
continue; }
748 shared_local.
SetSize(num_shared);
749 group_shared.
MakeJ();
752 for (
int i = 0, j = 0; i < conf_group.
Size(); i++)
756 if (entity_geom && (*entity_geom)[i] !=
geom) {
continue; }
766 for (
int i = 0; i < group_shared.
Size(); i++)
768 int size = group_shared.
RowSize(i);
769 int *row = group_shared.
GetRow(i);
772 ref_row.
Sort([&](
const int a,
const int b)
780 if (lsi_a != lsi_b) {
return lsi_a < lsi_b; }
782 return (el_loc_a & 0xf) < (el_loc_b & 0xf);
792 for (
int ent = 0; ent <
Dim; ent++)
806 for (
unsigned i = 0; i <
groups.size(); i++)
808 if (
groups[i].size() > 1 || !i)
811 group_map[i] = int_groups.
Insert(iset);
818 for (
int ent = 0; ent < 3; ent++)
823 ecg = group_map[ecg];
857 for (
int i = 0; i < slt.
Size(); i++)
862 int v[4], e[4], eo[4];
869 for (
int i = 0; i < slq.
Size(); i++)
879 for (
int ent = 0; ent <
Dim; ent++)
902 for (
int i = 0; i < shared.
conforming.Size(); i++)
906 MFEM_ASSERT(face != NULL,
"");
908 MFEM_ASSERT(face->
elem[0] >= 0 && face->
elem[1] >= 0,
"");
911 if (e[0]->rank ==
MyRank) { std::swap(e[0], e[1]); }
912 MFEM_ASSERT(e[0]->rank !=
MyRank && e[1]->rank ==
MyRank,
"");
918 for (
int i = 0; i < shared.
masters.Size(); i++)
924 if (sf.
element < 0) {
continue; }
926 MFEM_ASSERT(mf.
element >= 0,
"");
931 if (loc0 == loc1) {
continue; }
932 if (loc0) { std::swap(e[0], e[1]); }
939 MFEM_ASSERT(fnbr.
Size() <= bound,
"oops, bad upper bound");
953 for (
int i = 0; i < fnbr.
Size(); i++)
972 std::map<int, int> vert_map;
973 for (
int i = 0; i < fnbr.
Size(); i++)
981 for (
int k = 0; k < gi.
nv; k++)
983 int &v = vert_map[elem->
node[k]];
984 if (!v) { v = vert_map.size(); }
988 if (!i || elem->
rank != fnbr[i-1]->rank)
1004 for (
auto it = vert_map.begin(); it != vert_map.end(); ++it)
1017 for (
int i = 0, last_rank = -1; i < send_elems.
Size(); i++)
1020 if (c.
from != last_rank)
1028 c.
from = send_elems[i-1].from;
1034 for (
int i = 0; i < shared.
conforming.Size(); i++)
1040 if (e[0]->rank ==
MyRank) { std::swap(e[0], e[1]); }
1057 if (shared.
slaves.Size())
1063 MFEM_ASSERT(pmesh.
faces_info.Size() == nfaces,
"");
1066 for (
int i = nfaces; i < pmesh.
faces_info.Size(); i++)
1077 for (
int i = 0; i < shared.
masters.Size(); i++)
1083 if (sf.
element < 0) {
continue; }
1085 MFEM_ASSERT(mf.
element >= 0,
"");
1091 if (sloc == mloc) {
continue; }
1120 if (!sloc &&
Dim == 3)
1126 std::swap((*pm2)(0,1), (*pm2)(0,3));
1127 std::swap((*pm2)(1,1), (*pm2)(1,3));
1143 else if (!sloc &&
Dim == 2)
1152 MFEM_ASSERT(fi.
NCFace < 0,
"");
1181 bool removeAll =
true;
1184 for (
int i = 0; i < 8; i++)
1187 if (el.
child[i] >= 0)
1190 if (!
remove[i]) { removeAll =
false; }
1195 if (removeAll) {
return true; }
1198 for (
int i = 0; i < 8; i++)
1218 MFEM_WARNING(
"Can't prune 3D aniso meshes yet.");
1254 for (
int i = 0; i < refinements.
Size(); i++)
1258 "anisotropic parallel refinement not supported yet in 3D.");
1260 MFEM_VERIFY(
Iso ||
Dim < 3,
1261 "parallel refinement of 3D aniso meshes not supported yet.");
1268 for (
int i = 0; i < neighbors.
Size(); i++)
1270 send_ref[neighbors[i]].SetNCMesh(
this);
1278 for (
int i = 0; i < refinements.
Size(); i++)
1284 for (
int j = 0; j < ranks.
Size(); j++)
1286 send_ref[ranks[j]].AddRefinement(elem, ref.
ref_type);
1294 for (
int i = 0; i < refinements.
Size(); i++)
1301 for (
int j = 0; j < neighbors.
Size(); j++)
1311 for (
int i = 0; i < msg.
Size(); i++)
1326 MFEM_VERIFY(max_nc_level >= 1,
"'max_nc_level' must be 1 or greater.");
1333 long size = refinements.
Size(), glob_size;
1334 MPI_Allreduce(&size, &glob_size, 1, MPI_LONG, MPI_SUM,
MyComm);
1336 if (!glob_size) {
break; }
1344 MFEM_VERIFY(
Dim < 3 ||
Iso,
1345 "derefinement of 3D anisotropic meshes not implemented yet.");
1369 for (
int i = 0; i < derefs.
Size(); i++)
1371 int row = derefs[i];
1373 "invalid derefinement number.");
1378 int coarse_rank = INT_MAX;
1379 for (
int j = 0; j < size; j++)
1382 coarse_rank = std::min(coarse_rank, fine_rank);
1384 for (
int j = 0; j < size; j++)
1386 new_ranks[fine[j]] = coarse_rank;
1390 int target_elements = 0;
1391 for (
int i = 0; i < new_ranks.Size(); i++)
1393 if (new_ranks[i] ==
MyRank) { target_elements++; }
1407 for (
int i = 0; i < neighbors.
Size(); i++)
1409 send_deref[neighbors[i]].SetNCMesh(
this);
1416 for (
int i = 0; i < derefs.
Size(); i++)
1419 int parent =
elements[old_elements[fine[0]]].parent;
1423 for (
int j = 0; j < ranks.
Size(); j++)
1425 send_deref[ranks[j]].AddDerefinement(parent, new_ranks[fine[0]]);
1435 for (
int i = 0; i < old_elements.
Size(); i++)
1437 elements[old_elements[i]].index = i;
1442 old_elements.
Copy(coarse);
1443 for (
int i = 0; i < derefs.
Size(); i++)
1446 int parent =
elements[old_elements[fine[0]]].parent;
1455 for (
int j = 0; j < neighbors.
Size(); j++)
1465 for (
int i = 0; i < msg.
Size(); i++)
1483 for (
int i = 0; i < coarse.
Size(); i++)
1489 index = -1 -
elements[coarse[i]].rank;
1499 for (
int i = 0; i < coarse.
Size(); i++)
1513 template<
typename Type>
1515 const Table &deref_table)
1528 for (
int i = 0; i < deref_table.
Size(); i++)
1530 const int* fine = deref_table.
GetRow(i);
1531 int size = deref_table.
RowSize(i);
1532 MFEM_ASSERT(size <= 8,
"");
1534 int ranks[8], min_rank = INT_MAX, max_rank = INT_MIN;
1535 for (
int j = 0; j < size; j++)
1538 min_rank = std::min(min_rank, ranks[j]);
1539 max_rank = std::max(max_rank, ranks[j]);
1543 if (min_rank != max_rank)
1546 for (
int j = 0; j < size; j++)
1553 for (
int j = 0; j < size; j++)
1555 Type *data = &elem_data[fine[j]];
1557 int rnk = ranks[j], len = 1;
1563 for (
int k = 0; k < neigh.
Size(); k++)
1565 MPI_Request* req =
new MPI_Request;
1566 MPI_Isend(data, len, datatype, neigh[k], 292,
MyComm, req);
1572 MPI_Request* req =
new MPI_Request;
1573 MPI_Irecv(data, len, datatype, rnk, 292,
MyComm, req);
1580 for (
int i = 0; i < requests.
Size(); i++)
1582 MPI_Wait(requests[i], MPI_STATUS_IGNORE);
1589 ParNCMesh::SynchronizeDerefinementData<int>(
Array<int> &,
const Table &);
1601 for (
int i = 0; i < deref_table.
Size(); i++)
1603 const int *fine = deref_table.
GetRow(i),
1604 size = deref_table.
RowSize(i);
1609 for (
int j = 0; j < size; j++)
1617 for (
int k = 0; k <
Dim; k++)
1620 splits[k] >= max_nc_level)
1622 leaf_ok[fine[j]] = 0;
break;
1634 for (
int i = 0; i < deref_table.
Size(); i++)
1636 const int* fine = deref_table.
GetRow(i),
1637 size = deref_table.
RowSize(i);
1639 for (
int j = 0; j < size; j++)
1641 if (!leaf_ok[fine[j]])
1643 level_ok[i] = 0;
break;
1660 if (!custom_partition)
1666 long local_elems =
NElements, total_elems = 0;
1667 MPI_Allreduce(&local_elems, &total_elems, 1, MPI_LONG, MPI_SUM,
MyComm);
1669 long first_elem_global = 0;
1670 MPI_Scan(&local_elems, &first_elem_global, 1, MPI_LONG, MPI_SUM,
MyComm);
1671 first_elem_global -= local_elems;
1677 new_ranks[i] =
Partition(first_elem_global + (j++), total_elems);
1690 "Size of the partition array must match the number "
1691 "of local mesh elements (ParMesh::GetNE()).");
1694 custom_partition->
Copy(new_ranks);
1703 for (
int i = 0; i < old_elements.
Size(); i++)
1716 bool sfc = (target_elements >= 0);
1733 int begin = 0, end = 0;
1753 for (
int i = 0; i < rank_neighbors.
Size(); i++)
1755 int elem = rank_neighbors[i];
1763 recv_ghost_ranks[rank].SetNCMesh(
this);
1772 NeighborElementRankMessage::Map::iterator it;
1773 for (it = recv_ghost_ranks.begin(); it != recv_ghost_ranks.end(); ++it)
1776 for (
int i = 0; i < msg.
Size(); i++)
1780 new_ranks[ghost_index] = msg.
values[i];
1784 recv_ghost_ranks.clear();
1795 int received_elements = 0;
1801 received_elements++;
1803 el.
rank = new_ranks[i];
1806 int nsent = 0, nrecv = 0;
1813 owned_elements.
Sort([&](
const int a,
const int b)
1822 int begin = 0, end = 0;
1826 int rank =
elements[owned_elements[begin]].rank;
1827 while (end < owned_elements.
Size() &&
1828 elements[owned_elements[end]].rank == rank) { end++; }
1833 rank_elems.
MakeRef(&owned_elements[begin], end - begin);
1844 for (
int i = 0; i < batch.
Size(); i++)
1846 int elem = batch[i];
1892 while (received_elements < target_elements)
1901 for (
int i = 0; i < msg.
Size(); i++)
1903 int elem_rank = msg.
values[i];
1906 if (elem_rank ==
MyRank) { received_elements++; }
1928 MPI_Request barrier = MPI_REQUEST_NULL;
1940 for (
int i = 0; i < msg.
Size(); i++)
1952 if (barrier != MPI_REQUEST_NULL)
1954 MPI_Test(&barrier, &done, MPI_STATUS_IGNORE);
1960 int err = MPI_Ibarrier(
MyComm, &barrier);
1962 MFEM_VERIFY(err == MPI_SUCCESS,
"");
1963 MFEM_VERIFY(barrier != MPI_REQUEST_NULL,
"");
1974 int glob_sent, glob_recv;
1975 MPI_Reduce(&nsent, &glob_sent, 1, MPI_INT, MPI_SUM, 0,
MyComm);
1976 MPI_Reduce(&nrecv, &glob_recv, 1, MPI_INT, MPI_SUM, 0,
MyComm);
1980 MFEM_ASSERT(glob_sent == glob_recv,
1981 "(glob_sent, glob_recv) = ("
1982 << glob_sent <<
", " << glob_recv <<
")");
1989 const Table &old_element_dofs,
1990 long old_global_offset,
1997 RebalanceDofMessage::Map::iterator it;
2007 for (
int i = 0; i < ne; i++)
2028 RebalanceDofMessage::Map::iterator it;
2033 nd += msg.
dofs.size();
2044 for (
unsigned i = 0; i < msg.
elem_ids.size(); i++)
2048 for (
unsigned i = 0; i < msg.
dofs.size(); i++)
2061 : ncmesh(other.ncmesh), include_ref_types(other.include_ref_types)
2069 data.Append(value & 0xff);
2070 data.Append((value >> 8) & 0xff);
2071 data.Append((value >> 16) & 0xff);
2072 data.Append((value >> 24) & 0xff);
2078 return (
int) data[pos] +
2079 ((int) data[pos+1] << 8) +
2080 ((int) data[pos+2] << 16) +
2081 ((int) data[pos+3] << 24);
2086 for (
int i = 0; i < elements.
Size(); i++)
2088 int elem = elements[i];
2091 Element &el = ncmesh->elements[elem];
2092 if (el.
flag == flag) {
break; }
2101 Element &el = ncmesh->elements[elem];
2111 for (
int i = 0; i < 8; i++)
2113 if (el.
child[i] >= 0 && ncmesh->elements[el.
child[i]].flag)
2121 if (include_ref_types)
2126 for (
int i = 0; i < 8; i++)
2128 if (mask & (1 << i))
2130 EncodeTree(el.
child[i]);
2138 FlagElements(elements, 1);
2143 for (
int i = 0; i < ncmesh->root_state.Size(); i++)
2145 if (ncmesh->elements[i].flag)
2153 FlagElements(elements, 0);
2159 std::ostringstream oss;
2160 for (
int i = 0; i < ref_path.Size(); i++)
2162 oss <<
" elem " << ref_path[i] <<
" (";
2163 const Element &el = ncmesh->elements[ref_path[i]];
2164 for (
int j = 0; j <
GI[el.
Geom()].
nv; j++)
2166 if (j) { oss <<
", "; }
2167 oss << ncmesh->RetrieveNode(el, j);
2179 ref_path.Append(elem);
2181 int mask = data[pos++];
2188 Element &el = ncmesh->elements[elem];
2189 if (include_ref_types)
2191 int ref_type = data[pos++];
2194 ncmesh->RefineElement(elem, ref_type);
2196 else { MFEM_ASSERT(ref_type == el.
ref_type,
"") }
2200 MFEM_ASSERT(el.
ref_type != 0,
"Path not found:\n"
2201 << RefPath() <<
" mask = " << mask);
2204 for (
int i = 0; i < 8; i++)
2206 if (mask & (1 << i))
2208 DecodeTree(el.
child[i], pos, elements);
2213 ref_path.DeleteLast();
2220 while ((root = GetInt(pos)) >= 0)
2223 DecodeTree(root, pos, elements);
2229 write<int>(os, data.Size());
2230 os.write((
const char*) data.GetData(), data.Size());
2235 data.SetSize(read<int>(is));
2236 is.read((
char*) data.GetData(), data.Size());
2252 for (
unsigned i = 0; i <
groups.size(); i++)
2258 for (
int i = 0; i < ids[0].
Size(); i++)
2260 find_v[i].one = ids[0][i].index;
2274 for (
int j = 0; j < 2; j++)
2279 k = find_v[pos].two;
2290 for (
int i = 0; i < ids[1].
Size(); i++)
2292 find_e[i].one = ids[1][i].index;
2304 int v[4], e[4], eo[4], pos, k;
2306 for (
int j = 0; j < nfv; j++)
2310 k = find_v[pos].two;
2316 k = find_e[pos].two;
2331 for (
int i = 0; i < gi.
nv; i++)
2333 if (
nodes[el.
node[i]].vert_index ==
id.index)
2340 MFEM_ABORT(
"Vertex not found.");
2346 const int *ev =
GI[old.
Geom()].
edges[(int)
id.local];
2348 MFEM_ASSERT(node != NULL,
"Edge not found.");
2354 for (
int i = 0; i < gi.
ne; i++)
2356 const int* ev = gi.
edges[i];
2357 if ((el.
node[ev[0]] == node->
p1 && el.
node[ev[1]] == node->
p2) ||
2358 (el.
node[ev[1]] == node->
p1 && el.
node[ev[0]] == node->
p2))
2366 MFEM_ABORT(
"Edge not found.");
2372 const MeshId &first = ids[find[pos].two];
2373 while (++pos < find.Size() && ids[find[pos].two].index == first.
index)
2375 MeshId &other = ids[find[pos].two];
2383 std::map<int, int> stream_id;
2389 for (
int type = 0; type < 3; type++)
2391 for (
int i = 0; i < ids[type].
Size(); i++)
2393 elements.
Append(ids[type][i].element);
2405 for (
int i = 0; i < decoded.
Size(); i++)
2407 stream_id[decoded[i]] = i;
2412 for (
int type = 0; type < 3; type++)
2414 write<int>(os, ids[type].
Size());
2415 for (
int i = 0; i < ids[type].
Size(); i++)
2417 const MeshId&
id = ids[type][i];
2418 write<int>(os, stream_id[
id.element]);
2419 write<char>(os,
id.local);
2434 for (
int type = 0; type < 3; type++)
2436 int ne = read<int>(is);
2439 for (
int i = 0; i < ne; i++)
2441 int el_num = read<int>(is);
2442 int elem = elems[el_num];
2445 MFEM_VERIFY(!el.
ref_type,
"not a leaf element: " << el_num);
2447 MeshId &
id = ids[type][i];
2449 id.local = read<char>(is);
2457 id.index =
nodes[el.
node[(int)
id.local]].vert_index;
2462 const int* ev = gi.edges[(int)
id.local];
2464 MFEM_ASSERT(node && node->
HasEdge(),
"edge not found.");
2470 const int* fv = gi.faces[(int)
id.local];
2473 MFEM_ASSERT(face,
"face not found.");
2474 id.index = face->
index;
2484 std::map<GroupId, GroupId> stream_id;
2485 for (
int i = 0; i < ids.
Size(); i++)
2487 if (i && ids[i] == ids[i-1]) {
continue; }
2488 unsigned size = stream_id.size();
2489 GroupId &sid = stream_id[ids[i]];
2490 if (size != stream_id.size()) { sid = size; }
2494 write<short>(os, stream_id.size());
2495 for (std::map<GroupId, GroupId>::iterator
2496 it = stream_id.begin(); it != stream_id.end(); ++it)
2498 write<GroupId>(os, it->second);
2502 write<short>(os, group.size());
2503 for (
unsigned i = 0; i < group.size(); i++)
2505 write<int>(os, group[i]);
2511 write<short>(os, -1);
2516 write<int>(os, ids.
Size());
2517 for (
int i = 0; i < ids.
Size(); i++)
2519 write<GroupId>(os, stream_id[ids[i]]);
2525 int ngroups = read<short>(is);
2531 for (
int i = 0; i < ngroups; i++)
2533 int id = read<GroupId>(is);
2534 int size = read<short>(is);
2538 for (
int i = 0; i < size; i++)
2540 ranks[i] = read<int>(is);
2552 for (
int i = 0; i < ids.
Size(); i++)
2554 ids[i] = groups[read<GroupId>(is)];
2561 template<
class ValueType,
bool RefTypes,
int Tag>
2564 std::ostringstream ostream;
2570 eset.
Encode(tmp_elements);
2578 std::map<int, int> element_index;
2579 for (
int i = 0; i < decoded.
Size(); i++)
2581 element_index[decoded[i]] = i;
2584 write<int>(ostream, values.size());
2585 MFEM_ASSERT(
elements.size() == values.size(),
"");
2587 for (
unsigned i = 0; i < values.size(); i++)
2589 write<int>(ostream, element_index[
elements[i]]);
2590 write<ValueType>(ostream, values[i]);
2593 ostream.str().swap(data);
2596 template<
class ValueType,
bool RefTypes,
int Tag>
2599 std::istringstream istream(data);
2605 eset.
Decode(tmp_elements);
2607 int* el = tmp_elements.
GetData();
2611 int count = read<int>(istream);
2612 for (
int i = 0; i < count; i++)
2614 int index = read<int>(istream);
2615 MFEM_ASSERT(index >= 0 && (
size_t) index < values.size(),
"");
2616 values[
index] = read<ValueType>(istream);
2626 eset.SetNCMesh(ncmesh);
2631 eset.Decode(decoded);
2633 elem_ids.resize(decoded.
Size());
2634 for (
int i = 0; i < decoded.
Size(); i++)
2636 elem_ids[i] = eset.GetNCMesh()->elements[decoded[i]].index;
2640 static void write_dofs(std::ostream &os,
const std::vector<int> &dofs)
2642 write<int>(os, dofs.size());
2644 os.write((
const char*) dofs.data(), dofs.size() *
sizeof(int));
2647 static void read_dofs(std::istream &is, std::vector<int> &dofs)
2649 dofs.resize(read<int>(is));
2650 is.read((
char*) dofs.data(), dofs.size() *
sizeof(int));
2655 std::ostringstream stream;
2658 write<long>(stream, dof_offset);
2659 write_dofs(stream, dofs);
2661 stream.str().swap(data);
2666 std::istringstream stream(data);
2669 dof_offset = read<long>(stream);
2670 read_dofs(stream, dofs);
2677 elem_ids.resize(elems.
Size());
2678 for (
int i = 0; i < elems.
Size(); i++)
2680 elem_ids[i] = eset.GetNCMesh()->elements[elems[i]].index;
2693 for (
int i = 0; i < cle.
Size(); i++)
2701 debug_mesh.
ncmesh = copy;
2712 for (
int i = 0; i < 3; i++)
2729 return (elem_ids.capacity() + dofs.capacity()) *
sizeof(
int);
2732 template<
typename K,
typename V>
2733 static long map_memory_usage(
const std::map<K, V> &map)
2736 for (
typename std::map<K, V>::const_iterator
2737 it = map.begin(); it != map.end(); ++it)
2739 result += it->second.MemoryUsage();
2740 result +=
sizeof(std::pair<K, V>) + 3*
sizeof(
void*) +
sizeof(bool);
2748 for (
unsigned i = 0; i <
groups.size(); i++)
2750 groups_size +=
groups[i].capacity() *
sizeof(int);
2752 const int approx_node_size =
2753 sizeof(std::pair<CommGroup, GroupId>) + 3*
sizeof(
void*) +
sizeof(bool);
2754 return groups_size +
group_id.size() * approx_node_size;
2757 template<
typename Type,
int Size>
2758 static long arrays_memory_usage(
const Array<Type> (&arrays)[Size])
2761 for (
int i = 0; i < Size; i++)
2763 total += arrays[i].MemoryUsage();
2770 long total_groups_owners = 0;
2771 for (
int i = 0; i < 3; i++)
2807 << arrays_memory_usage(
entity_owner) <<
" entity_owner\n"
2833 #endif // MFEM_USE_MPI
NCList face_list
lazy-initialized list of faces, see GetFaceList
NCList edge_list
lazy-initialized list of edges, see GetEdgeList
int Partition(long index, long total_elements) const
Return the processor number for a global element number.
ElementSet(NCMesh *ncmesh=NULL, bool include_ref_types=false)
void SendRebalanceDofs(int old_ndofs, const Table &old_element_dofs, long old_global_offset, FiniteElementSpace *space)
Use the communication pattern from last Rebalance() to send element DOFs.
void DecodeMeshIds(std::istream &is, Array< MeshId > ids[])
void Create(ListOfIntegerSets &groups, int mpitag)
Set up the group topology given the list of sets of shared entities.
Array< char > element_type
int Size() const
Return the logical size of the array.
std::vector< ValueType > values
long PartitionFirstIndex(int rank, long total_elements) const
Return the global index of the first element owned by processor 'rank'.
void GetConformingSharedStructures(class ParMesh &pmesh)
CoarseFineTransformations transforms
storage for data returned by Get[De]RefinementTransforms()
void Recreate(const int n, const int *p)
Create an integer set from C-array 'p' of 'n' integers. Overwrites any existing set data...
int elem[2]
up to 2 elements sharing the face
virtual void Trim()
Save memory by releasing all non-essential and cached data.
char ref_type
refinement XYZ bit mask (7 = full isotropic)
std::vector< int > CommGroup
void CalcFaceOrientations()
int GetFaceVerticesEdges(const MeshId &face_id, int vert_index[4], int edge_index[4], int edge_orientation[4]) const
void Unique()
Removes duplicities from a sorted array. This requires operator== to be defined for T...
T * end()
STL-like end. Returns pointer after the last element of the array.
void SetDerefMatrixCodes(int parent, Array< int > &fine_coarse)
Helper struct to convert a C++ type to an MPI type.
void MakeSharedList(const NCList &list, NCList &shared)
char flag
generic flag/marker, can be used by algorithms
Array< char > face_orient
This holds in one place the constants about the geometries we support.
void MakeI(int nrows)
Next 7 methods are used together with the default constructor.
virtual void GetVertices(Array< int > &v) const =0
Returns element's vertices.
std::vector< int > elem_ids
void ChangeRemainingMeshIds(Array< MeshId > &ids, int pos, const Array< Pair< int, int > > &find)
virtual void BuildFaceList()
void Dump(std::ostream &os) const
virtual void Trim()
Save memory by releasing all non-essential and cached data.
int index
element number in the Mesh, -1 if refined
Array< Element * > face_nbr_elements
const NCList & GetSharedVertices()
GroupId GetSingletonGroup(int rank)
void GetFaceNeighbors(class ParMesh &pmesh)
Lists all edges/faces in the nonconforming mesh.
std::string RefPath() const
void SetNCMesh(ParNCMesh *pncmesh)
Set pointer to ParNCMesh (needed to encode the message).
std::map< int, NeighborElementRankMessage > Map
void Copy(Array ©) const
Create a copy of the internal array to the provided copy.
virtual void ElementSharesVertex(int elem, int local, int vnode)
unsigned matrix
index into NCList::point_matrices[geom]
std::vector< int > elements
T * GetData()
Returns the data.
void Load(std::istream &is)
Table derefinements
possible derefinements, see GetDerefinementTable
Data type dense matrix using column-major storage.
NCList vertex_list
lazy-initialized list of vertices, see GetVertexList
void InitDerefTransforms()
void GetRow(int i, Array< int > &row) const
Return row i in array row (the Table must be finalized)
char geom
Geometry::Type of the element (char for storage only)
void SynchronizeDerefinementData(Array< Type > &elem_data, const Table &deref_table)
int GetInt(int pos) const
virtual void Derefine(const Array< int > &derefs)
void SetElements(const Array< int > &elems, NCMesh *ncmesh)
bool operator<(const Pair< A, B > &p, const Pair< A, B > &q)
Comparison operator for class Pair, based on the first element only.
Array< Vert3 > shared_trias
void CountSplits(int elem, int splits[3]) const
Geometry::Type Geom() const
void ChangeEdgeMeshIdElement(NCMesh::MeshId &id, int elem)
Array< int > face_nbr_group
static int GetQuadOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
void GetLimitRefinements(Array< Refinement > &refinements, int max_level)
bool GroupContains(GroupId id, int rank) const
Return true if group 'id' contains the given rank.
int PrintMemoryDetail() const
void EncodeTree(int elem)
int spaceDim
dimensions of the elements and the vertex coordinates
const NCList & GetSharedEdges()
void DeleteAll()
Delete the whole array.
static int find_element_edge(const Element &el, int vn0, int vn1, bool abort=true)
int index
Mesh element number.
int master
master number (in Mesh numbering)
void Recv(int rank, int size, MPI_Comm comm)
Post-probe receive from processor 'rank' of message size 'size'.
const NCList & GetFaceList()
Return the current list of conforming and nonconforming faces.
const NCList & GetSharedList(int entity)
Helper to get shared vertices/edges/faces ('entity' == 0/1/2 resp.).
Array< NCFaceInfo > nc_faces_info
A parallel extension of the NCMesh class.
static int find_local_face(int geom, int a, int b, int c)
RebalanceDofMessage::Map send_rebalance_dofs
int GetNumFaces() const
Return the number of faces (3D), edges (2D) or vertices (1D).
void FlagElements(const Array< int > &elements, char flag)
void Decode(Array< int > &elements) const
Face * GetFace(Element &elem, int face_no)
void AddElementRank(int elem, int rank)
Array< int > face_nbr_elements_offset
int index
face number in the Mesh
void MakeFromList(int nrows, const Array< Connection > &list)
Array< Element * > shared_edges
Array< DenseMatrix * > aux_pm_store
Stores modified point matrices created by GetFaceNeighbors.
virtual void SetAttributes()
int Append(const T &el)
Append element 'el' to array, resize if necessary.
void RedistributeElements(Array< int > &new_ranks, int target_elements, bool record_comm)
void GetSubArray(int offset, int sa_size, Array< T > &sa) const
Copy sub array starting from offset out to the provided sa.
void FindNeighbors(int elem, Array< int > &neighbors, const Array< int > *search_set=NULL)
void ElementNeighborProcessors(int elem, Array< int > &ranks)
void AddConnection(int r, int c)
Array< int > old_index_or_rank
void CreateGroups(int nentities, Array< Connection > &index_rank, Array< GroupId > &entity_group)
RebalanceDofMessage::Map recv_rebalance_dofs
int InitialPartition(int index) const
Helper to get the partitioning when the serial mesh gets split initially.
Array< DenseMatrix * > point_matrices[Geometry::NumGeom]
List of unique point matrices for each slave geometry.
int Insert(IntegerSet &s)
Check to see if set 's' is in the list. If not append it to the end of the list. Returns the index of...
void Reserve(int capacity)
Ensures that the allocated size is at least the given size.
virtual void BuildEdgeList()
Identifies a vertex/edge/face in both Mesh and NCMesh.
T * begin()
STL-like begin. Returns pointer to the first element of the array.
int parent
parent element, -1 if this is a root element, -2 if free'd
std::map< int, NeighborRefinementMessage > Map
long MemoryUsage() const
Returns the number of bytes allocated for the array including any reserve.
static void IsendAll(MapT &rank_msg, MPI_Comm comm)
Helper to send all messages in a rank-to-message map container.
signed char local
local number within 'element'
void Sort()
Sorts the array in ascending order. This requires operator< to be defined for T.
int MyRank
used in parallel, or when loading a parallel file in serial
int FindSorted(const T &el) const
Do bisection search for 'el' in a sorted array; return -1 if not found.
int GetVDim() const
Returns vector dimension.
int Size() const
Returns the number of TYPE I elements.
virtual void ElementSharesFace(int elem, int local, int face)
A class for non-conforming AMR. The class is not used directly by the user, rather it is an extension...
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
int element
NCMesh::Element containing this vertex/edge/face.
virtual void GetBoundaryClosure(const Array< int > &bdr_attr_is_ess, Array< int > &bdr_vertices, Array< int > &bdr_edges)
Array< Vert4 > shared_quads
void Rebalance(const Array< int > *custom_partition=NULL)
static GeomInfo GI[Geometry::NumGeom]
static bool TestAllSent(MapT &rank_msg)
Return true if all messages in the map container were sent, otherwise return false, without waiting.
int slaves_end
slave faces
void DecodeGroups(std::istream &is, Array< GroupId > &ids)
void AdjustMeshIds(Array< MeshId > ids[], int rank)
virtual void ElementSharesEdge(int elem, int local, int enode)
bool Iso
true if the mesh only contains isotropic refinements
Array< int > tmp_neighbors
std::map< int, NeighborDerefinementMessage > Map
Array< Connection > entity_index_rank[3]
virtual void BuildEdgeList()
Array< GroupId > entity_owner[3]
void Encode(const Array< int > &elements)
Array< char > tmp_shared_flag
Array< Vertex > face_nbr_vertices
virtual void BuildFaceList()
Nonconforming edge/face within a bigger edge/face.
void DerefineElement(int elem)
Class FiniteElementSpace - responsible for providing FEM view of the mesh, mainly managing the set of...
OutStream err(std::cerr)
Global stream used by the library for standard error output. Initially it uses the same std::streambu...
void GetDebugMesh(Mesh &debug_mesh) const
bool PruneTree(int elem)
Internal. Recursive part of Prune().
int Find(const T &el) const
Return the first index where 'el' is found; return -1 if not found.
void AddAColumnInRow(int r)
void SetSize(int nsize)
Change the logical size of the array, keep existing entries.
mfem::Element * NewMeshElement(int geom) const
Array< int > leaf_sfc_index
natural tree ordering of leaf elements
static void RecvAll(MapT &rank_msg, MPI_Comm comm)
Helper to receive all messages in a rank-to-message map container.
Helper struct for defining a connectivity table, see Table::MakeFromList.
void Issend(int rank, MPI_Comm comm)
Non-blocking synchronous send to processor 'rank'. Returns immediately. Completion (MPI_Wait/Test) me...
void RefineElement(int elem, char ref_type)
virtual void Refine(const Array< Refinement > &refinements)
Array< int > boundary_faces
subset of all faces, set by BuildFaceList
int child[8]
2-8 children (if ref_type != 0)
Array< double > coordinates
void RecvRebalanceDofs(Array< int > &elements, Array< long > &dofs)
Receive element DOFs sent by SendRebalanceDofs().
Array< int > leaf_elements
finest elements, in Mesh ordering (+ ghosts)
Table group_svert
Shared objects in each group.
Array< int > boundary_layer
list of type 3 elements
Array< MeshId > conforming
signed char geom
Geometry::Type (faces only) (char to save RAM)
void Isend(int rank, MPI_Comm comm)
Non-blocking send to processor 'rank'. Returns immediately. Completion (as tested by MPI_Wait/Test) d...
int index(int i, int j, int nx, int ny)
void AddElementRank(int elem, int rank)
void InitOwners(int num, Array< GroupId > &entity_owner)
void ChangeVertexMeshIdElement(NCMesh::MeshId &id, int elem)
void DecodeTree(int elem, int &pos, Array< int > &elements) const
Array< FaceInfo > faces_info
void SetAttribute(const int attr)
Set element's attribute.
int NGroups() const
Return the number of groups.
void EncodeGroups(std::ostream &os, const Array< GroupId > &ids)
virtual void BuildVertexList()
void CalculatePMatrixGroups()
NCMesh * ncmesh
Optional non-conforming mesh extension.
Array< int > entity_elem_local[3]
T & Last()
Return the last element in the array.
void NeighborExpand(const Array< int > &elems, Array< int > &expanded, const Array< int > *search_set=NULL)
BlockArray< Element > elements
virtual void CheckDerefinementNCLevel(const Table &deref_table, Array< int > &level_ok, int max_nc_level)
void FindSetNeighbors(const Array< char > &elem_set, Array< int > *neighbors, Array< char > *neighbor_set=NULL)
Array< int > ghost_layer
list of elements whose 'element_type' == 2.
void GetEdgeVertices(const MeshId &edge_id, int vert_index[2], bool oriented=true) const
Return Mesh vertex indices of an edge identified by 'edge_id'.
void MakeRef(T *, int)
Make this Array a reference to a pointer.
static void WaitAllSent(MapT &rank_msg)
Helper to wait for all messages in a map container to be sent.
void EncodeMeshIds(std::ostream &os, Array< MeshId > ids[])
void AddConnections(int entity, int index, const Array< int > &ranks)
long MemoryUsage() const
Return total number of bytes allocated.
void MakeSharedTable(int ngroups, int ent, Array< int > &shared_local, Table &group_shared, Array< char > *entity_geom=NULL, char geom=0)
static int GetTriOrientation(const int *base, const int *test)
Returns the orientation of "test" relative to "base".
Array< int > svert_lvert
Shared to local index mapping.
const NCList & GetSharedFaces()
std::map< int, RebalanceMessage > Map
static int get_face_orientation(Face &face, Element &e1, Element &e2, int local[2]=NULL)
virtual void BuildVertexList()
const NCList & GetEdgeList()
Return the current list of conforming and nonconforming edges.
Array< GroupId > entity_pmat_group[3]
void NeighborProcessors(Array< int > &neighbors)
long GroupsMemoryUsage() const
static bool IProbe(int &rank, int &size, MPI_Comm comm)
Non-blocking probe for incoming message of this type from any rank. If there is an incoming message...
Table send_face_nbr_elements
char ref_type
bit mask of X,Y,Z refinements (bits 0,1,2 respectively)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void InitFromNCMesh(const NCMesh &ncmesh)
Initialize vertices/elements/boundary/tables from a nonconforming mesh.
static int find_node(const Element &el, int node)
bool CheckElementType(int elem, int type)
Class for parallel meshes.
Abstract data type element.
GroupId GetGroupId(const CommGroup &group)
Data type line segment element.
const double * CalcVertexPos(int node) const
int rank
processor number (ParNCMesh), -1 if undefined/unknown
virtual void LimitNCLevel(int max_nc_level)
Parallel version of NCMesh::LimitNCLevel.
int node[8]
element corners (if ref_type == 0)
Array< unsigned char > data
encoded refinement (sub-)trees
static void Probe(int &rank, int &size, MPI_Comm comm)
Blocking probe for incoming message of this type from any rank. Returns the rank and message size...
virtual void Refine(const Array< Refinement > &refinements)
ParNCMesh(MPI_Comm comm, const NCMesh &ncmesh, int *part=NULL)
Construct by partitioning a serial NCMesh.
Array< GroupId > entity_conf_group[3]
Array< char > face_geom
face geometry by face index, set by OnMeshUpdated
void DofsToVDofs(Array< int > &dofs, int ndofs=-1) const