Mesh Oriented datABase  (version 5.6.0)
An array-based unstructured mesh library
moab::ParCommGraph Class Reference

#include <ParCommGraph.hpp>

+ Collaboration diagram for moab::ParCommGraph:

Public Types

enum  TypeGraph { INITIAL_MIGRATE , COVERAGE , DOF_BASED }
 

Public Member Functions

virtual ~ParCommGraph ()
 
 ParCommGraph (MPI_Comm joincomm, MPI_Group group1, MPI_Group group2, int coid1, int coid2)
 collective constructor, will be called on all sender tasks and receiver tasks More...
 
 ParCommGraph (const ParCommGraph &)
 copy constructor will copy only the senders, receivers, compid1, etc More...
 
ErrorCode compute_trivial_partition (std::vector< int > &numElemsPerTaskInGroup1)
 Based on the number of elements on each task in group 1, partition for group 2, trivially. More...
 
ErrorCode pack_receivers_graph (std::vector< int > &packed_recv_array)
 pack information about receivers view of the graph, for future sending to receiver root More...
 
bool is_root_sender ()
 
bool is_root_receiver ()
 
int sender (int index)
 
int receiver (int index)
 
int get_component_id1 ()
 
int get_component_id2 ()
 
int get_context_id ()
 
void set_context_id (int other_id)
 
EntityHandle get_cover_set ()
 
void set_cover_set (EntityHandle cover)
 
ErrorCode split_owned_range (int sender_rank, Range &owned)
 
ErrorCode split_owned_range (Range &owned)
 
ErrorCode send_graph (MPI_Comm jcomm)
 
ErrorCode send_graph_partition (ParallelComm *pco, MPI_Comm jcomm)
 
ErrorCode send_mesh_parts (MPI_Comm jcomm, ParallelComm *pco, Range &owned)
 
ErrorCode receive_comm_graph (MPI_Comm jcomm, ParallelComm *pco, std::vector< int > &pack_array)
 
ErrorCode receive_mesh (MPI_Comm jcomm, ParallelComm *pco, EntityHandle local_set, std::vector< int > &senders_local)
 
ErrorCode release_send_buffers ()
 
ErrorCode send_tag_values (MPI_Comm jcomm, ParallelComm *pco, Range &owned, std::vector< Tag > &tag_handles)
 
ErrorCode receive_tag_values (MPI_Comm jcomm, ParallelComm *pco, Range &owned, std::vector< Tag > &tag_handles)
 
const std::vector< int > & senders ()
 
const std::vector< int > & receivers ()
 
ErrorCode settle_send_graph (TupleList &TLcovIDs)
 
void SetReceivingAfterCoverage (std::map< int, std::set< int > > &idsFromProcs)
 
void settle_comm_by_ids (int comp, TupleList &TLBackToComp, std::vector< int > &valuesComp)
 
ErrorCode compute_partition (ParallelComm *pco, Range &owned, int met)
 
ErrorCode dump_comm_information (std::string prefix, int is_send)
 

Private Member Functions

void find_group_ranks (MPI_Group group, MPI_Comm join, std::vector< int > &ranks)
 find ranks of a group with respect to an encompassing communicator More...
 

Private Attributes

MPI_Comm comm
 
std::vector< int > senderTasks
 
std::vector< int > receiverTasks
 
bool rootSender
 
bool rootReceiver
 
int rankInGroup1
 
int rankInGroup2
 
int rankInJoin
 
int joinSize
 
int compid1
 
int compid2
 
int context_id
 
EntityHandle cover_set
 
std::map< int, std::vector< int > > recv_graph
 
std::map< int, std::vector< int > > recv_sizes
 
std::map< int, std::vector< int > > sender_graph
 
std::map< int, std::vector< int > > sender_sizes
 
std::vector< ParallelComm::Buffer * > localSendBuffs
 
int * comm_graph
 
std::map< int, Rangesplit_ranges
 
std::vector< MPI_Request > sendReqs
 
std::vector< int > corr_tasks
 
std::vector< int > corr_sizes
 
TypeGraph graph_type
 
std::map< int, std::vector< int > > involved_IDs_map
 
std::map< int, std::vector< int > > map_index
 
std::map< int, std::vector< int > > map_ptr
 

Detailed Description

Definition at line 55 of file ParCommGraph.hpp.

Member Enumeration Documentation

◆ TypeGraph

Enumerator
INITIAL_MIGRATE 
COVERAGE 
DOF_BASED 

Definition at line 58 of file ParCommGraph.hpp.

59  {
61  COVERAGE,
62  DOF_BASED
63  };

Constructor & Destructor Documentation

◆ ~ParCommGraph()

moab::ParCommGraph::~ParCommGraph ( )
virtual

Definition at line 73 of file ParCommGraph.cpp.

74 {
75  // TODO Auto-generated destructor stub
76 }

◆ ParCommGraph() [1/2]

moab::ParCommGraph::ParCommGraph ( MPI_Comm  joincomm,
MPI_Group  group1,
MPI_Group  group2,
int  coid1,
int  coid2 
)

collective constructor, will be called on all sender tasks and receiver tasks

Parameters
[in]joincommjoint MPI communicator that covers both sender and receiver MPI groups
[in]group1MPI group formed with sender tasks; (sender usually loads the mesh in a migrate scenario)
[in]group2MPI group formed with receiver tasks; (receiver usually receives the mesh in a migrate scenario)
[in]coid1sender component unique identifier in a coupled application (in climate simulations could be the Atmosphere Comp id = 5 or Ocean Comp ID , 17)
[in]coid2receiver component unique identifier in a coupled application (it is usually the coupler, 2, in E3SM)

this graph will be formed on sender and receiver tasks, and, in principle, will hold info about how the local entities are distributed on the other side

Its major role is to help in migration of data, from component to the coupler and vice-versa; Each local entity has a corresponding task (or tasks) on the other side, to where the data needs to be sent

important data stored in ParCommGraph, immediately it is created

  • all sender and receiver tasks ids, with respect to the joint communicator
  • local rank in sender and receiver group (-1 if not part of the respective group)
  • rank in the joint communicator (from 0)

Definition at line 20 of file ParCommGraph.cpp.

21  : comm( joincomm ), compid1( coid1 ), compid2( coid2 )
22 {
23  // find out the tasks from each group, in the joint communicator
24  find_group_ranks( group1, comm, senderTasks );
26 
27  rootSender = rootReceiver = false;
28  rankInGroup1 = rankInGroup2 = rankInJoin = -1; // not initialized, or not part of the group
29 
30  int mpierr = MPI_Group_rank( group1, &rankInGroup1 );
31  if( MPI_SUCCESS != mpierr || rankInGroup1 == MPI_UNDEFINED ) rankInGroup1 = -1;
32 
33  mpierr = MPI_Group_rank( group2, &rankInGroup2 );
34  if( MPI_SUCCESS != mpierr || rankInGroup2 == MPI_UNDEFINED ) rankInGroup2 = -1;
35 
36  mpierr = MPI_Comm_rank( comm, &rankInJoin );
37  if( MPI_SUCCESS != mpierr ) // it should be a fatal error
38  rankInJoin = -1;
39 
40  mpierr = MPI_Comm_size( comm, &joinSize );
41  if( MPI_SUCCESS != mpierr ) // it should be a fatal error
42  joinSize = -1;
43 
44  if( 0 == rankInGroup1 ) rootSender = true;
45  if( 0 == rankInGroup2 ) rootReceiver = true;
47  comm_graph = NULL;
48  context_id = -1;
49  cover_set = 0; // refers to nothing yet
50 }

References comm, comm_graph, context_id, cover_set, find_group_ranks(), graph_type, INITIAL_MIGRATE, joinSize, rankInGroup1, rankInGroup2, rankInJoin, receiverTasks, rootReceiver, rootSender, and senderTasks.

◆ ParCommGraph() [2/2]

moab::ParCommGraph::ParCommGraph ( const ParCommGraph src)

copy constructor will copy only the senders, receivers, compid1, etc

Definition at line 53 of file ParCommGraph.cpp.

54 {
55  comm = src.comm;
56  senderTasks = src.senderTasks; // these are the sender tasks in joint comm
57  receiverTasks = src.receiverTasks; // these are all the receiver tasks in joint comm
58  rootSender = src.rootSender;
59  rootReceiver = src.rootReceiver;
60  rankInGroup1 = src.rankInGroup1;
61  rankInGroup2 = src.rankInGroup2; // group 1 is sender, 2 is receiver
62  rankInJoin = src.rankInJoin;
63  joinSize = src.joinSize;
64  compid1 = src.compid1;
65  compid2 = src.compid2;
66  comm_graph = NULL;
67  graph_type = src.graph_type;
68  context_id = src.context_id;
69  cover_set = src.cover_set;
70  return;
71 }

References comm, comm_graph, compid1, compid2, context_id, cover_set, graph_type, joinSize, rankInGroup1, rankInGroup2, rankInJoin, receiverTasks, rootReceiver, rootSender, and senderTasks.

Member Function Documentation

◆ compute_partition()

ErrorCode moab::ParCommGraph::compute_partition ( ParallelComm pco,
Range owned,
int  met 
)

Definition at line 1134 of file ParCommGraph.cpp.

