]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
Fix of bug TC6.2.0: SIGSEGV at presentation creation for Iso Surfaces, Deformed Shape...
authorouv <ouv@opencascade.com>
Fri, 17 Dec 2010 15:31:41 +0000 (15:31 +0000)
committerouv <ouv@opencascade.com>
Fri, 17 Dec 2010 15:31:41 +0000 (15:31 +0000)
15 files changed:
src/PIPELINE/Makefile.am
src/PIPELINE/VISU_CellDataToPointData.cxx [new file with mode: 0644]
src/PIPELINE/VISU_CellDataToPointData.hxx [new file with mode: 0644]
src/PIPELINE/VISU_DeformationPL.cxx
src/PIPELINE/VISU_DeformationPL.hxx
src/PIPELINE/VISU_DeformedShapeAndScalarMapPL.cxx
src/PIPELINE/VISU_DeformedShapeAndScalarMapPL.hxx
src/PIPELINE/VISU_DeformedShapePL.cxx
src/PIPELINE/VISU_DeformedShapePL.hxx
src/PIPELINE/VISU_IsoSurfacesPL.cxx
src/PIPELINE/VISU_IsoSurfacesPL.hxx
src/PIPELINE/VISU_OptionalDeformationPL.cxx
src/PIPELINE/VISU_PipeLineUtils.hxx
src/PIPELINE/VISU_Plot3DPL.cxx
src/PIPELINE/VISU_Plot3DPL.hxx

index ee6aa5b63a84766b328bac8efacd5d2645a5a39b..6aec80a165e40b29d5fd245e4c998ff0b1760328 100644 (file)
@@ -71,7 +71,8 @@ salomeinclude_HEADERS= \
        VISU_ElnoAssembleFilter.hxx \
        VISU_DeformationPL.hxx \
        VISU_OptionalDeformationPL.hxx \
-       VISU_XYPlotActor.hxx
+       VISU_XYPlotActor.hxx \
+       VISU_CellDataToPointData.hxx
 
 dist_libVisuPipeLine_la_SOURCES= \
        VISU_MapperHolder.cxx \
@@ -118,7 +119,8 @@ dist_libVisuPipeLine_la_SOURCES= \
        VISU_ElnoAssembleFilter.cxx \
        VISU_DeformationPL.cxx \
        VISU_OptionalDeformationPL.cxx\
-       VISU_XYPlotActor.cxx
+       VISU_XYPlotActor.cxx \
+       VISU_CellDataToPointData.cxx
 
 libVisuPipeLine_la_CPPFLAGS= \
        $(VTK_INCLUDES) \
