]> SALOME platform Git repositories - tools/medcoupling.git/commitdiff
Salome HOME
[Partial load] First implementation of a partial parallel load (of mesh and field...
authorAnida Khizar <anida.khizar@cea.fr>
Wed, 11 Jan 2023 14:48:19 +0000 (15:48 +0100)
committerAnida Khizar <anida.khizar@cea.fr>
Mon, 30 Jan 2023 15:48:05 +0000 (16:48 +0100)
12 files changed:
src/MEDLoader/MEDFileField1TS.cxx
src/MEDLoader/MEDFileField1TS.hxx
src/MEDLoader/MEDFileFieldInternal.cxx
src/MEDLoader/MEDFileMesh.cxx
src/MEDLoader/MEDFileMesh.hxx
src/MEDLoader/MEDFileMeshElt.cxx
src/MEDLoader/MEDFileMeshElt.hxx
src/MEDLoader/MEDFileMeshLL.cxx
src/MEDLoader/MEDFileMeshLL.hxx
src/MEDLoader/MEDFilterEntity.hxx [new file with mode: 0644]
src/ParaMEDLoader/ParaMEDFileMesh.cxx
src/ParaMEDLoader/ParaMEDFileMesh.hxx

index bd3db05258066545f31384b7a9c7a8c386968020..1c390258296021832bb80fce28bb234a770ec39b 100644 (file)
@@ -1976,6 +1976,7 @@ MEDFileAnyTypeField1TS *MEDFileAnyTypeField1TS::New(med_idt fid, bool loadAll)
   return ret.retn();
 }
 
