]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
IPAL15171: IOLS. Error SIGSEGV on creating Stream Lines on maill.3.med.
authorjfa <jfa@opencascade.com>
Tue, 15 May 2007 11:59:39 +0000 (11:59 +0000)
committerjfa <jfa@opencascade.com>
Tue, 15 May 2007 11:59:39 +0000 (11:59 +0000)
src/PIPELINE/VISUPipeline.hxx
src/PIPELINE/VISU_UsedPointsFilter.cxx
src/PIPELINE/VISU_UsedPointsFilter.hxx

index 897c68a250813437cbdd9f5618daf2ff889f8b89..334195fa1a0439f023892915839929d4c0382a95 100755 (executable)
@@ -44,4 +44,4 @@
  #define VISU_PIPELINE_EXPORT
 #endif
 
-#endif
\ No newline at end of file
+#endif
index 6fb0ffda5ea9222a61fe1bd3c01fb699f660468e..b88d18eda8d3cfd700e0fe6ffb531a8607c93387 100644 (file)
@@ -17,7 +17,7 @@
 //  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
+//  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 //
 // File:    VISU_StreamLinesPL.cxx
@@ -27,6 +27,8 @@
 
 #include "VISU_UsedPointsFilter.hxx"
 
+#include "VISU_ConvertorUtils.hxx"
+
 #include <vtkObjectFactory.h>
 #include <vtkUnstructuredGrid.h>
 #include <vtkPointData.h>
 #include <vtkPoints.h>
 #include <vtkIdList.h>
 
+#include <vtkInformation.h>
+#include <vtkInformationVector.h>
+
 
 //----------------------------------------------------------------------------
 vtkStandardNewMacro(VISU_UsedPointsFilter);
 
 
 //----------------------------------------------------------------------------
-VISU_UsedPointsFilter
-::VISU_UsedPointsFilter()
+VISU_UsedPointsFilter::VISU_UsedPointsFilter()
 {
   PercentsOfUsedPoints = 1.0;
 }
 
-
 //----------------------------------------------------------------------------
-VISU_UsedPointsFilter
-::~VISU_UsedPointsFilter()
+VISU_UsedPointsFilter::~VISU_UsedPointsFilter()
 {}
 
 
 //----------------------------------------------------------------------------
-void
-VISU_UsedPointsFilter
-::Execute()
+int VISU_UsedPointsFilter::RequestData (vtkInformation *vtkNotUsed(request),
+                                        vtkInformationVector **inputVector,
+                                        vtkInformationVector *outputVector)
 {
-  vtkDataSet *anInput = this->GetInput();
-  vtkPointSet *anOutput = this->GetOutput();
-  anOutput->GetPointData()->CopyAllOff();
-  anOutput->GetCellData()->CopyAllOff();
-  anOutput->CopyStructure(anInput);
+  // get the info objects
+  vtkInformation *inInfo = inputVector[0]->GetInformationObject(0);
+  vtkInformation *outInfo = outputVector->GetInformationObject(0);
+
+  // get the input and ouptut
+  vtkDataSet *input = vtkDataSet::SafeDownCast(inInfo->Get(vtkDataObject::DATA_OBJECT()));
+  vtkPointSet *output = vtkPointSet::SafeDownCast(outInfo->Get(vtkDataObject::DATA_OBJECT()));
+
+  //////////
+  output->GetPointData()->CopyAllOff();
+  output->GetCellData()->CopyAllOff();
+  output->CopyStructure(input);
 
   vtkPoints* aPoints = vtkPoints::New();
   vtkIdList *anIdList = vtkIdList::New();
-  vtkIdType iEnd = anInput->GetNumberOfPoints();
-  for(vtkIdType i = 0; i < iEnd; i++){
-    anInput->GetPointCells(i,anIdList);
-    if(anIdList->GetNumberOfIds() > 0)
-      aPoints->InsertNextPoint(anInput->GetPoint(i));
+  vtkIdType iEnd = input->GetNumberOfPoints();
+  for (vtkIdType i = 0; i < iEnd; i++) {
+    input->GetPointCells(i, anIdList);
+    if (anIdList->GetNumberOfIds() > 0)
+      aPoints->InsertNextPoint(input->GetPoint(i));
   }
   vtkPoints* aNewPoints = vtkPoints::New();
   iEnd = aPoints->GetNumberOfPoints();
-  if (PercentsOfUsedPoints > 0){
+  if (PercentsOfUsedPoints > 0) {
     vtkIdType anOffset = vtkIdType(1.0/PercentsOfUsedPoints);
-    if(anOffset < 1) anOffset = 1;
-    for(vtkIdType i = 0; i < iEnd; i += anOffset)
+    if (anOffset < 1) anOffset = 1;
+    for (vtkIdType i = 0; i < iEnd; i += anOffset)
       aNewPoints->InsertNextPoint(aPoints->GetPoint(i));
   }
-  anOutput->SetPoints(aNewPoints);
+  output->SetPoints(aNewPoints);
   aNewPoints->Delete();
   aPoints->Delete();
-}
-
 
-//----------------------------------------------------------------------------
+  return 1;
+}
index 3795da85a2f450a25a64fd54864f0f5ce1f64736..380e0fc9de1472b1697b548e702e97cd3d0472a1 100644 (file)
@@ -17,7 +17,7 @@
 //  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
+//  See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
 //
 //
 // File:    VISU_UsedPointsFilter.hxx
 #ifndef VISU_UsedPointsFilter_HeaderFile
 #define VISU_UsedPointsFilter_HeaderFile
 
-#include <vtkDataSetToUnstructuredGridFilter.h>
+#include <vtkPointSetAlgorithm.h>
 
 
 //----------------------------------------------------------------------------
-class VISU_UsedPointsFilter : public vtkDataSetToUnstructuredGridFilter
+class VISU_UsedPointsFilter : public vtkPointSetAlgorithm
 {
 public:
-  vtkTypeMacro(VISU_UsedPointsFilter, vtkDataSetToUnstructuredGridFilter);
+  vtkTypeMacro(VISU_UsedPointsFilter, vtkPointSetAlgorithm);
 
-  static 
-  VISU_UsedPointsFilter* 
+  static
+  VISU_UsedPointsFilter*
   New();
 
   vtkSetMacro(PercentsOfUsedPoints,float);
@@ -50,9 +50,7 @@ protected:
   virtual
   ~VISU_UsedPointsFilter();
 
-  virtual
-  void
-  Execute();
+  int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *);
 
   float PercentsOfUsedPoints;
 };