]> SALOME platform Git repositories - modules/paravis.git/commitdiff
Salome HOME
[EDF26799] : Addition of GetRidOffDebugArrays property to remove numbering output...
authorAnthony Geay <anthony.geay@edf.fr>
Fri, 27 Jan 2023 08:27:15 +0000 (09:27 +0100)
committerAnthony Geay <anthony.geay@edf.fr>
Fri, 27 Jan 2023 08:27:15 +0000 (09:27 +0100)
src/Insitu/VisualizationLibrary/visu.cxx
src/Plugins/MEDReader/plugin/MEDLoaderForPV/MEDFileFieldRepresentationTree.cxx
src/Plugins/MEDReader/plugin/MEDLoaderForPV/MEDFileFieldRepresentationTree.hxx
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkFileSeriesGroupReader.cxx
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkMEDReader.cxx
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkMEDReader.h
src/Plugins/MEDReader/plugin/ParaViewPlugin/Resources/MEDReaderServer.xml

index 7d3f3b0d6366fadebddd5e9205dbb20aa84579da..a8c6957e1b1867fad4eca5982973a84909edae37 100644 (file)
@@ -147,7 +147,7 @@ void Visualization::ConvertToVTK(MEDCoupling::MEDCouplingFieldDouble* field, vtk
   TimeKeeper TK(0);
   int tmp1,tmp2;
   double tmp3(f->getTime(tmp1,tmp2));
-  vtkDataSet *ret(Tree.buildVTKInstance(false,tmp3,meshName,TK));
+  vtkDataSet *ret(Tree.buildVTKInstance(false,tmp3,meshName,true,TK));
   VTKGrid = ret;
 }
 
index cea876abb748a10779dd98d4aadf48c0c722b2a8..1b93ae4af6d74f7b292b9f6f71b147c334cb2365 100644 (file)
@@ -851,8 +851,8 @@ vtkStructuredGrid *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInter
   da->Delete();
   return ret;
 }
-vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDFileMeshes *meshes, ExportedTinyInfo *internalInfo) const
+
+vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDFileMeshes *meshes, bool debugArrays, ExportedTinyInfo *internalInfo) const
 {
   vtkDataSet *ret(0);
   //_fsp->isDataSetSupportEqualToThePreviousOne(i,globs);
@@ -904,7 +904,7 @@ vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolatio
       famCells->decrRef();
     }
   ptMML2->retrieveNumberIdsOnCells(numCells,noCpyNumCells);
-  if(numCells)
+  if(numCells && debugArrays)
     {
       vtkMCIdTypeArray *vtkTab(vtkMCIdTypeArray::New());
       vtkTab->SetNumberOfComponents(1);
@@ -929,7 +929,7 @@ vtkDataSet *MEDFileFieldRepresentationLeaves::buildVTKInstanceNoTimeInterpolatio
       famNodes->decrRef();
     }
   ptMML2->retrieveNumberIdsOnNodes(numNodes,noCpyNumNodes);
-  if(numNodes)
+  if(numNodes && debugArrays)
     {
       vtkMCIdTypeArray *vtkTab(vtkMCIdTypeArray::New());
       vtkTab->SetNumberOfComponents(1);
@@ -1376,7 +1376,7 @@ std::vector<double> MEDFileFieldRepresentationTree::getTimeSteps(int& lev0, cons
   return leaf.getTimeSteps(tk);
 }
 
-vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk, ExportedTinyInfo *internalInfo) const
+vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk, bool debugArrays, ExportedTinyInfo *internalInfo) const
 {
   int lev0,lev1,lev2;
   const MEDFileFieldRepresentationLeaves& leaf(getTheSingleActivated(lev0,lev1,lev2));
@@ -1432,7 +1432,7 @@ vtkDataSet *MEDFileFieldRepresentationTree::buildVTKInstance(bool isStdOrMode, d
     tr=new MEDStdTimeReq((int)zeTimeId);
   else
     tr=new MEDModeTimeReq(tk.getTheVectOfBool(),tk.getPostProcessedTime());
-  vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms,internalInfo));
+  vtkDataSet *ret(leaf.buildVTKInstanceNoTimeInterpolation(tr,_fields,_ms,debugArrays,internalInfo));
   delete tr;
   return ret;
 }