+
 MEDFileAnyTypeField1TS *MEDFileAnyTypeField1TS::New(const std::string& fileName, const std::string& fieldName, bool loadAll)
 {
   MEDFileUtilities::AutoFid fid(OpenMEDFileForRead(fileName));
@@ -2018,6 +2019,20 @@ MEDFileAnyTypeField1TS *MEDFileAnyTypeField1TS::NewAdv(med_idt fid, const std::s
   return ret.retn();
 }
 
+MEDFileAnyTypeField1TS *MEDFileAnyTypeField1TS::NewAdv(const std::string& fileName, const std::string& fieldName, int iteration, int order, const MEDFileMeshes *ms)
+{
+  MEDFileUtilities::AutoFid fid(OpenMEDFileForRead(fileName));
+  return NewAdv(fid,fieldName,iteration,order,ms);
+}
+
+MEDFileAnyTypeField1TS *MEDFileAnyTypeField1TS::NewAdv(med_idt fid, const std::string& fieldName, int iteration, int order, const MEDFileMeshes *ms)
+{
+  MCAuto<MEDFileAnyTypeField1TSWithoutSDA> c(BuildContentFrom(fid,fieldName,iteration,order,true,ms,0));
+  MCAuto<MEDFileAnyTypeField1TS> ret(BuildNewInstanceFromContent(c,fid));
+  ret->loadGlobals(fid);
+  return ret.retn();
+}
+
 MEDFileAnyTypeField1TSWithoutSDA *MEDFileAnyTypeField1TS::BuildContentFrom(med_idt fid, const std::string& fieldName, int iteration, int order, bool loadAll, const MEDFileMeshes *ms, const MEDFileEntities *entities)
 {
   med_field_type typcha;
index 55864e3f702ceca2695956c5f5fa79db22814fea..752e85658265441e8236794ab97a8df794319eb2 100644 (file)
@@ -289,6 +289,9 @@ namespace MEDCoupling
     MEDLOADER_EXPORT static MEDFileAnyTypeField1TS *New(med_idt fid, const std::string& fieldName, int iteration, int order, bool loadAll=true);
     MEDLOADER_EXPORT static MEDFileAnyTypeField1TS *NewAdv(const std::string& fileName, const std::string& fieldName, int iteration, int order, bool loadAll, const MEDFileEntities *entities);
     MEDLOADER_EXPORT static MEDFileAnyTypeField1TS *NewAdv(med_idt fid, const std::string& fieldName, int iteration, int order, bool loadAll, const MEDFileEntities *entities);
+    MEDLOADER_EXPORT static MEDFileAnyTypeField1TS *NewAdv(const std::string& fileName, const std::string& fieldName, int iteration, int order, const MEDFileMeshes *ms);
+    MEDLOADER_EXPORT static MEDFileAnyTypeField1TS *NewAdv(med_idt fid, const std::string& fieldName, int iteration, int order, const MEDFileMeshes *ms);
+
     MEDLOADER_EXPORT int getDimension() const;
     MEDLOADER_EXPORT int getIteration() const;
     MEDLOADER_EXPORT int getOrder() const;
@@ -374,6 +377,7 @@ namespace MEDCoupling
     MEDLOADER_EXPORT static typename MLFieldTraits<T>::F1TSType *New(const std::string& fileName, const std::string& fieldName, int iteration, int order, bool loadAll=true);
     MEDLOADER_EXPORT static typename MLFieldTraits<T>::F1TSType *New(med_idt fid, const std::string& fieldName, int iteration, int order, bool loadAll=true);
     MEDLOADER_EXPORT static typename MLFieldTraits<T>::F1TSType *New(const typename MLFieldTraits<T>::F1TSWSDAType& other, bool shallowCopyOfContent);
+
   public:
     MEDLOADER_EXPORT static typename Traits<T>::ArrayType *ReturnSafelyTypedDataArray(MCAuto<DataArray>& arr);
     MEDLOADER_EXPORT typename Traits<T>::ArrayType *getFieldWithProfile(TypeOfField type, int meshDimRelToMax, const MEDFileMesh *mesh, DataArrayIdType *&pfl) const;
index 640b423938f7f6580e677d69f3d8e4ae9d3448f0..80e6f83a1c07e7b6dd7ac63a1dc5c4cc670ea563 100644 (file)
@@ -30,6 +30,8 @@
 #include "MEDCouplingFieldTemplate.hxx"
 #include "MEDCouplingFieldDouble.hxx"
 
+#include "MEDFilterEntity.hxx"
+
 #include "CellModel.hxx"
 
 // From MEDLOader.cxx TU
@@ -591,47 +593,17 @@ void MEDFileFieldPerMeshPerTypePerDisc::goReadZeValuesInFile(med_idt fid, const
       INTERP_KERNEL::AutoPtr<char> pflname(MEDLoaderBase::buildEmptyString(MED_NAME_SIZE)),locname(MEDLoaderBase::buildEmptyString(MED_NAME_SIZE));
       med_int profilesize,nbi;
       med_int overallNval(MEDfieldnValueWithProfile(fid,fieldName.c_str(),iteration,order,menti,mgeoti,FromIdType<int>(_profile_it+1),MED_COMPACT_PFLMODE,pflname,&profilesize,locname,&nbi));
-      const SlicePartDefinition *spd(dynamic_cast<const SlicePartDefinition *>(pd));
-      if(spd)
-        {
-          mcIdType start,stop,step;
-          spd->getSlice(start,stop,step);
-          mcIdType nbOfEltsToLoad(DataArray::GetNumberOfItemGivenBES(start,stop,step,"MEDFileFieldPerMeshPerTypePerDisc::goReadZeValuesInFile"));
-          med_filter filter=MED_FILTER_INIT;
-          MEDFILESAFECALLERRD0(MEDfilterBlockOfEntityCr,(fid,/*nentity*/overallNval,/*nvaluesperentity*/nbi,/*nconstituentpervalue*/nbOfCompo,
-                                                         MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
-                                                         /*start*/ToMedInt(start+1),/*stride*/ToMedInt(step),/*count*/1,/*blocksize*/ToMedInt(nbOfEltsToLoad),
-                                                         /*lastblocksize=useless because count=1*/0,&filter));
-          MEDFILESAFECALLERRD0(MEDfieldValueAdvancedRd,(fid,fieldName.c_str(),iteration,order,menti,mgeoti,&filter,startFeedingPtr));
-          MEDfilterClose(&filter);
-          return ;
-        }
-      const DataArrayPartDefinition *dpd(dynamic_cast<const DataArrayPartDefinition *>(pd));
-      if(dpd)
-        {
-          dpd->checkConsistencyLight();
-          MCAuto<DataArrayIdType> myIds(dpd->toDAI());
-          mcIdType a(myIds->getMinValueInArray()),b(myIds->getMaxValueInArray());
-          myIds=myIds->deepCopy();// WARNING deep copy here because _pd is modified by applyLin !!!
-          myIds->applyLin(1,-a);
-          mcIdType nbOfEltsToLoad(b-a+1);
-          med_filter filter=MED_FILTER_INIT;
-          {//TODO : manage int32 !
-            MCAuto<DataArrayDouble> tmp(DataArrayDouble::New());
-            tmp->alloc(nbOfEltsToLoad,nbOfCompo);
-            MEDFILESAFECALLERRD0(MEDfilterBlockOfEntityCr,(fid,/*nentity*/overallNval,/*nvaluesperentity*/nbi,/*nconstituentpervalue*/nbOfCompo,
-                                                           MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
-                                                           /*start*/ToMedInt(a+1),/*stride*/1,/*count*/1,/*blocksize*/ToMedInt(nbOfEltsToLoad),
-                                                           /*lastblocksize=useless because count=1*/0,&filter));
-            MEDFILESAFECALLERRD0(MEDfieldValueAdvancedRd,(fid,fieldName.c_str(),iteration,order,menti,mgeoti,&filter,reinterpret_cast<unsigned char *>(tmp->getPointer())));
-            MCAuto<DataArrayDouble> feeder(DataArrayDouble::New());
-            feeder->useExternalArrayWithRWAccess(reinterpret_cast<double *>(startFeedingPtr),_nval,nbOfCompo);
-            feeder->setContigPartOfSelectedValues(0,tmp,myIds);
-          }
-          MEDfilterClose(&filter);
-        }
-      else
-        throw INTERP_KERNEL::Exception("Not implemented yet for not slices!");
+
+      {//TODO : manage int32 !
+        pd->checkConsistencyLight();
+        MEDFilterEntity filter;
+        filter.init();
+        filter.fill(fid,/*nentity*/overallNval,/*nvaluesperentity*/nbi,/*nconstituentpervalue*/nbOfCompo,
+                  MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
+                  pd);
+        MEDFILESAFECALLERRD0(MEDfieldValueAdvancedRd,(fid,fieldName.c_str(),iteration,order,menti,mgeoti,filter.getPtr(),startFeedingPtr));
+        filter.close();
+      }
     }
 }
 
@@ -748,6 +720,7 @@ void MEDFileFieldPerMeshPerTypePerDisc::loadBigArray(med_idt fid, const MEDFileF
   throw INTERP_KERNEL::Exception("Error on array reading ! Unrecognized type of field ! Should be in FLOAT64 FLOAT32 INT32 or INT64 !");
 }
 
+
 /*!
  * Set a \c this->_start **and** \c this->_end keeping the same delta between the two.
  */
@@ -3162,7 +3135,10 @@ MEDFileFieldPerMesh::MEDFileFieldPerMesh(med_idt fid, MEDFileAnyTypeField1TSWith
         {
           const PartDefinition *pd(0);
           if(mmu)
-            pd=mmu->getPartDefAtLevel(mmu->getRelativeLevOnGeoType(typmai2[iter0->current()]),typmai2[iter0->current()]);
+            {
+              pd=mmu->getPartDefAtLevel(mmu->getRelativeLevOnGeoType(typmai2[iter0->current()]),typmai2[iter0->current()]);
+            }
+
           _field_pm_pt.push_back(MEDFileFieldPerMeshPerType::NewOnRead(fid,this,ON_CELLS,typmai2[iter0->current()],nasc,pd));
           if(nbProfile>0)
             setMeshName(name0);
index 8bed8693955aff1625fc29849116cf1f26c2c437..ba45f1c30c317e89ed86bf89c3bacc3af015081a 100644 (file)
@@ -2502,6 +2502,17 @@ MEDFileUMesh *MEDFileUMesh::LoadPartOf(med_idt fid, const std::string& mName, co
   return ret.retn();
 }
 
+/*!
+ * This method loads from file with name \a fileName the mesh called \a mName as New does. The difference is that
+ * here only a part of cells contained in the file (those in the \a distrib vector) will be loaded. Those selected cells must have the same geometrical type.
+ */
+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)
+{
+  MCAuto<MEDFileUMesh> ret(MEDFileUMesh::New());
+  ret->loadPartUMeshFromFileFromUserDistrib(fid,mName,type,distrib,dt,it,mrs);
+  return ret.retn();
+}
+
 /*!
  * This method is an helper to load only consecutive nodes chunk of data of MED file pointed by \a fileName.
  * Consecutive chunk is specified classicaly by start (included) stop (excluded) format with \a startNodeId and \a stopNodeId respectively.
@@ -2889,6 +2900,30 @@ void MEDFileUMesh::loadPartUMeshFromFile(med_idt fid, const std::string& mName,
   dispatchLoadedPart(fid,loaderl2,mName,mrs);
 }
 
+/*!
+ * This method loads only a part of specified cells in the \a distrib vector
+ * See MEDFileUMesh::LoadPartOf for detailed description.
+ *
+ * \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)
+{
+  MEDFileUMeshL2 loaderl2;
+  MEDCoupling::MEDCouplingMeshType meshType;
+  int dummy0,dummy1;
+  std::string dummy2;
+  MEDCoupling::MEDCouplingAxisType dummy3;
+  INTERP_KERNEL::AutoCppPtr<MeshOrStructMeshCls> mid(MEDFileUMeshL2::GetMeshIdFromName(fid,mName,meshType,dummy3,dummy0,dummy1,dummy2));
+  if(meshType!=UNSTRUCTURED)
+    {
+      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);
+  dispatchLoadedPart(fid,loaderl2,mName,mrs);
+}
+
+
 /*!
  * \brief Write joints in a file
  */
index 171f6aab2024e1630059d6312d47c181e3e14fb9..9b8820f5ff9100d2e50ac52901ccf9c60dc33ef0 100644 (file)
@@ -271,6 +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 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; }
@@ -382,6 +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 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;
@@ -392,6 +394,7 @@ MCAuto<DataArrayDouble>& coords, MCAuto<PartDefinition>& partCoords, MCAuto<Data
     void changeFamilyIdArr(mcIdType oldId, mcIdType newId);
     std::list< MCAuto<DataArrayIdType> > getAllNonNullFamilyIds() const;
     MCAuto<MEDFileUMeshSplitL1>& checkAndGiveEntryInSplitL1(int meshDimRelToMax, MEDCouplingPointSet *m);
+
   private:
     static const char SPE_FAM_STR_EXTRUDED_MESH[];
   private:
index c61e94d0a2fb608e6a30f39b893a18e612126f17..313b859fa943f0cc6562735ea4bf8753f59e0f77 100644 (file)
 #include "InterpKernelAutoPtr.hxx"
 #include "CellModel.hxx"
 
+#include "MEDFilterEntity.hxx"
+
 #include <iostream>
 
+
 // From MEDLOader.cxx TU
 extern med_geometry_type typmai3[INTERP_KERNEL::NORM_MAXTYPE];
 
@@ -120,6 +123,20 @@ MEDFileUMeshPerType *MEDFileUMeshPerType::NewPart(med_idt fid, const char *mName
   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)
+{
+  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);
+  return ret.retn();
+}
+
 std::size_t MEDFileUMeshPerType::getHeapMemorySizeWithoutChildren() const
 {
   return MEDFileUMeshPerTypeCommon::getHeapMemorySizeWithoutChildren()+0;
@@ -180,6 +197,23 @@ MEDFileUMeshPerType::MEDFileUMeshPerType(med_idt fid, const char *mName, int dt,
   loadPolyh(fid,mName,dt,it,mdim,curNbOfElem,geoElt,entity,mrs);
 }
 
+void MEDFileUMeshPerType::loadPart(med_idt fid, const char *mName, int dt, int it, int mdim, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
+                                   med_entity_type entity, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs)
+{
+  med_bool changement,transformation;
+  mcIdType curNbOfElem(MEDmeshnEntity(fid,mName,dt,it,entity,geoElt,MED_CONNECTIVITY,MED_NODAL,&changement,&transformation));
+  const INTERP_KERNEL::CellModel& cm(INTERP_KERNEL::CellModel::GetCellModel(type));
+  MCAuto<DataArrayIdType> listOfIds=DataArrayIdType::New();
+  listOfIds->useArray(distrib.data(),false,DeallocType::C_DEALLOC,distrib.size(),1);
+  _pd=PartDefinition::New(listOfIds);
+  if(!cm.isDynamic())
+    {
+      loadPartStaticType(fid,mName,dt,it,mdim,curNbOfElem,geoElt,type,entity,mrs);
+    }
+  else
+    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::loadPart : not implemented yet for the dynamic type !");
+}
+
 void MEDFileUMeshPerType::loadPart(med_idt fid, const char *mName, int dt, int it, int mdim, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
                                    med_entity_type entity, mcIdType strt, mcIdType end, mcIdType step, MEDFileMeshReadSelector *mrs)
 {
@@ -189,12 +223,13 @@ void MEDFileUMeshPerType::loadPart(med_idt fid, const char *mName, int dt, int i
   _pd=PartDefinition::New(strt,end,step);
   if(!cm.isDynamic())
     {
-      loadPartStaticType(fid,mName,dt,it,mdim,curNbOfElem,geoElt,type,entity,strt,end,step,mrs);
+      loadPartStaticType(fid,mName,dt,it,mdim,curNbOfElem,geoElt,type,entity,mrs);
     }
   else
     throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::loadPart : not implemented yet for the dynamic type !");
 }
 
+
 void MEDFileUMeshPerType::loadFromStaticType(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
                                              med_entity_type entity, MEDFileMeshReadSelector *mrs)
 {
@@ -210,50 +245,51 @@ void MEDFileUMeshPerType::loadFromStaticType(med_idt fid, const char *mName, int
 }
 
 void MEDFileUMeshPerType::loadPartStaticType(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
-                                             med_entity_type entity, mcIdType strt, mcIdType end, mcIdType step, MEDFileMeshReadSelector *mrs)
+                                             med_entity_type entity, MEDFileMeshReadSelector *mrs)
 {
-  if(strt<0)
-    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::loadPartStaticType : start pos is negative !");
-  if(end>curNbOfElem)
-    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::loadPartStaticType : end is after the authorized range !");
-  mcIdType nbOfEltsToLoad(DataArray::GetNumberOfItemGivenBES(strt,end,step,"MEDFileUMeshPerType::loadPartStaticType"));
   _m=MEDCoupling1SGTUMesh::New(mName,type);
   MEDCoupling1SGTUMesh *mc(dynamic_cast<MEDCoupling1SGTUMesh *>((MEDCoupling1GTUMesh *)_m));
   MCAuto<DataArrayMedInt> conn(DataArrayMedInt::New());
   mcIdType nbOfNodesPerCell(mc->getNumberOfNodesPerCell());
+  if(!_pd)
+    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::loadPartStaticType : no part definition !");
+  mcIdType nbOfEltsToLoad(_pd->getNumberOfElems());
   conn->alloc(nbOfNodesPerCell*nbOfEltsToLoad,1);
-  med_filter filter=MED_FILTER_INIT;
-  MEDfilterBlockOfEntityCr(fid,/*nentity*/ToMedInt(curNbOfElem),/*nvaluesperentity*/1,/*nconstituentpervalue*/ToMedInt(nbOfNodesPerCell),
-                           MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
-                           /*start*/ToMedInt(strt+1),/*stride*/ToMedInt(step),/*count*/1,/*blocksize*/ToMedInt(nbOfEltsToLoad),
-                           /*lastblocksize=useless because count=1*/0,&filter);
-  MEDFILESAFECALLERRD0(MEDmeshElementConnectivityAdvancedRd,(fid,mName,dt,it,entity,geoElt,MED_NODAL,&filter,conn->getPointer()));
-  MEDfilterClose(&filter);
+  MEDFilterEntity filter;
+  filter.init();
+  filter.fill(fid,/*nentity*/curNbOfElem,/*nvaluesperentity*/1,/*nconstituentpervalue*/nbOfNodesPerCell,
+      MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
+      _pd);
+  MEDFILESAFECALLERRD0(MEDmeshElementConnectivityAdvancedRd,(fid,mName,dt,it,entity,geoElt,MED_NODAL,filter.getPtr(),conn->getPointer()));
+  filter.close();
   std::transform(conn->begin(),conn->end(),conn->getPointer(),std::bind(std::plus<med_int>(),std::placeholders::_1,-1));
   mc->setNodalConnectivity(FromMedIntArray<mcIdType>(conn));
-  loadPartOfCellCommonPart(fid,mName,strt,end,step,dt,it,mdim,curNbOfElem,geoElt,entity,mrs);
+  loadPartOfCellCommonPart(fid,mName,dt,it,mdim,curNbOfElem,geoElt,entity,mrs);
 }
 
-void MEDFileUMeshPerType::loadPartOfCellCommonPart(med_idt fid, const char *mName, mcIdType strt, mcIdType stp, mcIdType step, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, med_entity_type entity, MEDFileMeshReadSelector *mrs)
+
+void MEDFileUMeshPerType::loadPartOfCellCommonPart(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, med_entity_type entity, MEDFileMeshReadSelector *mrs)
 {
   med_bool changement,transformation;
+  if(!_pd)
+    throw INTERP_KERNEL::Exception("MEDFileUMeshPerType::loadPartOfCellCommonPart : no part definition !");
+  mcIdType nbOfEltsToLoad(_pd->getNumberOfElems());
   _fam=0;
-  mcIdType nbOfEltsToLoad(DataArray::GetNumberOfItemGivenBES(strt,stp,step,"MEDFileUMeshPerType::loadPartOfCellCommonPart"));
   if(MEDmeshnEntity(fid,mName,dt,it,entity,geoElt,MED_FAMILY_NUMBER,MED_NODAL,&changement,&transformation)>0)
     {
       if(!mrs || mrs->isCellFamilyFieldReading())
         {
           MCAuto<DataArrayMedInt> miFam(DataArrayMedInt::New());
           miFam->alloc(nbOfEltsToLoad,1);
-          med_filter filter=MED_FILTER_INIT;
-          MEDfilterBlockOfEntityCr(fid,/*nentity*/ToMedInt(curNbOfElem),/*nvaluesperentity*/1,/*nconstituentpervalue*/1,
-                                   MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
-                                   /*start*/ToMedInt(strt+1),/*stride*/ToMedInt(step),/*count*/1,/*blocksize*/ToMedInt(nbOfEltsToLoad),
-                                   /*lastblocksize=useless because count=1*/0,&filter);
-          if(MEDmeshEntityAttributeAdvancedRd(fid,mName,MED_FAMILY_NUMBER,dt,it,entity,geoElt,&filter,miFam->getPointer())!=0)
+          MEDFilterEntity filter;
+          filter.init();
+          filter.fill(fid,/*nentity*/curNbOfElem,/*nvaluesperentity*/1,/*nconstituentpervalue*/1,
+              MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
+              _pd);
+          if(MEDmeshEntityAttributeAdvancedRd(fid,mName,MED_FAMILY_NUMBER,dt,it,entity,geoElt,filter.getPtr(),miFam->getPointer())!=0)
             miFam->fillWithZero();
           _fam=FromMedIntArray<mcIdType>(miFam);
-          MEDfilterClose(&filter);
+          filter.close();
         }
     }
   _num=0;
@@ -263,15 +299,15 @@ void MEDFileUMeshPerType::loadPartOfCellCommonPart(med_idt fid, const char *mNam
         {
           MCAuto<DataArrayMedInt> miNum(DataArrayMedInt::New());
           miNum->alloc(nbOfEltsToLoad,1);
-          med_filter filter=MED_FILTER_INIT;
-          MEDfilterBlockOfEntityCr(fid,/*nentity*/ToMedInt(curNbOfElem),/*nvaluesperentity*/1,/*nconstituentpervalue*/1,
-                                   MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
-                                   /*start*/ToMedInt(strt+1),/*stride*/ToMedInt(step),/*count*/1,/*blocksize*/ToMedInt(nbOfEltsToLoad),
-                                   /*lastblocksize=useless because count=1*/0,&filter);
-          if(MEDmeshEntityAttributeAdvancedRd(fid,mName,MED_NUMBER,dt,it,entity,geoElt,&filter,miNum->getPointer())!=0)
+          MEDFilterEntity filter;
+          filter.init();
+          filter.fill(fid,/*nentity*/curNbOfElem,/*nvaluesperentity*/1,/*nconstituentpervalue*/1,
+              MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
+              _pd);
+          if(MEDmeshEntityAttributeAdvancedRd(fid,mName,MED_NUMBER,dt,it,entity,geoElt,filter.getPtr(),miNum->getPointer())!=0)
             miNum->fillWithZero();
           _num=FromMedIntArray<mcIdType>(miNum);
-          MEDfilterClose(&filter);
+          filter.close();
         }
     }
   _names=0;
@@ -281,16 +317,16 @@ void MEDFileUMeshPerType::loadPartOfCellCommonPart(med_idt fid, const char *mNam
         {
           _names=DataArrayAsciiChar::New();
           _names->alloc(nbOfEltsToLoad+1,MED_SNAME_SIZE);//not a bug to avoid the memory corruption due to last \0 at the end
-          med_filter filter=MED_FILTER_INIT;
-          MEDfilterBlockOfEntityCr(fid,/*nentity*/ToMedInt(curNbOfElem),/*nvaluesperentity*/1,/*nconstituentpervalue*/1,
-                                   MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
-                                   /*start*/ToMedInt(strt+1),/*stride*/ToMedInt(step),/*count*/1,/*blocksize*/ToMedInt(nbOfEltsToLoad),
-                                   /*lastblocksize=useless because count=1*/0,&filter);
-          if(MEDmeshEntityAttributeAdvancedRd(fid,mName,MED_NAME,dt,it,entity,geoElt,&filter,_names->getPointer())!=0)
+          MEDFilterEntity filter;
+          filter.init();
+          filter.fill(fid,/*nentity*/curNbOfElem,/*nvaluesperentity*/1,/*nconstituentpervalue*/1,
+              MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
+              _pd);
+          if(MEDmeshEntityAttributeAdvancedRd(fid,mName,MED_NAME,dt,it,entity,geoElt,filter.getPtr(),_names->getPointer())!=0)
             _names=0;
           else
             _names->reAlloc(nbOfEltsToLoad);//not a bug to avoid the memory corruption due to last \0 at the end
-          MEDfilterClose(&filter);
+          filter.close();
         }
     }
 }
index 746a0a9b8ac4c15d0a98f8880b07ec293fa5ab87..cc65f166db1924ee08c66352b4cd7f5ea900881c 100644 (file)
@@ -58,6 +58,8 @@ 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);
+
     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);
     std::size_t getHeapMemorySizeWithoutChildren() const;
