]> SALOME platform Git repositories - modules/paravis.git/commitdiff
Salome HOME
Into MEDReader plugin NewGroupsNames filter linked with GroupAsMultiBlock filter
authorAnthony Geay <anthony.geay@edf.fr>
Wed, 12 Aug 2020 06:33:13 +0000 (08:33 +0200)
committerAnthony Geay <anthony.geay@edf.fr>
Wed, 12 Aug 2020 06:33:13 +0000 (08:33 +0200)
src/Plugins/MEDReader/plugin/MEDReaderIO/CMakeLists.txt
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkGroupsNames.cxx [new file with mode: 0644]
src/Plugins/MEDReader/plugin/MEDReaderIO/vtkGroupsNames.h [new file with mode: 0644]
src/Plugins/MEDReader/plugin/ParaViewPlugin/Resources/MEDReaderServer.xml
src/Plugins/MEDReader/plugin/Test/testMEDReader22.py

index 6d671f2bd06df78e8a3fd6d8bb52f8926153619e..3d9bea0b72dd8e20b70a1c9a3d2658809a4236b3 100644 (file)
@@ -27,6 +27,7 @@ set(classes
   vtkMEDReader
   vtkPVMetaDataInformation
   vtkGroupAsMultiBlock
+  vtkGroupsNames
   vtkUgSelectCellIds
 )
 
diff --git a/src/Plugins/MEDReader/plugin/MEDReaderIO/vtkGroupsNames.cxx b/src/Plugins/MEDReader/plugin/MEDReaderIO/vtkGroupsNames.cxx
new file mode 100644 (file)
index 0000000..3ba7530
--- /dev/null
@@ -0,0 +1,147 @@
+// Copyright (C) 2020  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 (EDF R&D)
+
+#include "vtkGroupsNames.h"
+#include "ExtractGroupHelper.h"
+#include "vtkMEDReader.h"
+#include "vtkUgSelectCellIds.h"
+#include "MEDFileFieldRepresentationTree.hxx"
+#include "VTKMEDTraits.hxx"
+
+#include <vtkCellArray.h>
+#include <vtkCellData.h>
+#include "vtkMutableDirectedGraph.h"
+#include "vtkInformationDataObjectMetaDataKey.h"
+#include <vtkCompositeDataToUnstructuredGridFilter.h>
+#include <vtkMultiBlockDataGroupFilter.h>
+#include <vtkDoubleArray.h>
+#include <vtkInformation.h>
+#include <vtkInformationVector.h>
+#include <vtkNew.h>
+#include <vtkPVGlyphFilter.h>
+#include <vtkPointData.h>
+#include <vtkPolyData.h>
+#include <vtkUnstructuredGrid.h>
+#include <vtkMultiBlockDataSet.h>
+#include <vtkCellCenters.h>
+#include <vtkGlyphSource2D.h>
+#include "vtkTable.h"
+#include "vtkStringArray.h"
+
+class vtkGroupsNamesInternal : public ExtractGroupInternal
+{
+};
+
+vtkStandardNewMacro(vtkGroupsNames);
+
+vtkGroupsNames::vtkGroupsNames():Internal(new ExtractGroupInternal)
+{
+}
+
+vtkGroupsNames::~vtkGroupsNames()
+{
+  delete this->Internal;
+}
+
+int vtkGroupsNames::RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector)
+{
+  vtkInformation *outInfo(outputVector->GetInformationObject(0));
+  vtkInformation *inputInfo(inputVector[0]->GetInformationObject(0));
+  if(!ExtractGroupInternal::IndependantIsInformationOK(vtkMEDReader::META_DATA(),inputInfo))
+  {
+    vtkErrorMacro("No SIL Data available ! The source of this filter must be MEDReader !");
+    return 0;
+  }
+  vtkMutableDirectedGraph *mdg(vtkMutableDirectedGraph::SafeDownCast(inputInfo->Get(vtkMEDReader::META_DATA())));
+  this->Internal->loadFrom(mdg);
+  return 1;
+}
+
+int vtkGroupsNames::RequestData(vtkInformation *vtkNotUsed(request), vtkInformationVector **inputVector, vtkInformationVector *outputVector)
+{
+  vtkInformation* inputInfo=inputVector[0]->GetInformationObject(0);
+  vtkDataSet *input(nullptr);
+  vtkMultiBlockDataSet *inputMB(vtkMultiBlockDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT())));
+  if(inputMB)
+  {
+    if(inputMB->GetNumberOfBlocks()!=1)
+    {
+      vtkErrorMacro("vtkGroupsNames::RequestData : input has not the right number of parts ! Expected 1 !");
+      return 0;
+    }
+    input = vtkDataSet::SafeDownCast(inputMB->GetBlock(0));
+  }
+  else
+  {
+    input = vtkDataSet::SafeDownCast(inputInfo->Get(vtkDataObject::DATA_OBJECT()));
+  }
+  if(!input)
+  {
+    vtkErrorMacro("vtkGroupsNames::RequestData : input is neither a DataSet nor a MultiblockDataSet !");
+    return 0;
+  }
+  vtkInformation *info(input->GetInformation());
+  vtkInformation *outInfo(outputVector->GetInformationObject(0));
+  vtkTable *output(vtkTable::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT())));
+  if(!output)
+  {
+    vtkErrorMacro("vtkGroupsNames::RequestData : something wrong !");
+    return 0; 
+  }
+  //
+  vtkUnstructuredGrid *inputc(vtkUnstructuredGrid::SafeDownCast(input));
+  if(!inputc)
+  {
+    vtkErrorMacro("vtkGroupsNames::RequestData : Only implemented for vtkUnstructuredGrid !");
+    return 0; 
+  }
+  //
+  vtkDataArray *famIdsGen(inputc->GetCellData()->GetScalars(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME));
+  if(!famIdsGen)
+  {
+    vtkErrorMacro(<< "vtkGroupsNames::RequestData : Fail to locate " << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME << " array !");
+    return 0;
+  }
+  using vtkMCIdTypeArray = MEDFileVTKTraits<mcIdType>::VtkType;
+  vtkMCIdTypeArray *famIdsArr(vtkMCIdTypeArray::SafeDownCast(famIdsGen));
+  if(!famIdsArr)
+  {
+    vtkErrorMacro(<< "vtkGroupsNames::RequestData : cell array " << MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME << " is exepected to be IdType !");
+    return 0; 
+  }
+  // Let's go !
+  std::vector< std::pair<std::string,std::vector<int> > > allGroups(this->Internal->getAllGroups());
+  vtkNew<vtkIntArray> blockIdArr;
+  blockIdArr->SetName("Block ID");
+  blockIdArr->SetNumberOfTuples(allGroups.size());
+  vtkNew<vtkStringArray> groupNameArr;
+  groupNameArr->SetName("Group Name");
+  groupNameArr->SetNumberOfTuples(allGroups.size());
+  int blockId(0);
+  for(const auto& grp : allGroups)
+  {
+    blockIdArr->SetValue(blockId,blockId);
+    groupNameArr->SetValue(blockId,vtkStdString(grp.first));
+    blockId++;
+  }
+  output->AddColumn(blockIdArr);
+  output->AddColumn(groupNameArr);
+  return 1;
+}
diff --git a/src/Plugins/MEDReader/plugin/MEDReaderIO/vtkGroupsNames.h b/src/Plugins/MEDReader/plugin/MEDReaderIO/vtkGroupsNames.h
new file mode 100644 (file)
index 0000000..14af366
--- /dev/null
@@ -0,0 +1,39 @@
+// Copyright (C) 2020  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 (EDF R&D)
+
+#pragma once
+
+#include "vtkTableAlgorithm.h"
+
+class ExtractGroupInternal;
+
+class vtkGroupsNames : public vtkTableAlgorithm
+{
+public:
+    static vtkGroupsNames *New();
+    vtkTypeMacro(vtkGroupsNames, vtkTableAlgorithm);
+protected:
+    vtkGroupsNames();
+    ~vtkGroupsNames();
+    int RequestInformation(vtkInformation *request, vtkInformationVector **inputVector, vtkInformationVector *outputVector) override;
+    int RequestData(vtkInformation*, vtkInformationVector**, vtkInformationVector*) override;
+private:
+    ExtractGroupInternal *Internal;
+};
index 44cdcba47ffa44d79278e513fb5924f0b7d115c5..9c4e031f3d469c03f646b89494162f782ce3a8da 100644 (file)
      </InputProperty>
     </SourceProxy>
 
