]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
[ParaMEDFileMesh] parallel partial loading is now enabled even if some procs have... akr/paramedfile
authorAnida Khizar <anida.khizar@cea.fr>
Wed, 9 Oct 2024 08:54:41 +0000 (10:54 +0200)
committerAnida Khizar <anida.khizar@cea.fr>
Fri, 11 Oct 2024 14:06:11 +0000 (16:06 +0200)
src/MEDLoader/MEDFileMeshElt.cxx
src/MEDLoader/MEDFileMeshElt.hxx
src/MEDLoader/MEDFileMeshLL.cxx
src/ParaMEDLoader/ParaMEDFileMesh.cxx
src/ParaMEDMEMTest/ParaMEDMEMTest.hxx
src/ParaMEDMEMTest/ParaMEDMEMTest_MEDLoader.cxx

index 7995b599c5867c305b9ee86028c3595d4b95d844..b3f59b7cd96edf0ad1e6af65fe4e0de8520518e4 100644 (file)
@@ -112,29 +112,20 @@ MEDFileUMeshPerType *MEDFileUMeshPerType::New(med_idt fid, const char *mName, in
 
 MEDFileUMeshPerType *MEDFileUMeshPerType::NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, mcIdType strt, mcIdType stp, mcIdType step, MEDFileMeshReadSelector *mrs)
 {
-  int geoElt2i((int)geoElt2);
-  if(geoElt2i<0 || geoElt2i>=INTERP_KERNEL::NORM_MAXTYPE)
-    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : Not recognized MEDCoupling/MEDLoader geometric type !");
-  med_geometry_type geoElt(typmai3[geoElt2]);
+  med_geometry_type geoElt;
   med_entity_type whichEntity;
-  if(!isExisting(fid,mName,dt,it,geoElt,whichEntity))
-    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : The specified geo type is not present in the specified mesh !");
+  getMedTypes(fid, mName, dt, it, geoElt2, geoElt, whichEntity);
+  
   MCAuto<MEDFileUMeshPerType> ret(new MEDFileUMeshPerType);
   ret->loadPart(fid,mName,dt,it,mdim,geoElt,geoElt2,whichEntity,strt,stp,step,mrs);
   return ret.retn();
 }
 
-MEDFileUMeshPerType *MEDFileUMeshPerType::NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs)
+
+MEDFileUMeshPerType *MEDFileUMeshPerType::NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, med_geometry_type geoElt, med_entity_type whichEntity, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs)
 {
-  int geoElt2i((int)geoElt2);
-  if(geoElt2i<0 || geoElt2i>=INTERP_KERNEL::NORM_MAXTYPE)
-    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : Not recognized MEDCoupling/MEDLoader geometric type !");
-  med_geometry_type geoElt(typmai3[geoElt2]);
-  med_entity_type whichEntity;
-  if(!isExisting(fid,mName,dt,it,geoElt,whichEntity))
-    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::NewPart : The specified geo type is not present in the specified mesh !");
-  MCAuto<MEDFileUMeshPerType> ret(new MEDFileUMeshPerType);
-  ret->loadPart(fid,mName,dt,it,mdim,geoElt,geoElt2,whichEntity,distrib,mrs);
+   MCAuto<MEDFileUMeshPerType> ret(new MEDFileUMeshPerType);
+  ret->loadPart(fid,mName,dt,it,mdim,geoElt,geoElt2,whichEntity,distrib,mrs);  
   return ret.retn();
 }
 
@@ -150,6 +141,16 @@ std::vector<const BigMemoryObject *> MEDFileUMeshPerType::getDirectChildrenWithN
   return ret;
 }
 