index d2aec99a2d81919d595c8e426d3d3ffe0c13a93f..ceed4fc9bdd5bc8e378abb5839f05f1f7f0ad5ed 100644 (file)
@@ -122,7 +122,7 @@ public:
   std::vector<double> getTimeSteps(const TimeKeeper& tk) const;
   std::vector< std::pair<int,int> > getTimeStepsInCoarseMEDFileFormat(std::vector<double>& ts) const;
   std::string getHumanReadableOverviewOfTS() const;
-  vtkDataSet *buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDFileMeshes *meshes, ExportedTinyInfo *internalInfo=0) const;
+  vtkDataSet *buildVTKInstanceNoTimeInterpolation(const MEDTimeReq *tr, const MEDCoupling::MEDFileFieldGlobsReal *globs, const MEDCoupling::MEDFileMeshes *meshes, bool debugArrays, ExportedTinyInfo *internalInfo=0) const;
 private:
   vtkUnstructuredGrid *buildVTKInstanceNoTimeInterpolationUnstructured(MEDCoupling::MEDUMeshMultiLev *mm) const;
   vtkRectilinearGrid *buildVTKInstanceNoTimeInterpolationCartesian(MEDCoupling::MEDCMeshMultiLev *mm) const;
@@ -154,7 +154,7 @@ public:
   //
   std::string getDftMeshName() const;
   std::vector<double> getTimeSteps(int& lev0, const TimeKeeper& tk) const;
-  vtkDataSet *buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk, ExportedTinyInfo *internalInfo=0) const;
+  vtkDataSet *buildVTKInstance(bool isStdOrMode, double timeReq, std::string& meshName, const TimeKeeper& tk, bool debugArrays, ExportedTinyInfo *internalInfo=0) const;
   void printMySelf(std::ostream& os) const;
   std::map<std::string,bool> dumpState() const;
   //non const methods
index d666796ddd367e4f23a8ab7eedc7d9c38d82d5e3..ee58aaefd8ecede0b500e0fe6b1ecd5031397b4e 100644 (file)
@@ -206,6 +206,7 @@ int vtkFileSeriesGroupReader::RequestData(vtkInformation* vtkNotUsed(request),
       localReader->GenerateVectors(exposedReader->GetGenerateVect());
       localReader->ChangeMode(exposedReader->GetIsStdOrMode());
       localReader->GhostCellGeneratorCallForPara(exposedReader->GetGCGCP());
+      localReader->GetRidOffDebugArrays(exposedReader->GetRemoveDebugArrays());
 
       // Configure the localReader for usage with the files
       localReader->SetFileName(this->GetFileName(i + offFile));
index 7737fcbbca9a65f0edca75132bbc202dcd0420ac..cb9d7c20bec254d7bd6ae1d07dc18528c10354e1 100755 (executable)
@@ -139,6 +139,7 @@ void vtkMEDReader::Reload()
   this->IsStdOrMode = false;
   this->GenerateVect = false;
   this->GCGCP = true;
+  this->RemoveDebugArrays = false;
   this->FieldSelection->RemoveAllArrays();
   this->TimeFlagSelection->RemoveAllArrays();
   this->Modified();
@@ -185,6 +186,19 @@ void vtkMEDReader::GhostCellGeneratorCallForPara(int gcgcp)
     }
 }
 
