]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
Implementation of ParaUMesh::getCellIdsLyingOnNodes V9_5_0a2
authorAnthony Geay <anthony.geay@edf.fr>
Mon, 2 Mar 2020 15:47:26 +0000 (16:47 +0100)
committerAnthony Geay <anthony.geay@edf.fr>
Mon, 30 Mar 2020 08:39:04 +0000 (10:39 +0200)
src/MEDCoupling/MCAuto.hxx
src/MEDCoupling/MEDCouplingMemArray.cxx
src/MEDCoupling/MEDCouplingMemArray.hxx
src/ParaMEDMEM/CMakeLists.txt
src/ParaMEDMEM/CommInterface.hxx
src/ParaMEDMEM/ParaMESH.cxx
src/ParaMEDMEM/ParaMESH.hxx
src/ParaMEDMEM/ParaUMesh.cxx [new file with mode: 0644]
src/ParaMEDMEM/ParaUMesh.hxx [new file with mode: 0644]
src/ParaMEDMEM_Swig/ParaMEDMEMCommon.i

index ff32e25a130849b8a06fd381ee4c3ad1eaab35d1..d1a8a2892c44d9a43bd207bf47367ed3ff3d6606 100644 (file)
@@ -18,8 +18,7 @@
 //
 // Author : Anthony Geay (CEA/DEN)
 
-#ifndef __PARAMEDMEM_MEDCOUPLINGAUTOREFCOUNTOBJECTPTR_HXX__
-#define __PARAMEDMEM_MEDCOUPLINGAUTOREFCOUNTOBJECTPTR_HXX__
+#pragma once
 
 #include "MEDCouplingRefCountObject.hxx"
 #include "InterpKernelException.hxx"
@@ -35,6 +34,7 @@ namespace MEDCoupling
     MCAuto(const MCAuto& other):_ptr(0) { referPtr(other._ptr); }
     MCAuto(T *ptr=0):_ptr(ptr) { }
     ~MCAuto() { destroyPtr(); }
+    void checkNotNull() const { if(!_ptr) throw INTERP_KERNEL::Exception("Pointer is nullptr !"); }
     bool isNull() const { return _ptr==0; }
     bool isNotNull() const { return !isNull(); }
     void nullify() { destroyPtr(); _ptr=0; }
@@ -50,6 +50,7 @@ namespace MEDCoupling
     operator T *() { return _ptr; }
     operator const T *() const { return _ptr; }
     T *retn() { if(_ptr) _ptr->incrRef(); return _ptr; }
+    T *retnConstCast() const { if(_ptr) _ptr->incrRef(); return _ptr; }
     T *iAmATrollConstCast() const { return _ptr; }
   private:
     void referPtr(T *ptr) { _ptr=ptr; if(_ptr) _ptr->incrRef(); }
@@ -121,5 +122,3 @@ namespace MEDCoupling
     const T *_ptr;
   };
 }
-
-#endif
index 309f3a376995e6d5e2972fd3f01e357ee0012144..baa8dbd292b300f3a4fae94592de66549ec7c953 100755 (executable)
@@ -60,6 +60,13 @@ template class MEDCOUPLING_EXPORT MEDCoupling::DataArrayTuple<mcIdType>;
 template class MEDCOUPLING_EXPORT MEDCoupling::DataArrayTuple<double>;
 template class MEDCOUPLING_EXPORT MEDCoupling::DataArrayTuple<float>;
 
+void MEDCoupling::DACheckNbOfTuplesAndComp(const DataArray *da, mcIdType nbOfTuples, std::size_t nbOfCompo, const std::string& msg)
+{
+  if(!da)
+    throw INTERP_KERNEL::Exception("DACheckNbOfTuplesAndComp : null input object !");
+  da->checkNbOfTuplesAndComp(nbOfTuples,nbOfCompo,msg);
+}
+
 template<mcIdType SPACEDIM>
 void DataArrayDouble::findCommonTuplesAlg(const double *bbox, mcIdType nbNodes, mcIdType limitNodeId, double prec, DataArrayIdType *c, DataArrayIdType *cI) const
 {
index e5f4a3873e28f03da14186d483d6e78216c8876a..030f2c013dfd31c52f1ea27d36d0d4f68e787ba3 100755 (executable)
@@ -134,8 +134,11 @@ namespace MEDCoupling
     static mcIdType GetPosOfItemGivenBESRelativeNoThrow(T value, T begin, T end, T step);
   };
 
+  class DataArray;
   class DataArrayByte;
 
+  MEDCOUPLING_EXPORT void DACheckNbOfTuplesAndComp(const DataArray *da, mcIdType nbOfTuples, std::size_t nbOfCompo, const std::string& msg);
+
   class MEDCOUPLING_EXPORT DataArray : public RefCountObject, public TimeLabel
   {
   public:
index 7324dd19d0a4f87f8e30eb529668d4df0dc6d240..5a7c28443763c81e3360b329572754eb6d5384be 100644 (file)
@@ -1,4 +1,4 @@
-# Copyright (C) 2012-2019  CEA/DEN, EDF R&D
+# Copyright (C) 2012-2020  CEA/DEN, EDF R&D
 #
 # This library is free software; you can redistribute it and/or
 # modify it under the terms of the GNU Lesser General Public
@@ -40,6 +40,7 @@ SET(paramedmem_SOURCES
   ProcessorGroup.cxx
   MPIProcessorGroup.cxx
   ParaMESH.cxx
+  ParaUMesh.cxx
   ComponentTopology.cxx
   MPIAccess.cxx
   InterpolationMatrix.cxx
index e76f22e6efe438b99b348fdd0938a17c960b6eb3..9d96cde93e15658800dc66a46b8707f9f6164c16 100644 (file)
 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 
-#ifndef __COMMINTERFACE_HXX__
-#define __COMMINTERFACE_HXX__
+#pragma once
+
+#include "ParaIdType.hxx"
+#include "MEDCouplingMemArray.hxx"
 
 #include <mpi.h>
 
-#include "ParaIdType.hxx"
+#include <memory>
 
 namespace MEDCoupling
 {
@@ -78,10 +80,12 @@ namespace MEDCoupling
     int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
                   void* recvbuf, int recvcount, MPI_Datatype recvtype,
                   MPI_Comm comm) const { return MPI_Allgather(sendbuf,sendcount, sendtype, recvbuf, recvcount, recvtype, comm); }
+    int allGatherV(const void *sendbuf, int sendcount, MPI_Datatype sendtype, void *recvbuf, const int recvcounts[],
+                   const int displs[], MPI_Datatype recvtype, MPI_Comm comm) { return MPI_Allgatherv(sendbuf,sendcount,sendtype,recvbuf,recvcounts,displs,recvtype,comm); }
     int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype,
                  void* recvbuf, int recvcount, MPI_Datatype recvtype,
                  MPI_Comm comm) const { return MPI_Alltoall(sendbuf, sendcount, sendtype, recvbuf, recvcount, recvtype, comm); }
-    int allToAllV(void* sendbuf, int* sendcounts, int* senddispls,
+    int allToAllV(const void* sendbuf, int* sendcounts, int* senddispls,
                   MPI_Datatype sendtype, void* recvbuf, int* recvcounts,
                   int* recvdispls, MPI_Datatype recvtype, 
                   MPI_Comm comm) const { return MPI_Alltoallv(sendbuf, sendcounts, senddispls, sendtype, recvbuf, recvcounts, recvdispls, recvtype, comm); }
@@ -89,7 +93,46 @@ namespace MEDCoupling
     int reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype,
                MPI_Op op, int root, MPI_Comm comm) const { return MPI_Reduce(sendbuf, recvbuf, count, datatype, op, root, comm); }
     int allReduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) const { return MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm); }
+  public:
+
+    /*!
+    * \a counts is expected to be an array of array length. This method returns an array of split array.
+    */
+    static std::unique_ptr<mcIdType[]> SplitArrayOfLength(const std::unique_ptr<mcIdType[]>& counts, std::size_t countsSz, int rk, int size)
+    {
+      std::unique_ptr<mcIdType[]> ret(new mcIdType[countsSz]);
+      for(std::size_t i=0;i<countsSz;++i)
+      {
+        mcIdType a,b;
+        DataArray::GetSlice(0,counts[i],1,rk,size,a,b);
+        ret[i] = b-a;
+      }
+      return ret;
+    }
+
+    /*!
+    * Helper of alltoallv and allgatherv
+    */
+    template<class T>
+    static std::unique_ptr<int []> ToIntArray(const std::unique_ptr<T []>& arr, std::size_t size)
+    {
+      std::unique_ptr<int []> ret(new int[size]);
+      std::copy(arr.get(),arr.get()+size,ret.get());
+      return ret;
+    }
+    
+    /*!
+    * Helper of alltoallv and allgatherv
+    */
+    static std::unique_ptr<int []> ComputeOffset(const std::unique_ptr<int []>& counts, std::size_t sizeOfCounts)
+    {
+      std::unique_ptr<int []> ret(new int[sizeOfCounts]);
+      ret[0] = 0;
+      for(std::size_t i = 1 ; i < sizeOfCounts ; ++i)
+      {
+        ret[i] = ret[i-1] + counts[i-1];
+      }
+      return ret;
+    }
   };
 }
-
-#endif /*COMMINTERFACE_HXX_*/
index 1fe4ddf9cb6481cc4b6b6b1e9edd7cc8ac175075..f0ed51e7591a1f65c05ac7cc686b4d7d927d7a1c 100644 (file)
@@ -36,37 +36,22 @@ namespace MEDCoupling
   ParaMESH::ParaMESH( MEDCouplingPointSet *subdomain_mesh, MEDCouplingPointSet *subdomain_face,
             DataArrayIdType *CorrespElt_local2global, DataArrayIdType *CorrespFace_local2global,
             DataArrayIdType *CorrespNod_local2global, const ProcessorGroup& proc_group ):
-    _cell_mesh(subdomain_mesh),
-    _face_mesh(subdomain_face),
     _my_domain_id(proc_group.myRank()),
-    _block_topology (new BlockTopology(proc_group, subdomain_mesh->getNumberOfCells())),
-    _explicit_topology(0),
-    _node_global(CorrespNod_local2global),
-    _face_global(CorrespFace_local2global),
-    _cell_global(CorrespElt_local2global)
+    _block_topology(new BlockTopology(proc_group, subdomain_mesh->getNumberOfCells())),
+    _explicit_topology(nullptr)
   {
-    if(_cell_mesh)
-      _cell_mesh->incrRef();
-    if(_face_mesh)
-      _face_mesh->incrRef();
-    if(CorrespElt_local2global)
-      CorrespElt_local2global->incrRef();
-    if(CorrespFace_local2global)
-      CorrespFace_local2global->incrRef();
-    if(CorrespNod_local2global)
-      CorrespNod_local2global->incrRef();
+    _cell_mesh.takeRef(subdomain_mesh);
+    _face_mesh.takeRef(subdomain_face);
+    _node_global.takeRef(CorrespNod_local2global);
+    _face_global.takeRef(CorrespFace_local2global);
+    _cell_global.takeRef(CorrespElt_local2global);
   }
 
   ParaMESH::ParaMESH( MEDCouplingPointSet *mesh, const ProcessorGroup& proc_group, const std::string& name):
-    _cell_mesh(mesh),
-    _face_mesh(0),
     _my_domain_id(proc_group.myRank()),
-    _block_topology (new BlockTopology(proc_group, mesh->getNumberOfCells())),
-    _node_global(0),
-    _face_global(0)
+    _block_topology(new BlockTopology(proc_group, mesh->getNumberOfCells()))
   {
-    if(_cell_mesh)
-      _cell_mesh->incrRef();
+    _cell_mesh.takeRef(mesh);
     mcIdType nb_elem=mesh->getNumberOfCells();
     _explicit_topology=new BlockTopology(proc_group,nb_elem);
     mcIdType nbOfCells=mesh->getNumberOfCells();
@@ -82,41 +67,17 @@ namespace MEDCoupling
 
   void ParaMESH::setNodeGlobal(DataArrayIdType *nodeGlobal)
   {
-    if(nodeGlobal!=_node_global)
-      {
-        if(_node_global)
-          _node_global->decrRef();
-        _node_global=nodeGlobal;
-        if(_node_global)
-          _node_global->incrRef();
-      }
+    _node_global.takeRef(nodeGlobal);
   }
 
   void ParaMESH::setCellGlobal(DataArrayIdType *cellGlobal)
   {
-    if(cellGlobal!=_cell_global)
-      {
-        if(_cell_global)
-          _cell_global->decrRef();
-        _cell_global=cellGlobal;
-        if(_cell_global)
-          _cell_global->incrRef();
-      }
+    _cell_global.takeRef(cellGlobal);
   }
 
   ParaMESH::~ParaMESH()
   {
-    if(_cell_mesh)
-      _cell_mesh->decrRef();
-    if(_face_mesh)
-      _face_mesh->decrRef();
     delete _block_topology;
-    if(_node_global)
-      _node_global->decrRef();
-    if(_cell_global)
-      _cell_global->decrRef();
-    if(_face_global)
-      _face_global->decrRef();
     delete _explicit_topology;
   }
 
index eee91dfe93d14322303495316431e2c91e2ed478..f642a95fda8c46246cc13d918f91f93a2441dc08 100644 (file)
@@ -17,8 +17,7 @@
 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 
-#ifndef __PARAMESH_HXX__
-#define __PARAMESH_HXX__
+#pragma once
 
 #include "MEDCouplingPointSet.hxx"
 #include "ProcessorGroup.hxx"
@@ -61,21 +60,21 @@ namespace MEDCoupling
     void setCellGlobal(DataArrayIdType *cellGlobal);
     Topology* getTopology() const { return _explicit_topology; }
     bool isStructured() const { return _cell_mesh->isStructured(); }
-    MEDCouplingPointSet *getCellMesh() const { return _cell_mesh; }
-    MEDCouplingPointSet *getFaceMesh() const { return _face_mesh; }
+    MEDCouplingPointSet *getCellMesh() const { return _cell_mesh.iAmATrollConstCast(); }
+    MEDCouplingPointSet *getFaceMesh() const { return _face_mesh.iAmATrollConstCast(); }
     BlockTopology* getBlockTopology() const { return _block_topology; }
 
-    DataArrayIdType* getGlobalNumberingNodeDA() const { if(_node_global) _node_global->incrRef(); return _node_global; }
-    DataArrayIdType* getGlobalNumberingFaceDA() const { if(_face_global) _face_global->incrRef(); return _face_global; }
-    DataArrayIdType* getGlobalNumberingCellDA() const { if(_cell_global) _cell_global->incrRef(); return _cell_global; }
-    const mcIdType* getGlobalNumberingNode() const { if(_node_global) return _node_global->getConstPointer(); return 0; }
-    const mcIdType* getGlobalNumberingFace() const { if(_face_global) return _face_global->getConstPointer(); return 0; }
-    const mcIdType* getGlobalNumberingCell() const { if(_cell_global) return _cell_global->getConstPointer(); return 0; }
+    DataArrayIdType* getGlobalNumberingNodeDA() const { return _node_global.retnConstCast(); }
+    DataArrayIdType* getGlobalNumberingFaceDA() const { return _face_global.retnConstCast(); }
+    DataArrayIdType* getGlobalNumberingCellDA() const { return _cell_global.retnConstCast(); }
+    const mcIdType* getGlobalNumberingNode() const { if(_node_global) return _node_global->getConstPointer(); return nullptr; }
+    const mcIdType* getGlobalNumberingFace() const { if(_face_global) return _face_global->getConstPointer(); return nullptr; }
+    const mcIdType* getGlobalNumberingCell() const { if(_cell_global) return _cell_global->getConstPointer(); return nullptr; }
 
   private:
     //mesh object underlying the ParaMESH object
-    MEDCouplingPointSet *_cell_mesh ;
-    MEDCouplingPointSet *_face_mesh ;
+    MCAuto<MEDCouplingPointSet> _cell_mesh ;
+    MCAuto<MEDCouplingPointSet> _face_mesh ;
 
     //id of the local grid
     int _my_domain_id;
@@ -84,10 +83,8 @@ namespace MEDCoupling
     MEDCoupling::BlockTopology* _block_topology;
     Topology*  _explicit_topology;
     // pointers to global numberings
-    DataArrayIdType* _node_global;
-    DataArrayIdType* _face_global;
-    DataArrayIdType* _cell_global;
+    MCAuto<DataArrayIdType> _node_global;
+    MCAuto<DataArrayIdType> _face_global;
+    MCAuto<DataArrayIdType> _cell_global;
   };
 }
-
-#endif
diff --git a/src/ParaMEDMEM/ParaUMesh.cxx b/src/ParaMEDMEM/ParaUMesh.cxx
new file mode 100644 (file)
index 0000000..0b939fd
--- /dev/null
@@ -0,0 +1,151 @@
+//
+// Copyright (C) 2020  CEA/DEN, EDF R&D
+//
+// This library is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+//
+// This library is distributed in the hope that it will be useful,
+// but WITHOUT ANY WARRANTY; without even the implied warranty of
+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+// Lesser General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License along with this library; if not, write to the Free Software
+// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
+//
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+// Author : Anthony Geay (EDF R&D)
+
+#include "ParaUMesh.hxx"
+#include "ProcessorGroup.hxx"
+#include "MPIProcessorGroup.hxx"
+#include "Topology.hxx"
+#include "BlockTopology.hxx"
+#include "CommInterface.hxx"
+#include "MEDCouplingMemArray.hxx"
+
+#include "mpi.h"
+
+#include <fstream>
+#include <sstream>
+#include <numeric>
+#include <memory>
+#include <vector>
+
+using namespace std;
+using namespace MEDCoupling;
+
+ParaUMesh::ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds)
+{
+  _mesh.takeRef(mesh);
+  _cell_global.takeRef(globalCellIds);
+  _node_global.takeRef(globalNodeIds);
+  _mesh.checkNotNull();
+  _cell_global.checkNotNull();
+  _node_global.checkNotNull();
+  _mesh->checkConsistencyLight();
+  if(_mesh->getNumberOfNodes() != _node_global->getNumberOfTuples())
+    throw INTERP_KERNEL::Exception("ParaUMesh constructor : mismatch between # nodes and len of global # nodes.");
+  if(_mesh->getNumberOfCells() != _cell_global->getNumberOfTuples())
+    throw INTERP_KERNEL::Exception("ParaUMesh constructor : mismatch between # cells and len of global # cells.");
+}
+
+/*!
+* This method computes the cells part of distributed mesh lying on \a globalNodeIds nodes.
+* The input \a globalNodeIds are not supposed to reside on the current process.
+*/
+MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
+{
+  if(fullyIn)
+    throw INTERP_KERNEL::Exception("ParaUMesh::getCellIdsLyingOnNodes : not implemented yet for fullyIn == True !");
+  MPI_Comm comm(MPI_COMM_WORLD);
+  CommInterface ci;
+  int size;
+  ci.commSize(comm,&size);
+  std::unique_ptr<mcIdType[]> nbOfElems(new mcIdType[size]),nbOfElems2(new mcIdType[size]),nbOfElems3(new mcIdType[size]);
+  mcIdType nbOfNodeIdsLoc(globalNodeIds->getNumberOfTuples());
+  ci.allGather(&nbOfNodeIdsLoc,1,MPI_ID_TYPE,nbOfElems.get(),1,MPI_ID_TYPE,comm);
+  std::vector< MCAuto<DataArrayIdType> > tabs(size);
+  // loop to avoid to all procs to have all the nodes per proc
+  for(int subDiv = 0 ; subDiv < size ; ++subDiv)
+  {
+    std::unique_ptr<mcIdType[]> nbOfElemsSp(CommInterface::SplitArrayOfLength(nbOfElems,size,subDiv,size));
+    mcIdType nbOfNodeIdsSum(std::accumulate(nbOfElemsSp.get(),nbOfElemsSp.get()+size,0));
+    std::unique_ptr<mcIdType[]> allGlobalNodeIds(new mcIdType[nbOfNodeIdsSum]);
+    std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElemsSp,size) );
+    std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) );
+    mcIdType startGlobalNodeIds,endGlobalNodeIds;
+    DataArray::GetSlice(0,globalNodeIds->getNumberOfTuples(),1,subDiv,size,startGlobalNodeIds,endGlobalNodeIds);
+    ci.allGatherV(globalNodeIds->begin()+startGlobalNodeIds,endGlobalNodeIds-startGlobalNodeIds,MPI_ID_TYPE,allGlobalNodeIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm);
+    mcIdType offset(0);
+    for(int curRk = 0 ; curRk < size ; ++curRk)
+    {
+      MCAuto<DataArrayIdType> globalNodeIdsOfCurProc(DataArrayIdType::New());
+      globalNodeIdsOfCurProc->useArray(allGlobalNodeIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElemsSp[curRk],1);
+      offset += nbOfElemsSp[curRk];
+      MCAuto<DataArrayIdType> globalNodeIdsCaptured(_node_global->buildIntersection(globalNodeIdsOfCurProc));
+      MCAuto<DataArrayIdType> localNodeIdsToLocate(_node_global->findIdForEach(globalNodeIdsCaptured->begin(),globalNodeIdsCaptured->end()));
+      MCAuto<DataArrayIdType> localCellCaptured(_mesh->getCellIdsLyingOnNodes(localNodeIdsToLocate->begin(),localNodeIdsToLocate->end(),false));
+      MCAuto<DataArrayIdType> localCellCapturedGlob(_cell_global->selectByTupleIdSafe(localCellCaptured->begin(),localCellCaptured->end()));
+      if(tabs[curRk].isNull())
+        tabs[curRk] = localCellCapturedGlob;
+      else
+        tabs[curRk]->insertAtTheEnd(localCellCapturedGlob->begin(),localCellCapturedGlob->end());
+    }
+  }
+  for(int curRk = 0 ; curRk < size ; ++curRk)
+  {
+    tabs[curRk] = tabs[curRk]->buildUniqueNotSorted();
+    nbOfElems3[curRk] = tabs[curRk]->getNumberOfTuples();
+  }
+  std::vector<const DataArrayIdType *> tabss(tabs.begin(),tabs.end());
+  MCAuto<DataArrayIdType> cells(DataArrayIdType::Aggregate(tabss));
+  ci.allToAll(nbOfElems3.get(),1,MPI_ID_TYPE,nbOfElems2.get(),1,MPI_ID_TYPE,comm);
+  mcIdType nbOfCellIdsSum(std::accumulate(nbOfElems2.get(),nbOfElems2.get()+size,0));
+  MCAuto<DataArrayIdType> cellIdsFromProcs(DataArrayIdType::New());
+  cellIdsFromProcs->alloc(nbOfCellIdsSum,1);
+  {
+    std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElems3,size) ),nbOfElemsOutInt( CommInterface::ToIntArray<mcIdType>(nbOfElems2,size) );
+    std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ), offsetsOut( CommInterface::ComputeOffset(nbOfElemsOutInt,size) );
+    ci.allToAllV(cells->begin(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,
+                 cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_ID_TYPE,comm);
+  }
+  cellIdsFromProcs->sort();
+  return cellIdsFromProcs;
+}
+
+/*!
+ */
+MCAuto<ParaUMesh> ParaUMesh::redistributeCells(const DataArrayIdType *globalCellIds) const
+{
+  MPI_Comm comm(MPI_COMM_WORLD);
+  CommInterface ci;
+  int size;
+  ci.commSize(comm,&size);
+  std::unique_ptr<mcIdType[]> nbOfElems(new mcIdType[size]),nbOfElems2(new mcIdType[size]);
+  mcIdType nbOfNodeIdsLoc(globalCellIds->getNumberOfTuples());
+  ci.allGather(&nbOfNodeIdsLoc,1,MPI_ID_TYPE,nbOfElems.get(),1,MPI_ID_TYPE,comm);
+  mcIdType nbOfCellIdsSum(std::accumulate(nbOfElems.get(),nbOfElems.get()+size,0));
+  std::unique_ptr<mcIdType[]> allGlobalCellIds(new mcIdType[nbOfCellIdsSum]);
+  std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElems,size) );
+  std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) );
+  ci.allGatherV(globalCellIds->begin(),nbOfNodeIdsLoc,MPI_ID_TYPE,allGlobalCellIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm);
+  mcIdType offset(0);
+  for(int curRk = 0 ; curRk < size ; ++curRk)
+  {
+    MCAuto<DataArrayIdType> globalCellIdsOfCurProc(DataArrayIdType::New());
+    globalCellIdsOfCurProc->useArray(allGlobalCellIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElems[curRk],1);
+    offset += nbOfElems[curRk];
+    // prepare all2all session
+    MCAuto<DataArrayIdType> globalCellIdsCaptured(_cell_global->buildIntersection(globalCellIdsOfCurProc));// OK for the global cellIds
+    MCAuto<DataArrayIdType> localCellIdsCaptured(_node_global->findIdForEach(globalCellIdsCaptured->begin(),globalCellIdsCaptured->end()));
+    MCAuto<MEDCouplingUMesh> meshPart(_mesh->buildPartOfMySelf(localCellIdsCaptured->begin(),localCellIdsCaptured->end(),true));
+    MCAuto<DataArrayIdType> o2n(meshPart->zipCoordsTraducer());// OK for the mesh
+    MCAuto<DataArrayIdType> n2o(o2n->invertArrayO2N2N2O(meshPart->getNumberOfNodes()));
+    MCAuto<DataArrayIdType> globalNodeIdsPart(_node_global->selectByTupleIdSafe(n2o->begin(),n2o->end())); // OK for the global nodeIds
+  }
+  
+}
diff --git a/src/ParaMEDMEM/ParaUMesh.hxx b/src/ParaMEDMEM/ParaUMesh.hxx
new file mode 100644 (file)
index 0000000..3d0bff3
--- /dev/null
@@ -0,0 +1,49 @@
+// Copyright (C) 2020  CEA/DEN, EDF R&D
+//
+// This library is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+//
+// This library is distributed in the hope that it will be useful,
+// but WITHOUT ANY WARRANTY; without even the implied warranty of
+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+// Lesser General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License along with this library; if not, write to the Free Software
+// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
+//
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+// Author : Anthony Geay (EDF R&D)
+
+#pragma once
+
+#include "MEDCouplingUMesh.hxx"
+#include "ProcessorGroup.hxx"
+#include "MEDCouplingMemArray.hxx"
+
+#include <string>
+#include <vector>
+
+namespace MEDCoupling
+{
+  /*!
+   * Parallel representation of an unstructured mesh.
+   *
+   * This class is very specific to the requirement of parallel code computations.
+   */
+  class ParaUMesh
+  {
+  public:
+    ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds);
+    MCAuto<DataArrayIdType> getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const;
+    MCAuto<ParaUMesh> redistributeCells(const DataArrayIdType *globalCellIds) const;
+    virtual ~ParaUMesh() { }
+  private:
+    MCAuto<MEDCouplingUMesh> _mesh;
+    MCAuto<DataArrayIdType> _cell_global;
+    MCAuto<DataArrayIdType> _node_global;
+  };
+}
index a75ae3074f2dbe7aaf21a0a72ac871eb5e137440..054d680b936fca4f2022b4cbd6cb2f9464b4e60e 100644 (file)
@@ -34,6 +34,7 @@
 #include "ParaFIELD.hxx"
 #include "ICoCoMEDField.hxx"
 #include "ComponentTopology.hxx"
+#include "ParaUMesh.hxx"
 
 using namespace INTERP_KERNEL;
 using namespace MEDCoupling;
@@ -41,7 +42,6 @@ using namespace ICoCo;
 %}
 
 %include "InterpolationOptions.hxx"
-%include "CommInterface.hxx"
 %include "ProcessorGroup.hxx"
 %include "DECOptions.hxx"
 %include "ParaMESH.hxx"
@@ -57,8 +57,87 @@ using namespace ICoCo;
 %rename(ICoCoMEDField) ICoCo::MEDField;
 %include "ICoCoMEDField.hxx"
 
+%newobject MEDCoupling::ParaUMesh::getCellIdsLyingOnNodes;
+
 %nodefaultctor;
 
+namespace MEDCoupling
+{
+  class CommInterface
+  {
+  public:
+    CommInterface();
+    virtual ~CommInterface();
+    int worldSize() const;
+    int commSize(MPI_Comm comm, int* size) const;
+    int commRank(MPI_Comm comm, int* rank) const;
+    int commGroup(MPI_Comm comm, MPI_Group* group) const;
+    int groupIncl(MPI_Group group, int size, int* ranks, MPI_Group* group_output) const;
+    int commCreate(MPI_Comm comm, MPI_Group group, MPI_Comm* comm_output) const;
+    int groupFree(MPI_Group* group) const;
+    int commFree(MPI_Comm* comm) const;
+
+    int send(void* buffer, int count, MPI_Datatype datatype, int target, int tag, MPI_Comm comm) const;
+    int recv(void* buffer, int count, MPI_Datatype datatype, int source, int tag, MPI_Comm comm, MPI_Status* status) const;
+    int sendRecv(void* sendbuf, int sendcount, MPI_Datatype sendtype, 
+                 int dest, int sendtag, void* recvbuf, int recvcount, 
+                 MPI_Datatype recvtype, int source, int recvtag, MPI_Comm comm,
+                 MPI_Status* status);
+
+    int Isend(void* buffer, int count, MPI_Datatype datatype, int target,
+              int tag, MPI_Comm comm, MPI_Request *request) const;
+    int Irecv(void* buffer, int count, MPI_Datatype datatype, int source,
+              int tag, MPI_Comm comm, MPI_Request* request) const;
+
+    int wait(MPI_Request *request, MPI_Status *status) const;
+    int test(MPI_Request *request, int *flag, MPI_Status *status) const;
+    int requestFree(MPI_Request *request) const;
+    int waitany(int count, MPI_Request *array_of_requests, int *index, MPI_Status *status) const;
+    int testany(int count, MPI_Request *array_of_requests, int *index, int *flag, MPI_Status *status) const;
+    int waitall(int count, MPI_Request *array_of_requests, MPI_Status *array_of_status) const { return MPI_Waitall(count, array_of_requests, array_of_status); }
+    int testall(int count, MPI_Request *array_of_requests, int *flag, MPI_Status *array_of_status) const;
+    int waitsome(int incount, MPI_Request *array_of_requests,int *outcount, int *array_of_indices, MPI_Status *array_of_status) const;
+    int testsome(int incount, MPI_Request *array_of_requests, int *outcount,
+                 int *array_of_indices, MPI_Status *array_of_status) const;
+    int probe(int source, int tag, MPI_Comm comm, MPI_Status *status) const;
+    int Iprobe(int source, int tag, MPI_Comm comm, int *flag, MPI_Status *status) const;
+    int cancel(MPI_Request *request) const;
+    int testCancelled(MPI_Status *status, int *flag) const;
+    int barrier(MPI_Comm comm) const;
+    int errorString(int errorcode, char *string, int *resultlen) const;
+    int getCount(MPI_Status *status, MPI_Datatype datatype, int *count) const;
+
+    int broadcast(void* buffer, int count, MPI_Datatype datatype, int root, MPI_Comm comm) const;
+    int allGather(void* sendbuf, int sendcount, MPI_Datatype sendtype,
+                  void* recvbuf, int recvcount, MPI_Datatype recvtype,
+                  MPI_Comm comm) const;
+    int allToAll(void* sendbuf, int sendcount, MPI_Datatype sendtype,
+                 void* recvbuf, int recvcount, MPI_Datatype recvtype,
+                 MPI_Comm comm) const;
+    int allToAllV(void* sendbuf, int* sendcounts, int* senddispls,
+                  MPI_Datatype sendtype, void* recvbuf, int* recvcounts,
+                  int* recvdispls, MPI_Datatype recvtype, 
+                  MPI_Comm comm) const;
+
+    int reduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, int root, MPI_Comm comm) const;
+    int allReduce(void* sendbuf, void* recvbuf, int count, MPI_Datatype datatype, MPI_Op op, MPI_Comm comm) const;
+  };
+
+  class ParaUMesh
+  {
+  public:
+    ParaUMesh(MEDCouplingUMesh *mesh, DataArrayIdType *globalCellIds, DataArrayIdType *globalNodeIds);
+    %extend
+    {
+      DataArrayIdType *getCellIdsLyingOnNodes(const DataArrayIdType *globalNodeIds, bool fullyIn) const
+      { 
+        MCAuto<DataArrayIdType> ret(self->getCellIdsLyingOnNodes(globalNodeIds,fullyIn));
+        return ret.retn();
+      }
+    }
+  };
+}
+
 /* This object can be used only if MED_ENABLE_FVM is defined*/
 #ifdef MED_ENABLE_FVM
 class NonCoincidentDEC : public DEC