+void MEDFileUMeshPerType::getMedTypes(med_idt fid, const char *mName, int dt, int it, INTERP_KERNEL::NormalizedCellType MCgeoElt, med_geometry_type& geoElt, med_entity_type& whichEntity)
+{
+  int geoElt2i((int)MCgeoElt);
+  if(geoElt2i<0 || geoElt2i>=INTERP_KERNEL::NORM_MAXTYPE)
+    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::getMedTypes : Not recognized MEDCoupling/MEDLoader geometric type !");
+  geoElt = typmai3[MCgeoElt];
+  if(!isExisting(fid,mName,dt,it,geoElt,whichEntity))
+    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::getMedTypes : The specified geo type is not present in the specified mesh !");
+}
+
 bool MEDFileUMeshPerType::isExisting(med_idt fid, const char *mName, int dt, int it, med_geometry_type geoElt, med_entity_type& whichEntity)
 {
   static const med_entity_type entities[3]={MED_CELL,MED_DESCENDING_FACE,MED_DESCENDING_EDGE};
index 7975553b009ab09355149a6f5c998714cdc65d9f..ad5d46b71aa28d555c3a83365ab705c49cbc6ba0 100644 (file)
@@ -58,10 +58,11 @@ namespace MEDCoupling
   public:
     static MEDFileUMeshPerType *New(med_idt fid, const char *mName, int dt, int it, int mdim, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType geoElt2, MEDFileMeshReadSelector *mrs);
     static MEDFileUMeshPerType *NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, mcIdType strt, mcIdType stp, mcIdType step, MEDFileMeshReadSelector *mrs);
-    static MEDFileUMeshPerType *NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs);
+    static MEDFileUMeshPerType *NewPart(med_idt fid, const char *mName, int dt, int it, int mdim, INTERP_KERNEL::NormalizedCellType geoElt2, med_geometry_type geoElt, med_entity_type whichEntity, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs);
 
     std::string getClassName() const override { return std::string("MEDFileUMeshPerType"); }
     static bool isExisting(med_idt fid, const char *mName, int dt, int it, med_geometry_type geoElt, med_entity_type& whichEntity);
+    static void getMedTypes(med_idt fid, const char *mName, int dt, int it, INTERP_KERNEL::NormalizedCellType MCgeoElt, med_geometry_type& geoElt, med_entity_type& whichEntity);
     std::size_t getHeapMemorySizeWithoutChildren() const;
     std::vector<const BigMemoryObject *> getDirectChildrenWithNull() const;
     int getDim() const;