1135 {
1136  // we are on a task on sender, and need to compute a new partition;
1137  // primary cells need to be distributed to nb receivers tasks
1138  // first, we will use graph partitioner, with zoltan;
1139  // in the graph that we need to build, the first layer of ghosts is needed;
1140  // can we avoid that ? For example, we can find out from each boundary edge/face what is the
1141  // other cell (on the other side), then form the global graph, and call zoltan in parallel met 1
1142  // would be a geometric partitioner, and met 2 would be a graph partitioner for method 1 we do
1143  // not need any ghost exchange
1144 
1145  // find first edges that are shared
1146  if( owned.empty() ) return MB_SUCCESS; // nothing to do? empty partition is not allowed, maybe we should return
1147  // error?
1148  Core* mb = (Core*)pco->get_moab();
1149 
1150  double t1, t2, t3;
1151  t1 = MPI_Wtime();
1152  int primaryDim = mb->dimension_from_handle( *owned.rbegin() );
1153  int interfaceDim = primaryDim - 1; // should be 1 or 2
1154  Range sharedEdges;
1155 
1156  std::vector< int > shprocs( MAX_SHARING_PROCS );
1157  std::vector< EntityHandle > shhandles( MAX_SHARING_PROCS );
1158 
1159  Tag gidTag = mb->globalId_tag();
1160  int np;
1161  unsigned char pstatus;
1162 
1163  std::multimap< int, int > extraGraphEdges;
1164  // std::map<int, int> adjCellsId;
1165  std::map< int, int > extraCellsProc;
1166  // if method is 2, no need to do the exchange for adjacent cells across partition boundary
1167  // these maps above will be empty for method 2 (geometry)
1168  if( 1 == met )
1169  {
1170  MB_CHK_ERR( pco->get_shared_entities( /*int other_proc*/ -1, sharedEdges, interfaceDim,
1171  /*const bool iface*/ true ) );
1172 
1173 #ifdef VERBOSE
1174  std::cout << " on sender task " << pco->rank() << " number of shared interface cells " << sharedEdges.size()
1175  << "\n";
1176 #endif
1177  // find to what processors we need to send the ghost info about the edge
1178  // first determine the local graph; what elements are adjacent to each cell in owned range
1179  // cells that are sharing a partition interface edge, are identified first, and form a map
1180  TupleList TLe; // tuple list for cells
1181  TLe.initialize( 2, 0, 1, 0, sharedEdges.size() ); // send to, id of adj cell, remote edge
1182  TLe.enableWriteAccess();
1183 
1184  std::map< EntityHandle, int > edgeToCell; // from local boundary edge to adjacent cell id
1185  // will be changed after
1186  for( Range::iterator eit = sharedEdges.begin(); eit != sharedEdges.end(); eit++ )
1187  {
1188  EntityHandle edge = *eit;
1189  // get the adjacent cell
1190  Range adjEnts;
1191  MB_CHK_ERR( mb->get_adjacencies( &edge, 1, primaryDim, false, adjEnts ) );
1192  if( adjEnts.size() > 0 )
1193  {
1194  EntityHandle adjCell = adjEnts[0];
1195  int gid;
1196  MB_CHK_ERR( mb->tag_get_data( gidTag, &adjCell, 1, &gid ) );
1197  MB_CHK_ERR( pco->get_sharing_data( edge, shprocs.data(), shhandles.data(), pstatus, np ) );
1198  int n = TLe.get_n();
1199  TLe.vi_wr[2 * n] = shprocs[0];
1200  TLe.vi_wr[2 * n + 1] = gid;
1201  TLe.vul_wr[n] = shhandles[0]; // the remote edge corresponding to shared edge
1202  edgeToCell[edge] = gid; // store the map between edge and local id of adj cell
1203  TLe.inc_n();
1204  }
1205  }
1206 
1207 #ifdef VERBOSE
1208  std::stringstream ff2;
1209  ff2 << "TLe_" << pco->rank() << ".txt";
1210  TLe.print_to_file( ff2.str().c_str() );
1211 #endif
1212  // send the data to the other processors:
1213  ( pco->proc_config().crystal_router() )->gs_transfer( 1, TLe, 0 );
1214  // on receiver side, each local edge will have the remote cell adjacent to it!
1215 
1216  int ne = TLe.get_n();
1217  for( int i = 0; i < ne; i++ )
1218  {
1219  int sharedProc = TLe.vi_rd[2 * i]; // this info is coming from here, originally
1220  int remoteCellID = TLe.vi_rd[2 * i + 1]; // this is the id of the remote cell, on sharedProc
1221  EntityHandle localCell = TLe.vul_rd[i]; // this is now local edge/face on this proc
1222  int localCellId = edgeToCell[localCell]; // this is the local cell adjacent to edge/face
1223  // now, we will need to add to the graph the pair <localCellId, remoteCellID>
1224  std::pair< int, int > extraAdj = std::make_pair( localCellId, remoteCellID );
1225  extraGraphEdges.insert( extraAdj );
1226  // adjCellsId [edgeToCell[localCell]] = remoteCellID;
1227  extraCellsProc[remoteCellID] = sharedProc;
1228 #ifdef VERBOSE
1229  std::cout << "local ID " << edgeToCell[localCell] << " remote cell ID: " << remoteCellID << "\n";
1230 #endif
1231  }
1232  }
1233  t2 = MPI_Wtime();
1234  if( rootSender ) std::cout << " time preparing the input for Zoltan:" << t2 - t1 << " seconds. \n";
1235  // so adj cells ids; need to call zoltan for parallel partition
1236 #ifdef MOAB_HAVE_ZOLTAN
1237  ZoltanPartitioner* mbZTool = new ZoltanPartitioner( mb, pco );
1238  if( 1 <= met ) // partition in zoltan, either graph or geometric partitioner
1239  {
1240  std::map< int, Range > distribution; // how to distribute owned elements by processors in receiving groups
1241  // in how many tasks do we want to be distributed?
1242  int numNewPartitions = (int)receiverTasks.size();
1243  Range primaryCells = owned.subset_by_dimension( primaryDim );
1244  MB_CHK_ERR( mbZTool->partition_owned_cells( primaryCells, extraGraphEdges, extraCellsProc, numNewPartitions,
1245  distribution, met ) );
1246  for( std::map< int, Range >::iterator mit = distribution.begin(); mit != distribution.end(); mit++ )
1247  {
1248  int part_index = mit->first;
1249  assert( part_index < numNewPartitions );
1250  split_ranges[receiverTasks[part_index]] = mit->second;
1251  }
1252  }
1253  // delete now the partitioner
1254  delete mbZTool;
1255 #endif
1256  t3 = MPI_Wtime();
1257  if( rootSender ) std::cout << " time spent by Zoltan " << t3 - t2 << " seconds. \n";
1258  return MB_SUCCESS;
1259 }

References moab::Range::begin(), moab::ProcConfig::crystal_router(), moab::Range::empty(), moab::TupleList::enableWriteAccess(), moab::Range::end(), moab::ParallelComm::get_moab(), moab::TupleList::get_n(), moab::ParallelComm::get_shared_entities(), moab::ParallelComm::get_sharing_data(), moab::TupleList::inc_n(), moab::TupleList::initialize(), MAX_SHARING_PROCS, mb, MB_CHK_ERR, MB_SUCCESS, ZoltanPartitioner::partition_owned_cells(), moab::TupleList::print_to_file(), moab::ParallelComm::proc_config(), moab::ParallelComm::rank(), moab::Range::rbegin(), receiverTasks, rootSender, moab::Range::size(), split_ranges, moab::Range::subset_by_dimension(), moab::TupleList::vi_rd, moab::TupleList::vi_wr, moab::TupleList::vul_rd, and moab::TupleList::vul_wr.

Referenced by iMOAB_SendMesh().

◆ compute_trivial_partition()

ErrorCode moab::ParCommGraph::compute_trivial_partition ( std::vector< int > &  numElemsPerTaskInGroup1)

Based on the number of elements on each task in group 1, partition for group 2, trivially.

Operations: it is called on every receiver task; decides how are all elements distributed

Note: establish how many elements are sent from each task in group 1 to tasks in group 2 This call is usually made on a root / master process, and will construct local maps that are member data, which contain the communication graph, in both directions Also, number of elements migrated/exchanged between each sender/receiver

Parameters
[in]numElemsPerTaskInGroup1(std::vector<int> &) number of elements on each sender task

Definition at line 100 of file ParCommGraph.cpp.

101 {
102 
103  recv_graph.clear();
104  recv_sizes.clear();
105  sender_graph.clear();
106  sender_sizes.clear();
107 
108  if( numElemsPerTaskInGroup1.size() != senderTasks.size() )
109  return MB_FAILURE; // each sender has a number of elements that it owns
110 
111  // first find out total number of elements to be sent from all senders
112  int total_elems = 0;
113  std::vector< int > accum;
114  accum.push_back( 0 );
115 
116  int num_senders = (int)senderTasks.size();
117 
118  for( size_t k = 0; k < numElemsPerTaskInGroup1.size(); k++ )
119  {
120  total_elems += numElemsPerTaskInGroup1[k];
121  accum.push_back( total_elems );
122  }
123 
124  int num_recv = ( (int)receiverTasks.size() );
125  // in trivial partition, every receiver should get about total_elems/num_receivers elements
126  int num_per_receiver = (int)( total_elems / num_recv );
127  int leftover = total_elems - num_per_receiver * num_recv;
128 
129  // so receiver k will receive [starts[k], starts[k+1] ) interval
130  std::vector< int > starts;
131  starts.resize( num_recv + 1 );
132  starts[0] = 0;
133  for( int k = 0; k < num_recv; k++ )
134  {
135  starts[k + 1] = starts[k] + num_per_receiver;
136  if( k < leftover ) starts[k + 1]++;
137  }
138 
139  // each sender will send to a number of receivers, based on how the
140  // arrays starts[0:num_recv] and accum[0:sendr] overlap
141  int lastUsedReceiverRank = 0; // first receiver was not treated yet
142  for( int j = 0; j < num_senders; j++ )
143  {
144  // we could start the receiver loop with the latest receiver that received from previous
145  // sender
146  for( int k = lastUsedReceiverRank; k < num_recv; k++ )
147  {
148  // if overlap:
149  if( starts[k] < accum[j + 1] && starts[k + 1] > accum[j] )
150  {
151  recv_graph[receiverTasks[k]].push_back( senderTasks[j] );
152  sender_graph[senderTasks[j]].push_back( receiverTasks[k] );
153 
154  // we still need to decide what is the overlap
155  int sizeOverlap = 1; // at least 1, for sure
156  // 1
157  if( starts[k] >= accum[j] ) // one end is starts[k]
158  {
159  if( starts[k + 1] >= accum[j + 1] ) // the other end is accum[j+1]
160  sizeOverlap = accum[j + 1] - starts[k];
161  else //
162  sizeOverlap = starts[k + 1] - starts[k];
163  }
164  else // one end is accum[j]
165  {
166  if( starts[k + 1] >= accum[j + 1] ) // the other end is accum[j+1]
167  sizeOverlap = accum[j + 1] - accum[j];
168  else
169  sizeOverlap = starts[k + 1] - accum[j];
170  }
171  recv_sizes[receiverTasks[k]].push_back( sizeOverlap ); // basically, task k will receive from
172  // sender j, sizeOverlap elems
173  sender_sizes[senderTasks[j]].push_back( sizeOverlap );
174  if( starts[k] > accum[j + 1] )
175  {
176  lastUsedReceiverRank = k - 1; // so next k loop will start a little higher, we
177  // probably finished with first few receivers (up
178  // to receiver lastUsedReceiverRank)
179  break; // break the k loop, we distributed all elements from sender j to some
180  // receivers
181  }
182  }
183  }
184  }
185 
186  return MB_SUCCESS;
187 }

References MB_SUCCESS, receiverTasks, recv_graph, recv_sizes, sender_graph, sender_sizes, and senderTasks.

Referenced by iMOAB_SendMesh().

◆ dump_comm_information()

