]> SALOME platform Git repositories - modules/paravis.git/commitdiff
Salome HOME
MEDReader is now working in Parallel for *.med files.
authorAnthony Geay <anthony.geay@edf.fr>
Fri, 3 Oct 2014 09:15:01 +0000 (11:15 +0200)
committerAnthony Geay <anthony.geay@edf.fr>
Fri, 3 Oct 2014 09:15:01 +0000 (11:15 +0200)
src/Plugins/MEDReader/IO/CMakeLists.txt
src/Plugins/MEDReader/IO/MEDFileFieldRepresentationTree.cxx
src/Plugins/MEDReader/IO/MEDFileFieldRepresentationTree.hxx
src/Plugins/MEDReader/IO/vtkMEDReader.cxx

index 091b533542b15ac9b88633df02c7abe59eba797b..ce8030ae4f4a960a6a8fda636b98f97a89923641 100644 (file)
@@ -23,6 +23,10 @@ INCLUDE_DIRECTORIES(
   ${MED_ROOT_DIR}/include/salome
   )
 
+IF(PARAVIEW_USE_MPI)
+  ADD_DEFINITIONS("-DMEDREADER_USE_MPI")
+ENDIF(PARAVIEW_USE_MPI)
+
 SET(MEDReader_CLASSES vtkMEDReader vtkExtractGroup vtkELNOMeshFilter vtkELNOSurfaceFilter vtkELNOFilter vtkExtractCellType)
 
 SET(MEDReader_SRCS)
index 26a9c6cd08d1faa0c63bdaba0f03a269303d7417..11bb8031f8641e14b08d79f7dabe03d9702fabf3 100644 (file)
 #include "MEDFileData.hxx"
 #include "SauvReader.hxx"
 
+#ifdef MEDREADER_USE_MPI
+  #include "ParaMEDFileMesh.hxx"
+#endif
+
 #include "vtkXMLUnstructuredGridWriter.h"//
 
 #include "vtkUnstructuredGrid.h"
@@ -1087,12 +1091,17 @@ int MEDFileFieldRepresentationTree::getMaxNumberOfTimeSteps() const
 /*!
  * 
  */
-void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv)
+void MEDFileFieldRepresentationTree::loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv, int iPart, int nbOfParts)
 {
   if(isMEDOrSauv)
     {
+#ifdef MEDREADER_USE_MPI
+      _ms=ParaMEDFileMeshes::New(iPart,nbOfParts,fileName);
+      _fields=MEDFileFields::LoadPartOf(fileName,false,_ms);//false is important to not read the values
+#else
       _ms=MEDFileMeshes::New(fileName);
       _fields=MEDFileFields::New(fileName,false);//false is important to not read the values
+#endif
     }
   else
     {
index 329ff6cb00db0a3d6f88d2b24fa016a474b5942b..d184aa5791c53e016431d25d055ea4790508272b 100644 (file)
@@ -149,7 +149,7 @@ public:
   vtkDataSet *buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk) const;
   void printMySelf(std::ostream& os) const;
   //non const methods
-  void loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv);
+  void loadMainStructureOfFile(const char *fileName, bool isMEDOrSauv, int iPart, int nbOfParts);
   void removeEmptyLeaves();
   // static methods
   static bool IsFieldMeshRegardingInfo(const std::vector<std::string>& compInfos);
index 769e87a4b802a452a251d4acdc0daf0a1e05cd73..7b526f4d898f8b8d7545510d36cecafaa17f7be5 100644 (file)
@@ -257,7 +257,12 @@ void vtkMEDReader::SetFileName(const char *fname)
         }
       if(this->Internal->Tree.getNumberOfLeavesArrays()==0)
         {
-          this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),this->Internal->IsMEDOrSauv);
+          int iPart(-1),nbOfParts(-1);
+#ifdef MEDREADER_USE_MPI
+          MPI_Comm_rank(MPI_COMM_WORLD,&iPart);
+          MPI_Comm_size(MPI_COMM_WORLD,&nbOfParts);
+#endif
+          this->Internal->Tree.loadMainStructureOfFile(this->Internal->FileName.c_str(),this->Internal->IsMEDOrSauv,iPart,nbOfParts);
           if(!this->Internal->PK.arePropertiesOnTreeToSetAfter())
             this->Internal->Tree.activateTheFirst();//This line manually initialize the status of server (this) with the remote client.
           this->Internal->TK.setMaxNumberOfTimeSteps(this->Internal->Tree.getMaxNumberOfTimeSteps());