]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
With memory optimization
authorAnthony Geay <anthony.geay@edf.fr>
Tue, 10 Mar 2020 21:36:24 +0000 (22:36 +0100)
committerAnthony Geay <anthony.geay@edf.fr>
Tue, 10 Mar 2020 21:36:24 +0000 (22:36 +0100)
src/ParaMEDMEM/ParaUMesh.cxx

index fed10a5f87cc3fae3c80c822e3069dc9282f463f..ad6803606f2c3748203a87168ee9f2ef04aca501 100644 (file)
@@ -65,41 +65,53 @@ MCAuto<DataArrayIdType> ParaUMesh::getCellIdsLyingOnNodes(DataArrayIdType *globa
   CommInterface ci;
   int size;
   ci.commSize(comm,&size);
-  std::unique_ptr<mcIdType[]> nbOfElems(new mcIdType[size]),nbOfElems2(new mcIdType[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);
-  mcIdType nbOfNodeIdsSum(std::accumulate(nbOfElems.get(),nbOfElems.get()+size,0));
   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>(nbOfElems,size) );
+    std::unique_ptr<int[]> nbOfElemsInt( CommInterface::ToIntArray<mcIdType>(nbOfElemsSp,size) );
     std::unique_ptr<int[]> offsetsIn( CommInterface::ComputeOffset(nbOfElemsInt,size) );
-    ci.allGatherV(globalNodeIds->begin(),globalNodeIds->getNumberOfTuples(),MPI_ID_TYPE,allGlobalNodeIds.get(),nbOfElemsInt.get(),offsetsIn.get(),MPI_ID_TYPE,comm);
+    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,nbOfElems[curRk],1);
-      offset += nbOfElems[curRk];
+      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()));
-      nbOfElems[curRk] = localCellCapturedGlob->getNumberOfTuples();
-      tabs[curRk] = localCellCapturedGlob;
+      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(nbOfElems.get(),1,MPI_ID_TYPE,nbOfElems2.get(),1,MPI_ID_TYPE,comm);
+  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>(nbOfElems,size) ),nbOfElemsOutInt( CommInterface::ToIntArray<mcIdType>(nbOfElems2,size) );
+    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->getPointer(),nbOfElemsOutInt.get(),offsetsOut.get(),MPI_ID_TYPE,comm);
   }
   cellIdsFromProcs->sort();
   return cellIdsFromProcs;