+    <SourceProxy name="GroupsNames" class="vtkGroupsNames" label="Groups Names">
+      <InputProperty name="Input" command="SetInputConnection">
+       <ProxyGroupDomain name="groups">
+         <Group name="sources"/>
+         <Group name="filters"/>
+       </ProxyGroupDomain>
+       <DataTypeDomain name="input_type">
+         <DataType value="vtkDataSet"/>
+       </DataTypeDomain>
+       <Documentation>
+         This property specifies the input to the Level Scalars filter.
+       </Documentation>
+     </InputProperty>
+    </SourceProxy>
+
 
   </ProxyGroup>
 
index 92ee70f6d7192fcb91844a00959ccdb30ac93196..343c5dc5b53df92b23f96719ac10fa52cff7c7f6 100644 (file)
@@ -70,6 +70,16 @@ def test():
     reader.AllArrays = ['TS0/mesh/ComSup0/field@@][@@P0','TS0/mesh/ComSup0/field2@@][@@P1','TS0/mesh/ComSup0/mesh@@][@@P0']
     reader.AllTimeSteps = ['0000']
 
+    groupsNames = GroupsNames(Input=reader)
+    groupsNames.UpdatePipeline()
+    MyAssert( servermanager.Fetch(groupsNames).GetNumberOfBlocks() == 1 )
+    gn = servermanager.Fetch(groupsNames).GetBlock(0)
+    blockIDS = numpy_support.vtk_to_numpy(gn.GetColumnByName("Block ID"))
+    MyAssert(np.all(blockIDS==np.array([0,1,2,3],dtype=np.int32)))
+    grpNames = gn.GetColumnByName("Group Name")
+    MyAssert(grpNames.GetNumberOfTuples()==4)
+    MyAssert([grpNames.GetValue(i) for i in range(4)] == ['BottomLeft','BottomRight','TopLeft','TopRight'])
+    
     groupsAsMultiBlocks = GroupsAsMultiBlocks(Input=reader)
     groupsAsMultiBlocks.UpdatePipeline()
     blocks = servermanager.Fetch(groupsAsMultiBlocks)