@@ -72,19 +74,21 @@ namespace MEDCoupling
                         med_entity_type entity, MEDFileMeshReadSelector *mrs);
     void loadPart(med_idt fid, const char *mName, int dt, int it, int mdim, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
                   med_entity_type entity, mcIdType strt, mcIdType end, mcIdType step, MEDFileMeshReadSelector *mrs);
+    void loadPart(med_idt fid, const char *mName, int dt, int it, int mdim, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
+                  med_entity_type entity, const std::vector<mcIdType>& distrib, MEDFileMeshReadSelector *mrs);
     void loadFromStaticType(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
                             med_entity_type entity, MEDFileMeshReadSelector *mrs);
     void loadPartStaticType(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, INTERP_KERNEL::NormalizedCellType type,
-                            med_entity_type entity, mcIdType strt, mcIdType end, mcIdType step, MEDFileMeshReadSelector *mrs);
+                            med_entity_type entity, MEDFileMeshReadSelector *mrs);
     void loadPolyg(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType arraySize, med_geometry_type geoElt,
                    med_entity_type entity, MEDFileMeshReadSelector *mrs);
     void loadPolyh(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType connFaceLgth, med_geometry_type geoElt,
                    med_entity_type entity, MEDFileMeshReadSelector *mrs);
-    void loadPartOfCellCommonPart(med_idt fid, const char *mName, mcIdType strt, mcIdType stp, mcIdType step, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, med_entity_type entity, MEDFileMeshReadSelector *mrs);
+    void loadPartOfCellCommonPart(med_idt fid, const char *mName, int dt, int it, int mdim, mcIdType curNbOfElem, med_geometry_type geoElt, med_entity_type entity, MEDFileMeshReadSelector *mrs);
+
   private:
     MCAuto<MEDCoupling1GTUMesh> _m;
     MCAuto<PartDefinition> _pd;
   };
 }