ErrorCode moab::ParCommGraph::dump_comm_information ( std::string  prefix,
int  is_send 
)

Definition at line 1364 of file ParCommGraph.cpp.

1365 {
1366  //
1367  if( -1 != rankInGroup1 && 1 == is_send ) // it is a sender task
1368  {
1369  std::ofstream dbfile;
1370  std::stringstream outf;
1371  outf << prefix << "_sender_" << rankInGroup1 << "_joint" << rankInJoin << "_type_" << (int)graph_type << ".txt";
1372  dbfile.open( outf.str().c_str() );
1373 
1374  if( graph_type == COVERAGE )
1375  {
1376  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
1377  mit != involved_IDs_map.end(); mit++ )
1378  {
1379  int receiver_proc = mit->first;
1380  std::vector< int >& eids = mit->second;
1381  dbfile << "receiver: " << receiver_proc << " size:" << eids.size() << "\n";
1382  }
1383  }
1384  else if( graph_type == INITIAL_MIGRATE ) // just after migration
1385  {
1386  for( std::map< int, Range >::iterator mit = split_ranges.begin(); mit != split_ranges.end(); mit++ )
1387  {
1388  int receiver_proc = mit->first;
1389  Range& eids = mit->second;
1390  dbfile << "receiver: " << receiver_proc << " size:" << eids.size() << "\n";
1391  }
1392  }
1393  else if( graph_type == DOF_BASED ) // just after migration, or from computeGraph
1394  {
1395  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
1396  mit != involved_IDs_map.end(); mit++ )
1397  {
1398  int receiver_proc = mit->first;
1399  dbfile << "receiver: " << receiver_proc << " size:" << mit->second.size() << "\n";
1400  }
1401  }
1402  dbfile.close();
1403  }
1404  if( -1 != rankInGroup2 && 0 == is_send ) // it is a receiver task
1405  {
1406  std::ofstream dbfile;
1407  std::stringstream outf;
1408  outf << prefix << "_receiver_" << rankInGroup2 << "_joint" << rankInJoin << "_type_" << (int)graph_type
1409  << ".txt";
1410  dbfile.open( outf.str().c_str() );
1411 
1412  if( graph_type == COVERAGE )
1413  {
1414  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
1415  mit != involved_IDs_map.end(); mit++ )
1416  {
1417  int sender_proc = mit->first;
1418  std::vector< int >& eids = mit->second;
1419  dbfile << "sender: " << sender_proc << " size:" << eids.size() << "\n";
1420  }
1421  }
1422  else if( graph_type == INITIAL_MIGRATE ) // just after migration
1423  {
1424  for( std::map< int, Range >::iterator mit = split_ranges.begin(); mit != split_ranges.end(); mit++ )
1425  {
1426  int sender_proc = mit->first;
1427  Range& eids = mit->second;
1428  dbfile << "sender: " << sender_proc << " size:" << eids.size() << "\n";
1429  }
1430  }
1431  else if( graph_type == DOF_BASED ) // just after migration
1432  {
1433  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
1434  mit != involved_IDs_map.end(); mit++ )
1435  {
1436  int sender_proc = mit->first;
1437  dbfile << "receiver: " << sender_proc << " size:" << mit->second.size() << "\n";
1438  }
1439  }
1440  dbfile.close();
1441  }
1442  return MB_SUCCESS;
1443 }

References COVERAGE, DOF_BASED, graph_type, INITIAL_MIGRATE, involved_IDs_map, MB_SUCCESS, rankInGroup1, rankInGroup2, rankInJoin, moab::Range::size(), and split_ranges.

◆ find_group_ranks()

void moab::ParCommGraph::find_group_ranks ( MPI_Group  group,
MPI_Comm  join,
std::vector< int > &  ranks 
)
private

find ranks of a group with respect to an encompassing communicator

Operations: Local, usually called on root process of the group

Parameters
[in]joincomm(MPI_Comm)
[in]group(MPI_Group)
[out]ranks( std::vector<int>) ranks with respect to the joint communicator

Definition at line 81 of file ParCommGraph.cpp.

82 {
83  MPI_Group global_grp;
84  MPI_Comm_group( joincomm, &global_grp );
85 
86  int grp_size;
87 
88  MPI_Group_size( group, &grp_size );
89  std::vector< int > rks( grp_size );
90  ranks.resize( grp_size );
91 
92  for( int i = 0; i < grp_size; i++ )
93  rks[i] = i;
94 
95  MPI_Group_translate_ranks( group, grp_size, rks.data(), global_grp, ranks.data() );
96  MPI_Group_free( &global_grp );
97  return;
98 }

Referenced by ParCommGraph().

◆ get_component_id1()

int moab::ParCommGraph::get_component_id1 ( )
inline

Definition at line 153 of file ParCommGraph.hpp.

154  {
155  return compid1;
156  }

References compid1.

◆ get_component_id2()

int moab::ParCommGraph::get_component_id2 ( )
inline

Definition at line 157 of file ParCommGraph.hpp.

158  {
159  return compid2;
160  }

References compid2.

◆ get_context_id()

int moab::ParCommGraph::get_context_id ( )
inline

Definition at line 162 of file ParCommGraph.hpp.

163  {
164  return context_id;
165  }

References context_id.

◆ get_cover_set()

EntityHandle moab::ParCommGraph::get_cover_set ( )
inline

Definition at line 171 of file ParCommGraph.hpp.

172  {
173  return cover_set;
174  }

References cover_set.

Referenced by iMOAB_ReceiveElementTag(), and iMOAB_SendElementTag().

◆ is_root_receiver()

bool moab::ParCommGraph::is_root_receiver ( )
inline

Definition at line 138 of file ParCommGraph.hpp.

139  {
140  return rootReceiver;
141  }

References rootReceiver.

◆ is_root_sender()

bool moab::ParCommGraph::is_root_sender ( )
inline

Definition at line 133 of file ParCommGraph.hpp.

134  {
135  return rootSender;
136  }

References rootSender.

Referenced by send_graph(), send_graph_partition(), and send_mesh_parts().

◆ pack_receivers_graph()

ErrorCode moab::ParCommGraph::pack_receivers_graph ( std::vector< int > &  packed_recv_array)

pack information about receivers view of the graph, for future sending to receiver root

Operations: Local, called on root process of the senders group

Parameters
[out]packed_recv_arraypacked data will be sent to the root of receivers, and distributed from there, and will have this information, for each receiver, concatenated receiver 1 task, number of senders for receiver 1, then sender tasks for receiver 1, receiver 2 task, number of senders for receiver 2, sender tasks for receiver 2, etc Note: only the root of senders will compute this, and send it over to the receiver root, which will distribute it over each receiver; We do not pack the sizes of data to be sent, only the senders for each of the receivers (could be of size O(n^2) , where n is the number of tasks ; but in general, it should be more like O(n) ). Each sender sends to a "finite" number of receivers, and each receiver receives from a finite number of senders). We need this info to decide how to set up the send/receive waiting game for non-blocking communication )

Definition at line 189 of file ParCommGraph.cpp.

190 {
191  // it will basically look at local data, to pack communication graph, each receiver task will
192  // have to post receives for each sender task that will send data to it; the array will be
193  // communicated to root receiver, and eventually distributed to receiver tasks
194 
195  /*
196  * packed_array will have receiver, number of senders, then senders, etc
197  */
198  if( recv_graph.size() < receiverTasks.size() )
199  {
200  // big problem, we have empty partitions in receive
201  std::cout << " WARNING: empty partitions, some receiver tasks will receive nothing.\n";
202  }
203  for( std::map< int, std::vector< int > >::iterator it = recv_graph.begin(); it != recv_graph.end(); it++ )
204  {
205  int recv = it->first;
206  std::vector< int >& senders = it->second;
207  packed_recv_array.push_back( recv );
208  packed_recv_array.push_back( (int)senders.size() );
209 
210  for( int k = 0; k < (int)senders.size(); k++ )
211  packed_recv_array.push_back( senders[k] );
212  }
213 
214  return MB_SUCCESS;
215 }

References MB_SUCCESS, receiverTasks, recv_graph, and senders().

Referenced by send_graph().

◆ receive_comm_graph()

ErrorCode moab::ParCommGraph::receive_comm_graph ( MPI_Comm  jcomm,
ParallelComm pco,
std::vector< int > &  pack_array 
)

Definition at line 358 of file ParCommGraph.cpp.

359 {
360  // first, receive from sender_rank 0, the communication graph (matrix), so each receiver
361  // knows what data to expect
362  MPI_Comm receive = pco->comm();
363  int size_pack_array, ierr;
364  MPI_Status status;
365  if( rootReceiver )
366  {
367  /*
368  * MPI_Probe(
369  int source,
370  int tag,
371  MPI_Comm comm,
372  MPI_Status* status)
373  *
374  */
375  int mtag = compid2;
376  ierr = MPI_Probe( sender( 0 ), mtag, jcomm, &status );
377  if( 0 != ierr )
378  {
379  std::cout << " MPI_Probe failure: " << ierr << "\n";
380  return MB_FAILURE;
381  }
382  // get the count of data received from the MPI_Status structure
383  ierr = MPI_Get_count( &status, MPI_INT, &size_pack_array );
384  if( 0 != ierr )
385  {
386  std::cout << " MPI_Get_count failure: " << ierr << "\n";
387  return MB_FAILURE;
388  }
389 #ifdef VERBOSE
390  std::cout << " receive comm graph size: " << size_pack_array << "\n";
391 #endif
392  pack_array.resize( size_pack_array );
393  ierr = MPI_Recv( pack_array.data(), size_pack_array, MPI_INT, sender( 0 ), mtag, jcomm, &status );
394  if( 0 != ierr ) return MB_FAILURE;
395 #ifdef VERBOSE
396  std::cout << " receive comm graph ";
397  for( int k = 0; k < (int)pack_array.size(); k++ )
398  std::cout << " " << pack_array[k];
399  std::cout << "\n";
400 #endif
401  }
402 
403  // now broadcast this whole array to all receivers, so they know what to expect
404  ierr = MPI_Bcast( &size_pack_array, 1, MPI_INT, 0, receive );
405  if( 0 != ierr ) return MB_FAILURE;
406  pack_array.resize( size_pack_array );
407  ierr = MPI_Bcast( pack_array.data(), size_pack_array, MPI_INT, 0, receive );
408  if( 0 != ierr ) return MB_FAILURE;
409  return MB_SUCCESS;
410 }

References moab::ParallelComm::comm(), compid2, MB_SUCCESS, rootReceiver, and sender().

Referenced by iMOAB_ReceiveMesh().

◆ receive_mesh()

