]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
Some cleanup agy/para_mesh_opt_2
authorAnthony Geay <anthony.geay@edf.fr>
Sun, 26 Apr 2020 22:28:57 +0000 (00:28 +0200)
committerAnthony Geay <anthony.geay@edf.fr>
Sun, 26 Apr 2020 22:28:57 +0000 (00:28 +0200)
src/ParaMEDMEM/CommInterface.cxx
src/ParaMEDMEM/CommInterface.hxx
src/ParaMEDMEM/ParaUMesh.cxx

index 725b76fb1c10101c11051a9096152de806bcb834..6fda70ddd5fd54a9371daca96e436c68ae7a8823 100644 (file)
@@ -105,6 +105,13 @@ namespace MEDCoupling
     this->allToAllArraysT<double>(comm,arrays,arraysOut);
   }
 
+  void CommInterface::allToAllArrays(MPI_Comm comm, const std::vector< MCAuto<DataArrayDouble> >& arrays, MCAuto<DataArrayDouble>& arraysOut) const
+  {
+    std::unique_ptr<mcIdType[]> notUsed1;
+    mcIdType notUsed2;
+    this->allToAllArraysT2<double>(comm,arrays,arraysOut,notUsed1,notUsed2);
+  }
+
   /*!
   * Generalized AllToAll collective communication.
   */
@@ -112,4 +119,5 @@ namespace MEDCoupling
   {
     this->allToAllArraysT<mcIdType>(comm,arrays,arraysOut);
   }
+
 }
index d795185cce53faba26442249a124f4ea148761dc..e969f93940947d35b5faffcee367bb671a2ba13d 100644 (file)
@@ -122,8 +122,9 @@ namespace MEDCoupling
     int allGatherArrays(MPI_Comm comm, const DataArrayIdType *array, std::unique_ptr<mcIdType[]>& result, std::unique_ptr<mcIdType[]>& resultIndex) const;
     void allToAllArrays(MPI_Comm comm, const std::vector< MCAuto<DataArrayIdType> >& arrays, std::vector< MCAuto<DataArrayIdType> >& arraysOut) const;
     void allToAllArrays(MPI_Comm comm, const std::vector< MCAuto<DataArrayDouble> >& arrays, std::vector< MCAuto<DataArrayDouble> >& arraysOut) const;
+    void allToAllArrays(MPI_Comm comm, const std::vector< MCAuto<DataArrayDouble> >& arrays, MCAuto<DataArrayDouble>& arraysOut) const;
     template<class T>
-    void allToAllArraysT(MPI_Comm comm, const std::vector< MCAuto<typename Traits<T>::ArrayType> >& arrays, std::vector< MCAuto<typename Traits<T>::ArrayType> >& arraysOut) const
+    int allToAllArraysT2(MPI_Comm comm, const std::vector< MCAuto<typename Traits<T>::ArrayType> >& arrays, MCAuto<typename Traits<T>::ArrayType>& arrayOut, std::unique_ptr<mcIdType[]>& nbOfElems2, mcIdType& nbOfComponents) const
     {
       using DataArrayT = typename Traits<T>::ArrayType;
       int size;
@@ -132,8 +133,9 @@ namespace MEDCoupling
         throw INTERP_KERNEL::Exception("AllToAllArrays : internal error ! Invalid size of input array.");
         
       std::vector< const DataArrayT *> arraysBis(FromVecAutoToVecOfConst<DataArrayT>(arrays));
-      std::unique_ptr<mcIdType[]> nbOfElems2(new mcIdType[size]),nbOfElems3(new mcIdType[size]);
-      mcIdType nbOfComponents(std::numeric_limits<mcIdType>::max());
+      std::unique_ptr<mcIdType[]> nbOfElems3(new mcIdType[size]);
+      nbOfElems2.reset(new mcIdType[size]);
+      nbOfComponents = std::numeric_limits<mcIdType>::max();
       for(int curRk = 0 ; curRk < size ; ++curRk)
       {
         mcIdType curNbOfCompo( ToIdType( arrays[curRk]->getNumberOfComponents() ) );
@@ -150,15 +152,26 @@ namespace MEDCoupling
       }
       this->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<DataArrayT> cellIdsFromProcs(DataArrayT::New());
-      cellIdsFromProcs->alloc(nbOfCellIdsSum/nbOfComponents,nbOfComponents);
+      arrayOut = DataArrayT::New();
+      arrayOut->alloc(nbOfCellIdsSum/nbOfComponents,nbOfComponents);
       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) );
       {
         MCAuto<DataArrayT> arraysAcc(DataArrayT::Aggregate(arraysBis));
         this->allToAllV(arraysAcc->begin(),nbOfElemsInt.get(),offsetsIn.get(),ParaTraits<T>::MPIDataType,
-                        cellIdsFromProcs->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),ParaTraits<T>::MPIDataType,comm);
+                        arrayOut->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),ParaTraits<T>::MPIDataType,comm);
       }
