]> SALOME platform Git repositories - modules/paravis.git/commitdiff
Salome HOME
Implementation of the '22884: EDF 9622 MED: Filter pipeline Gauss Points + ExtractGro... V7_6_0b1
authorrnv <rnv@opencascade.com>
Tue, 14 Apr 2015 12:49:55 +0000 (15:49 +0300)
committerrnv <rnv@opencascade.com>
Tue, 14 Apr 2015 12:49:55 +0000 (15:49 +0300)
src/Plugins/MEDReader/IO/CMakeLists.txt
src/Plugins/MEDReader/IO/vtkMEDQuadraturePointsGenerator.cxx [new file with mode: 0644]
src/Plugins/MEDReader/IO/vtkMEDQuadraturePointsGenerator.h [new file with mode: 0644]
src/Plugins/MEDReader/ParaViewPlugin/Resources/MEDReaderServer.xml

index 0ecbbbb9343f445d5515518110578c1c829b35c7..48865a678c22a193871d911d97db609b0b35063a 100644 (file)
@@ -27,7 +27,7 @@ IF(HDF5_IS_PARALLEL)
   ADD_DEFINITIONS("-DMEDREADER_USE_MPI")
 ENDIF(HDF5_IS_PARALLEL)
 
-SET(MEDReader_CLASSES vtkMEDReader vtkExtractGroup vtkELNOMeshFilter vtkELNOSurfaceFilter vtkELNOFilter vtkExtractCellType)
+SET(MEDReader_CLASSES vtkMEDReader vtkExtractGroup vtkELNOMeshFilter vtkELNOSurfaceFilter vtkELNOFilter vtkExtractCellType vtkMEDQuadraturePointsGenerator)
 
 SET(MEDReader_SRCS)
 SET(MEDReader_HDRS)