ErrorCode moab::ParCommGraph::receive_mesh ( MPI_Comm  jcomm,
ParallelComm pco,
EntityHandle  local_set,
std::vector< int > &  senders_local 
)

Definition at line 412 of file ParCommGraph.cpp.

416 {
417  ErrorCode rval;
418  int ierr;
419  MPI_Status status;
420  // we also need to fill corresponding mesh info on the other side
421  corr_tasks = senders_local;
422  Range newEnts;
423 
424  Tag orgSendProcTag; // this will be a tag set on the received mesh, with info about from what
425  // task / PE the
426  // primary element came from, in the joint communicator ; this will be forwarded by coverage
427  // mesh
428  int defaultInt = -1; // no processor, so it was not migrated from somewhere else
429  MB_CHK_SET_ERR( pco->get_moab()->tag_get_handle( "orig_sending_processor", 1, MB_TYPE_INTEGER, orgSendProcTag,
430  MB_TAG_DENSE | MB_TAG_CREAT, &defaultInt ),
431  "can't create original sending processor tag" );
432  int mtag = compid2;
433  if( !senders_local.empty() )
434  {
435  for( size_t k = 0; k < senders_local.size(); k++ )
436  {
437  int sender1 = senders_local[k];
438  // first receive the size of the buffer using probe
439  /*
440  * MPI_Probe(
441  int source,
442  int tag,
443  MPI_Comm comm,
444  MPI_Status* status)
445  *
446  */
447  ierr = MPI_Probe( sender1, mtag, jcomm, &status );
448  if( 0 != ierr )
449  {
450  std::cout << " MPI_Probe failure in ParCommGraph::receive_mesh " << ierr << "\n";
451  return MB_FAILURE;
452  }
453  // get the count of data received from the MPI_Status structure
454  int size_pack;
455  ierr = MPI_Get_count( &status, MPI_CHAR, &size_pack );
456  if( 0 != ierr )
457  {
458  std::cout << " MPI_Get_count failure in ParCommGraph::receive_mesh " << ierr << "\n";
459  return MB_FAILURE;
460  }
461 
462  /* ierr = MPI_Recv (&size_pack, 1, MPI_INT, sender1, 1, jcomm, &status);
463  if (0!=ierr) return MB_FAILURE;*/
464  // now resize the buffer, then receive it
465  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( size_pack );
466  // buffer->reserve(size_pack);
467 
468  ierr = MPI_Recv( buffer->mem_ptr, size_pack, MPI_UNSIGNED_CHAR, sender1, mtag, jcomm, &status );
469  if( 0 != ierr )
470  {
471  std::cout << " MPI_Recv failure in ParCommGraph::receive_mesh " << ierr << "\n";
472  return MB_FAILURE;
473  }
474  // now unpack the buffer we just received
475  Range entities;
476  std::vector< std::vector< EntityHandle > > L1hloc, L1hrem;
477  std::vector< std::vector< int > > L1p;
478  std::vector< EntityHandle > L2hloc, L2hrem;
479  std::vector< unsigned int > L2p;
480 
481  buffer->reset_ptr( sizeof( int ) );
482  std::vector< EntityHandle > entities_vec( entities.size() );
483  std::copy( entities.begin(), entities.end(), entities_vec.begin() );
484  rval = pco->unpack_buffer( buffer->buff_ptr, false, -1, -1, L1hloc, L1hrem, L1p, L2hloc, L2hrem, L2p,
485  entities_vec );
486  delete buffer;
487  if( MB_SUCCESS != rval ) return rval;
488 
489  std::copy( entities_vec.begin(), entities_vec.end(), range_inserter( entities ) );
490  // we have to add them to the local set
491  rval = pco->get_moab()->add_entities( local_set, entities );
492  if( MB_SUCCESS != rval ) return rval;
493  // corr_sizes is the size of primary entities received
494  Range verts = entities.subset_by_dimension( 0 );
495  Range local_primary_ents = subtract( entities, verts );
496  if( local_primary_ents.empty() )
497  {
498  // it is possible that all ents sent were vertices (point cloud)
499  // then consider primary entities the vertices
500  local_primary_ents = verts;
501  }
502  else
503  {
504  // set a tag with the original sender for the primary entity
505  // will be used later for coverage mesh
506  std::vector< int > orig_senders( local_primary_ents.size(), sender1 );
507  rval = pco->get_moab()->tag_set_data( orgSendProcTag, local_primary_ents, orig_senders.data() );
508  }
509  corr_sizes.push_back( (int)local_primary_ents.size() );
510 
511  newEnts.merge( entities );
512  // make these in split ranges
513  split_ranges[sender1] = local_primary_ents;
514 
515 #ifdef VERBOSE
516  std::ostringstream partial_outFile;
517 
518  partial_outFile << "part_send_" << sender1 << "."
519  << "recv" << rankInJoin << ".vtk";
520 
521  // the mesh contains ghosts too, but they are not part of mat/neumann set
522  // write in serial the file, to see what tags are missing
523  std::cout << " writing from receiver " << rankInJoin << " from sender " << sender1
524  << " entities: " << entities.size() << std::endl;
525  rval = pco->get_moab()->write_file( partial_outFile.str().c_str(), 0, 0, &local_set,
526  1 ); // everything on local set received
527  if( MB_SUCCESS != rval ) return rval;
528 #endif
529  }
530  }
531  // make sure adjacencies are updated on the new elements
532 
533  if( newEnts.empty() )
534  {
535  std::cout << " WARNING: this task did not receive any entities \n";
536  }
537  // in order for the merging to work, we need to be sure that the adjacencies are updated
538  // (created)
539  Range local_verts = newEnts.subset_by_type( MBVERTEX );
540  newEnts = subtract( newEnts, local_verts );
541  Core* mb = (Core*)pco->get_moab();
542  AEntityFactory* adj_fact = mb->a_entity_factory();
543  if( !adj_fact->vert_elem_adjacencies() )
544  adj_fact->create_vert_elem_adjacencies();
545  else
546  {
547  for( Range::iterator it = newEnts.begin(); it != newEnts.end(); ++it )
548  {
549  EntityHandle eh = *it;
550  const EntityHandle* conn = NULL;
551  int num_nodes = 0;
552  MB_CHK_ERR( mb->get_connectivity( eh, conn, num_nodes ) );
553  adj_fact->notify_create_entity( eh, conn, num_nodes );
554  }
555  }
556 
557  return MB_SUCCESS;
558 }

References moab::Interface::add_entities(), moab::Range::begin(), buffer, compid2, corr_sizes, corr_tasks, moab::AEntityFactory::create_vert_elem_adjacencies(), moab::Range::empty(), moab::Range::end(), ErrorCode, moab::ParallelComm::get_moab(), mb, MB_CHK_ERR, MB_CHK_SET_ERR, MB_SUCCESS, MB_TAG_CREAT, MB_TAG_DENSE, MB_TYPE_INTEGER, MBVERTEX, moab::Range::merge(), moab::AEntityFactory::notify_create_entity(), rankInJoin, moab::Range::size(), split_ranges, moab::Range::subset_by_dimension(), moab::Range::subset_by_type(), moab::subtract(), moab::Interface::tag_get_handle(), moab::Interface::tag_set_data(), moab::ParallelComm::unpack_buffer(), moab::AEntityFactory::vert_elem_adjacencies(), and moab::Interface::write_file().

Referenced by iMOAB_ReceiveMesh().

◆ receive_tag_values()

ErrorCode moab::ParCommGraph::receive_tag_values ( MPI_Comm  jcomm,
ParallelComm pco,
Range owned,
std::vector< Tag > &  tag_handles 
)

Get the size of the specified tag in bytes

Definition at line 787 of file ParCommGraph.cpp.