index 5baacd48a2b313fa1aa33be200ae8174379564c6..4d81121d78f605c0f15149e272fd9aa1add44584 100644 (file)
@@ -636,14 +636,16 @@ void MEDFileUMeshL2::loadPartFromUserDistrib(med_idt fid, const MeshOrStructMesh
       MCAuto<DataArrayIdType> fni(DataArrayIdType::BuildListOfSwitchedOn(fetchedNodeIds));
       MCAuto< MapKeyVal<mcIdType, mcIdType> > o2n(fni->invertArrayN2O2O2NOptimized());
       for(std::vector< std::vector< MCAuto<MEDFileUMeshPerType> > >::const_iterator it0=_per_type_mesh.begin();it0!=_per_type_mesh.end();it0++)
-        for(std::vector< MCAuto<MEDFileUMeshPerType> >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
-          (*it1)->getMesh()->renumberNodesInConn(o2n->data());
+          for(std::vector< MCAuto<MEDFileUMeshPerType> >::const_iterator it1=(*it0).begin();it1!=(*it0).end();it1++)
+            (*it1)->getMesh()->renumberNodesInConn(o2n->data());
+
 
       // loading coordinates of fetched nodes
       std::vector<mcIdType> distribNodes;
       for(std::map<mcIdType,mcIdType>::const_iterator mapIter = o2n->data().begin(); mapIter != o2n->data().end(); ++mapIter)
-          distribNodes.push_back(mapIter->first);
-      this->loadPartCoords(fid,infosOnComp,mName,dt,it,distribNodes);
+          distribNodes.push_back(mapIter->first);      
+      if(distribNodes.size() != 0)
+         this->loadPartCoords(fid,infosOnComp,mName,dt,it,distribNodes);
     }
   else
     throw INTERP_KERNEL::Exception("MEDFileUMeshL2::loadPartFromUserDistrib: multiple load sessions not handled!");
@@ -692,8 +694,15 @@ void MEDFileUMeshL2::loadPartOfConnectivityFromUserDistrib(med_idt fid, int mdim
   std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>::const_iterator iter;
   for (iter = distrib.begin(); iter != distrib.end(); iter++)
     {
-        MCAuto<MEDFileUMeshPerType> tmp(MEDFileUMeshPerType::NewPart(fid,mName.c_str(),dt,it,mdim,iter->first/*type*/,iter->second/*distrib over the current type*/,mrs));
-        _per_type_mesh[0].push_back(tmp);
+      INTERP_KERNEL::NormalizedCellType MCGeoType = iter->first;
+      med_geometry_type geoElt;
+      med_entity_type whichEntity;
+      MEDFileUMeshPerType::getMedTypes(fid, mName.c_str(), dt, it, MCGeoType, geoElt, whichEntity);
+      if(iter->second.size() !=0)
+       {
+         MCAuto<MEDFileUMeshPerType> tmp(MEDFileUMeshPerType::NewPart(fid,mName.c_str(),dt,it,mdim,MCGeoType,geoElt,whichEntity,iter->second/*distrib over the current type*/,mrs));
+         _per_type_mesh[0].push_back(tmp);
+       }         
     }
   sortTypes();
 }
index 0f442055c59784b7b49aaddefbeb53f3de0f9149..8b5cd40b27f0499fc38f340585ce77e77a3ee200 100644 (file)
@@ -35,25 +35,6 @@ extern med_geometry_type typmai3[INTERP_KERNEL::NORM_MAXTYPE];
 
 using namespace MEDCoupling;
 
-void checkDistribution(const MPI_Comm& com, mcIdType totalNumberOfElements, const std::vector<mcIdType>& distrib)
-{
-  mcIdType nbEltsInDistribLoc = distrib.size();
-  mcIdType nbEltsInDistribTot = -1;
-#ifdef HAVE_MPI
-  MPI_Allreduce(&nbEltsInDistribLoc, &nbEltsInDistribTot, 1, MPI_LONG, MPI_SUM, com);
-#else
-      throw INTERP_KERNEL::Exception("not(HAVE_MPI) incompatible with MPI_World_Size>1");
-#endif
-  if(nbEltsInDistribTot != totalNumberOfElements)
-    {
-      if(nbEltsInDistribTot > totalNumberOfElements)
-        throw INTERP_KERNEL::Exception("ParaMEDFileMesh : Some of your partitions overlap each other ! Each element in your distribution vector must appear only once ! ");
-      else
-        throw INTERP_KERNEL::Exception("ParaMEDFileMesh : The distribution does not cover the whole mesh ! Each element of the mesh must appear once in your distribution vector ");
-    }
-}
-
-
 MEDFileMesh *ParaMEDFileMesh::New(int iPart, int nbOfParts, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   MEDFileUtilities::CheckFileForRead(fileName);
@@ -169,13 +150,6 @@ MEDFileUMesh *ParaMEDFileUMesh::ParaNew(int iPart, int nbOfParts, const MPI_Comm
 MEDFileUMesh *ParaMEDFileUMesh::NewPrivate(med_idt fid, const MPI_Comm& com, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   MCAuto<MEDFileUMesh> ret;
-  for(std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>::const_iterator iter=distrib.begin(); iter!= distrib.end(); iter++)
-    {
-        med_geometry_type geoMedType(typmai3[iter->first /*current geometric type*/]);
-        med_bool changement,transformation;
-        med_int totalNumberOfElements(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_CELL,geoMedType,MED_CONNECTIVITY,MED_NODAL,&changement,&transformation));
-        checkDistribution(com,totalNumberOfElements,iter->second /*distrib over this geometric type*/);
-    }
   ret=MEDFileUMesh::LoadPartOfFromUserDistrib(fid,mName,distrib,dt,it,mrs);
   return ret.retn();
 }
index 9f12434ddb4ed682355b74a3502a13cae0181199..d6c23d0009eed993e1fe9f44c1a2fcfe0cfb0b3b 100644 (file)
@@ -92,6 +92,7 @@ class ParaMEDMEMTest : public CppUnit::TestFixture
   CPPUNIT_TEST(testParallelLoad3);   // 2 procs
   CPPUNIT_TEST(testParallelLoad4);   // 2 procs
   CPPUNIT_TEST(testParallelLoad5);   // 2 procs
+  CPPUNIT_TEST(testParallelLoad6);   // 3 procs
   CPPUNIT_TEST_SUITE_END();
 
 public:
@@ -163,6 +164,7 @@ public:
   void testParallelLoad3();
   void testParallelLoad4();
   void testParallelLoad5();
+  void testParallelLoad6();
 
   std::string getTmpDirectory();
   std::string makeTmpFile( const std::string&, const std::string& = "" );
index c992562d0e1a5943177e6bf797b73ce3fe6a4676..a6e90417c1b7b461ac6f2532f2d6d5da87d4a83a 100644 (file)
@@ -131,6 +131,67 @@ MEDCouplingUMesh* genLocMeshMultipleTypes3()
   return ret.retn();
 }
 
+/*
+ * Generate a 2D mesh that is supposed to match the part that will be loaded by proc0 in testParallelLoad6
+ */
+MEDCouplingUMesh* genPartialLocMeshMultipleTypes1()
+{
+  MCAuto<MEDCouplingUMesh> ret= MEDCouplingUMesh::New("mesh",2);
+  double coords[8] = {0.,2.,  1.,2.,  0.,3.,   1.,3.};
+  DataArrayDouble *myCoords=DataArrayDouble::New();
+  myCoords->alloc(4,2);
+  std::copy(coords,coords+8,myCoords->getPointer());
+  ret->setCoords(myCoords);
+  myCoords->decrRef();
+  mcIdType conn[4]={0,1,3,2};
+  ret->allocateCells(1);
+  ret->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn);
+  ret->finishInsertingCells();
+  return ret.retn();
+}
+
+/*
+ * Generate a 2D mesh that is supposed to match the part that will be loaded by proc1 in testParallelLoad6
+ */
+MEDCouplingUMesh* genPartialLocMeshMultipleTypes2()
+{
+  MCAuto<MEDCouplingUMesh> ret= MEDCouplingUMesh::New("mesh",2);
+  double coords[6] = {0.,1.,  1.,1.,  1.,2.};
+  DataArrayDouble *myCoords=DataArrayDouble::New();
+  myCoords->alloc(3,2);
+  std::copy(coords,coords+6,myCoords->getPointer());
+  ret->setCoords(myCoords);
+  myCoords->decrRef();
+  mcIdType conn[3]={0,1,2};
+  ret->allocateCells(1);
+  ret->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
+  ret->finishInsertingCells();
+  return ret.retn();
+}
+
+
+/*
+ * Generate a 2D mesh that is supposed to match the part that will be loaded by proc2 in testParallelLoad6
+ */
+MEDCouplingUMesh* genPartialLocMeshMultipleTypes3()
+{
+  MCAuto<MEDCouplingUMesh> ret= MEDCouplingUMesh::New("mesh",2);
+  double coords[12] = {1.,0.,  2.,0.,  1.,1.,  2.,1.,  1.,2.,  2.,2. };
+  DataArrayDouble *myCoords=DataArrayDouble::New();
+  myCoords->alloc(6,2);
+  std::copy(coords,coords+12,myCoords->getPointer());
+  ret->setCoords(myCoords);
+  myCoords->decrRef();
+  mcIdType conn[7]={0,1,3,  2,3,5,4};
+  ret->allocateCells(2);
+  ret->insertNextCell(INTERP_KERNEL::NORM_TRI3,3,conn);
+  ret->insertNextCell(INTERP_KERNEL::NORM_QUAD4,4,conn+3);
+  ret->finishInsertingCells();
+  return ret.retn();
+}
+
+
+
 /*
  * Generate a 2D field that is supposed to match the local field loaded by each proc in testParallelLoad4
  */