diff --git a/src/PIPELINE/VISU_CellDataToPointData.cxx b/src/PIPELINE/VISU_CellDataToPointData.cxx
new file mode 100644 (file)
index 0000000..46e22d3
--- /dev/null
@@ -0,0 +1,137 @@
+/*=========================================================================
+
+  Program:   Visualization Toolkit
+  Module:    $RCSfile$
+
+  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
+  All rights reserved.
+  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
+
+     This software is distributed WITHOUT ANY WARRANTY; without even
+     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
+     PURPOSE.  See the above copyright notice for more information.
+
+ Note from SALOME:
+ This file is a part of VTK library
+ It has been renamed and modified for SALOME project
+
+=========================================================================*/
+
+#include "VISU_CellDataToPointData.hxx"
+
+#include <vtkCellData.h>
+#include <vtkDataSet.h>
+#include <vtkIdList.h>
+#include <vtkInformation.h>
+#include <vtkInformationVector.h>
+#include <vtkObjectFactory.h>
+#include <vtkPointData.h>
+
+vtkCxxRevisionMacro(VISU_CellDataToPointData, "$Revision$");
+vtkStandardNewMacro(VISU_CellDataToPointData);
+
+//----------------------------------------------------------------------------
+// Instantiate object so that cell data is not passed to output.
+VISU_CellDataToPointData::VISU_CellDataToPointData()
+{
+  this->PassCellData = 0;
+}
+
+#define VTK_MAX_CELLS_PER_POINT 4096
+
+//----------------------------------------------------------------------------
+int VISU_CellDataToPointData::RequestData(
+  vtkInformation*,
+  vtkInformationVector** inputVector,
+  vtkInformationVector* outputVector)
+{
+  vtkInformation* info = outputVector->GetInformationObject(0);
+  vtkDataSet *output = vtkDataSet::SafeDownCast(
+    info->Get(vtkDataObject::DATA_OBJECT()));
+
+  vtkInformation* inInfo = inputVector[0]->GetInformationObject(0);
+  vtkDataSet *input = vtkDataSet::SafeDownCast(
+    inInfo->Get(vtkDataObject::DATA_OBJECT()));
+
+  vtkIdType cellId, ptId;
+  vtkIdType numCells, numPts;
+  vtkCellData *inPD=input->GetCellData();
+  vtkPointData *outPD=output->GetPointData();
+  vtkIdList *cellIds;
+  double weight;
+  double *weights;
+
+  vtkDebugMacro(<<"Mapping cell data to point data");
+
+  // First, copy the input to the output as a starting point
+  output->CopyStructure( input );
+
+  cellIds = vtkIdList::New();
+  cellIds->Allocate(VTK_MAX_CELLS_PER_POINT);
+
+  if ( (numPts=input->GetNumberOfPoints()) < 1 )
+    {
+    vtkDebugMacro(<<"No input point data!");
+    cellIds->Delete();
+    return 1;
+    }
+  weights = new double[VTK_MAX_CELLS_PER_POINT];
+  
+  // Pass the point data first. The fields and attributes
+  // which also exist in the cell data of the input will
+  // be over-written during CopyAllocate
+  output->GetPointData()->CopyGlobalIdsOff();
+  output->GetPointData()->PassData(input->GetPointData());
+  output->GetPointData()->CopyFieldOff("vtkGhostLevels");
+
+  // notice that inPD and outPD are vtkCellData and vtkPointData; respectively.
+  // It's weird, but it works.
+  outPD->InterpolateAllocate(inPD,numPts);
+
+  int abort=0;
+  vtkIdType progressInterval=numPts/20 + 1;
+  for (ptId=0; ptId < numPts && !abort; ptId++)
+    {
+    if ( !(ptId % progressInterval) )
+      {
+      this->UpdateProgress(static_cast<double>(ptId)/numPts);
+      abort = GetAbortExecute();
+      }
+
+    input->GetPointCells(ptId, cellIds);
+    numCells = cellIds->GetNumberOfIds();
+    if ( numCells > 0 )
+      {
+      weight = 1.0 / numCells;
+      for (cellId=0; cellId < numCells; cellId++)
+        {
+        weights[cellId] = weight;
+        }
+      outPD->InterpolatePoint(inPD, ptId, cellIds, weights);
+      }
+    else
+      {
+      outPD->NullPoint(ptId);
+      }
+    }
+
+  if ( !this->PassCellData )
+    {
+    output->GetCellData()->CopyAllOff();
+    output->GetCellData()->CopyFieldOn("vtkGhostLevels");
+    }
+  output->GetCellData()->PassData(input->GetCellData());
+
+  cellIds->Delete();
+  delete [] weights;
+  
+  return 1;
+}
+
+//----------------------------------------------------------------------------
+void VISU_CellDataToPointData::PrintSelf(ostream& os, vtkIndent indent)
+{
+  this->Superclass::PrintSelf(os,indent);
+
+  os << indent << "Pass Cell Data: " << (this->PassCellData ? "On\n" : "Off\n");
+}
diff --git a/src/PIPELINE/VISU_CellDataToPointData.hxx b/src/PIPELINE/VISU_CellDataToPointData.hxx
new file mode 100644 (file)
index 0000000..ed4fb38
--- /dev/null
@@ -0,0 +1,74 @@
+/*=========================================================================
+
+  Program:   Visualization Toolkit
+  Module:    $RCSfile$
+
+  Copyright (c) Ken Martin, Will Schroeder, Bill Lorensen
+  All rights reserved.
+  See Copyright.txt or http://www.kitware.com/Copyright.htm for details.
+
+     This software is distributed WITHOUT ANY WARRANTY; without even
+     the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR
+     PURPOSE.  See the above copyright notice for more information.
+
+=========================================================================*/
+// .NAME VISU_CellDataToPointData - map cell data to point data
+// .SECTION Description
+// VISU_CellDataToPointData is a filter that transforms cell data (i.e., data
+// specified per cell) into point data (i.e., data specified at cell
+// points). The method of transformation is based on averaging the data
+// values of all cells using a particular point. Optionally, the input cell
+// data can be passed through to the output as well.
+
+// .SECTION Caveats
+// This filter is an abstract filter, that is, the output is an abstract type
+// (i.e., vtkDataSet). Use the convenience methods (e.g.,
+// GetPolyDataOutput(), GetStructuredPointsOutput(), etc.) to get the type
+// of output you want.
+
+// .SECTION See Also
+// vtkPointData vtkCellData vtkPointDataToCellData
+
+// note from SALOME:
+// This file is a part of VTK library
+// It has been renamed and modified for SALOME project
+
+#ifndef __VISU_CellDataToPointData_h
+#define __VISU_CellDataToPointData_h
+
+#include <vtkDataSetAlgorithm.h>
+
+class vtkDataSet;
+
+class VTK_GRAPHICS_EXPORT VISU_CellDataToPointData : public vtkDataSetAlgorithm
+{
+public:
+  static VISU_CellDataToPointData *New();
+  vtkTypeRevisionMacro(VISU_CellDataToPointData,vtkDataSetAlgorithm);
+  void PrintSelf(ostream& os, vtkIndent indent);
+
+  // Description:
+  // Control whether the input cell data is to be passed to the output. If
+  // on, then the input cell data is passed through to the output; otherwise,
+  // only generated point data is placed into the output.
+  vtkSetMacro(PassCellData,int);
+  vtkGetMacro(PassCellData,int);
+  vtkBooleanMacro(PassCellData,int);
+
+protected:
+  VISU_CellDataToPointData();
+  ~VISU_CellDataToPointData() {};
+
+  virtual int RequestData(vtkInformation* request,
+                          vtkInformationVector** inputVector,
+                          vtkInformationVector* outputVector);
+
+  int PassCellData;
+private:
+  VISU_CellDataToPointData(const VISU_CellDataToPointData&);  // Not implemented.
+  void operator=(const VISU_CellDataToPointData&);  // Not implemented.
+};
+
+#endif
+
+
index 0fbaa6c00501f24378376fef955df105fde77254..f55e795f702b3fddbc245fbb3ae35ba681e91c99 100755 (executable)
@@ -31,7 +31,6 @@
 #include <vtkDataSet.h>
 #include <vtkPassThroughFilter.h>
 #include <vtkWarpVector.h>