-
 #endif
index 656955a1c153ef01249fc0824ef4918bbfd33550..5535647f6313c4d881b3f0b5f162897046887266 100644 (file)
@@ -31,6 +31,7 @@
 #include "InterpKernelAutoPtr.hxx"
 #include "CellModel.hxx"
 
+#include "MEDFilterEntity.hxx"
 #include <set>
 #include <iomanip>
 
@@ -582,6 +583,39 @@ void MEDFileUMeshL2::loadPart(med_idt fid, const MeshOrStructMeshCls *mId, const
   }
 }
 
+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)
+{
+  int Mdim;
+  std::vector<std::string> infosOnComp(loadCommonPart(fid,mId,mName,dt,it,Mdim));
+  if(Mdim==-4)
+    return ;
+  loadPartOfConnectivityFromUserDistrib(fid,Mdim,mName,type,distrib,dt,it,mrs);
+  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();
+  mesh->computeNodeIdsAlg(fetchedNodeIds);
+
+  if(!mrs || mrs->getNumberOfCoordsLoadSessions()==1)
+    {
+      std::vector<mcIdType> distribNodes;
+      for(mcIdType i=0; i<nCoords; i++)
+        if(fetchedNodeIds[i])
+          distribNodes.push_back(i);
+
+      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;
+      mesh->renumberNodesInConn(mapO2N);
+      this->loadPartCoords(fid,infosOnComp,mName,dt,it,distribNodes);
+    }
+  else
+    throw INTERP_KERNEL::Exception("MEDFileUMeshL2::loadPartFromUserDistrib: multiple load sessions not handled!");
+}
+
 void MEDFileUMeshL2::loadConnectivity(med_idt fid, int mdim, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   _per_type_mesh.resize(1);
@@ -614,6 +648,14 @@ void MEDFileUMeshL2::loadPartOfConnectivity(med_idt fid, int mdim, const std::st
   sortTypes();
 }
 
+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)
+{
+  _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);
+}
+
 void MEDFileUMeshL2::loadCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it)
 {
   int spaceDim((int)infosOnComp.size());
@@ -662,29 +704,67 @@ void MEDFileUMeshL2::loadCoords(med_idt fid, const std::vector<std::string>& inf
     _coords->setInfoOnComponent(i,infosOnComp[i]);
 }
 
+void MEDFileUMeshL2::LoadPartCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, const std::vector<mcIdType>& distribNodes,
+MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords)
+{
+  med_int spaceDim((int)infosOnComp.size());
+  allocCoordsPartCoords(spaceDim,distribNodes,_coords,_part_coords);
+  _coords->setInfoOnComponents(infosOnComp);
+  fillPartCoords(fid,spaceDim,mName,dt,it,_coords,_part_coords,_fam_coords,_num_coords,_name_coords);
+}
+
 void MEDFileUMeshL2::LoadPartCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, mcIdType nMin, mcIdType nMax,
 MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords)
 {
-  med_bool changement,transformation;
-  med_int spaceDim((int)infosOnComp.size()),nCoords(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_NODE,MED_NONE,MED_COORDINATE,MED_NO_CMODE,&changement,&transformation));
+  med_int spaceDim((int)infosOnComp.size());
+  allocCoordsPartCoords(spaceDim,nMin,nMax,_coords,_part_coords);
+  _coords->setInfoOnComponents(infosOnComp);
+  fillPartCoords(fid,spaceDim,mName,dt,it,_coords,_part_coords,_fam_coords,_num_coords,_name_coords);
+}
+
+void MEDFileUMeshL2::allocCoordsPartCoords(mcIdType spaceDim, const std::vector<mcIdType>& nodeIds, MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords)
+{
+  mcIdType nbNodesToLoad(nodeIds.size());
+  _coords=DataArrayDouble::New();
+  _coords->alloc(nbNodesToLoad,spaceDim);
+
+  MCAuto<DataArrayIdType> nodeIdsArray=DataArrayIdType::New();
+  nodeIdsArray->useArray(nodeIds.data(),false,DeallocType::C_DEALLOC,nbNodesToLoad,1);
+  _part_coords=PartDefinition::New(nodeIdsArray);
+}
+
+void MEDFileUMeshL2::allocCoordsPartCoords(mcIdType spaceDim, mcIdType nMin, mcIdType nMax, MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords)
+{
   _coords=DataArrayDouble::New();
   mcIdType nbNodesToLoad(nMax-nMin);
   _coords->alloc(nbNodesToLoad,spaceDim);
-  med_filter filter=MED_FILTER_INIT,filter2=MED_FILTER_INIT;
-  MEDfilterBlockOfEntityCr(fid,/*nentity*/nCoords,/*nvaluesperentity*/1,/*nconstituentpervalue*/spaceDim,
-                           MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
-                           /*start*/ToMedInt(nMin+1),/*stride*/1,/*count*/1,/*blocksize*/ToMedInt(nbNodesToLoad),
-                           /*lastblocksize=useless because count=1*/0,&filter);
-  MEDFILESAFECALLERRD0(MEDmeshNodeCoordinateAdvancedRd,(fid,mName.c_str(),dt,it,&filter,_coords->getPointer()));
+
   _part_coords=PartDefinition::New(nMin,nMax,1);
-  MEDfilterClose(&filter);
-  MEDfilterBlockOfEntityCr(fid,nCoords,1,1,MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,
-                           MED_NO_PROFILE,ToMedInt(nMin+1),1,1,ToMedInt(nbNodesToLoad),0,&filter2);
+}
+
+void MEDFileUMeshL2::fillPartCoords(med_idt fid, mcIdType spaceDim, const std::string& mName, int dt, int it,
+                                    MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords)
+{
+  med_bool changement,transformation;
+  med_int nCoords(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_NODE,MED_NONE,MED_COORDINATE,MED_NO_CMODE,&changement,&transformation));
+  mcIdType nbNodesToLoad = _part_coords->getNumberOfElems();
+
+  MEDFilterEntity filter1, filter2;
+  filter1.init();
+  filter1.fill(fid,/*nentity*/nCoords,/*nvaluesperentity*/1,/*nconstituentpervalue*/spaceDim,
+                           MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE,
+                           _part_coords);
+  MEDFILESAFECALLERRD0(MEDmeshNodeCoordinateAdvancedRd,(fid,mName.c_str(),dt,it,filter1.getPtr(),_coords->getPointer()));
+  filter1.close();
+
+  filter2.init();
+  filter2.fill(fid,nCoords,1,1,
+               MED_ALL_CONSTITUENT,MED_FULL_INTERLACE,MED_COMPACT_STMODE,MED_NO_PROFILE, _part_coords);
   if(MEDmeshnEntity(fid,mName.c_str(),dt,it,MED_NODE,MED_NO_GEOTYPE,MED_FAMILY_NUMBER,MED_NODAL,&changement,&transformation)>0)
     {
       MCAuto<DataArrayMedInt> miFamCoord=DataArrayMedInt::New();
       miFamCoord->alloc(nbNodesToLoad,1);
-      MEDFILESAFECALLERRD0(MEDmeshEntityAttributeAdvancedRd,(fid,mName.c_str(),MED_FAMILY_NUMBER,dt,it,MED_NODE,MED_NO_GEOTYPE,&filter2,miFamCoord->getPointer()));
+      MEDFILESAFECALLERRD0(MEDmeshEntityAttributeAdvancedRd,(fid,mName.c_str(),MED_FAMILY_NUMBER,dt,it,MED_NODE,MED_NO_GEOTYPE,filter2.getPtr(),miFamCoord->getPointer()));
       _fam_coords=FromMedIntArray<mcIdType>(miFamCoord);
     }
   else