791 {
792  // opposite to sending, we will use blocking receives
793  int ierr;
794  MPI_Status status;
795  // basically, owned.size() needs to be equal to sum(corr_sizes)
796  // get info about the tag size, type, etc
797  Core* mb = (Core*)pco->get_moab();
798  // get info about the tag
799  //! Get the size of the specified tag in bytes
800  ErrorCode rval;
801  int total_bytes_per_entity = 0;
802  std::vector< int > vect_bytes_per_tag;
803 #ifdef VERBOSE
804  std::vector< int > tag_sizes;
805 #endif
806  for( size_t i = 0; i < tag_handles.size(); i++ )
807  {
808  int bytes_per_tag;
809  MB_CHK_ERR( mb->tag_get_bytes( tag_handles[i], bytes_per_tag ) );
810  total_bytes_per_entity += bytes_per_tag;
811  vect_bytes_per_tag.push_back( bytes_per_tag );
812 #ifdef VERBOSE
813  int tag_size;
814  MB_CHK_ERR( mb->tag_get_length( tag_handles[i], tag_size ) );
815  tag_sizes.push_back( tag_size );
816 #endif
817  }
818 
819  int mtag = compid1 + compid2;
820 
821  if( graph_type == INITIAL_MIGRATE )
822  {
823  // std::map<int, Range> split_ranges;
824  // rval = split_owned_range ( owned);MB_CHK_ERR ( rval );
825 
826  // use the buffers data structure to allocate memory for receiving the tags
827  for( std::map< int, Range >::iterator it = split_ranges.begin(); it != split_ranges.end(); it++ )
828  {
829  int sender_proc = it->first;
830  Range ents = it->second; // primary entities, with the tag data, we will receive
831  int size_buffer = 4 + total_bytes_per_entity *
832  (int)ents.size(); // hopefully, below 2B; if more, we have a big problem ...
833  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( size_buffer );
834 
835  buffer->reset_ptr( sizeof( int ) );
836 
837  *( (int*)buffer->mem_ptr ) = size_buffer;
838  // int size_pack = buffer->get_current_size(); // debug check
839 
840  ierr = MPI_Recv( buffer->mem_ptr, size_buffer, MPI_UNSIGNED_CHAR, sender_proc, mtag, jcomm, &status );
841  if( ierr != 0 ) return MB_FAILURE;
842  // now set the tag
843  // copy to tag
844 
845  for( size_t i = 0; i < tag_handles.size(); i++ )
846  {
847  rval = mb->tag_set_data( tag_handles[i], ents, (void*)( buffer->buff_ptr ) );
848  buffer->buff_ptr += vect_bytes_per_tag[i] * ents.size();
849  }
850  delete buffer; // no need for it afterwards
851  MB_CHK_ERR( rval );
852  }
853  }
854  else if( graph_type == COVERAGE ) // receive buffer, then extract tag data, in a loop
855  {
856  // we know that we will need to receive some tag data in a specific order (by ids stored)
857  // first, get the ids of the local elements, from owned Range; unpack the buffer in order
858  Tag gidTag = mb->globalId_tag();
859  std::vector< int > gids;
860  gids.resize( owned.size() );
861  MB_CHK_ERR( mb->tag_get_data( gidTag, owned, gids.data() ) );
862  std::map< int, EntityHandle > gidToHandle;
863  size_t i = 0;
864  for( Range::iterator it = owned.begin(); it != owned.end(); it++ )
865  {
866  EntityHandle eh = *it;
867  gidToHandle[gids[i++]] = eh;
868  }
869  //
870  // now, unpack the data and set it to the tag
871  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
872  mit != involved_IDs_map.end(); mit++ )
873  {
874  int sender_proc = mit->first;
875  std::vector< int >& eids = mit->second;
876  int size_buffer = 4 + total_bytes_per_entity *
877  (int)eids.size(); // hopefully, below 2B; if more, we have a big problem ...
878  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( size_buffer );
879  buffer->reset_ptr( sizeof( int ) );
880  *( (int*)buffer->mem_ptr ) = size_buffer; // this is really not necessary, it should receive this too
881 
882  // receive the buffer
883  ierr = MPI_Recv( buffer->mem_ptr, size_buffer, MPI_UNSIGNED_CHAR, sender_proc, mtag, jcomm, &status );
884  if( ierr != 0 ) return MB_FAILURE;
885 // start copy
886 #ifdef VERBOSE
887  std::ofstream dbfile;
888  std::stringstream outf;
889  outf << "recvFrom_" << sender_proc << "_on_proc_" << rankInJoin << ".txt";
890  dbfile.open( outf.str().c_str() );
891  dbfile << "recvFrom_" << sender_proc << " on proc " << rankInJoin << "\n";
892 #endif
893 
894  // copy tag data from buffer->buff_ptr
895  // data is arranged by tag , and repeat the loop for each entity ()
896  // maybe it should be arranged by entity now, not by tag (so one loop for entities,
897  // outside)
898 
899  for( std::vector< int >::iterator it = eids.begin(); it != eids.end(); ++it )
900  {
901  int eID = *it;
902  std::map< int, EntityHandle >::iterator mit2 = gidToHandle.find( eID );
903  if( mit2 == gidToHandle.end() )
904  {
905  std::cout << " on rank: " << rankInJoin << " cannot find entity handle with global ID " << eID
906  << "\n";
907  return MB_FAILURE;
908  }
909  EntityHandle eh = mit2->second;
910  for( i = 0; i < tag_handles.size(); i++ )
911  {
912  MB_CHK_ERR( mb->tag_set_data( tag_handles[i], &eh, 1, (void*)( buffer->buff_ptr ) ) );
913 #ifdef VERBOSE
914  dbfile << "global ID " << eID << " local handle " << mb->id_from_handle( eh ) << " vals: ";
915  double* vals = (double*)( buffer->buff_ptr );
916  for( int kk = 0; kk < tag_sizes[i]; kk++ )
917  {
918  dbfile << " " << *vals;
919  vals++;
920  }
921  dbfile << "\n";
922 #endif
923  buffer->buff_ptr += vect_bytes_per_tag[i];
924  }
925  }
926 
927  // delete receive buffer
928  delete buffer;
929 #ifdef VERBOSE
930  dbfile.close();
931 #endif
932  }
933  }
934  else if( graph_type == DOF_BASED )
935  {
936  // need to fill up the values for each tag, in the order desired, from the buffer received
937  //
938  // get all the tags, for all owned entities, and pack the buffers accordingly
939  // we do not want to get the tags by entity, it may be too expensive
940  std::vector< std::vector< double > > valuesTags;
941  valuesTags.resize( tag_handles.size() );
942  for( size_t i = 0; i < tag_handles.size(); i++ )
943  {
944  int bytes_per_tag;
945  MB_CHK_ERR( mb->tag_get_bytes( tag_handles[i], bytes_per_tag ) );
946  valuesTags[i].resize( owned.size() * bytes_per_tag / sizeof( double ) );
947  // fill the whole array, we will pick up from here
948  // we will fill this array, using data from received buffer
949  // rval = mb->tag_get_data(owned, (void*)( valuesTags[i].data() ) );MB_CHK_ERR ( rval );
950  }
951  // now, unpack the data and set the tags
952  sendReqs.resize( involved_IDs_map.size() );
953  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
954  mit != involved_IDs_map.end(); ++mit )
955  {
956  int sender_proc = mit->first;
957  std::vector< int >& eids = mit->second;
958  std::vector< int >& index_in_values = map_index[sender_proc];
959  std::vector< int >& index_ptr = map_ptr[sender_proc]; // this is eids.size()+1;
960  int size_buffer = 4 + total_bytes_per_entity *
961  (int)eids.size(); // hopefully, below 2B; if more, we have a big problem ...
962  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( size_buffer );
963  buffer->reset_ptr( sizeof( int ) );
964 
965  // receive the buffer
966  ierr = MPI_Recv( buffer->mem_ptr, size_buffer, MPI_UNSIGNED_CHAR, sender_proc, mtag, jcomm, &status );
967  if( ierr != 0 ) return MB_FAILURE;
968  // use the values in buffer to populate valuesTag arrays, fill it up!
969  int j = 0;
970  for( std::vector< int >::iterator it = eids.begin(); it != eids.end(); ++it, ++j )
971  {
972  for( size_t i = 0; i < tag_handles.size(); i++ )
973  {
974  // right now, move just doubles; but it could be any type of tag
975  double val = *( (double*)( buffer->buff_ptr ) );
976  buffer->buff_ptr += 8; // we know we are working with doubles only !!!
977  for( int k = index_ptr[j]; k < index_ptr[j + 1]; k++ )
978  valuesTags[i][index_in_values[k]] = val;
979  }
980  }
981  // we are done with the buffer in which we received tags, release / delete it
982  delete buffer;
983  }
984  // now we populated the values for all tags; set now the tags!
985  for( size_t i = 0; i < tag_handles.size(); i++ )
986  {
987  // we will fill this array, using data from received buffer
988  MB_CHK_ERR( mb->tag_set_data( tag_handles[i], owned, (void*)( valuesTags[i].data() ) ) );
989  }
990  }
991  return MB_SUCCESS;
992 }

References moab::Range::begin(), buffer, compid1, compid2, COVERAGE, DOF_BASED, moab::Range::end(), ErrorCode, moab::ParallelComm::get_moab(), graph_type, INITIAL_MIGRATE, involved_IDs_map, map_index, map_ptr, mb, MB_CHK_ERR, MB_SUCCESS, rankInJoin, sendReqs, moab::Range::size(), and split_ranges.

Referenced by iMOAB_ReceiveElementTag().

◆ receiver()

int moab::ParCommGraph::receiver ( int  index)
inline

Definition at line 148 of file ParCommGraph.hpp.

149  {
150  return receiverTasks[index];
151  }

References moab::index, and receiverTasks.

Referenced by iMOAB_ReceiveMesh(), and send_graph().

◆ receivers()

const std::vector< int >& moab::ParCommGraph::receivers ( )
inline

Definition at line 210 of file ParCommGraph.hpp.

211  {
212  return receiverTasks;
213  }

References receiverTasks.

Referenced by split_owned_range().

◆ release_send_buffers()

ErrorCode moab::ParCommGraph::release_send_buffers ( )

Definition at line 561 of file ParCommGraph.cpp.

562 {
563  int ierr, nsize = (int)sendReqs.size();
564  std::vector< MPI_Status > mult_status;
565  mult_status.resize( sendReqs.size() );
566  ierr = MPI_Waitall( nsize, sendReqs.data(), mult_status.data() );
567 
568  if( ierr != 0 ) return MB_FAILURE;
569  // now we can free all buffers
570  delete[] comm_graph;
571  comm_graph = NULL;
572  std::vector< ParallelComm::Buffer* >::iterator vit;
573  for( vit = localSendBuffs.begin(); vit != localSendBuffs.end(); ++vit )
574  delete( *vit );
575  localSendBuffs.clear();
576  return MB_SUCCESS;
577 }

References comm_graph, localSendBuffs, MB_SUCCESS, and sendReqs.

◆ send_graph()

ErrorCode moab::ParCommGraph::send_graph ( MPI_Comm  jcomm)

use tag 10 to send size and tag 20 to send the packed array

Definition at line 265 of file ParCommGraph.cpp.

266 {
267  if( is_root_sender() )
268  {
269  int ierr;
270  // will need to build a communication graph, because each sender knows now to which receiver
271  // to send data the receivers need to post receives for each sender that will send data to
272  // them will need to gather on rank 0 on the sender comm, global ranks of sender with
273  // receivers to send build communication matrix, each receiver will receive from what sender
274 
275  std::vector< int > packed_recv_array;
276  ErrorCode rval = pack_receivers_graph( packed_recv_array );
277  if( MB_SUCCESS != rval ) return rval;
278 
279  int size_pack_array = (int)packed_recv_array.size();
280  comm_graph = new int[size_pack_array + 1]; // this should be at least size 2
281  comm_graph[0] = size_pack_array;
282  for( int k = 0; k < size_pack_array; k++ )
283  comm_graph[k + 1] = packed_recv_array[k];
284  // will add 2 requests
285  /// use tag 10 to send size and tag 20 to send the packed array
286  sendReqs.resize( 1 );
287  // do not send the size in advance, because we use probe now
288  /*ierr = MPI_Isend( comm_graph.data(), 1, MPI_INT, receiver(0), 10, jcomm, sendReqs.data()); // we
289  have to use global communicator if (ierr!=0) return MB_FAILURE;*/
290  int mtag = compid2;
291  ierr = MPI_Isend( &comm_graph[1], size_pack_array, MPI_INT, receiver( 0 ), mtag, jcomm,
292  sendReqs.data() ); // we have to use global communicator
293  if( ierr != 0 ) return MB_FAILURE;
294  }
295  return MB_SUCCESS;
296 }

References comm_graph, compid2, ErrorCode, is_root_sender(), MB_SUCCESS, pack_receivers_graph(), receiver(), and sendReqs.

Referenced by iMOAB_SendMesh(), and send_graph_partition().

◆ send_graph_partition()

ErrorCode moab::ParCommGraph::send_graph_partition ( ParallelComm pco,
MPI_Comm  jcomm 
)

Definition at line 1263 of file ParCommGraph.cpp.