+      return size;
+    }
+
+    template<class T>
+    void allToAllArraysT(MPI_Comm comm, const std::vector< MCAuto<typename Traits<T>::ArrayType> >& arrays, std::vector< MCAuto<typename Traits<T>::ArrayType> >& arraysOut) const
+    {
+      using DataArrayT = typename Traits<T>::ArrayType;
+      MCAuto<DataArrayT> cellIdsFromProcs;
+      std::unique_ptr<mcIdType[]> nbOfElems2;
+      mcIdType nbOfComponents(0);
+      int size(this->allToAllArraysT2<T>(comm,arrays,cellIdsFromProcs,nbOfElems2,nbOfComponents));
       std::unique_ptr<mcIdType[]> offsetsOutIdType( CommInterface::ComputeOffset(nbOfElems2,size) );
       // build output arraysOut by spliting cellIdsFromProcs into parts
       arraysOut.resize(size);
index 0ceb35565b7cb5b3b4625f98ea52c2e66f4bce99..60d44a3231f14e5bc286f32b4cf9e43ba752c2cc 100644 (file)
@@ -132,30 +132,22 @@ MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodes(const DataArrayIdType
 }
 
 /*!
+ * Return part of \a this mesh split over COMM_WORLD. Part is defined by global cell ids array \a globaCellIds.
  */
 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]);
-  mcIdType nbOfCellsRequested(globalCellIds->getNumberOfTuples());
-  ci.allGather(&nbOfCellsRequested,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(),nbOfCellsRequested,MPI_ID_TYPE,allGlobalCellIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm);
-  mcIdType offset(0);
+  std::unique_ptr<mcIdType[]> allGlobalCellIds,allGlobalCellIdsIndex;
+  int size(ci.allGatherArrays(comm,globalCellIds,allGlobalCellIds,allGlobalCellIdsIndex));
   // Prepare ParaUMesh parts to be sent : compute for each proc the contribution of current rank.
   std::vector< MCAuto<DataArrayIdType> > globalCellIdsToBeSent(size),globalNodeIdsToBeSent(size);
   std::vector< MCAuto<MEDCouplingUMesh> > meshPartsToBeSent(size);
   for(int curRk = 0 ; curRk < size ; ++curRk)
   {
+    mcIdType offset(allGlobalCellIdsIndex[curRk]);
     MCAuto<DataArrayIdType> globalCellIdsOfCurProc(DataArrayIdType::New());
-    globalCellIdsOfCurProc->useArray(allGlobalCellIds.get()+offset,false,DeallocType::CPP_DEALLOC,nbOfElems[curRk],1);
-    offset += nbOfElems[curRk];
+    globalCellIdsOfCurProc->useArray(allGlobalCellIds.get()+offset,false,DeallocType::CPP_DEALLOC,allGlobalCellIdsIndex[curRk+1]-offset,1);
     // the key call is here : compute for rank curRk the cells to be sent
     MCAuto<DataArrayIdType> globalCellIdsCaptured(_cell_global->buildIntersection(globalCellIdsOfCurProc));// OK for the global cellIds
     MCAuto<DataArrayIdType> localCellIdsCaptured(_cell_global->findIdForEach(globalCellIdsCaptured->begin(),globalCellIdsCaptured->end()));
@@ -185,14 +177,13 @@ ParaUMesh *ParaUMesh::redistributeCells(const DataArrayIdType *globalCellIds) co
     ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayIdType>(connectivityToBeSent),connectivityReceived);
   }
   //coordinates
-  std::vector< MCAuto<DataArrayDouble> > coordsReceived;
+  MCAuto<DataArrayDouble> coords;
   {
     std::vector<const DataArrayDouble *> coordsToBeSent(UMeshCoordsIterator(0,&meshPartsToBeSent2),UMeshCoordsIterator(meshPartsToBeSent2.size(),&meshPartsToBeSent2));
-    ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayDouble>(coordsToBeSent),coordsReceived);
+    ci.allToAllArrays(comm,FromVecConstToVecAuto<DataArrayDouble>(coordsToBeSent),coords);
   }
   /////// Sort it all !
   // firstly deal with nodes.
-  MCAuto<DataArrayDouble> coords( DataArrayDouble::Aggregate(FromVecAutoToVecOfConst<DataArrayDouble>(coordsReceived)) );
   MCAuto<DataArrayIdType> aggregatedNodeIds( DataArrayIdType::Aggregate(FromVecAutoToVecOfConst<DataArrayIdType>(globalNodeIdsReceived)) );
   MCAuto<DataArrayIdType> aggregatedNodeIdsSorted(aggregatedNodeIds->copySorted());
   MCAuto<DataArrayIdType> nodeIdsIntoAggregatedIds(DataArrayIdType::FindPermutationFromFirstToSecondDuplicate(aggregatedNodeIdsSorted,aggregatedNodeIds));