Page MenuHomec4science

communicator_unity.hh
No OneTemporary

File Metadata

Created
Fri, Jan 10, 15:33

communicator_unity.hh

/**
* @file communicator_unity.hh
*
* @author Guillaume Anciaux <guillaume.anciaux@epfl.ch>
*
* @date Mon Oct 28 19:23:14 2013
*
* @brief Simple LM communicator used in the sequential case
*
* @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_UNITY_HH__
#define __LIBMULTISCALE_COMMUNICATOR_UNITY_HH__
/* -------------------------------------------------------------------------- */
#include "communicator.hh"
/* -------------------------------------------------------------------------- */
__BEGIN_LIBMULTISCALE__
/* -------------------------------------------------------------------------- */
class UnityMPI : public Communicator {
/* ------------------------------------------------------------------------ */
/* Constructors/Destructors */
/* ------------------------------------------------------------------------ */
public:
UnityMPI(){
DUMP("creating unity communicateur", DBG_INFO);
free_procs = 1;
};
~UnityMPI(){}
/* ------------------------------------------------------------------------ */
/* Methods */
/* ------------------------------------------------------------------------ */
UInt addGroup(UInt nb_procs){
if (free_procs == 0 && nb_procs > 1)
LM_FATAL("not enough free processors to make group");
--free_procs;
return 0;
};
bool amIinGroup(CommGroup i){return(lm_my_proc_id == 0);}
MPI_Comm getMpiGroup(CommGroup i){return MPI_COMM_WORLD;};
UInt getNBprocsOnGroup(CommGroup i){return 1;};
void synchronize(CommGroup group_index){};
UInt realRank( UInt i, CommGroup group){return 0;};
UInt groupRank( UInt i, CommGroup group){return 0;};
bool isInGroup(UInt i , CommGroup group) {return true;};
UInt getNBGroups() {return 1;};
//operation de communications
//si un seul proc pas de comm
void sendLocalGeometriesToGroup(Geometry & geom,CommGroup group){};
void receiveLocalGeometriesFromGroup(Geometry ** geom,CommGroup group){};
void sendCommunicationTable(std::vector<UInt> & com_with,CommGroup destgroup){};
void receiveCommunicationTable(std::vector<UInt> & com_with,CommGroup fromgroup){};
void sendUInts(CommBuffer<UInt> & d,UInt nb,UInt dest,CommGroup group,
const std::string & buf){};
void receiveUInts(CommBuffer<UInt> & d,UInt nb,UInt from,CommGroup group,const std::string & buf){};
void sendReals(CommBuffer<Real> & d,UInt nb,UInt dest,CommGroup group,
const std::string & buf){};
void receiveReals(CommBuffer<Real> & d,UInt nb,UInt from,CommGroup group,const std::string & buf){};
void reduceUInt(CommBuffer<UInt> & contrib,UInt nb,CommGroup group,
const std::string & comment,Operator op){};
void reduceReal(CommBuffer<Real> & contrib,UInt nb,CommGroup group,
const std::string & comment, Operator op){};
void allReduceUInt(CommBuffer<UInt> & contrib,UInt nb,CommGroup group,
const std::string & comment,Operator op){};
void allReduceReal(CommBuffer<Real> & contrib,UInt nb,CommGroup group,
const std::string & comment, Operator op){};
void waitForPendingComs(){};
/* ------------------------------------------------------------------------ */
/* Class Members */
/* ------------------------------------------------------------------------ */
protected:
UInt free_procs;
};
/* -------------------------------------------------------------------------- */
__END_LIBMULTISCALE__
#endif /* __LIBMULTISCALE_COMMUNICATOR_UNITY_HH__ */

Event Timeline