-#include <vtkCellDataToPointData.h>
 #include <vtkUnstructuredGrid.h>
 #ifdef _DEBUG_
 static int MYDEBUG = 0;
@@ -53,7 +52,7 @@ VISU_DeformationPL::VISU_DeformationPL():
   myVectorMergeFilter->SetMergingInputs(true);
   myInputPassFilter = vtkPassThroughFilter::New();
   myOutputPassFiler = vtkPassThroughFilter::New();
-  myCellDataToPointData = vtkCellDataToPointData::New();
+  myCellDataToPointData = VISU_CellDataToPointData::New();
   myCellDataToPointData->PassCellDataOn();
 
   myInputPassFilter->SetInput(vtkUnstructuredGrid::New());
index e8d665761af46cb324db741089f61da4844383dc..69c75c905a4b4c2bcc5469c2ece18ff80e2084b2 100755 (executable)
@@ -32,7 +32,7 @@ class vtkDataSet;
 class VISU_MergeFilter;
 class vtkPassThroughFilter;
 class vtkWarpVector;
-class vtkCellDataToPointData;
+class VISU_CellDataToPointData;
 
 class VISU_PIPELINE_EXPORT VISU_DeformationPL {
   
@@ -74,7 +74,7 @@ protected:
   vtkSmartPointer<VISU_MergeFilter> myVectorMergeFilter;
   vtkPassThroughFilter *myInputPassFilter;
   vtkPassThroughFilter *myOutputPassFiler;
-  vtkCellDataToPointData *myCellDataToPointData;
+  VISU_CellDataToPointData *myCellDataToPointData;
 
 private:
   vtkFloatingPointType myScaleFactor;
index 4f2891c9d971617339dd842f36303742235228a5..f533608dd41de8d9cceaefa5508e86bc8ea128e2 100644 (file)
@@ -39,7 +39,6 @@
 #include <vtkImplicitBoolean.h>
 #include <vtkImplicitFunction.h>
 #include <vtkUnstructuredGrid.h>
-#include <vtkCellDataToPointData.h>
 #include <vtkPointDataToCellData.h>
 #include <vtkImplicitFunctionCollection.h>
 
@@ -74,7 +73,7 @@ VISU_DeformedShapeAndScalarMapPL
 
   myScalarsFieldTransform = VISU_FieldTransform::New();
 
-  myCellDataToPointData = vtkCellDataToPointData::New();
+  myCellDataToPointData = VISU_CellDataToPointData::New();
   myScalarsElnoDisassembleFilter = VISU_ElnoDisassembleFilter::New();
 
   vtkImplicitBoolean* anImplicitBoolean = vtkImplicitBoolean::New();
index a59794db59bdda9cda4dfa7075b5c4735e5aacfd..0817c6a9f84c29d8da96e1f2b3cbd79285240835 100644 (file)
@@ -31,7 +31,7 @@
 class VISU_MergeFilter;
 class vtkWarpVector;
 class vtkUnstructuredGrid;
-class vtkCellDataToPointData;
+class VISU_CellDataToPointData;
 class vtkPointDataToCellData;
 class VISU_ElnoDisassembleFilter;
 class SALOME_ExtractGeometry;
@@ -168,7 +168,7 @@ private:
   vtkWarpVector  *myWarpVector;
   VISU_MergeFilter *myScalarsMergeFilter;
   vtkSmartPointer<vtkUnstructuredGrid> myScalars;
-  vtkCellDataToPointData* myCellDataToPointData;
+  VISU_CellDataToPointData* myCellDataToPointData;
   VISU_FieldTransform* myScalarsFieldTransform;
   VISU_Extractor* myScalarsExtractor;
   VISU_ElnoDisassembleFilter* myScalarsElnoDisassembleFilter;
index 3312723f25b0631385ff6ba0c227db855f70cd9c..f33dc1321c395c30cfb2fd85bbb230b0575e0f94 100644 (file)
@@ -46,7 +46,7 @@ VISU_DeformedShapePL
   SetIsFeatureEdgesAllowed(true);
 
   myWarpVector = vtkWarpVector::New();
-  myCellDataToPointData = vtkCellDataToPointData::New();
+  myCellDataToPointData = VISU_CellDataToPointData::New();
 }
 
 