@@ -693,7 +773,7 @@ MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<D
     {
       MCAuto<DataArrayMedInt> miNumCoord=DataArrayMedInt::New();
       miNumCoord->alloc(nbNodesToLoad,1);
-      MEDFILESAFECALLERRD0(MEDmeshEntityAttributeAdvancedRd,(fid,mName.c_str(),MED_NUMBER,dt,it,MED_NODE,MED_NO_GEOTYPE,&filter2,miNumCoord->getPointer()));
+      MEDFILESAFECALLERRD0(MEDmeshEntityAttributeAdvancedRd,(fid,mName.c_str(),MED_NUMBER,dt,it,MED_NODE,MED_NO_GEOTYPE,filter2.getPtr(),miNumCoord->getPointer()));
       _num_coords=FromMedIntArray<mcIdType>(miNumCoord);
     }
   else
@@ -702,15 +782,15 @@ MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<D
     {
       _name_coords=DataArrayAsciiChar::New();
       _name_coords->alloc(nbNodesToLoad+1,MED_SNAME_SIZE);//not a bug to avoid the memory corruption due to last \0 at the end
-      MEDFILESAFECALLERRD0(MEDmeshEntityAttributeAdvancedRd,(fid,mName.c_str(),MED_NAME,dt,it,MED_NODE,MED_NO_GEOTYPE,&filter2,_name_coords->getPointer()));
+      MEDFILESAFECALLERRD0(MEDmeshEntityAttributeAdvancedRd,(fid,mName.c_str(),MED_NAME,dt,it,MED_NODE,MED_NO_GEOTYPE,filter2.getPtr(),_name_coords->getPointer()));
       _name_coords->reAlloc(nbNodesToLoad);//not a bug to avoid the memory corruption due to last \0 at the end
     }
   else
     _name_coords=nullptr;
