Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F101569339
communicator.hh
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Tue, Feb 11, 16:11
Size
16 KB
Mime Type
text/x-c
Expires
Thu, Feb 13, 16:11 (2 d)
Engine
blob
Format
Raw Data
Handle
24187770
Attached To
rLIBMULTISCALE LibMultiScale
communicator.hh
View Options
/**
* @file communicator.hh
*
* @author Guillaume Anciaux <guillaume.anciaux@epfl.ch>
*
* @date Fri Nov 08 12:44:53 2013
*
* @brief Mother class of LM communicators
*
* @section LICENSE
*
* Copyright INRIA and CEA
*
* The LibMultiScale is a C++ parallel framework for the multiscale
* coupling methods dedicated to material simulations. This framework
* provides an API which makes it possible to program coupled simulations
* and integration of already existing codes.
*
* This Project was initiated in a collaboration between INRIA Futurs Bordeaux
* within ScAlApplix team and CEA/DPTA Ile de France.
* The project is now continued at the Ecole Polytechnique Fédérale de Lausanne
* within the LSMS/ENAC laboratory.
*
* This software is governed by the CeCILL-C license under French law and
* abiding by the rules of distribution of free software. You can use,
* modify and/ or redistribute the software under the terms of the CeCILL-C
* license as circulated by CEA, CNRS and INRIA at the following URL
* "http://www.cecill.info".
*
* As a counterpart to the access to the source code and rights to copy,
* modify and redistribute granted by the license, users are provided only
* with a limited warranty and the software's author, the holder of the
* economic rights, and the successive licensors have only limited
* liability.
*
* In this respect, the user's attention is drawn to the risks associated
* with loading, using, modifying and/or developing or reproducing the
* software by the user in light of its specific status of free software,
* that may mean that it is complicated to manipulate, and that also
* therefore means that it is reserved for developers and experienced
* professionals having in-depth computer knowledge. Users are therefore
* encouraged to load and test the software's suitability as regards their
* requirements in conditions enabling the security of their systems and/or
* data to be ensured and, more generally, to use and operate it in the
* same conditions as regards security.
*
* The fact that you are presently reading this means that you have had
* knowledge of the CeCILL-C license and that you accept its terms.
*
*/
#ifndef __LIBMULTISCALE_COMMUNICATOR_HH__
#define __LIBMULTISCALE_COMMUNICATOR_HH__
/* -------------------------------------------------------------------------- */
#include "comm_buffer.hh"
#include <mpi.h>
#include "tensor.hh"
/* -------------------------------------------------------------------------- */
__BEGIN_LIBMULTISCALE__
/* -------------------------------------------------------------------------- */
class
Geometry
;
/* -------------------------------------------------------------------------- */
class
Communicator
{
public
:
/* -------------------------------------------------------------------------- */
// Constructor/desctructor
/* -------------------------------------------------------------------------- */
virtual
~
Communicator
(){};
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
//! return the communicator used
static
Communicator
&
getCommunicator
(){
if
(
!
communicator
)
LM_THROW
(
"communicator not already created"
);
return
*
communicator
;
};
//! set the communicator
template
<
typename
T
>
static
void
createCommunicator
(){
if
(
communicator
)
LM_THROW
(
"communicator was already created"
);
communicator
=
new
T
();
};
/* -------------------------------------------------------------------------- */
//Virtual functions
/* -------------------------------------------------------------------------- */
//! test if am i in the group specified by ID
virtual
bool
amIinGroup
(
CommGroup
group_index
)
=
0
;
//! add a group of nb_proc processors
virtual
UInt
addGroup
(
UInt
nb_procs
)
=
0
;
//! return the equivalent MPI group from one group ID
virtual
MPI_Comm
getMpiGroup
(
CommGroup
group_index
)
=
0
;
//! barrier for a given group (no param is a barrier on all processors)
virtual
void
synchronize
(
CommGroup
group_index
=
group_all
)
=
0
;
//! return number of processors present in group specified
virtual
UInt
getNBprocsOnGroup
(
CommGroup
i
)
=
0
;
//! return my real rank from a group rank and group ID
virtual
UInt
realRank
(
UInt
i
,
CommGroup
group
)
=
0
;
//! return my group rank from a real rank and group ID
virtual
UInt
groupRank
(
UInt
i
,
CommGroup
group
)
=
0
;
virtual
bool
isInGroup
(
UInt
i
,
CommGroup
group
)
=
0
;
virtual
UInt
getNBGroups
()
=
0
;
//! wait for unsynchronisated comunications
virtual
void
waitForPendingComs
()
=
0
;
//! global comunication over a group to send geometries structures
virtual
void
sendLocalGeometriesToGroup
(
Geometry
&
geom
,
CommGroup
group
)
=
0
;
//! global comunication over a group to receive geometries structures
virtual
void
receiveLocalGeometriesFromGroup
(
Geometry
**
geom
,
CommGroup
group
)
=
0
;
//! global comunication Over a group to send boolean tables
virtual
void
sendCommunicationTable
(
std
::
vector
<
UInt
>
&
com_with
,
CommGroup
destgroup
)
=
0
;
//! global comunication over a group to receive boolean tables
virtual
void
receiveCommunicationTable
(
std
::
vector
<
UInt
>
&
com_with
,
CommGroup
fromgroup
)
=
0
;
//! receive array from a processor on by its ID in a given group
virtual
void
receiveUInts
(
CommBuffer
<
UInt
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
)
=
0
;
//! send array to a processor on by its ID in a given group
virtual
void
sendUInts
(
CommBuffer
<
UInt
>
&
i
,
UInt
nb
,
UInt
dest
,
CommGroup
group
,
const
std
::
string
&
comment
)
=
0
;
//! receive array from a processor on by its ID in a given group
virtual
void
receiveReals
(
CommBuffer
<
Real
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
)
=
0
;
//! send array to a processor on by its ID in a given group
virtual
void
sendReals
(
CommBuffer
<
Real
>
&
i
,
UInt
nb
,
UInt
dest
,
CommGroup
group
,
const
std
::
string
&
comment
)
=
0
;
//! reduction operation over a given group
virtual
void
reduceUInt
(
CommBuffer
<
UInt
>
&
contrib
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
)
=
0
;
//! reduction operation over a given group
virtual
void
reduceReal
(
CommBuffer
<
Real
>
&
contrib
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
)
=
0
;
//! reduction operation over a given group
virtual
void
allReduceUInt
(
CommBuffer
<
UInt
>
&
contrib
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
)
=
0
;
//! reduction operation over a given group
virtual
void
allReduceReal
(
CommBuffer
<
Real
>
&
contrib
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
)
=
0
;
//! printself function
virtual
void
printself
(
std
::
ostream
&
stream
){};
/* -------------------------------------------------------------------------- */
// Generic functions
/* -------------------------------------------------------------------------- */
//receive functions
//! receive an nb-data array from a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
receive
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! receive array from a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
receive
(
std
::
vector
<
T
>
&
d
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! receive an nb-data array from a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
receive
(
CommBuffer
<
T
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! receive array from a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
receive
(
CommBuffer
<
T
>
&
d
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! receive a C-style array from a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
receive
(
T
*
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
);
//sending functions
//! send an nb-data array to a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
send
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! send array to a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
send
(
std
::
vector
<
T
>
&
d
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! send an nb-data array to a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
send
(
CommBuffer
<
T
>
&
d
,
UInt
nb
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! send array to a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
send
(
CommBuffer
<
T
>
&
d
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
);
//! send a C-style array to a processor specified by its ID in a given group
template
<
typename
T
>
inline
void
send
(
T
*
d
,
UInt
nb
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
);
//reducing functions
//! reduce by sum from a commbuffer
template
<
typename
T
>
inline
void
reduce
(
CommBuffer
<
T
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
);
//! reduce from a vector
template
<
typename
T
>
inline
void
reduce
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
);
//! reduce a C-style array (for compatibility : involve additional copy)
template
<
typename
T
>
inline
void
reduce
(
T
*
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
);
//! reduce by sum from a commbuffer
template
<
typename
T
>
inline
void
allReduce
(
CommBuffer
<
T
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
);
//! reduce from a vector
template
<
typename
T
>
inline
void
allReduce
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
);
//! reduce a C-style array (for compatibility : involve additional copy)
template
<
typename
T
>
inline
void
allReduce
(
T
*
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
);
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
static
Communicator
*
communicator
;
};
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
receive
(
std
::
vector
<
T
>
&
d
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
CommBuffer
<
T
>
buf
(
d
);
receive
(
buf
,
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
receive
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
CommBuffer
<
T
>
buf
(
d
);
receive
(
buf
,
nb
,
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
send
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
CommBuffer
<
T
>
buf
(
d
);
send
(
buf
,
nb
,
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
send
(
std
::
vector
<
T
>
&
d
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
send
(
d
,
d
.
size
(),
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
receive
(
CommBuffer
<
T
>
&
d
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
receive
(
d
,
UINT_MAX
,
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
receive
(
CommBuffer
<
Real
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
receiveReals
(
d
,
nb
,
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
send
(
CommBuffer
<
Real
>
&
d
,
UInt
nb
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
){
sendReals
(
d
,
nb
,
to
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
receive
(
CommBuffer
<
UInt
>
&
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
receiveUInts
(
d
,
nb
,
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
send
(
CommBuffer
<
UInt
>
&
d
,
UInt
nb
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
){
sendUInts
(
d
,
nb
,
to
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
send
(
CommBuffer
<
T
>
&
d
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
send
(
d
,
d
.
size
(),
from
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
reduce
(
CommBuffer
<
Real
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
reduceReal
(
d
,
nb
,
group
,
comment
,
op
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
reduce
(
CommBuffer
<
UInt
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
reduceUInt
(
d
,
nb
,
group
,
comment
,
op
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
reduce
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
CommBuffer
<
T
>
buf
(
d
);
reduce
(
buf
,
nb
,
group
,
comment
,
op
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
reduce
(
T
*
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
CommBuffer
<
T
>
buf
;
buf
.
resize
(
nb
);
std
::
copy
(
d
,
d
+
nb
,
buf
.
buffer
());
reduce
(
buf
,
nb
,
group
,
comment
,
op
);
std
::
copy
(
buf
.
buffer
(),
buf
.
buffer
()
+
nb
,
d
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
allReduce
(
CommBuffer
<
Real
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
allReduceReal
(
d
,
nb
,
group
,
comment
,
op
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
allReduce
(
CommBuffer
<
UInt
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
allReduceUInt
(
d
,
nb
,
group
,
comment
,
op
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
allReduce
(
std
::
vector
<
T
>
&
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
CommBuffer
<
T
>
buf
(
d
);
allReduce
(
buf
,
nb
,
group
,
comment
,
op
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
allReduce
(
T
*
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
CommBuffer
<
T
>
buf
;
buf
.
resize
(
nb
);
std
::
copy
(
d
,
d
+
nb
,
buf
.
buffer
());
allReduce
(
buf
,
nb
,
group
,
comment
,
op
);
std
::
copy
(
buf
.
buffer
(),
buf
.
buffer
()
+
nb
,
d
);
}
/* -------------------------------------------------------------------------- */
template
<>
inline
void
Communicator
::
allReduce
(
Tensor
*
d
,
UInt
nb
,
CommGroup
group
,
const
std
::
string
&
comment
,
Operator
op
){
allReduce
(
d
->
getVec
(),
9
,
group
,
comment
,
op
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
send
(
T
*
d
,
UInt
nb
,
UInt
to
,
CommGroup
group
,
const
std
::
string
&
comment
){
CommBuffer
<
T
>
buf
;
buf
.
resize
(
nb
);
std
::
copy
(
d
,
d
+
nb
,
buf
.
buffer
());
send
(
buf
,
nb
,
to
,
group
,
comment
);
}
/* -------------------------------------------------------------------------- */
template
<
typename
T
>
inline
void
Communicator
::
receive
(
T
*
d
,
UInt
nb
,
UInt
from
,
CommGroup
group
,
const
std
::
string
&
comment
){
CommBuffer
<
T
>
buf
;
buf
.
resize
(
nb
);
receive
(
buf
,
nb
,
from
,
group
,
comment
);
std
::
copy
(
buf
.
buffer
(),
buf
.
buffer
()
+
nb
,
d
);
}
/* -------------------------------------------------------------------------- */
inline
std
::
ostream
&
operator
<<
(
std
::
ostream
&
os
,
Communicator
&
c
){
c
.
printself
(
os
);
return
os
;
}
/* -------------------------------------------------------------------------- */
__END_LIBMULTISCALE__
//#include "communicator_world.hh"
//__BEGIN_LIBMULTISCALE__
//static WorldCommunicator WORLD_COMMUNICATOR;
///* -------------------------------------------------------------------------- */
//__END_LIBMULTISCALE__
#endif
/* __LIBMULTISCALE_COMMUNICATOR_HH__ */
Event Timeline
Log In to Comment