1264 {
1265  // first, accumulate the info to root of sender; use gatherv
1266  // first, accumulate number of receivers from each sender, to the root receiver
1267  int numberReceivers =
1268  (int)split_ranges.size(); // these are ranges of elements to be sent to each receiver, from this task
1269  int nSenders = (int)senderTasks.size(); //
1270  // this sender will have to send to this many receivers
1271  std::vector< int > displs( 1 ); // displacements for gatherv
1272  std::vector< int > counts( 1 );
1273  if( is_root_sender() )
1274  {
1275  displs.resize( nSenders + 1 );
1276  counts.resize( nSenders );
1277  }
1278 
1279  int ierr = MPI_Gather( &numberReceivers, 1, MPI_INT, counts.data(), 1, MPI_INT, 0, pco->comm() );
1280  if( ierr != MPI_SUCCESS ) return MB_FAILURE;
1281  // compute now displacements
1282  if( is_root_sender() )
1283  {
1284  displs[0] = 0;
1285  for( int k = 0; k < nSenders; k++ )
1286  {
1287  displs[k + 1] = displs[k] + counts[k];
1288  }
1289  }
1290  std::vector< int > buffer;
1291  if( is_root_sender() ) buffer.resize( displs[nSenders] ); // the last one will have the total count now
1292 
1293  std::vector< int > recvs;
1294  for( std::map< int, Range >::iterator mit = split_ranges.begin(); mit != split_ranges.end(); mit++ )
1295  {
1296  recvs.push_back( mit->first );
1297  }
1298  ierr = MPI_Gatherv( recvs.data(), numberReceivers, MPI_INT, buffer.data(), counts.data(), displs.data(), MPI_INT, 0,
1299  pco->comm() );
1300  if( ierr != MPI_SUCCESS ) return MB_FAILURE;
1301 
1302  // now form recv_graph map; points from the
1303  // now form the graph to be sent to the other side; only on root
1304  if( is_root_sender() )
1305  {
1306 #ifdef GRAPH_INFO
1307  std::ofstream dbfileSender;
1308  std::stringstream outf;
1309  outf << "S_" << compid1 << "_R_" << compid2 << "_SenderGraph.txt";
1310  dbfileSender.open( outf.str().c_str() );
1311  dbfileSender << " number senders: " << nSenders << "\n";
1312  dbfileSender << " senderRank \treceivers \n";
1313  for( int k = 0; k < nSenders; k++ )
1314  {
1315  int indexInBuff = displs[k];
1316  int senderTask = senderTasks[k];
1317  dbfileSender << senderTask << "\t\t";
1318  for( int j = 0; j < counts[k]; j++ )
1319  {
1320  int recvTask = buffer[indexInBuff + j];
1321  dbfileSender << recvTask << " ";
1322  }
1323  dbfileSender << "\n";
1324  }
1325  dbfileSender.close();
1326 #endif
1327  for( int k = 0; k < nSenders; k++ )
1328  {
1329  int indexInBuff = displs[k];
1330  int senderTask = senderTasks[k];
1331  for( int j = 0; j < counts[k]; j++ )
1332  {
1333  int recvTask = buffer[indexInBuff + j];
1334  recv_graph[recvTask].push_back( senderTask ); // this will be packed and sent to root receiver, with
1335  // nonblocking send
1336  }
1337  }
1338 #ifdef GRAPH_INFO
1339  std::ofstream dbfile;
1340  std::stringstream outf2;
1341  outf2 << "S_" << compid1 << "_R_" << compid2 << "_RecvGraph.txt";
1342  dbfile.open( outf2.str().c_str() );
1343  dbfile << " number receivers: " << recv_graph.size() << "\n";
1344  dbfile << " receiverRank \tsenders \n";
1345  for( std::map< int, std::vector< int > >::iterator mit = recv_graph.begin(); mit != recv_graph.end(); mit++ )
1346  {
1347  int recvTask = mit->first;
1348  std::vector< int >& senders = mit->second;
1349  dbfile << recvTask << "\t\t";
1350  for( std::vector< int >::iterator vit = senders.begin(); vit != senders.end(); vit++ )
1351  dbfile << *vit << " ";
1352  dbfile << "\n";
1353  }
1354  dbfile.close();
1355 #endif
1356  // this is the same as trivial partition
1357  MB_CHK_ERR( send_graph( jcomm ) );
1358  }
1359 
1360  return MB_SUCCESS;
1361 }

References buffer, moab::ParallelComm::comm(), compid1, compid2, is_root_sender(), MB_CHK_ERR, MB_SUCCESS, recv_graph, send_graph(), senders(), senderTasks, and split_ranges.

Referenced by iMOAB_SendMesh().

◆ send_mesh_parts()

ErrorCode moab::ParCommGraph::send_mesh_parts ( MPI_Comm  jcomm,
ParallelComm pco,
Range owned 
)

not anymore !

Definition at line 300 of file ParCommGraph.cpp.

301 {
302 
303  ErrorCode rval;
304  if( split_ranges.empty() ) // in trivial partition
305  {
306  rval = split_owned_range( rankInGroup1, owned );
307  if( rval != MB_SUCCESS ) return rval;
308  // we know this on the sender side:
310  corr_sizes = sender_sizes[senderTasks[rankInGroup1]]; // another copy
311  }
312  int mtag = compid2;
313  int indexReq = 0;
314  int ierr; // MPI error
315  if( is_root_sender() ) indexReq = 1; // for sendReqs
316  sendReqs.resize( indexReq + split_ranges.size() );
317  for( std::map< int, Range >::iterator it = split_ranges.begin(); it != split_ranges.end(); it++ )
318  {
319  int receiver_proc = it->first;
320  Range ents = it->second;
321 
322  // add necessary vertices too
323  Range verts;
324  rval = pco->get_moab()->get_adjacencies( ents, 0, false, verts, Interface::UNION );
325  if( rval != MB_SUCCESS )
326  {
327  std::cout << " can't get adjacencies. for entities to send\n";
328  return rval;
329  }
330  ents.merge( verts );
331  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( ParallelComm::INITIAL_BUFF_SIZE );
332  buffer->reset_ptr( sizeof( int ) );
333  rval = pco->pack_buffer( ents, false, true, false, -1, buffer );
334  if( rval != MB_SUCCESS )
335  {
336  std::cout << " can't pack buffer for entities to send\n";
337  return rval;
338  }
339  int size_pack = buffer->get_current_size();
340 
341  // TODO there could be an issue with endian things; check !!!!!
342  // we are sending the size of the buffer first as an int!!!
343  /// not anymore !
344  /* ierr = MPI_Isend(buffer->mem_ptr, 1, MPI_INT, receiver_proc, 1, jcomm,
345  &sendReqs[indexReq]); // we have to use global communicator if (ierr!=0) return MB_FAILURE;
346  indexReq++;*/
347 
348  ierr = MPI_Isend( buffer->mem_ptr, size_pack, MPI_UNSIGNED_CHAR, receiver_proc, mtag, jcomm,
349  &sendReqs[indexReq] ); // we have to use global communicator
350  if( ierr != 0 ) return MB_FAILURE;
351  indexReq++;
352  localSendBuffs.push_back( buffer );
353  }
354  return MB_SUCCESS;
355 }

References buffer, compid2, corr_sizes, corr_tasks, ErrorCode, moab::Interface::get_adjacencies(), moab::ParallelComm::get_moab(), moab::ParallelComm::INITIAL_BUFF_SIZE, is_root_sender(), localSendBuffs, MB_SUCCESS, moab::Range::merge(), moab::ParallelComm::pack_buffer(), rankInGroup1, sender_graph, sender_sizes, senderTasks, sendReqs, split_owned_range(), split_ranges, and moab::Interface::UNION.

Referenced by iMOAB_SendMesh().

◆ send_tag_values()

ErrorCode moab::ParCommGraph::send_tag_values ( MPI_Comm  jcomm,
ParallelComm pco,
Range owned,
std::vector< Tag > &  tag_handles 
)

Get the size of the specified tag in bytes

Definition at line 581 of file ParCommGraph.cpp.