index ccbeb5823fc31ddd44218590ffa556dcb2f44e54..1176b5541036f53f2d3eb63e18cf09dd5e3541d1 100644 (file)
@@ -31,7 +31,7 @@
 #include "VISUPipeline.hxx"
 #include "VISU_ScalarMapPL.hxx"
 
-class vtkCellDataToPointData;
+class VISU_CellDataToPointData;
 class SALOME_Transform;
 class vtkWarpVector;
 
@@ -102,7 +102,7 @@ protected:
   vtkFloatingPointType myScaleFactor;
   vtkFloatingPointType myMapScaleFactor;
   vtkWarpVector *myWarpVector;
-  vtkCellDataToPointData* myCellDataToPointData;
+  VISU_CellDataToPointData* myCellDataToPointData;
 
 private:
   VISU_DeformedShapePL(const VISU_DeformedShapePL&);  // Not implemented.
index ae875a055fcb461c2450734ab06ffacc7472e0d5..99ee33f086614c320d77a84d8ddb8fa838c57cb7 100644 (file)
@@ -51,7 +51,7 @@ VISU_IsoSurfacesPL
 
   myContourFilter = vtkContourFilter::New();
 
-  myCellDataToPointData = vtkCellDataToPointData::New();
+  myCellDataToPointData = VISU_CellDataToPointData::New();
 }
 
 