-  MEDfilterClose(&filter2);
-  _coords->setInfoOnComponents(infosOnComp);
+  filter2.close();
 }
 
+
 /*!
  * For performance reasons LoadPartCoordsArray method calls LoadPartCoords
  */
@@ -745,6 +825,12 @@ void MEDFileUMeshL2::loadPartCoords(med_idt fid, const std::vector<std::string>&
   LoadPartCoords(fid,infosOnComp,mName,dt,it,nMin,nMax,_coords,_part_coords,_fam_coords,_num_coords,_name_coords);
 }
 
+void MEDFileUMeshL2::loadPartCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, const std::vector<mcIdType>& distribNodes)
+{
+  LoadPartCoords(fid,infosOnComp,mName,dt,it,distribNodes,_coords,_part_coords,_fam_coords,_num_coords,_name_coords);
+}
+
+
 void MEDFileUMeshL2::loadPartCoordsSlice(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, const DataArrayIdType *nodeIds, mcIdType nbOfCoordLS)
 {
   nodeIds->checkAllocated();
index c79c371c8cc9c015e34b26557d2b0496db81bb01..c04126080fec9fb6ed6ed72327df01b07472beb6 100644 (file)
@@ -121,10 +121,15 @@ 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 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 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);
+    void loadPartCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, const std::vector<mcIdType>& distribNodes);
     void loadPartCoordsSlice(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, const DataArrayIdType *nodeIds, mcIdType nbOfCoordLS);
     int getNumberOfLevels() const { return (int)_per_type_mesh.size(); }
     bool emptyLev(int levId) const { return _per_type_mesh[levId].empty(); }
@@ -138,10 +143,15 @@ namespace MEDCoupling
     MCAuto<DataArrayIdType> getCoordsGlobalNum() const { return _global_num_coords; }
     MCAuto<DataArrayAsciiChar> getCoordsName() const { return _name_coords; }
     static void WriteCoords(med_idt fid, const std::string& mname, int dt, int it, double time, const DataArrayDouble *coords, const DataArrayIdType *famCoords, const DataArrayIdType *numCoords, const DataArrayAsciiChar *nameCoords, const DataArrayIdType *globalNumCoords);