585 {
586  // basically, owned.size() needs to be equal to sum(corr_sizes)
587  // get info about the tag size, type, etc
588  int ierr;
589  Core* mb = (Core*)pco->get_moab();
590  // get info about the tag
591  //! Get the size of the specified tag in bytes
592  int total_bytes_per_entity = 0; // we need to know, to allocate buffers
593  ErrorCode rval;
594  std::vector< int > vect_bytes_per_tag;
595 #ifdef VERBOSE
596  std::vector< int > tag_sizes;
597 #endif
598  for( size_t i = 0; i < tag_handles.size(); i++ )
599  {
600  int bytes_per_tag;
601  MB_CHK_ERR( mb->tag_get_bytes( tag_handles[i], bytes_per_tag ) );
602  int tag_size1; // length
603  MB_CHK_ERR( mb->tag_get_length( tag_handles[i], tag_size1 ) );
604  if( graph_type == DOF_BASED )
605  bytes_per_tag = bytes_per_tag / tag_size1; // we know we have one double per tag , per ID sent;
606  // could be 8 for double, 4 for int, etc
607  total_bytes_per_entity += bytes_per_tag;
608  vect_bytes_per_tag.push_back( bytes_per_tag );
609 #ifdef VERBOSE
610  int tag_size;
611  MB_CHK_ERR( mb->tag_get_length( tag_handles[i], tag_size ) );
612  tag_sizes.push_back( tag_size );
613 #endif
614  }
615 
616  int mtag = compid1 + compid2; // used as mpi tag to differentiate a little the messages
617  int indexReq = 0;
618  if( graph_type == INITIAL_MIGRATE ) // original send
619  {
620  // use the buffers data structure to allocate memory for sending the tags
621  sendReqs.resize( split_ranges.size() );
622 
623  for( std::map< int, Range >::iterator it = split_ranges.begin(); it != split_ranges.end(); it++ )
624  {
625  int receiver_proc = it->first;
626  Range ents = it->second; // primary entities, with the tag data
627  int size_buffer = 4 + total_bytes_per_entity *
628  (int)ents.size(); // hopefully, below 2B; if more, we have a big problem ...
629  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( size_buffer );
630 
631  buffer->reset_ptr( sizeof( int ) );
632  for( size_t i = 0; i < tag_handles.size(); i++ )
633  {
634  // copy tag data to buffer->buff_ptr, and send the buffer (we could have used
635  // regular char arrays)
636  MB_CHK_ERR( mb->tag_get_data( tag_handles[i], ents, (void*)( buffer->buff_ptr ) ) );
637  // advance the butter
638  buffer->buff_ptr += vect_bytes_per_tag[i] * ents.size();
639  }
640  *( (int*)buffer->mem_ptr ) = size_buffer;
641  // int size_pack = buffer->get_current_size(); // debug check
642 
643  ierr = MPI_Isend( buffer->mem_ptr, size_buffer, MPI_UNSIGNED_CHAR, receiver_proc, mtag, jcomm,
644  &sendReqs[indexReq] ); // we have to use global communicator
645  if( ierr != 0 ) return MB_FAILURE;
646  indexReq++;
647  localSendBuffs.push_back( buffer ); // we will release them after nonblocking sends are completed
648  }
649  }
650  else if( graph_type == COVERAGE )
651  {
652  // we know that we will need to send some tag data in a specific order
653  // first, get the ids of the local elements, from owned Range; arrange the buffer in order
654  // of increasing global id
655  Tag gidTag = mb->globalId_tag();
656  std::vector< int > gids;
657  gids.resize( owned.size() );
658  MB_CHK_ERR( mb->tag_get_data( gidTag, owned, gids.data() ) );
659  std::map< int, EntityHandle > gidToHandle;
660  size_t i = 0;
661  for( Range::iterator it = owned.begin(); it != owned.end(); it++ )
662  {
663  EntityHandle eh = *it;
664  gidToHandle[gids[i++]] = eh;
665  }
666  // now, pack the data and send it
667  sendReqs.resize( involved_IDs_map.size() );
668  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
669  mit != involved_IDs_map.end(); mit++ )
670  {
671  int receiver_proc = mit->first;
672  std::vector< int >& eids = mit->second;
673  int size_buffer = 4 + total_bytes_per_entity *
674  (int)eids.size(); // hopefully, below 2B; if more, we have a big problem ...
675  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( size_buffer );
676  buffer->reset_ptr( sizeof( int ) );
677 #ifdef VERBOSE
678  std::ofstream dbfile;
679  std::stringstream outf;
680  outf << "from_" << rankInJoin << "_send_to_" << receiver_proc << ".txt";
681  dbfile.open( outf.str().c_str() );
682  dbfile << "from " << rankInJoin << " send to " << receiver_proc << "\n";
683 #endif
684  // copy tag data to buffer->buff_ptr, and send the buffer (we could have used regular
685  // char arrays) pack data by tag, to be consistent with above, even though we loop
686  // through the entities for each tag
687 
688  for( std::vector< int >::iterator it = eids.begin(); it != eids.end(); it++ )
689  {
690  int eID = *it;
691  EntityHandle eh = gidToHandle[eID];
692  for( i = 0; i < tag_handles.size(); i++ )
693  {
694  rval = mb->tag_get_data( tag_handles[i], &eh, 1, (void*)( buffer->buff_ptr ) );
695  if( rval != MB_SUCCESS )
696  {
697  delete buffer; // free parallel comm buffer first
698 
699  MB_SET_ERR( rval, "Tag get data failed" );
700  }
701 #ifdef VERBOSE
702  dbfile << "global ID " << eID << " local handle " << mb->id_from_handle( eh ) << " vals: ";
703  double* vals = (double*)( buffer->buff_ptr );
704  for( int kk = 0; kk < tag_sizes[i]; kk++ )
705  {
706  dbfile << " " << *vals;
707  vals++;
708  }
709  dbfile << "\n";
710 #endif
711  buffer->buff_ptr += vect_bytes_per_tag[i];
712  }
713  }
714 
715 #ifdef VERBOSE
716  dbfile.close();
717 #endif
718  *( (int*)buffer->mem_ptr ) = size_buffer;
719  // int size_pack = buffer->get_current_size(); // debug check
720  ierr = MPI_Isend( buffer->mem_ptr, size_buffer, MPI_UNSIGNED_CHAR, receiver_proc, mtag, jcomm,
721  &sendReqs[indexReq] ); // we have to use global communicator
722  if( ierr != 0 ) return MB_FAILURE;
723  indexReq++;
724  localSendBuffs.push_back( buffer ); // we will release them after nonblocking sends are completed
725  }
726  }
727  else if( graph_type == DOF_BASED )
728  {
729  // need to fill up the buffer, in the order desired, send it
730  // get all the tags, for all owned entities, and pack the buffers accordingly
731  // we do not want to get the tags by entity, it may be too expensive
732  std::vector< std::vector< double > > valuesTags;
733  valuesTags.resize( tag_handles.size() );
734  for( size_t i = 0; i < tag_handles.size(); i++ )
735  {
736  int bytes_per_tag;
737  MB_CHK_ERR( mb->tag_get_bytes( tag_handles[i], bytes_per_tag ) );
738  valuesTags[i].resize( owned.size() * bytes_per_tag / sizeof( double ) );
739  // fill the whole array, we will pick up from here
740  MB_CHK_ERR( mb->tag_get_data( tag_handles[i], owned, (void*)( valuesTags[i].data() ) ) );
741  }
742  // now, pack the data and send it
743  sendReqs.resize( involved_IDs_map.size() );
744  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin();
745  mit != involved_IDs_map.end(); ++mit )
746  {
747  int receiver_proc = mit->first;
748  std::vector< int >& eids = mit->second;
749  std::vector< int >& index_in_values = map_index[receiver_proc];
750  std::vector< int >& index_ptr = map_ptr[receiver_proc]; // this is eids.size()+1;
751  int size_buffer = 4 + total_bytes_per_entity *
752  (int)eids.size(); // hopefully, below 2B; if more, we have a big problem ...
753  ParallelComm::Buffer* buffer = new ParallelComm::Buffer( size_buffer );
754  buffer->reset_ptr( sizeof( int ) );
755 #ifdef VERBOSE
756  std::ofstream dbfile;
757  std::stringstream outf;
758  outf << "from_" << rankInJoin << "_send_to_" << receiver_proc << ".txt";
759  dbfile.open( outf.str().c_str() );
760  dbfile << "from " << rankInJoin << " send to " << receiver_proc << "\n";
761 #endif
762  // copy tag data to buffer->buff_ptr, and send the buffer
763  // pack data by tag, to be consistent with above
764  int j = 0;
765  for( std::vector< int >::iterator it = eids.begin(); it != eids.end(); it++, j++ )
766  {
767  int index_in_v = index_in_values[index_ptr[j]];
768  for( size_t i = 0; i < tag_handles.size(); i++ )
769  {
770  // right now, move just doubles; but it could be any type of tag
771  *( (double*)( buffer->buff_ptr ) ) = valuesTags[i][index_in_v];
772  buffer->buff_ptr += 8; // we know we are working with doubles only !!!
773  }
774  };
775  *( (int*)buffer->mem_ptr ) = size_buffer;
776  // int size_pack = buffer->get_current_size(); // debug check
777  ierr = MPI_Isend( buffer->mem_ptr, size_buffer, MPI_UNSIGNED_CHAR, receiver_proc, mtag, jcomm,
778  &sendReqs[indexReq] ); // we have to use global communicator
779  if( ierr != 0 ) return MB_FAILURE;
780  indexReq++;
781  localSendBuffs.push_back( buffer ); // we will release them after nonblocking sends are completed
782  }
783  }
784  return MB_SUCCESS;
785 }

References moab::Range::begin(), buffer, compid1, compid2, COVERAGE, DOF_BASED, moab::Range::end(), ErrorCode, moab::ParallelComm::get_moab(), graph_type, INITIAL_MIGRATE, involved_IDs_map, localSendBuffs, map_index, map_ptr, mb, MB_CHK_ERR, MB_SET_ERR, MB_SUCCESS, rankInJoin, sendReqs, moab::Range::size(), and split_ranges.

Referenced by iMOAB_SendElementTag().

◆ sender()

int moab::ParCommGraph::sender ( int  index)
inline

Definition at line 143 of file ParCommGraph.hpp.

144  {
145  return senderTasks[index];
146  }

References moab::index, and senderTasks.

Referenced by receive_comm_graph().

◆ senders()

const std::vector< int >& moab::ParCommGraph::senders ( )
inline

Definition at line 206 of file ParCommGraph.hpp.

207  {
208  return senderTasks;
209  } // reference copy; refers to sender tasks in joint comm

References senderTasks.

Referenced by pack_receivers_graph(), and send_graph_partition().

◆ set_context_id()

void moab::ParCommGraph::set_context_id ( int  other_id)
inline

Definition at line 166 of file ParCommGraph.hpp.

167  {
168  context_id = other_id;
169  }

References context_id.

◆ set_cover_set()

void moab::ParCommGraph::set_cover_set ( EntityHandle  cover)
inline

Definition at line 175 of file ParCommGraph.hpp.

176  {
177  cover_set = cover;
178  }

References cover_set.

◆ SetReceivingAfterCoverage()

void moab::ParCommGraph::SetReceivingAfterCoverage ( std::map< int, std::set< int > > &  idsFromProcs)

Definition at line 1024 of file ParCommGraph.cpp.

1026 {
1027  for( auto mt = idsFromProcs.begin(); mt != idsFromProcs.end(); ++mt )
1028  {
1029  int fromProc = mt->first;
1030  std::set< int >& setIds = mt->second;
1031  involved_IDs_map[fromProc].resize( setIds.size() );
1032  std::vector< int >& listIDs = involved_IDs_map[fromProc];
1033  size_t indx = 0;
1034  for( std::set< int >::iterator st = setIds.begin(); st != setIds.end(); st++ )
1035  {
1036  int valueID = *st;
1037  listIDs[indx++] = valueID;
1038  }
1039  }
1040  graph_type = COVERAGE;
1041  return;
1042 }

References COVERAGE, graph_type, and involved_IDs_map.

◆ settle_comm_by_ids()

void moab::ParCommGraph::settle_comm_by_ids ( int  comp,
TupleList TLBackToComp,
std::vector< int > &  valuesComp 
)

Definition at line 1044 of file ParCommGraph.cpp.

