]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
[Partial load] Allowing partial load of meshes with multiple geometrical types
authorAnida Khizar <anida.khizar@cea.fr>
Thu, 16 Mar 2023 09:59:41 +0000 (10:59 +0100)
committerAnida Khizar <anida.khizar@cea.fr>
Mon, 27 Mar 2023 12:17:11 +0000 (14:17 +0200)
src/MEDLoader/MEDFileMesh.cxx
src/MEDLoader/MEDFileMesh.hxx
src/MEDLoader/MEDFileMeshLL.cxx
src/MEDLoader/MEDFileMeshLL.hxx
src/ParaMEDLoader/ParaMEDFileMesh.cxx
src/ParaMEDLoader/ParaMEDFileMesh.hxx
src/ParaMEDMEMTest/ParaMEDMEMTest_MEDLoader.cxx

index d6c47707ff9ad935288fe2632f1d0a63d0732ab3..2da4f8c74ed46877c8091cd8d161235c61eba5e3 100644 (file)
@@ -2515,10 +2515,10 @@ MEDFileUMesh *MEDFileUMesh::LoadPartOf(med_idt fid, const std::string& mName, co
  * \param [in] mrs - the request for what to be loaded.
  * \return MEDFileUMesh * - a new instance of MEDFileUMesh. The caller is to delete this mesh using decrRef() as it is no more needed.
  */
-MEDFileUMesh *MEDFileUMesh::LoadPartOfFromUserDistrib(med_idt fid, const std::string& mName,  INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
+MEDFileUMesh *MEDFileUMesh::LoadPartOfFromUserDistrib(med_idt fid, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   MCAuto<MEDFileUMesh> ret(MEDFileUMesh::New());
-  ret->loadPartUMeshFromFileFromUserDistrib(fid,mName,type,distrib,dt,it,mrs);
+  ret->loadPartUMeshFromFileFromUserDistrib(fid,mName,distrib,dt,it,mrs);
   return ret.retn();
 }
 
@@ -2915,7 +2915,7 @@ void MEDFileUMesh::loadPartUMeshFromFile(med_idt fid, const std::string& mName,
  *
  * \sa loadLL
  */
-void MEDFileUMesh::loadPartUMeshFromFileFromUserDistrib(med_idt fid, const std::string& mName, INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
+void MEDFileUMesh::loadPartUMeshFromFileFromUserDistrib(med_idt fid, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   MEDFileUMeshL2 loaderl2;
   MEDCoupling::MEDCouplingMeshType meshType;
@@ -2928,7 +2928,7 @@ void MEDFileUMesh::loadPartUMeshFromFileFromUserDistrib(med_idt fid, const std::
       std::ostringstream oss; oss << "loadPartUMeshFromFileFromUserDistrib : Trying to load as unstructured an existing mesh with name '" << mName << "' !";
       throw INTERP_KERNEL::Exception(oss.str().c_str());
     }
-  loaderl2.loadPartFromUserDistrib(fid,mid,mName,type,distrib,dt,it,mrs);
+  loaderl2.loadPartFromUserDistrib(fid,mid,mName,distrib,dt,it,mrs);
   dispatchLoadedPart(fid,loaderl2,mName,mrs);
 }
 
index 9b8820f5ff9100d2e50ac52901ccf9c60dc33ef0..24b3e0459c68ec64efea08aa4c2bf19d847917e6 100644 (file)
@@ -271,7 +271,7 @@ namespace MEDCoupling
     MEDLOADER_EXPORT std::string getClassName() const override { return std::string("MEDFileUMesh"); }
     MEDLOADER_EXPORT static MEDFileUMesh *LoadPartOf(const std::string& fileName, const std::string& mName, const std::vector<INTERP_KERNEL::NormalizedCellType>& types, const std::vector<mcIdType>& slicPerTyp, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
     MEDLOADER_EXPORT static MEDFileUMesh *LoadPartOf(med_idt fid, const std::string& mName, const std::vector<INTERP_KERNEL::NormalizedCellType>& types, const std::vector<mcIdType>& slicPerTyp, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
-    MEDLOADER_EXPORT static MEDFileUMesh *LoadPartOfFromUserDistrib(med_idt fid, const std::string& mName, INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
+    MEDLOADER_EXPORT static MEDFileUMesh *LoadPartOfFromUserDistrib(med_idt fid, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
     MEDLOADER_EXPORT static void LoadPartCoords(const std::string& fileName, const std::string& mName, int dt, int it, const std::vector<std::string>& infosOnComp, mcIdType startNodeId, mcIdType stopNodeId,
 MCAuto<DataArrayDouble>& coords, MCAuto<PartDefinition>& partCoords, MCAuto<DataArrayIdType>& famCoords, MCAuto<DataArrayIdType>& numCoords, MCAuto<DataArrayAsciiChar>& nameCoords);
     MEDLOADER_EXPORT static const char *GetSpeStr4ExtMesh() { return SPE_FAM_STR_EXTRUDED_MESH; }
@@ -383,7 +383,7 @@ MCAuto<DataArrayDouble>& coords, MCAuto<PartDefinition>& partCoords, MCAuto<Data
     MEDFileUMesh();
     MEDFileUMesh(med_idt fid, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs);
     void loadPartUMeshFromFile(med_idt fid, const std::string& mName, const std::vector<INTERP_KERNEL::NormalizedCellType>& types, const std::vector<mcIdType>& slicPerTyp, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
-    void loadPartUMeshFromFileFromUserDistrib(med_idt fid, const std::string& mName, INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
+    void loadPartUMeshFromFileFromUserDistrib(med_idt fid, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
     void loadLL(med_idt fid, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs);
     void dispatchLoadedPart(med_idt fid, const MEDFileUMeshL2& loaderl2, const std::string& mName, MEDFileMeshReadSelector *mrs);
     const MEDFileUMeshSplitL1 *getMeshAtLevSafe(int meshDimRelToMaxExt) const;
index 3bc9116a83ada1929d0608edb76838b8b18cad34..b2ea8df34908ed6bf86bb081716f0185f6683fe9 100644 (file)
@@ -588,7 +588,7 @@ void MEDFileUMeshL2::loadPart(med_idt fid, const MeshOrStructMeshCls *mId, const
  * First, we load the connectivity of nodes. Second, we load coordinates of nodes lying in the specified cells.
  * \throw exception if multiple load sessions are requested
  */
-void MEDFileUMeshL2::loadPartFromUserDistrib(med_idt fid, const MeshOrStructMeshCls *mId, const std::string& mName, INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
+void MEDFileUMeshL2::loadPartFromUserDistrib(med_idt fid, const MeshOrStructMeshCls *mId, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   int Mdim;
   std::vector<std::string> infosOnComp(loadCommonPart(fid,mId,mName,dt,it,Mdim));
@@ -596,36 +596,32 @@ void MEDFileUMeshL2::loadPartFromUserDistrib(med_idt fid, const MeshOrStructMesh
     return ;
 
   /* First step : loading connectivity of nodes, ie building a new mesh of one geometrical type with only the specified cells in distrib */
-  loadPartOfConnectivityFromUserDistrib(fid,Mdim,mName,type,distrib,dt,it,mrs);
+  loadPartOfConnectivityFromUserDistrib(fid,Mdim,mName,distrib,dt,it,mrs);
 
   /* Second step : loading nodes */
   med_bool changement,transformation;
   mcIdType nCoords(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_NODE,MED_NONE,MED_COORDINATE,MED_NO_CMODE,&changement,&transformation));
   std::vector<bool> fetchedNodeIds(nCoords,false);
-  MEDCoupling1GTUMesh *mesh = _per_type_mesh[0][0/*only one type in our case*/]->getMesh();
-  // for each node in the original mesh, which ones are laying on the partial mesh defined above
-  mesh->computeNodeIdsAlg(fetchedNodeIds);
+  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()->computeNodeIdsAlg(fetchedNodeIds);    // for each node in the original mesh, which ones are laying on the current single geometrical type partial mesh
 
   if(!mrs || mrs->getNumberOfCoordsLoadSessions()==1)
     {
-      std::vector<mcIdType> distribNodes;
-      for(mcIdType i=0; i<nCoords; i++)
-        if(fetchedNodeIds[i])
-          distribNodes.push_back(i);
-
-      // map defining the renumbering of each node
-      mcIdType nbOfNodesToLoad = distribNodes.size();
-      MCAuto< MapKeyVal<mcIdType, mcIdType> > o2n(MapKeyVal<mcIdType, mcIdType>::New());
-      std::map<mcIdType, mcIdType>& mapO2N(o2n->data());
-      for(mcIdType node=0; node<nbOfNodesToLoad; node++)
-        mapO2N[distribNodes[node]] = node;
-
       // renumbering nodes inside the connectivity of the partial mesh:
       // until now, the numbering used in the connectivity of the cells was that of the integral mesh,
       // so it might be sparsed as some original nodes are missing in the partial mesh,
       // thus we want each node to be renumbered so that the sequence of their numbers form a range
-      mesh->renumberNodesInConn(mapO2N);
+      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());
+
       // 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);
     }
   else
@@ -668,12 +664,17 @@ void MEDFileUMeshL2::loadPartOfConnectivity(med_idt fid, int mdim, const std::st
  * This method builds a new mesh of single geometrical type based on the partition of cells \a distrib, from mesh \a mName in file \a fid.
  * This distribution is not necessarily a slice.
  */
-void MEDFileUMeshL2::loadPartOfConnectivityFromUserDistrib(med_idt fid, int mdim, const std::string& mName, INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
+void MEDFileUMeshL2::loadPartOfConnectivityFromUserDistrib(med_idt fid, int mdim, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   _per_type_mesh.resize(1);
   _per_type_mesh[0].clear();
-  MCAuto<MEDFileUMeshPerType> tmp(MEDFileUMeshPerType::NewPart(fid,mName.c_str(),dt,it,mdim,type,distrib,mrs));
-  _per_type_mesh[0].push_back(tmp);
+  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);
+    }
+  sortTypes();
 }
 
 void MEDFileUMeshL2::loadCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it)
index 9efa9f3b50c881f3a271bcde671c6694ed3d217a..95246ff2d657660a8648da210ac11bb4acc14189 100644 (file)
@@ -121,11 +121,11 @@ namespace MEDCoupling
     std::vector<std::string> loadCommonPart(med_idt fid, const MeshOrStructMeshCls *mId, const std::string& mName, int dt, int it, int& Mdim);
     void loadAll(med_idt fid, const MeshOrStructMeshCls *mId, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs);
     void loadPart(med_idt fid, const MeshOrStructMeshCls *mId, const std::string& mName, const std::vector<INTERP_KERNEL::NormalizedCellType>& types, const std::vector<mcIdType>& slicPerTyp, int dt, int it, MEDFileMeshReadSelector *mrs);
-    void loadPartFromUserDistrib(med_idt fid, const MeshOrStructMeshCls *mId, const std::string& mName, INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs);
+    void loadPartFromUserDistrib(med_idt fid, const MeshOrStructMeshCls *mId, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs);
 
     void loadConnectivity(med_idt fid, int mdim, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs);
     void loadPartOfConnectivity(med_idt fid, int mdim, const std::string& mName, const std::vector<INTERP_KERNEL::NormalizedCellType>& types, const std::vector<mcIdType>& slicPerTyp, int dt, int it, MEDFileMeshReadSelector *mrs);
-    void loadPartOfConnectivityFromUserDistrib(med_idt fid, int mdim, const std::string& mName,  INTERP_KERNEL::NormalizedCellType type, const std::vector<mcIdType>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs);
+    void loadPartOfConnectivityFromUserDistrib(med_idt fid, int mdim, const std::string& mName, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>& distrib, int dt, int it, MEDFileMeshReadSelector *mrs);
 
     void loadCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it);
     void loadPartCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, mcIdType nMin, mcIdType nMax);
index 7ff676952cc47adbc3d5a45559ca1352e27e1177..91972cccdeb7c13babd497bb5a29435d43785b38 100644 (file)
@@ -132,7 +132,7 @@ MEDFileUMesh *ParaMEDFileUMesh::New(int iPart, int nbOfParts, const std::string&
  * \throw exception if the mesh contains multiple types of cells
  * \throw exception if the partition of the mesh cells defined by \a distrib does not cover the whole mesh
  */
-MEDFileUMesh *ParaMEDFileUMesh::ParaNew(const std::vector<mcIdType>& distrib, const MPI_Comm& com, const MPI_Info& nfo, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
+MEDFileUMesh *ParaMEDFileUMesh::ParaNew(const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>> &distrib, const MPI_Comm& com, const MPI_Info& nfo, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   MEDFileUtilities::CheckFileForRead(fileName);
 #ifdef HDF5_IS_PARALLEL
@@ -176,16 +176,17 @@ MEDFileUMesh *ParaMEDFileUMesh::ParaNew(int iPart, int nbOfParts, const MPI_Comm
  * Loads mesh \a mName in parallel using a custom partition of the mesh cells among the processes.
  * See ParaMEDFileUMesh::ParaNew for detailed description.
  */
-MEDFileUMesh *ParaMEDFileUMesh::NewPrivate(med_idt fid, const MPI_Comm& com, const std::vector<mcIdType>& distrib, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
+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;
-  INTERP_KERNEL::NormalizedCellType geoType;
-  getSingleGeometricType(fileName,mName,geoType);
-  med_geometry_type geoMedType(typmai3[geoType]);
-  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,distrib);
-  ret=MEDFileUMesh::LoadPartOfFromUserDistrib(fid,mName,geoType,distrib,dt,it,mrs);
+  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();
 }
 
@@ -250,12 +251,12 @@ MEDFileField1TS *ParaMEDFileField1TS::NewPrivate(med_idt fid, const MPI_Comm& co
 {
   INTERP_KERNEL::NormalizedCellType geoType;
   getSingleGeometricType(MEDFileWritable::FileNameFromFID(fid),mName,geoType);
-  med_geometry_type geoMedType(typmai3[geoType]);
-  med_bool changement,transformation;
   if(loc==ON_CELLS) //if distribution is on nodes, no fast way to check it (as a node can be shared by multiple processors)
     {
-        med_int totalNumberOfElements(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_NODE,geoMedType,MED_FAMILY_NUMBER,MED_NODAL,&changement,&transformation));
-        checkDistribution(com,totalNumberOfElements,distrib);
+      med_geometry_type geoMedType(typmai3[geoType]);
+      med_bool changement,transformation;
+      med_int totalNumberOfElements(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_CELL,geoMedType,MED_FAMILY_NUMBER,MED_NODAL,&changement,&transformation));
+      checkDistribution(com,totalNumberOfElements,distrib);
     }
   std::vector<std::pair<TypeOfField,INTERP_KERNEL::NormalizedCellType>> tmp={ {loc, geoType} };
   INTERP_KERNEL::AutoCppPtr<MEDFileEntities> entities(MEDFileEntities::BuildFrom(&tmp));
index 8c3d83974a710460d1b303897e7b0e2fb572b2f8..18a3dc8616bf52eed23fa4b8bac0971250c02d3b 100644 (file)
 
 #include <string>
 #include <vector>
+#include <map>
 #include "MCIdType.hxx"
 #include "MEDCouplingRefCountObject.hxx"
+#include "NormalizedGeometricTypes"
 
 namespace MEDCoupling
 {
@@ -50,11 +52,11 @@ namespace MEDCoupling
   public:
     static MEDFileUMesh *New(int iPart, int nbOfParts, const std::string& fileName, const std::string& mName, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
     static MEDFileUMesh *ParaNew(int iPart, int nbOfParts, const MPI_Comm& com, const MPI_Info& nfo, const std::string& fileName, const std::string& mName, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
-    static MEDFileUMesh *ParaNew(const std::vector<mcIdType>& distrib, const MPI_Comm& com, const MPI_Info& nfo, const std::string& fileName, const std::string& mName, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
+    static MEDFileUMesh *ParaNew(const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>&, const MPI_Comm& com, const MPI_Info& nfo, const std::string& fileName, const std::string& mName, int dt=-1, int it=-1, MEDFileMeshReadSelector *mrs=0);
 
   private:
     static MEDFileUMesh *NewPrivate(med_idt fid, int iPart, int nbOfParts, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs);
-    static MEDFileUMesh *NewPrivate(med_idt fid, const MPI_Comm& com, const std::vector<mcIdType>& distrib, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs);
+    static MEDFileUMesh *NewPrivate(med_idt fid, const MPI_Comm& com, const std::map<INTERP_KERNEL::NormalizedCellType,std::vector<mcIdType>>&, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs);
   };
 
 
index bb779bdcb60edaf272908c23b73cbcd8e75eb7e5..2420a4fa0b59c539c718d30f3f9226809469c94a 100644 (file)
@@ -96,18 +96,17 @@ void ParaMEDMEMTest::testParallelLoad1()
   if(size!=2)
     return ;
 
-  std::vector<mcIdType> distrib;
+  std::map<INTERP_KERNEL::NormalizedCellType, std::vector<mcIdType>> distrib;
   if (rank == 0)
-      distrib = {0,1,4,5};    //c++ type of indexing: index starts from zero!
+      distrib = { {INTERP_KERNEL::NORM_QUAD4,{0,1,4,5}/*c++ type of indexing: index starts from zero!*/} };
   else
-      distrib = {2,3,6,7};
+      distrib = { {INTERP_KERNEL::NORM_QUAD4,{2,3,6,7}} };
 
   std::string filename=INTERP_TEST::getResourceFile("SimpleTest2D.med");
   MCAuto<MEDFileUMesh> mu = ParaMEDFileUMesh::ParaNew(distrib, MPI_COMM_WORLD, MPI_INFO_NULL, filename, "mesh");
   MCAuto<MEDCouplingUMesh> meshRef = genLocMesh(rank);
   CPPUNIT_ASSERT(mu->getMeshAtLevel(0)->isEqual(meshRef,1e-12));
   MPI_Barrier(MPI_COMM_WORLD);
-
 }