@@ -228,8 +289,10 @@ void ParaMEDMEMTest::testParallelLoad2()
    distrib= { {INTERP_KERNEL::NORM_TRI3,{0,1}} , {INTERP_KERNEL::NORM_QUAD4,{1,3}} };
 
   std::string filename=INTERP_TEST::getResourceFile("Test2DMultiGeoType.med");
+  //partial loading
   MCAuto<MEDFileUMesh> mu = ParaMEDFileUMesh::ParaNew(distrib, MPI_COMM_WORLD, MPI_INFO_NULL, filename, "mesh");
   MCAuto<MEDCouplingUMesh> mesh = mu->getMeshAtLevel(0);
+
   MEDCouplingUMesh *meshRef;
   if(rank==0)
       meshRef=genLocMeshMultipleTypes1();
@@ -242,8 +305,8 @@ void ParaMEDMEMTest::testParallelLoad2()
   int allEqual = -1;
   MPI_Allreduce(&equal, &allEqual, 1, MPI_INT,MPI_SUM,MPI_COMM_WORLD);
   CPPUNIT_ASSERT(allEqual==3);
-
   meshRef->decrRef();
+
   MPI_Barrier(MPI_COMM_WORLD);
 }
 
@@ -381,4 +444,48 @@ void ParaMEDMEMTest::testParallelLoad5()
 }
 
 
