1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459 460 461 462 463
|
/////////////////////////////////////////////////////////////
// //
// Copyright (c) 2003-2014 by The University of Queensland //
// Centre for Geoscience Computing //
// http://earth.uq.edu.au/centre-geoscience-computing //
// //
// Primary Business: Brisbane, Queensland, Australia //
// Licensed under the Open Software License version 3.0 //
// http://www.apache.org/licenses/LICENSE-2.0 //
// //
/////////////////////////////////////////////////////////////
//--- point to point communication primitives for TML_Comm ---
#include "tml/message/packed_message.h"
//--- STL ---
#include <map>
using std::map;
/*!
send a C-array of data with known dimensions
\param data the data to be sent
\param ndata the size of the array (nr. of elements)
\param dest the rank of the destination process the data is sent to
\param tag the message tag
\warning no checks
*/
template <typename T>
void TML_Comm::send_array(T* data,int ndata,int dest,int tag)
{
MPI_Send(data,ndata,GetType(*data),dest,tag,m_comm);
}
/*!
receive a C-array of data with known dimensions
\param data the data to be received
\param ndata the number of integers to be received
\param source the rank of the destination process the data comes from
\param tag the message tag
\warning no checks
*/
template <typename T>
void TML_Comm::receive_array(T* data,int ndata,int source,int tag)
{
MPI_Recv(data,ndata,GetType(*data),source,tag,m_comm,&m_status);
}
/*!
send and receive a C-array of data with known dimensions
\param send_data the data to be sent
\param send_count the number of integers to be sent
\param recv_data the data to be received
\param recv_count the number of integers to be received
\param dest the rank of the destination process the data is sent to
\param source the rank of the destination process the data comes from
\param tag the message tag
\warning no checks
*/
template <typename T,typename P>
void TML_Comm::sendrecv_array(T* send_data,int send_count,P* recv_data,int recv_count,int dest,int source,int tag)
{
MPI_Sendrecv(send_data,send_count,GetType(*send_data),dest,tag,recv_data,recv_count,GetType(*recv_data),source,tag,m_comm,&m_status);
}
/*!
send single data
\param data the data to be sent
\param dest the rank of the destination process the data is sent to
\param tag the message tag
*/
template <typename T>
void TML_Comm::send(T data,int dest,int tag)
{
MPI_Send(&data,1,GetType(data),dest,tag,m_comm);
}
/*!
receive single data
\param data the data to be received
\param source the rank of the destination process the data comes from
\param tag the message tag
\warning no checks
*/
template <typename T>
void TML_Comm::receive(T& data,int source ,int tag)
{
MPI_Recv(&data,1,GetType(data),source,tag,m_comm,&m_status);
}
/*!
send and receive single data
\param send_data the data to be sent
\param recv_data the data to be received
\param dest the rank of the destination process the data is sent to
\param source the rank of the destination process the data comes from
\param tag the message tag
\warning no checks
*/template <typename T,typename P>
void TML_Comm::sendrecv(T send_data,P& recv_data,int dest,int source,int tag)
{
MPI_Sendrecv(&send_data,1,GetType(send_data),dest,tag,&recv_data,1,GetType(recv_data),source,tag,m_comm,&m_status);
}
/*!
Send an STL container or anything that has iterators, begin() and end(). Uses 2 send operations for size and data.
\param data the data to be sent
\param dest the rank of the destination process the data is sent to
\param tag the message tag
*/
template <typename T>
void TML_Comm::send_cont(const T& data, int dest, int tag)
{
int data_size=data.size();
// setup buffer
typename T::value_type *buffer=new typename T::value_type[data_size];
// put data into buffer
int count=0;
for(typename T::const_iterator iter=data.begin();
iter!=data.end();
iter++){
void* buf=reinterpret_cast<void*>(&(buffer[count])); // get memory adress for buffer element
new(buf)(typename T::value_type)(*iter); // initialize object at this adress
// the placement new stuff (see Stroustrup, p.255ff) is necessary
// because assignment buffer[count]=(*iter) doesn't work if the
// T::value_type contains constant member data. Initialization by
// calling the copy constructor directly doesn't work for builtin
// types
count++;
}
//send size
send(data_size,dest,tag);
//send data
send_array(buffer,data_size,dest,tag+1024);
// clean up
delete [] buffer;
}
/*!
Receive an STL container or anything that has iterators, begin() and end(). Single item
insert (a.insert(p,t)) is used instead of range insert to be compatible with both
sequence and associative containers;
\param data the data to be received
\param source the rank of the destination process the data is coming from
\param tag the message tag
*/
template <typename T>
void TML_Comm::receive_cont(T& data,int source,int tag)
{
int data_size;
//get size
receive(data_size,source,tag);
// setup recv buffer
typename T::value_type *buffer=new typename T::value_type[data_size];
//get data
receive_array(buffer,data_size,source,tag+1024);
// insert into container
for(int i=0;i<data_size;i++){
data.insert(data.end(),buffer[i]);
}
delete [] buffer;
}
/*!
send and receive an STL container or anything that has iterators, begin() and end(). Uses 2 MPI_Sendrecv operations for size and data
\param send_data the data to be sent
\param recv_data the data to be received
\param dest the rank of the destination process the data is sent to
\param source the rank of the destination process the data comes from
\param tag the message tag
*/
template <typename T,typename P>
void TML_Comm::sendrecv_cont(T send_data,P& recv_data,int dest,int source,int tag)
{
int send_count=send_data.size();
int recv_count;
// setup send buffer
typename T::value_type *send_buffer=new typename T::value_type[send_count];
// put data into send buffer
int count=0;
for(typename T::const_iterator iter=send_data.begin();
iter!=send_data.end();
iter++){
void* buf=reinterpret_cast<void*>(&(send_buffer[count])); // get memory adress for buffer element
new(buf)(typename T::value_type)(*iter); // initialize object at this adress (see send_cont)
count++;
}
//send/receive size
sendrecv(send_count,recv_count,dest,source,tag);
// check for recv from Null process -> set recv_count to 0 if so
if(source==MPI_PROC_NULL){
recv_count=0;
}
// setup recv buffer
typename T::value_type *recv_buffer=new typename T::value_type[recv_count];
//send/receive data
sendrecv_array(send_buffer,send_count,recv_buffer,recv_count,dest,source,tag+1024);
// insert into container
for(int i=0;i<recv_count;i++){
recv_data.insert(recv_data.end(),recv_buffer[i]);
}
delete [] send_buffer;
delete [] recv_buffer;
}
/*!
send and receive an STL container or anything that has iterators, begin() and end(), erase() and insert(). The input data is replaced (via erase/insert) by the output data. Uses separate send/receive buffers and MPI_Sendrecv instead of a single buffer and MPI_Sendrecv_replace to accomodate different size send and received messages. Uses 2 MPI_Sendrecv operations for size and data.
\param data the data to be sent and received
\param dest the rank of the destination process the data is sent to
\param source the rank of the destination process the data comes from
\param tag the message tag
*/
template <typename T>
void TML_Comm::sendrecv_cont_replace(T& data,int dest,int source,int tag)
{
int send_count=data.size();
int recv_count;
// setup send buffer
typename T::value_type *send_buffer=new typename T::value_type[send_count];
// put data into send buffer
int count=0;
for(typename T::const_iterator iter=data.begin();
iter!=data.end();
iter++){
void* buf=reinterpret_cast<void*>(&(send_buffer[count])); // get memory adress for buffer element
new(buf)(typename T::value_type)(*iter); // initialize object at this adress
count++;
}
//send/receive size
sendrecv(send_count,recv_count,dest,source,tag);
// check for recv from Null process -> set recv_count to 0 if so
if(source==MPI_PROC_NULL){
recv_count=0;
}
// setup recv buffer
typename T::value_type *recv_buffer=new typename T::value_type[recv_count];
//send/receive data
sendrecv_array(send_buffer,send_count,recv_buffer,recv_count,dest,source,tag+1024);
// replace data
data.erase(data.begin(),data.end());
// insert into container
for(int i=0;i<recv_count;i++){
data.insert(data.end(),recv_buffer[i]);
}
delete [] send_buffer;
delete [] recv_buffer;
}
/*!
Send an STL container or anything that has iterators, begin() and end(). Uses 2 send operations for size and data. If the "checked" option is set, a checked message buffer is used to detect type mismatches between pack and unpack operations. The TML_Packed_Message has to grow dynamically. It is not possible to predetermine the size of the message because while T always has a size() function, it is not guaranteed that sizeof() or size() work for the valuetype of T.
\param data the data to be sent
\param dest the rank of the destination process the data is sent to
\param checked use checked message buffer if true
\param tag the message tag
\warning checked stuff not yet implemented
*/
template <typename T>
void TML_Comm::send_cont_packed(T data, int dest ,bool checked ,int tag)
{
TML_Packed_Message* msg=new TML_Packed_Message(m_comm);
int nb_data=data.size();
msg->pack(nb_data); // pack number of items first
// pack data
for(typename T::const_iterator iter=data.begin();
iter!=data.end();
iter++){
msg->pack(*iter);
}
//send size
send(msg->size(),dest,tag);
//send data
send_array(msg->buffer(),msg->size(),dest,tag+1024);
delete msg;
}
/*!
Receive an STL container or anything that has iterators and insert(). Uses 2 receive operations for size and data. If the "checked" option is set, a checked message buffer is used to detect type mismatches between pack and unpack operations.
\param data the data to be received
\param source the rank of the destination process the data is coming from
\param checked use checked message buffer if true
\param tag the message tag
\warning checked stuff not yet implemented
*/
template <typename T>
void TML_Comm::receive_cont_packed(T& data,int source,bool checked,int tag)
{
int msg_size; // total size of the message
int nb_data; // number of data (i.e. the expected data.size())
//get size
receive(msg_size,source,tag);
TML_Packed_Message* msg=new TML_Packed_Message(m_comm,msg_size);
//get data
receive_array(msg->buffer(),msg_size,source,tag+1024);
// extract nuber of items
nb_data=msg->pop_int();
// unpack data
for(int i=0;i<nb_data;i++){
typename T::value_type tv;
msg->unpack(tv);
data.insert(data.end(),tv);
}
delete msg;
}
/*!
send and receive an STL container or anything that has iterators, begin() and end() and insert(). Uses 2 MPI_Sendrecv operations for size and data. If the "checked" option is set, a checked message buffer is used to detect type mismatches between pack and unpack operations.
\param send_data the data to be sent
\param recv_data the data to be received
\param dest the rank of the destination process the data is sent to
\param source the rank of the destination process the data comes from
\param tag the message tag
\warning checked stuff not yet implemented
*/
template <typename T,typename P>
void TML_Comm::sendrecv_cont_packed(T send_data, P& recv_data, int dest, int source, bool checked,int tag)
{
TML_Packed_Message* send_msg=new TML_Packed_Message(m_comm);
int send_nb_data=send_data.size();
int send_msg_size;
int recv_nb_data;
int recv_msg_size;
send_msg->pack(send_nb_data); // pack number of items first
// pack data
for(typename T::const_iterator iter=send_data.begin();
iter!=send_data.end();
iter++){
send_msg->pack(*iter);
}
send_msg_size=send_msg->size();
//send/receive size
sendrecv(send_msg_size,recv_msg_size,dest,source,tag);
// check for recv from Null process -> set recv_count to 0 if so
if(source==MPI_PROC_NULL){
recv_msg_size=0;
}
// setup receive message
TML_Packed_Message* recv_msg=new TML_Packed_Message(m_comm,recv_msg_size);
//send/receive data
sendrecv_array(send_msg->buffer(),send_msg_size,recv_msg->buffer(),recv_msg_size,dest,source,tag+1024);
if(source!=MPI_PROC_NULL){ // if the source exists
// extract nuber of items
recv_nb_data=recv_msg->pop_int();
// unpack data
typename T::value_type tv;
for(int i=0;i<recv_nb_data;i++){
recv_msg->unpack(tv);
recv_data.insert(recv_data.end(),tv);
}
}
delete send_msg;
delete recv_msg;
}
/*!
send and receive an STL container or anything that has iterators, begin() and end() and insert(). In input data is replaced by the received data. Uses 2 MPI_Sendrecv operations for size and data. If the "checked" option is set, a checked message buffer is used to detect type mismatches between pack and unpack operations.
\param data the data to be sent and received
\param dest the rank of the destination process the data is sent to
\param source the rank of the destination process the data comes from
\param tag the message tag
\warning checked stuff not yet implemented
*/
template <typename T>
void TML_Comm::sendrecv_cont_packed_replace(T& data, int dest, int source, bool checked,int tag)
{
TML_Packed_Message* send_msg=new TML_Packed_Message(m_comm);
int send_nb_data=data.size();
int send_msg_size;
int recv_nb_data;
int recv_msg_size;
send_msg->pack(send_nb_data); // pack number of items first
// pack data
for(typename T::const_iterator iter=data.begin();
iter!=data.end();
iter++){
send_msg->pack(*iter);
}
send_msg_size=send_msg->size();
//send/receive size
sendrecv(send_msg_size,recv_msg_size,dest,source,tag);
// check for recv from Null process -> set recv_count to 0 if so
if(source==MPI_PROC_NULL){
recv_msg_size=0;
}
// setup receive message
TML_Packed_Message* recv_msg=new TML_Packed_Message(m_comm,recv_msg_size);
//send/receive data
sendrecv_array(send_msg->buffer(),send_msg_size,recv_msg->buffer(),recv_msg_size,dest,source,tag+1024);
// extract nuber of items
recv_nb_data=recv_msg->pop_int();
// replace data
data.erase(data.begin(),data.end());
for(int i=0;i<recv_nb_data;i++){
typename T::value_type tv;
recv_msg->unpack(tv);
data.insert(data.end(),tv);
}
delete send_msg;
delete recv_msg;
}
|