+    static void LoadPartCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, const std::vector<mcIdType>& distribNodes,
+    MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords);
     static void LoadPartCoords(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, mcIdType nMin, mcIdType nMax,
-MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords);
+    MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords);
     static void LoadPartCoordsArray(med_idt fid, const std::vector<std::string>& infosOnComp, const std::string& mName, int dt, int it, const DataArrayIdType *nodeIds,
 MCAuto<DataArrayDouble>& _coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords);
+    static void allocCoordsPartCoords(mcIdType spaceDim, mcIdType nMin, mcIdType nMax, MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords);
+    static void allocCoordsPartCoords(mcIdType spaceDim, const std::vector<mcIdType>& nodeIds, MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords);
+    static void fillPartCoords(med_idt fid, mcIdType spaceDim, const std::string& mName, int dt, int it,MCAuto<DataArrayDouble>& _coords, MCAuto<PartDefinition>& _part_coords, MCAuto<DataArrayIdType>& _fam_coords, MCAuto<DataArrayIdType>& _num_coords, MCAuto<DataArrayAsciiChar>& _name_coords);
   private:
     void sortTypes();
   private:
diff --git a/src/MEDLoader/MEDFilterEntity.hxx b/src/MEDLoader/MEDFilterEntity.hxx
new file mode 100644 (file)
index 0000000..4a61d6e
--- /dev/null
@@ -0,0 +1,89 @@
+// Copyright (C) 2007-2022  CEA/DEN, EDF R&D
+//
+// This library is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+//
+// This library is distributed in the hope that it will be useful,
+// but WITHOUT ANY WARRANTY; without even the implied warranty of
+// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+// Lesser General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License along with this library; if not, write to the Free Software
+// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
+//
+// See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
+//
+// Author : Anthony Geay (CEA/DEN)
+
+#ifndef __MEDFILTERENTITY_HXX__
+#define __MEDFILTERENTITY_HXX__
+
+#include "MEDCouplingPartDefinition.hxx"
+#include "med.h"
+
+namespace MEDCoupling
+{
+  class MEDFilterEntity
+  {
+  public:
+    MEDFilterEntity() { _filter = new med_filter; }
+    ~MEDFilterEntity() { delete _filter; }
+    void init() { *_filter = MED_FILTER_INIT; }
+    inline void fill(med_idt fid, mcIdType nbOfEntity, mcIdType nbOfValuesPerEntity, mcIdType nbOfConstituentPerValue,
+                     const med_int constituentSelect, const med_switch_mode switchMode, const med_storage_mode storageMode, const char * const profileName,
+                     const PartDefinition* pd);
+    const med_filter *getPtr() const { return _filter; }
+    void close() { MEDfilterClose(_filter); }
+  private:
+    med_filter *_filter;
+  };
+}
+
+void MEDCoupling::MEDFilterEntity::fill(med_idt fid, mcIdType nbOfEntity, mcIdType nbOfValuesPerEntity, mcIdType nbOfConstituentPerValue,
+                                          const med_int constituentSelect, const med_switch_mode switchMode, const med_storage_mode storageMode, const char * const profileName,
+                                          const PartDefinition* pd)
+{
+  const SlicePartDefinition *spd(dynamic_cast<const SlicePartDefinition *>(pd));
+  if(spd)
+    {
+      //Here, pd contains a slice, so it's more efficient to define a filter of block
+      //(which will load contiguous values)
+      mcIdType nbOfEltsToLoad = spd->getNumberOfElems();
+      mcIdType strt,end,step;
+      spd->getSlice(strt,end,step);
+      if(strt<0)
+        throw INTERP_KERNEL::Exception("MEDFilterEntity::fill : start pos is negative !");
+      if(end>nbOfEntity)
+        throw INTERP_KERNEL::Exception("MEDFilterEntity::fill : end is after the authorized range !");
+      MEDfilterBlockOfEntityCr(fid,ToMedInt(nbOfEntity),ToMedInt(nbOfValuesPerEntity),ToMedInt(nbOfConstituentPerValue),
+                               constituentSelect,switchMode,storageMode,profileName,
+                               /*start*/ToMedInt(strt+1),/*stride*/ToMedInt(step),/*count*/1,/*blocksize*/ToMedInt(nbOfEltsToLoad),
+                               /*lastblocksize=useless because count=1*/0,_filter);
+      return;
+    }
+  const DataArrayPartDefinition *dpd(dynamic_cast<const DataArrayPartDefinition *>(pd));
+  if(dpd)
+    {
+      mcIdType nbOfEltsToLoad = dpd->getNumberOfElems();
+
+      //convert to fortran indexing
+      std::vector<mcIdType> dpdPlus1;
+      std::copy(pd->toDAI()->begin(), pd->toDAI()->end(), std::back_inserter(dpdPlus1));
+      std::for_each(dpdPlus1.begin(), dpdPlus1.end(), [](mcIdType &node){ node+=1; });
+
+      //Here, pd contains a random selection of non-contiguous values:
+      //we need to use a more generic filter (less efficient)
+      MEDfilterEntityCr(fid,ToMedInt(nbOfEntity),ToMedInt(nbOfValuesPerEntity),ToMedInt(nbOfConstituentPerValue),
+                            constituentSelect,switchMode,storageMode,profileName,
+                            ToMedInt(nbOfEltsToLoad), dpdPlus1.data(),
+                            _filter);
+      return;
+    }
+  throw INTERP_KERNEL::Exception("MEDFilterEntity::fill : empty part definition !");
+}
+
+
+#endif
index c5b64689d761ee99367f755c97b8e7998c3d4fc4..364141bd0777a6bd784977a05ff47885cee44df9 100644 (file)
 #include "MEDFileMesh.hxx"
 #include "MEDFileMeshLL.hxx"
 #include "MEDLoader.hxx"
+#include <MEDFileField1TS.hxx>
+#include <iostream>
+#include <fstream>
+
+
+// From MEDLOader.cxx TU
+extern med_geometry_type typmai3[INTERP_KERNEL::NORM_MAXTYPE];
 
 using namespace MEDCoupling;
 
@@ -73,18 +80,98 @@ MEDFileUMesh *ParaMEDFileUMesh::New(int iPart, int nbOfParts, const std::string&
   return ParaMEDFileUMesh::NewPrivate(fid,iPart,nbOfParts,fileName,mName,dt,it,mrs);
 }
 