+/*!
+ * Test case to load a 2D mesh with multiple geometric types in parallel on 3 procs.
+ * Some procs may have empty partition to load for some types
+ */
+void ParaMEDMEMTest::testParallelLoad6()
+{
+         int size;
+         int rank;
+         MPI_Comm_size(MPI_COMM_WORLD,&size);
+         MPI_Comm_rank(MPI_COMM_WORLD,&rank);
+         //
+         if(size!=3)
+           return ;
+
+         std::map<INTERP_KERNEL::NormalizedCellType, std::vector<mcIdType>> distrib;
+         // independant numerotation for each geometric type!
+         if (rank == 0)
+          distrib = { {INTERP_KERNEL::NORM_TRI3,{}} , {INTERP_KERNEL::NORM_QUAD4,{2}} };
+        else if(rank == 1)
+          distrib = { {INTERP_KERNEL::NORM_TRI3,{2}}  , {INTERP_KERNEL::NORM_QUAD4,{}} };
+        else
+          distrib= { {INTERP_KERNEL::NORM_TRI3,{0}} , {INTERP_KERNEL::NORM_QUAD4,{1}} };
+
+         std::string filename=INTERP_TEST::getResourceFile("Test2DMultiGeoType.med");
+         MCAuto<MEDFileUMesh> mu = ParaMEDFileUMesh::ParaNew(distrib, MPI_COMM_WORLD, MPI_INFO_NULL, filename, "mesh");
+         MCAuto<MEDCouplingUMesh> mesh = mu->getMeshAtLevel(0);
+
+         MEDCouplingUMesh *meshRef;
+         if(rank==0)
+             meshRef=genPartialLocMeshMultipleTypes1();
+         else if(rank==1)
+             meshRef=genPartialLocMeshMultipleTypes2();
+         else
+             meshRef=genPartialLocMeshMultipleTypes3();
+         //checking that all 3 procs have correctly loaded their part
+         int equal = (int)mesh->isEqual(meshRef,1e-12);
+         int allEqual = -1;
+         MPI_Allreduce(&equal, &allEqual, 1, MPI_INT,MPI_SUM,MPI_COMM_WORLD);
+         CPPUNIT_ASSERT(allEqual==3);
+         meshRef->decrRef();
+
+         MPI_Barrier(MPI_COMM_WORLD);
+}
+