+void vtkMEDReader::GetRidOffDebugArrays(int rmda)
+{
+  if ( !this->Internal )
+    return;
+
+  bool newVal(rmda!=0);
+  if(newVal!=this->RemoveDebugArrays)
+    {
+      this->RemoveDebugArrays=newVal;
+      this->Modified();
+    }
+}
+
 const char *vtkMEDReader::GetSeparator()
 {
   return MEDFileFieldRepresentationLeavesArrays::ZE_SEP;
@@ -345,21 +359,21 @@ int vtkMEDReader::RequestData(vtkInformation *request, vtkInformationVector ** /
       this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
 #else
       if (this->DistributeWithMPI && this->GCGCP)
-       {
-         vtkSmartPointer<vtkGhostCellsGenerator> gcg(vtkSmartPointer<vtkGhostCellsGenerator>::New());
-         {
-           vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,&ti));
-           gcg->SetInputData(ret);
-           ret->Delete();
-         }
-          // To be checked
-         // gcg->SetUseGlobalPointIds(true);
-         gcg->SetBuildIfRequired(false);
-         gcg->Update();
-         output->SetBlock(0,gcg->GetOutput());
-       }
+      {
+        vtkSmartPointer<vtkGhostCellsGenerator> gcg(vtkSmartPointer<vtkGhostCellsGenerator>::New());
+        {
+          vtkDataSet *ret(RetrieveDataSetAtTime(reqTS,&ti));
+          gcg->SetInputData(ret);
+          ret->Delete();
+        }
+              // To be checked
+        // gcg->SetUseGlobalPointIds(true);
+        gcg->SetBuildIfRequired(false);
+        gcg->Update();
+        output->SetBlock(0,gcg->GetOutput());
+      }
       else
-       this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
+             this->FillMultiBlockDataSetInstance(output,reqTS,&ti);
 #endif
       if(!ti.empty())
         {
@@ -545,7 +559,7 @@ vtkDataSet *vtkMEDReader::RetrieveDataSetAtTime(double reqTS, ExportedTinyInfo *
   if( !this->Internal )
     return 0;
   std::string meshName;
-  vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->IsStdOrMode,reqTS,meshName,this->Internal->TK,internalInfo));
+  vtkDataSet *ret(this->Internal->Tree.buildVTKInstance(this->IsStdOrMode,reqTS,meshName,this->Internal->TK,!this->RemoveDebugArrays,internalInfo));
   if(this->GenerateVect)
     {
       vtkGenerateVectors::Operate(ret->GetPointData());
index 03b1cb479f12f3fc8ca3d9eb6499c6b7ef25c24e..d4d81a822e195b530a398cf26d29c21f934a0c57 100755 (executable)
@@ -100,6 +100,12 @@ class VTK_EXPORT vtkMEDReader : public vtkMultiBlockDataSetAlgorithm
   // Default is true
   void GhostCellGeneratorCallForPara(int);
   vtkGetMacro(GCGCP, bool);
+  
+  // Description
+  // Control if mesh debug arrays should be removed or not
+  // Default is true
+  void GetRidOffDebugArrays(int);
+  vtkGetMacro(RemoveDebugArrays, bool);
 
 
  protected:
@@ -128,6 +134,7 @@ class VTK_EXPORT vtkMEDReader : public vtkMultiBlockDataSetAlgorithm
   bool GenerateVect = false;
   bool GCGCP = true;
   bool DistributeWithMPI = true;
+  bool RemoveDebugArrays = false;
 };
 
 #endif //__vtkMEDReader_h_
index b1522a2dc076ab102337a0a4c552085b264bcf4d..df9bb0e495b937f159cf96ecc8fddaa1e5349c33 100644 (file)
         <BooleanDomain name="bool"/>
       </IntVectorProperty>
 
+      <IntVectorProperty name="GetRidOffDebugArrays"
+                        label="Get Rid Off Debug Arrays"
+                        command="GetRidOffDebugArrays"
+                        number_of_elements="1"
+                        default_values="0"
+                                         panel_visibility="advanced">
+        <Documentation>
+          This property tells if arrays of mesh ids and number ids should be removed or not from output dataset. By default these arrays are present.
+        </Documentation>
+        <BooleanDomain name="bool"/>
+      </IntVectorProperty>
+
    </SourceProxy>
   </ProxyGroup>
   <ProxyGroup name="sources">
           <Property name="TimesFlagsInfo" />
           <Property name="TimesFlagsStatus" />
           <Property name="GhostCellGeneratorCallForPara" />
+          <Property name="GetRidOffDebugArrays" />
         </ExposedProperties>
       </SubProxy>
       <StringVectorProperty animateable="0"