1045 {
1046  // settle comm graph on comp
1047  if( rootSender || rootReceiver ) std::cout << " settle comm graph by id on component " << comp << "\n";
1048  int n = TLBackToComp.get_n();
1049  // third_method = true; // do not rely only on involved_IDs_map.size(); this can be 0 in some
1050  // cases
1051  std::map< int, std::set< int > > uniqueIDs;
1052  for( int i = 0; i < n; i++ )
1053  {
1054  int to_proc = TLBackToComp.vi_wr[3 * i + 2];
1055  int globalId = TLBackToComp.vi_wr[3 * i + 1];
1056  uniqueIDs[to_proc].insert( globalId );
1057  }
1058 
1059  // Vector to store element
1060  // with respective present index
1061  std::vector< std::pair< int, int > > vp;
1062  vp.reserve( valuesComp.size() );
1063 
1064  // Inserting element in pair vector
1065  // to keep track of previous indexes in valuesComp
1066  for( size_t i = 0; i < valuesComp.size(); ++i )
1067  {
1068  vp.push_back( std::make_pair( valuesComp[i], i ) );
1069  }
1070  // Sorting pair vector
1071  sort( vp.begin(), vp.end() );
1072 
1073  // vp[i].first, second
1074 
1075  // count now how many times some value appears in ordered (so in valuesComp)
1076  for( auto it = uniqueIDs.begin(); it != uniqueIDs.end(); ++it )
1077  {
1078  int procId = it->first;
1079  std::set< int >& nums = it->second;
1080  std::vector< int >& indx = map_ptr[procId];
1081  std::vector< int >& indices = map_index[procId];
1082  indx.resize( nums.size() + 1 );
1083  int indexInVp = 0;
1084  int indexVal = 0;
1085  indx[0] = 0; // start from 0
1086  for( auto sst = nums.begin(); sst != nums.end(); ++sst, ++indexVal )
1087  {
1088  int val = *sst;
1089  involved_IDs_map[procId].push_back( val );
1090  indx[indexVal + 1] = indx[indexVal];
1091  while( ( indexInVp < (int)valuesComp.size() ) && ( vp[indexInVp].first <= val ) ) // should be equal !
1092  {
1093  if( vp[indexInVp].first == val )
1094  {
1095  indx[indexVal + 1]++;
1096  indices.push_back( vp[indexInVp].second );
1097  }
1098  indexInVp++;
1099  }
1100  }
1101  }
1102 #ifdef VERBOSE
1103  std::stringstream f1;
1104  std::ofstream dbfile;
1105  f1 << "Involve_" << comp << "_" << rankInJoin << ".txt";
1106  dbfile.open( f1.str().c_str() );
1107  for( auto mit = involved_IDs_map.begin(); mit != involved_IDs_map.end(); ++mit )
1108  {
1109  int corrTask = mit->first;
1110  std::vector< int >& corrIds = mit->second;
1111  std::vector< int >& indx = map_ptr[corrTask];
1112  std::vector< int >& indices = map_index[corrTask];
1113 
1114  dbfile << " towards proc " << corrTask << " \n";
1115  for( int i = 0; i < (int)corrIds.size(); i++ )
1116  {
1117  dbfile << corrIds[i] << " [" << indx[i] << "," << indx[i + 1] << ") : ";
1118  for( int j = indx[i]; j < indx[i + 1]; j++ )
1119  dbfile << indices[j] << " ";
1120  dbfile << "\n";
1121  }
1122  dbfile << " \n";
1123  }
1124  dbfile.close();
1125 #endif
1126 
1128  // now we need to fill back and forth information, needed to fill the arrays
1129  // for example, from spectral to involved_IDs_map, in case we want to send data from
1130  // spectral to phys
1131 }

References DOF_BASED, moab::GeomUtil::first(), moab::TupleList::get_n(), graph_type, involved_IDs_map, map_index, map_ptr, rankInJoin, rootReceiver, rootSender, and moab::TupleList::vi_wr.

Referenced by iMOAB_ComputeCommGraph().

◆ settle_send_graph()

ErrorCode moab::ParCommGraph::settle_send_graph ( TupleList TLcovIDs)

Definition at line 996 of file ParCommGraph.cpp.

997 {
998  // fill involved_IDs_map with data
999  // will have "receiving proc" and global id of element
1000  int n = TLcovIDs.get_n();
1001  graph_type = COVERAGE; // do not rely only on involved_IDs_map.size(); this can be 0 in some cases
1002  for( int i = 0; i < n; i++ )
1003  {
1004  int to_proc = TLcovIDs.vi_wr[2 * i];
1005  int globalIdElem = TLcovIDs.vi_wr[2 * i + 1];
1006  involved_IDs_map[to_proc].push_back( globalIdElem );
1007  }
1008 #ifdef VERBOSE
1009  for( std::map< int, std::vector< int > >::iterator mit = involved_IDs_map.begin(); mit != involved_IDs_map.end();
1010  ++mit )
1011  {
1012  std::cout << " towards task " << mit->first << " send: " << mit->second.size() << " cells " << std::endl;
1013  for( size_t i = 0; i < mit->second.size(); i++ )
1014  {
1015  std::cout << " " << mit->second[i];
1016  }
1017  std::cout << std::endl;
1018  }
1019 #endif
1020  return MB_SUCCESS;
1021 }

References COVERAGE, moab::TupleList::get_n(), graph_type, involved_IDs_map, MB_SUCCESS, and moab::TupleList::vi_wr.

◆ split_owned_range() [1/2]

ErrorCode moab::ParCommGraph::split_owned_range ( int  sender_rank,
Range owned 
)

Definition at line 217 of file ParCommGraph.cpp.

218 {
219  int senderTask = senderTasks[sender_rank];
220  std::vector< int >& distribution = sender_sizes[senderTask];
221  std::vector< int >& receivers = sender_graph[senderTask];
222  if( distribution.size() != receivers.size() ) //
223  return MB_FAILURE;
224 
225  Range current = owned; // get the full range first, then we will subtract stuff, for
226  // the following ranges
227 
228  Range rleftover = current;
229  for( size_t k = 0; k < receivers.size(); k++ )
230  {
231  Range newr;
232  newr.insert( current.begin(), current.begin() + distribution[k] );
233  split_ranges[receivers[k]] = newr;
234 
235  rleftover = subtract( current, newr );
236  current = rleftover;
237  }
238 
239  return MB_SUCCESS;
240 }

References moab::Range::begin(), moab::Range::insert(), MB_SUCCESS, receivers(), sender_graph, sender_sizes, senderTasks, split_ranges, and moab::subtract().

Referenced by send_mesh_parts().

◆ split_owned_range() [2/2]

ErrorCode moab::ParCommGraph::split_owned_range ( Range owned)

Definition at line 243 of file ParCommGraph.cpp.

244 {
245  if( corr_tasks.size() != corr_sizes.size() ) //
246  return MB_FAILURE;
247 
248  Range current = owned; // get the full range first, then we will subtract stuff, for
249  // the following ranges
250 
251  Range rleftover = current;
252  for( size_t k = 0; k < corr_tasks.size(); k++ )
253  {
254  Range newr;
255  newr.insert( current.begin(), current.begin() + corr_sizes[k] );
256  split_ranges[corr_tasks[k]] = newr;
257 
258  rleftover = subtract( current, newr );
259  current = rleftover;
260  }
261 
262  return MB_SUCCESS;
263 }

References moab::Range::begin(), corr_sizes, corr_tasks, moab::Range::insert(), MB_SUCCESS, split_ranges, and moab::subtract().

Member Data Documentation

◆ comm

MPI_Comm moab::ParCommGraph::comm
private

Definition at line 244 of file ParCommGraph.hpp.

Referenced by ParCommGraph().

◆ comm_graph

int* moab::ParCommGraph::comm_graph
private

Definition at line 269 of file ParCommGraph.hpp.

Referenced by ParCommGraph(), release_send_buffers(), and send_graph().

◆ compid1

int moab::ParCommGraph::compid1
private

◆ compid2

◆ context_id

int moab::ParCommGraph::context_id
private

Definition at line 252 of file ParCommGraph.hpp.

Referenced by get_context_id(), ParCommGraph(), and set_context_id().

◆ corr_sizes

std::vector< int > moab::ParCommGraph::corr_sizes
private

Definition at line 283 of file ParCommGraph.hpp.

Referenced by receive_mesh(), send_mesh_parts(), and split_owned_range().

◆ corr_tasks

std::vector< int > moab::ParCommGraph::corr_tasks
private

Definition at line 281 of file ParCommGraph.hpp.

Referenced by receive_mesh(), send_mesh_parts(), and split_owned_range().

◆ cover_set

EntityHandle moab::ParCommGraph::cover_set
private

Definition at line 253 of file ParCommGraph.hpp.

Referenced by get_cover_set(), ParCommGraph(), and set_cover_set().

◆ graph_type

◆ involved_IDs_map

std::map< int, std::vector< int > > moab::ParCommGraph::involved_IDs_map
private

◆ joinSize

int moab::ParCommGraph::joinSize
private

Definition at line 250 of file ParCommGraph.hpp.

Referenced by ParCommGraph().

◆ localSendBuffs

std::vector< ParallelComm::Buffer* > moab::ParCommGraph::localSendBuffs
private

Definition at line 266 of file ParCommGraph.hpp.

Referenced by release_send_buffers(), send_mesh_parts(), and send_tag_values().

◆ map_index

std::map< int, std::vector< int > > moab::ParCommGraph::map_index
private

Definition at line 293 of file ParCommGraph.hpp.

Referenced by receive_tag_values(), send_tag_values(), and settle_comm_by_ids().

◆ map_ptr

std::map< int, std::vector< int > > moab::ParCommGraph::map_ptr
private

Definition at line 294 of file ParCommGraph.hpp.

Referenced by receive_tag_values(), send_tag_values(), and settle_comm_by_ids().

◆ rankInGroup1

int moab::ParCommGraph::rankInGroup1
private

Definition at line 249 of file ParCommGraph.hpp.

Referenced by dump_comm_information(), ParCommGraph(), and send_mesh_parts().

◆ rankInGroup2

int moab::ParCommGraph::rankInGroup2
private

Definition at line 249 of file ParCommGraph.hpp.

Referenced by dump_comm_information(), and ParCommGraph().

◆ rankInJoin

int moab::ParCommGraph::rankInJoin
private

◆ receiverTasks

std::vector< int > moab::ParCommGraph::receiverTasks
private

◆ recv_graph

std::map< int, std::vector< int > > moab::ParCommGraph::recv_graph
private

◆ recv_sizes

std::map< int, std::vector< int > > moab::ParCommGraph::recv_sizes
private

Definition at line 260 of file ParCommGraph.hpp.

Referenced by compute_trivial_partition().

◆ rootReceiver

bool moab::ParCommGraph::rootReceiver
private

◆ rootSender

bool moab::ParCommGraph::rootSender
private

◆ sender_graph

std::map< int, std::vector< int > > moab::ParCommGraph::sender_graph
private

Definition at line 262 of file ParCommGraph.hpp.

Referenced by compute_trivial_partition(), send_mesh_parts(), and split_owned_range().

◆ sender_sizes

std::map< int, std::vector< int > > moab::ParCommGraph::sender_sizes
private

Definition at line 264 of file ParCommGraph.hpp.

Referenced by compute_trivial_partition(), send_mesh_parts(), and split_owned_range().

◆ senderTasks

std::vector< int > moab::ParCommGraph::senderTasks
private

◆ sendReqs

std::vector< MPI_Request > moab::ParCommGraph::sendReqs
private

◆ split_ranges

std::map< int, Range > moab::ParCommGraph::split_ranges
private

The documentation for this class was generated from the following files: