]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
WIP
authorAnthony Geay <anthony.geay@edf.fr>
Fri, 6 Mar 2020 22:23:34 +0000 (23:23 +0100)
committerAnthony Geay <anthony.geay@edf.fr>
Fri, 6 Mar 2020 22:23:34 +0000 (23:23 +0100)
src/ParaMEDMEM/CommInterface.hxx
src/ParaMEDMEM/ParaUMesh.cxx

index e76f22e6efe438b99b348fdd0938a17c960b6eb3..76692c86ec99f68d05402b39e1a23036db6d4405 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 <mpi.h>
 
-#include "ParaIdType.hxx"
+#include <memory>
 
 namespace MEDCoupling
 {
@@ -89,7 +90,21 @@ 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:
+    static std::unique_ptr<int []> ToIntArray(const std::unique_ptr<long []>& arr, std::size_t size)
+    {
+      std::unique_ptr<int []> ret(new int[size]);
+      std::copy(arr.get(),arr.get()+size,ret.get());
+    }
+    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 e1d873d8f8046a92fc3e86120dc14c406b44c596..1db9288dc643118341579f4c7b7dd611eb91b530 100644 (file)
@@ -24,6 +24,7 @@
 #include "MPIProcessorGroup.hxx"
 #include "Topology.hxx"
 #include "BlockTopology.hxx"
+#include "CommInterface.hxx"
 #include "MEDCouplingMemArray.hxx"
 
 #include "mpi.h"
@@ -56,25 +57,41 @@ MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodes(DataArrayIdType *globa
   int rk,size;
   MPI_Comm_rank(comm,&rk);
   MPI_Comm_size(comm,&size);
-  std::unique_ptr<long[]> nbOfElems(new long[size]);
+  std::unique_ptr<long[]> nbOfElems(new long[size]),nbOfElems2(new long[size]);
   long nbOfNodeIdsLoc(FromIdType<long>(globalNodeIds->getNumberOfTuples()));
   MPI_Allgather(&nbOfNodeIdsLoc,1,MPI_LONG,nbOfElems.get(),size,MPI_LONG,comm);
   long nbOfNodeIdsSum(std::accumulate(nbOfElems.get(),nbOfElems.get()+size,0L));
-  std::unique_ptr<int []> allGlobalNodeIds(new int[nbOfNodeIdsSum]);
-  MPI_Allgather(globalNodeIds->begin(),1,MPI_INT,allGlobalNodeIds.get(),size,MPI_INT,comm);
-  // compute for each proc the cell contribution among _mesh cells
-  mcIdType offset(0);
-  for(int curRk = 0 ; curRk < size ; ++curRk)
+  std::vector< MCAuto<DataArrayIdType> > tabs(size);
   {
-      MCAuto<DataArrayIdType> globalNodeIdsOfCurProc(DataArrayIdType::New());
-      globalNodeIdsOfCurProc->useArray(allGlobalNodeIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElems[curRk],1);
-      offset += nbOfElems[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()));
-      
+    std::unique_ptr<int []> allGlobalNodeIds(new int[nbOfNodeIdsSum]);
+    MPI_Allgather(globalNodeIds->begin(),1,MPI_INT,allGlobalNodeIds.get(),1,MPI_INT,comm);
+    // compute for each proc the cell contribution among _mesh cells
+    long offset(0);
+    for(int curRk = 0 ; curRk < size ; ++curRk)
+    {
+        MCAuto<DataArrayIdType> globalNodeIdsOfCurProc(DataArrayIdType::New());
+        globalNodeIdsOfCurProc->useArray(allGlobalNodeIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElems[curRk],1);
+        offset += nbOfElems[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()));
+        nbOfElems[curRk] = FromIdType<long>(localCellCapturedGlob->getNumberOfTuples());
+        tabs[curRk] = localCellCapturedGlob;
+    }
   }
-  //
-  MCAuto<DataArrayIdType> local;
+  std::vector<const DataArrayIdType *> tabss(tabs.begin(),tabs.end());
+  MCAuto<DataArrayIdType> cells(DataArrayIdType::Aggregate(tabss));
+  MPI_Alltoall(nbOfElems.get(),1,MPI_LONG,nbOfElems2.get(),1,MPI_LONG,comm);
+  long nbOfCellIdsSum(std::accumulate(nbOfElems2.get(),nbOfElems2.get()+size,0L));
+  MCAuto<DataArrayIdType> cellIdsFromProcs(DataArrayIdType::New());
+  cellIdsFromProcs->alloc(nbOfCellIdsSum,1);
+  {
+    std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray(nbOfElems,size) ),nbOfElemsOutInt( CommInterface::ToIntArray(nbOfElems2,size) );
+    std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) ), offsetsOut( CommInterface::ComputeOffset(CommInterface::ToIntArray(nbOfElems2,size),size) );
+    MPI_Alltoallv(cells->begin(),nbOfElemsInt.get(),offsetsIn.get(),MPI_INT,
+                  cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_INT,comm);
+  }
+  cellIdsFromProcs->sort();
+  return cellIdsFromProcs;
 }