index 81651dd09e8469753049cc6e5d38b02c7bc5314f..195a985acc413e3ae6686288aae32789c5edd98c 100644 (file)
@@ -32,7 +32,7 @@
 #include "VISU_ScalarMapPL.hxx"
 
 class vtkContourFilter;
-class vtkCellDataToPointData;
+class VISU_CellDataToPointData;
 
 
 //----------------------------------------------------------------------------
@@ -120,7 +120,7 @@ protected:
   int myNbParts;
   vtkFloatingPointType myRange[2];
   bool myIsRangeFixed;
-  vtkCellDataToPointData* myCellDataToPointData;
+  VISU_CellDataToPointData* myCellDataToPointData;
   vtkContourFilter *myContourFilter;
 
 private:
index ae8527c3eba98f2c9f05106688cac97213aff7ff..58ddb5e2acea215557adba2bf440a23ef398e472 100755 (executable)
@@ -29,7 +29,6 @@
 #include <vtkDataSet.h>
 #include <vtkPassThroughFilter.h>
 #include <vtkWarpVector.h>
-#include <vtkCellDataToPointData.h>
 #ifdef _DEBUG_
 static int MYDEBUG = 0;
 #else
index c3aef9e3d457df56cac1c753de599f48dbdec316..ca2e2d542a0fc036869d82aacff4f98180d8f9aa 100644 (file)
 
 #include "VISUPipeline.hxx"
 #include "VISU_ConvertorUtils.hxx"
+#include "VISU_CellDataToPointData.hxx"
 
 #include <vtkProperty.h>
 #include <vtkObjectFactory.h>
 #include <vtkDataSetMapper.h>
 #include <vtkUnstructuredGrid.h>
 
-#include <vtkCellDataToPointData.h>
 #include <vtkPointData.h>
 #include <vtkCellData.h>
 #include <vtkPolyData.h>
@@ -76,7 +76,7 @@ namespace VISU
   template<class TOutputFilter> 
   void
   CellDataToPoint(TOutputFilter* theOutputFilter, 
-                  vtkCellDataToPointData *theCellDataToPointData,
+                  VISU_CellDataToPointData *theCellDataToPointData,
                   vtkDataSet* theDataSet)
 
   {
index 08846958fc615b054ab8bc1e72cf707bb3e2c495..dbb39bf0da3bb7959b8259abead31d6965f5ac4e 100644 (file)
@@ -33,7 +33,6 @@
 #include <vtkCutter.h>
 #include <vtkPlane.h>
 
-#include <vtkCellDataToPointData.h>
 #include <vtkGeometryFilter.h>
 #include <vtkContourFilter.h>
 #include <vtkWarpScalar.h>
@@ -47,7 +46,7 @@ vtkStandardNewMacro(VISU_Plot3DPL);
 //----------------------------------------------------------------------------
 VISU_Plot3DPL
 ::VISU_Plot3DPL():
-  myCellDataToPointData(vtkCellDataToPointData::New()),
+  myCellDataToPointData(VISU_CellDataToPointData::New()),
   myAppendPolyData(vtkAppendPolyData::New()),
   myGeometryFilter(vtkGeometryFilter::New()),
   myContourFilter(vtkContourFilter::New()),
index ab3d259462a30733ab651992b1fbef8a8a654cb0..cb98e963b389845b9e59b7f3bc0bce1d005ef7eb 100644 (file)
@@ -35,7 +35,7 @@
 class vtkWarpScalar;
 class vtkContourFilter;
 class vtkGeometryFilter;
-class vtkCellDataToPointData;
+class VISU_CellDataToPointData;
 
 
 //----------------------------------------------------------------------------
@@ -150,7 +150,7 @@ protected:
   vtkFloatingPointType myPosition, myScaleFactor, myMapScaleFactor;
   VISU_CutPlanesPL::PlaneOrientation myOrientation;
 
-  vtkSmartPointer<vtkCellDataToPointData> myCellDataToPointData;
+  vtkSmartPointer<VISU_CellDataToPointData> myCellDataToPointData;
   vtkSmartPointer<vtkAppendPolyData> myAppendPolyData;
   vtkSmartPointer<vtkGeometryFilter> myGeometryFilter;
   vtkSmartPointer<vtkContourFilter> myContourFilter;