-// MPI_COMM_WORLD, MPI_INFO_NULL 
-MEDFileUMesh *ParaMEDFileUMesh::ParaNew(int iPart, int nbOfParts, const MPI_Comm& com, const MPI_Info& nfo, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
+/*!
+ * Opens the given file in parallel to load a specific part of it (the selection of cells is done via \a distrib)
+ * \param [in] distrib - vector of cells we want to load
+ * \param [in] com - group of MPI processes that will read the file
+ * \param [in] nfo- MPI info object (used to manage MPI routines)
+ * \param [in] filename - name of the file we want to read
+ * \param [in] mName - name of the mesh we want to read
+ * \param [in] dt - order at which to read the mesh
+ * \param [in] it - iteration at which to read the mesh
+ * \param [in] mrs - object used to read additional low-level information
+ * \return MEDFileUMesh* - a new instance of MEDFileUMesh. The
+ *         caller is to delete this mesh using decrRef() as it is no more needed.
+ */
+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)
 {
   MEDFileUtilities::CheckFileForRead(fileName);
 #ifdef HDF5_IS_PARALLEL
   MEDFileUtilities::AutoFid fid(MEDparFileOpen(fileName.c_str(),MED_ACC_RDONLY,com,nfo));
 #else
   MEDFileUtilities::AutoFid fid(MEDfileOpen(fileName.c_str(),MED_ACC_RDONLY));
+#endif
+  return ParaMEDFileUMesh::NewPrivate(fid,com,distrib,fileName,mName,dt,it,mrs);
+}
+
+
+/*!
+ * Opens the given file in parallel to partially load the file. The mesh will be equally and linearly distributed among all processes:
+ * the list of cells will be divided into \a nbOfParts slices and only slice \a iPart will be loaded.
+ * \param [in] iPart - part of the mesh that will be loaded
+ * \param [in] nbOfParts - total number of parts in which to divide the mesh
+ * \param [in] com - group of MPI processes that will read the file
+ * \param [in] nfo- MPI info object (used to manage MPI routines)
+ * \param [in] filename - name of the file we want to read
+ * \param [in] mName - name of the mesh we want to read
+ * \param [in] dt - Time order at which to read the mesh
+ * \param [in] it - Time iteration at which to read the mesh
+ * \param [in] mrs - object used to read additional low-level information
+ * \return MEDFileUMesh* - a new instance of MEDFileUMesh. The
+ *         caller is to delete this mesh using decrRef() as it is no more needed.
+ */
+MEDFileUMesh *ParaMEDFileUMesh::ParaNew(int iPart, int nbOfParts, 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
+  MEDFileUtilities::AutoFid fid(MEDparFileOpen(fileName.c_str(),MED_ACC_RDONLY,com,nfo)); // MPI_COMM_WORLD, MPI_INFO_NULL
+#else
+  MEDFileUtilities::AutoFid fid(MEDfileOpen(fileName.c_str(),MED_ACC_RDONLY));
 #endif
   return ParaMEDFileUMesh::NewPrivate(fid,iPart,nbOfParts,fileName,mName,dt,it,mrs);
 }
 
+/*!
+ * Partially load the file (using a custom distribution of the cells)
+ * See ParaMEDFileUMesh::ParaNew for detailed description.
+ * \throw exception if the mesh contains multiple types
+ */
+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)
+{
+  MCAuto<MEDFileUMesh> ret;
+  int meshDim, spaceDim;
+  mcIdType numberOfNodes;
+  std::vector< std::vector< std::pair<INTERP_KERNEL::NormalizedCellType,int> > > typesDistrib(GetUMeshGlobalInfo(fileName,mName,meshDim,spaceDim,numberOfNodes));
+  std::size_t numberOfTypes = typesDistrib.size();
+  if(numberOfTypes != 1)
+    throw INTERP_KERNEL::Exception("ParaMEDFileMesh::NewPrivate : only mesh with single geometrical type are supported with given distribution !");
+
+  INTERP_KERNEL::NormalizedCellType geoType = typesDistrib[0][0].first;
+  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));
+  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::NewPrivate : Some of your partitions overlap each other ! Each element in your distribution vector must appear only once ! ");
+      else
+        throw INTERP_KERNEL::Exception("ParaMEDFileMesh::NewPrivate : The distribution does not cover the whole mesh ! Each element of the mesh must appear once in your distribution vector ");
+    }
+  ret=MEDFileUMesh::LoadPartOfFromUserDistrib(fid,mName,geoType,distrib,dt,it,mrs);
+  return ret.retn();
+}
+
+/*!
+ * Partially load the mesh \a mName inside the file named \a fileName (using a linear distribution of the cells)
+ * See ParaMEDFileUMesh::ParaNew for detailed description.
+ */
 MEDFileUMesh *ParaMEDFileUMesh::NewPrivate(med_idt fid, int iPart, int nbOfParts, const std::string& fileName, const std::string& mName, int dt, int it, MEDFileMeshReadSelector *mrs)
 {
   MCAuto<MEDFileUMesh> ret;
@@ -106,6 +193,20 @@ MEDFileUMesh *ParaMEDFileUMesh::NewPrivate(med_idt fid, int iPart, int nbOfParts
   return ret.retn();
 }
 
+MEDFileField1TS *ParaMEDFileField1TS::ParaNew(MEDFileUMesh* partMesh, const std::string& fileName, const std::string& mName, const std::string& fName, int dt, int it)
+{
+  MCAuto<MEDFileMeshes> ms(MEDFileMeshes::New());
+  ms->pushMesh(partMesh);
+  MCAuto<MEDFileAnyTypeField1TS> partFile(MEDFileAnyTypeField1TS::NewAdv(fileName,fName,dt,it,ms));
+
+  MCAuto<MEDFileField1TS> ret(MEDCoupling::DynamicCast<MEDFileAnyTypeField1TS,MEDFileField1TS>(partFile));
+  if(ret.isNotNull())
+    return ret.retn();
+  else
+    throw INTERP_KERNEL::Exception("ParaMEDFileField1TS::ParaNew : only FLOAT64 field supported for the moment !");
+}
+
+
 MEDFileMeshes *ParaMEDFileMeshes::New(int iPart, int nbOfParts, const std::string& fileName)
 {
   std::vector<std::string> ms(GetMeshNames(fileName));
index da482e3067063fcaebf61c1775d0998b43e3a35e..22abf48dd596c7dc486c9226613812739d530c4f 100644 (file)
@@ -26,6 +26,8 @@
 #include "mpi.h"
 
 #include <string>
+#include <vector>
+#include "MCIdType.hxx"
 
 namespace MEDCoupling
 {
@@ -33,6 +35,7 @@ namespace MEDCoupling
   class MEDFileUMesh;
   class MEDFileMeshes;
   class MEDFileMeshReadSelector;
+  class MEDFileField1TS;
 
   class ParaMEDFileMesh
   {
@@ -46,16 +49,28 @@ 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);
+
   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);
   };
 
+
   class ParaMEDFileMeshes
   {
   public:
     static MEDFileMeshes *New(int iPart, int nbOfParts, const std::string& fileName);
     static MEDFileMeshes *ParaNew(int iPart, int nbOfParts, const MPI_Comm& com, const MPI_Info& nfo, const std::string& fileName);
   };
+
+  class ParaMEDFileField1TS
+  {
+  public:
+      static MEDFileField1TS *ParaNew(MEDFileUMesh* partMesh, const std::string& fileName, const std::string& mName, const std::string& fName, int dt=-1, int it=-1);
+
+  };
+
 }
 
 #endif