diff --git a/src/Plugins/MEDReader/IO/vtkMEDQuadraturePointsGenerator.cxx b/src/Plugins/MEDReader/IO/vtkMEDQuadraturePointsGenerator.cxx
new file mode 100644 (file)
index 0000000..452123b
--- /dev/null
@@ -0,0 +1,134 @@
+// Copyright (C) 2010-2015  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 : Roman NIKOLAEV
+
+//Local includes
+#include "vtkMEDQuadraturePointsGenerator.h"
+#include "MEDFileFieldRepresentationTree.hxx"
+
+//VTK includes
+#include <vtkObjectFactory.h>
+#include <vtkInformation.h>
+#include <vtkInformationVector.h>
+#include <vtkUnstructuredGrid.h>
+#include <vtkCellData.h>
+#include <vtkPointData.h>
+#include <vtkInformationQuadratureSchemeDefinitionVectorKey.h>
+#include <vtkQuadratureSchemeDefinition.h>
+
+//-----------------------------------------------------------------------------
+vtkStandardNewMacro(vtkMEDQuadraturePointsGenerator);
+
+//-----------------------------------------------------------------------------
+vtkMEDQuadraturePointsGenerator::vtkMEDQuadraturePointsGenerator()
+{
+}
+
+//-----------------------------------------------------------------------------
+vtkMEDQuadraturePointsGenerator::~vtkMEDQuadraturePointsGenerator()
+{}
+
+
+//-----------------------------------------------------------------------------
+int vtkMEDQuadraturePointsGenerator::RequestData(
+      vtkInformation* request,
+      vtkInformationVector **input,
+      vtkInformationVector *output)
+{
+  if (this->Superclass::RequestData(request, input, output) == 0 )
+    {
+      return 0;
+    }    
+
+  //Fill MED internal array
+  vtkDataObject *tmpDataObj;
+  
+  // Get the input.
+  tmpDataObj = input[0]->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
+  vtkUnstructuredGrid *usgIn = vtkUnstructuredGrid::SafeDownCast(tmpDataObj);
+  // Get the output.
+  tmpDataObj = output->GetInformationObject(0)->Get(vtkDataObject::DATA_OBJECT());
+  vtkPolyData *pdOut = vtkPolyData::SafeDownCast(tmpDataObj);
+  vtkDataArray* offsets = this->GetInputArrayToProcess(0, input);
+
+  if (usgIn == NULL || pdOut == NULL || offsets == NULL )
+    {
+    vtkErrorMacro("Filter data has not been configured correctly. Aborting.");
+    return 1;
+    }
+  
+  vtkInformation *info = offsets->GetInformation();
+  vtkInformationQuadratureSchemeDefinitionVectorKey *key 
+    = vtkQuadratureSchemeDefinition::DICTIONARY();
+  if (!key->Has(info))
+    {
+    vtkErrorMacro(
+    << "Dictionary is not present in array "
+    << offsets->GetName() << " " << offsets
+    << " Aborting.");
+    return 0;
+    }
+
+  vtkIdType nCells = usgIn->GetNumberOfCells();
+  int dictSize = key->Size(info);
+  vtkQuadratureSchemeDefinition **dict
+    = new vtkQuadratureSchemeDefinition *[dictSize];
+  key->GetRange(info, dict, 0, 0, dictSize);  
+
+  // Loop over all fields to map the internal MED cell array to the points array
+  int nCArrays = usgIn->GetCellData()->GetNumberOfArrays();
+  for (int i = 0; i<nCArrays; ++i)
+    {
+    vtkDataArray* array = usgIn->GetCellData()->GetArray(i);
+    if ( !array )
+      {
+      continue;
+      }
+    std::string arrName = array->GetName();
+    if ( arrName == MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_CELL_NAME ) 
+      {        
+        vtkDataArray *out_id_cells = array->NewInstance();
+        out_id_cells->SetName(MEDFileFieldRepresentationLeavesArrays::FAMILY_ID_NODE_NAME);
+        out_id_cells->SetNumberOfComponents(array->GetNumberOfComponents());
+        out_id_cells->CopyComponentNames( array );
+        for (int cellId = 0; cellId < nCells; cellId++)
+          {
+          int cellType = usgIn->GetCellType(cellId);
+
+          // a simple check to see if a scheme really exists for this cell type.
+          // should not happen if the cell type has not been modified.
+          if (dict[cellType] == NULL)
+            {
+            continue;
+            }
+
+          int np = dict[cellType]->GetNumberOfQuadraturePoints();
+          for (int id = 0; id < np; id++)
+            {
+            out_id_cells->InsertNextTuple(cellId, array);
+            }
+          }
+        out_id_cells->Squeeze();
+        pdOut->GetPointData()->AddArray(out_id_cells);
+        out_id_cells->Delete();
+      }          
+    }
+  delete[] dict;  
+  return 1;
+}
diff --git a/src/Plugins/MEDReader/IO/vtkMEDQuadraturePointsGenerator.h b/src/Plugins/MEDReader/IO/vtkMEDQuadraturePointsGenerator.h
new file mode 100644 (file)
index 0000000..9509831
--- /dev/null
@@ -0,0 +1,46 @@
+// Copyright (C) 2010-2015  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 : Roman NIKOLAEV
+
+#ifndef vtkMEDQuadraturePointsGenerator_h
+#define vtkMEDQuadraturePointsGenerator_h
+
+#include <vtkFiltersGeneralModule.h> // For export macro
+#include <vtkQuadraturePointsGenerator.h>
+
+class vtkInformation;
+class vtkInformationVector;
+
+class VTKFILTERSGENERAL_EXPORT vtkMEDQuadraturePointsGenerator : public vtkQuadraturePointsGenerator
+{
+public:
+  vtkTypeMacro(vtkMEDQuadraturePointsGenerator,vtkQuadraturePointsGenerator);
+  static vtkMEDQuadraturePointsGenerator *New();
+
+protected:
+
+  int RequestData(vtkInformation *req, vtkInformationVector **input, vtkInformationVector *output);
+  vtkMEDQuadraturePointsGenerator();
+  virtual ~vtkMEDQuadraturePointsGenerator();
+private:
+  vtkMEDQuadraturePointsGenerator(const vtkMEDQuadraturePointsGenerator &); // Not implemented
+  void operator=(const vtkMEDQuadraturePointsGenerator &); // Not implemented
+};
+
+#endif
index 80d730b3f62b3f4bea7755cdcdf025310246d73f..64b3169f29b4aedfa6abee0bf9527826f6e280d2 100644 (file)
       </Hints>
     </SourceProxy>
 
-    <SourceProxy name="GaussPoints" class="vtkQuadraturePointsGenerator" label="Gauss Points">
+    <SourceProxy name="GaussPoints" class="vtkMEDQuadraturePointsGenerator" label="Gauss Points">
       <Documentation
         long_help="Create a point set with data at Gauss points."
         short_help="Create a point set with data at Gauss points.">