From: apo Date: Fri, 2 May 2008 11:04:57 +0000 (+0000) Subject: To introduce ElnoAssembleFilter which can be used as terminal point for the ELNO... X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=e34d15600df70ecc635f771a43e6f6231502b8f3;p=modules%2Fvisu.git To introduce ElnoAssembleFilter which can be used as terminal point for the ELNO pipeline that is started by ElnoDisassembleFilter --- diff --git a/src/PIPELINE/Makefile.am b/src/PIPELINE/Makefile.am index 7b51bbe4..ef891d78 100644 --- a/src/PIPELINE/Makefile.am +++ b/src/PIPELINE/Makefile.am @@ -106,6 +106,7 @@ dist_libVisuPipeLine_la_SOURCES= \ VISU_OpenGLElnoMapper.cxx \ VISU_ElnoGeometryFilter.cxx \ VISU_ElnoDisassembleFilter.cxx \ + VISU_ElnoAssembleFilter.cxx \ VISU_ElnoExtractScalars.cxx \ VISU_ElnoMapperHolder.cxx \ VISU_ElnoScalarMapPL.cxx \ diff --git a/src/PIPELINE/VISUPipeLine.cxx b/src/PIPELINE/VISUPipeLine.cxx index 692c7e74..8ab6e28b 100644 --- a/src/PIPELINE/VISUPipeLine.cxx +++ b/src/PIPELINE/VISUPipeLine.cxx @@ -49,6 +49,7 @@ #include "VISU_OpenGLElnoMapper.hxx" #include "VISU_ElnoDisassembleFilter.hxx" +#include "VISU_ElnoAssembleFilter.hxx" #include "VISU_ElnoExtractScalars.hxx" #include @@ -253,6 +254,7 @@ main(int argc, char** argv) VISU_ElnoDisassembleFilter* aDisassembleFilter = VISU_ElnoDisassembleFilter::New(); aDisassembleFilter->SetInput( anUnstructuredGrid ); + aDisassembleFilter->SetShrinkFactor( 0.999 ); VISU::WriteToFile( aDisassembleFilter->GetOutput(), "/data/apo/elno_from_disassemble.vtk" ); VISU_ElnoExtractScalars* anExtractScalars = VISU_ElnoExtractScalars::New(); @@ -284,6 +286,10 @@ main(int argc, char** argv) aPlane->Delete(); VISU::WriteToFile( aCutter->GetOutput(), "/data/apo/elno_from_cutter.vtk" ); + VISU_ElnoAssembleFilter* anAssembleFilter = VISU_ElnoAssembleFilter::New(); + anAssembleFilter->SetInput( aCutter->GetOutput() ); + VISU::WriteToFile( anAssembleFilter->GetPolyDataOutput(), "/data/apo/elno_from_assemble.vtk" ); + //VISU_ElnoWarpVector* anElnoWarpVector = VISU_ElnoWarpVector::New(); //anElnoWarpVector->SetInput( anUnstructuredGridIDMapper->GetUnstructuredGridOutput() ); @@ -310,7 +316,8 @@ main(int argc, char** argv) aMapper->SetScalarModeToUsePointData(); //aMapper->SetInput( aGeometryFilter->GetOutput() ); - aMapper->SetInput( aCutter->GetOutput() ); + //aMapper->SetInput( aCutter->GetOutput() ); + aMapper->SetInput( anAssembleFilter->GetOutput() ); VISU_LookupTable* aMapperTable( VISU_LookupTable::New() ); aMapperTable->SetHueRange( 0.667, 0.0 ); diff --git a/src/PIPELINE/VISU_ElnoAssembleFilter.cxx b/src/PIPELINE/VISU_ElnoAssembleFilter.cxx new file mode 100644 index 00000000..6351d772 --- /dev/null +++ b/src/PIPELINE/VISU_ElnoAssembleFilter.cxx @@ -0,0 +1,154 @@ +// Copyright (C) 2003 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, +// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS +// +// 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. +// +// 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 + +#include "VISU_ElnoAssembleFilter.hxx" +#include "VISU_PipeLineUtils.hxx" +#include "VISU_ElnoMeshValue.hxx" + +#include +#include +#include +#include +#include +#include +#include +#include + + +//---------------------------------------------------------------------------- +vtkStandardNewMacro( VISU_ElnoAssembleFilter ); + + +//---------------------------------------------------------------------------- +VISU_ElnoAssembleFilter::VISU_ElnoAssembleFilter() +{ + this->SetInputArrayToProcess( 0, // idx + 0, // port + 0, // connection + vtkDataObject::FIELD_ASSOCIATION_POINTS, // field association + "ELNO_POINT_COORDS" ); // name +} + + +//---------------------------------------------------------------------------- +VISU_ElnoAssembleFilter::~VISU_ElnoAssembleFilter() +{} + + +//---------------------------------------------------------------------------- +namespace +{ + //---------------------------------------------------------------------------- + template < int points_type, int elno_type > + int Execute2( vtkPointSet *theInput, + vtkPointSet *theOutput, + vtkDataArray *theElnoPointCoords ) + { + theOutput->CopyStructure( theInput ); + + vtkCellData *aCellData = theOutput->GetCellData(); + aCellData->PassData( theInput->GetCellData() ); + + vtkPointData *aPointData = theOutput->GetPointData(); + aPointData->PassData( theInput->GetPointData() ); + + vtkPoints *anInputPoints = theInput->GetPoints(); + vtkPoints *aPoints = anInputPoints->New( elno_type ); + vtkIdType aNbPoints = theInput->GetNumberOfPoints(); + aPoints->SetNumberOfPoints( aNbPoints ); + + typedef typename VISU::TL::TEnum2VTKArrayType< elno_type >::TResult TPointsDataArray; + typedef typename VISU::TL::TEnum2VTKBasicType< elno_type >::TResult TPointsDataType; + TPointsDataArray* anOutputPointsArray = TPointsDataArray::SafeDownCast( aPoints->GetData() ); + + TPointsDataArray* anElnoPointCoords = TPointsDataArray::SafeDownCast( theElnoPointCoords ); + + for ( vtkIdType aPointId = 0; aPointId < aNbPoints; aPointId++ ) { + TPointsDataType aCoords[ 3 ]; + anElnoPointCoords->GetTupleValue( aPointId, aCoords ); + anOutputPointsArray->SetTupleValue( aPointId, aCoords ); + } + + theOutput->SetPoints( aPoints ); + + return 1; + } + + + //---------------------------------------------------------------------------- + template < int points_type > + int Execute( vtkPointSet *theInput, + vtkPointSet *theOutput, + vtkDataArray *theElnoPointCoords ) + { + switch( theElnoPointCoords->GetDataType() ){ + case VTK_DOUBLE: + return Execute2< points_type, VTK_DOUBLE >( theInput, theOutput, theElnoPointCoords ); + case VTK_FLOAT: + return Execute2< points_type, VTK_FLOAT >( theInput, theOutput, theElnoPointCoords ); + case VTK_INT: + return Execute2< points_type, VTK_INT >( theInput, theOutput, theElnoPointCoords ); + case VTK_LONG: + return Execute2< points_type, VTK_LONG >( theInput, theOutput, theElnoPointCoords ); + default: + break; + } + + return 0; + } + + + //---------------------------------------------------------------------------- +} + + +//---------------------------------------------------------------------------- +int VISU_ElnoAssembleFilter::RequestData( vtkInformation *vtkNotUsed(request), + vtkInformationVector **inputVector, + vtkInformationVector *outputVector ) +{ + // get the info objects + vtkInformation *inInfo = inputVector[0]->GetInformationObject(0); + vtkInformation *outInfo = outputVector->GetInformationObject(0); + + // get the input and ouptut + vtkPointSet *anInput = vtkPointSet::SafeDownCast( inInfo->Get( vtkDataObject::DATA_OBJECT() ) ); + vtkPointSet *anOutput = vtkPointSet::SafeDownCast( outInfo->Get( vtkDataObject::DATA_OBJECT() ) ); + + vtkDataArray *anElnoPointCoords = this->GetInputArrayToProcess( 0, inputVector ); + + vtkPoints *aPoints = anInput->GetPoints(); + switch( aPoints->GetDataType() ){ + case VTK_DOUBLE: + return ::Execute< VTK_DOUBLE >( anInput, anOutput, anElnoPointCoords ); + case VTK_FLOAT: + return ::Execute< VTK_FLOAT >( anInput, anOutput, anElnoPointCoords ); + case VTK_INT: + return ::Execute< VTK_INT >( anInput, anOutput, anElnoPointCoords ); + case VTK_LONG: + return ::Execute< VTK_LONG >( anInput, anOutput, anElnoPointCoords ); + default: + break; + } + + return 0; +} + + +//---------------------------------------------------------------------------- diff --git a/src/PIPELINE/VISU_ElnoAssembleFilter.hxx b/src/PIPELINE/VISU_ElnoAssembleFilter.hxx new file mode 100644 index 00000000..4aca5d5b --- /dev/null +++ b/src/PIPELINE/VISU_ElnoAssembleFilter.hxx @@ -0,0 +1,43 @@ +// Copyright (C) 2003 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, +// CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS +// +// 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. +// +// 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 + +#ifndef VISU_ElnoAssembleFilter_H +#define VISU_ElnoAssembleFilter_H + +#include + +class VISU_ElnoAssembleFilter : public vtkPointSetAlgorithm +{ +public: + typedef vtkPointSetAlgorithm Superclass; + + static VISU_ElnoAssembleFilter *New(); + +protected: + VISU_ElnoAssembleFilter(); + ~VISU_ElnoAssembleFilter(); + + int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *); + +private: + VISU_ElnoAssembleFilter(const VISU_ElnoAssembleFilter&); // Not implemented. + void operator=(const VISU_ElnoAssembleFilter&); // Not implemented. +}; + +#endif diff --git a/src/PIPELINE/VISU_ElnoDisassembleFilter.cxx b/src/PIPELINE/VISU_ElnoDisassembleFilter.cxx index e7ecfedf..b73b5214 100644 --- a/src/PIPELINE/VISU_ElnoDisassembleFilter.cxx +++ b/src/PIPELINE/VISU_ElnoDisassembleFilter.cxx @@ -49,6 +49,8 @@ VISU_ElnoDisassembleFilter::VISU_ElnoDisassembleFilter() 0, // connection vtkDataObject::FIELD_ASSOCIATION_CELLS, // field association "ELNO_COMPONENT_MAPPER" ); // name + + this->myShrinkFactor = -1.0; } @@ -57,58 +59,87 @@ VISU_ElnoDisassembleFilter::~VISU_ElnoDisassembleFilter() {} +//---------------------------------------------------------------------------- +void VISU_ElnoDisassembleFilter::SetShrinkFactor( vtkFloatingPointType theValue ) +{ + if ( VISU::CheckIsSameValue( theValue, myShrinkFactor ) ) + return; + + myShrinkFactor = theValue; + this->Modified(); +} + + +//---------------------------------------------------------------------------- +vtkFloatingPointType VISU_ElnoDisassembleFilter::GetShrinkFactor() +{ + return myShrinkFactor; +} + + //---------------------------------------------------------------------------- namespace { - template < int points_type, int elno_type > - int Execute2( vtkUnstructuredGrid *theInput, - vtkUnstructuredGrid *theOutput, - vtkDataArray *theElnoDataArray, - vtkDataArray *theElnoDataMapper ) + //---------------------------------------------------------------------------- + template < int points_type, class TPointsDataArray, int elno_type, class TElnoDataArray > + void SimpleExecute( vtkCellArray *theConnectivity, + vtkPointData *theInputPointData, + vtkPointData *theOutputPointData, + TPointsDataArray *theInputPointsArray, + TPointsDataArray *theOutputPointsArray, + VISU::TGetElnoNodeData< elno_type >& theGetElnoNodeData, + TElnoDataArray* theElnoPointDataArray, + TPointsDataArray *theElnoPointCoords ) { - VISU::TGetElnoNodeData< elno_type > aGetElnoNodeData( theElnoDataArray, theElnoDataMapper ); - - vtkCellArray *aConnectivity = vtkCellArray::New(); - aConnectivity->DeepCopy( theInput->GetCells() ); - - vtkPoints *anInputPoints = theInput->GetPoints(); - vtkPoints *aPoints = anInputPoints->New( anInputPoints->GetDataType() ); - vtkIdType aNbCells = aConnectivity->GetNumberOfCells(); - vtkIdType aNbPoints = aConnectivity->GetNumberOfConnectivityEntries() - aNbCells; - aPoints->Allocate( aNbPoints ); - - typedef typename VISU::TL::TEnum2VTKArrayType< points_type >::TResult TPointsDataArray; typedef typename VISU::TL::TEnum2VTKBasicType< points_type >::TResult TPointsDataType; - TPointsDataArray* anInputPointsArray = TPointsDataArray::SafeDownCast( anInputPoints->GetData() ); - TPointsDataArray* anOutputPointsArray = TPointsDataArray::SafeDownCast( aPoints->GetData() ); - - vtkPointData *anInputPointData = theInput->GetPointData(); - vtkPointData *aPointData = theOutput->GetPointData(); - aPointData->Allocate( aNbPoints ); - - typedef typename VISU::TL::TEnum2VTKArrayType< elno_type >::TResult TElnoDataArray; typedef typename VISU::TL::TEnum2VTKBasicType< elno_type >::TResult TElnoDataType; - TElnoDataArray* anElnoPointsDataArray = TElnoDataArray::New(); - anElnoPointsDataArray->SetNumberOfComponents( aGetElnoNodeData.getNbComp() ); - anElnoPointsDataArray->SetNumberOfTuples( aNbPoints ); - anElnoPointsDataArray->SetName( "ELNO_POINT_DATA" ); - TPointsDataArray* aPointCoords = TPointsDataArray::New(); - aPointCoords->SetNumberOfComponents( 3 ); - aPointCoords->SetNumberOfTuples( aNbPoints ); - aPointCoords->SetName( "ELNO_POINT_COORDS" ); + theConnectivity->InitTraversal(); + vtkIdType aNbPts = 0, *aPts = 0; + for ( vtkIdType aCellId = 0; theConnectivity->GetNextCell( aNbPts, aPts ); aCellId++ ) { + for ( vtkIdType aPntId = 0; aPntId < aNbPts; aPntId++ ) { + TPointsDataType aCoords[ 3 ]; + vtkIdType aCurrentPntId = aPts[ aPntId ]; + theInputPointsArray->GetTupleValue( aCurrentPntId, aCoords ); + + aPts[ aPntId ] = theOutputPointsArray->InsertNextTupleValue( aCoords ); + vtkIdType aNewPntId = aPts[ aPntId ]; - aConnectivity->InitTraversal(); + theElnoPointCoords->SetTupleValue( aNewPntId, aCoords ); + + theOutputPointData->CopyData( theInputPointData, aCurrentPntId, aNewPntId ); + + TElnoDataType* anElnoData = theGetElnoNodeData( aCellId, aPntId ); + theElnoPointDataArray->SetTupleValue( aNewPntId, anElnoData ); + } + } + } + + + //---------------------------------------------------------------------------- + template < int points_type, class TPointsDataArray, int elno_type, class TElnoDataArray > + void ShrunkExecute( vtkCellArray *theConnectivity, + vtkPointData *theInputPointData, + vtkPointData *theOutputPointData, + TPointsDataArray *theInputPointsArray, + TPointsDataArray *theOutputPointsArray, + VISU::TGetElnoNodeData< elno_type >& theGetElnoNodeData, + TElnoDataArray* theElnoPointDataArray, + TPointsDataArray *theElnoPointCoords, + vtkFloatingPointType theShrinkFactor ) + { + typedef typename VISU::TL::TEnum2VTKBasicType< points_type >::TResult TPointsDataType; + typedef typename VISU::TL::TEnum2VTKBasicType< elno_type >::TResult TElnoDataType; + + theConnectivity->InitTraversal(); vtkIdType aNbPts = 0, *aPts = 0; - for ( vtkIdType aCellId = 0; aConnectivity->GetNextCell( aNbPts, aPts ); aCellId++ ) { + for ( vtkIdType aCellId = 0; theConnectivity->GetNextCell( aNbPts, aPts ); aCellId++ ) { - TPointsDataType aCenter[ 3 ] = { TPointsDataType(), - TPointsDataType(), - TPointsDataType() }; + TPointsDataType aCenter[ 3 ] = { TPointsDataType(), TPointsDataType(), TPointsDataType() }; for ( vtkIdType aPntId = 0; aPntId < aNbPts; aPntId++ ) { TPointsDataType aCoords[ 3 ]; - anInputPointsArray->GetTupleValue( aPts[ aPntId ], aCoords ); + theInputPointsArray->GetTupleValue( aPts[ aPntId ], aCoords ); aCenter[ 0 ] += aCoords[ 0 ]; aCenter[ 1 ] += aCoords[ 1 ]; @@ -122,30 +153,92 @@ namespace for ( vtkIdType aPntId = 0; aPntId < aNbPts; aPntId++ ) { TPointsDataType aCoords[ 3 ]; vtkIdType aCurrentPntId = aPts[ aPntId ]; - anInputPointsArray->GetTupleValue( aCurrentPntId, aCoords ); + theInputPointsArray->GetTupleValue( aCurrentPntId, aCoords ); TPointsDataType aNewCoords[ 3 ]; - static vtkFloatingPointType SHRINK_FACTOR = 0.999; aNewCoords[ 0 ] = aCenter[ 0 ] + - TPointsDataType( SHRINK_FACTOR * ( aCoords[ 0 ] - aCenter[ 0 ] ) ); + TPointsDataType( theShrinkFactor * ( aCoords[ 0 ] - aCenter[ 0 ] ) ); aNewCoords[ 1 ] = aCenter[ 1 ] + - TPointsDataType( SHRINK_FACTOR * ( aCoords[ 1 ] - aCenter[ 1 ] ) ); + TPointsDataType( theShrinkFactor * ( aCoords[ 1 ] - aCenter[ 1 ] ) ); aNewCoords[ 2 ] = aCenter[ 2 ] + - TPointsDataType( SHRINK_FACTOR * ( aCoords[ 2 ] - aCenter[ 2 ] ) ); + TPointsDataType( theShrinkFactor * ( aCoords[ 2 ] - aCenter[ 2 ] ) ); - aPts[ aPntId ] = anOutputPointsArray->InsertNextTupleValue( aNewCoords ); + aPts[ aPntId ] = theOutputPointsArray->InsertNextTupleValue( aNewCoords ); vtkIdType aNewPntId = aPts[ aPntId ]; - aPointCoords->SetTupleValue( aNewPntId, aCoords ); + theElnoPointCoords->SetTupleValue( aNewPntId, aCoords ); - aPointData->CopyData( anInputPointData, aCurrentPntId, aNewPntId ); + theOutputPointData->CopyData( theInputPointData, aCurrentPntId, aNewPntId ); - TElnoDataType* anElnoData = aGetElnoNodeData( aCellId, aPntId ); - anElnoPointsDataArray->SetTupleValue( aNewPntId, anElnoData ); + TElnoDataType* anElnoData = theGetElnoNodeData( aCellId, aPntId ); + theElnoPointDataArray->SetTupleValue( aNewPntId, anElnoData ); } } + } + + + //---------------------------------------------------------------------------- + template < int points_type, int elno_type > + int Execute2( vtkUnstructuredGrid *theInput, + vtkUnstructuredGrid *theOutput, + vtkDataArray *theElnoDataArray, + vtkDataArray *theElnoDataMapper, + vtkFloatingPointType theShrinkFactor ) + { + VISU::TGetElnoNodeData< elno_type > aGetElnoNodeData( theElnoDataArray, theElnoDataMapper ); + + vtkCellArray *aConnectivity = vtkCellArray::New(); + aConnectivity->DeepCopy( theInput->GetCells() ); + + vtkPoints *anInputPoints = theInput->GetPoints(); + vtkPoints *aPoints = anInputPoints->New( anInputPoints->GetDataType() ); + vtkIdType aNbCells = aConnectivity->GetNumberOfCells(); + vtkIdType aNbPoints = aConnectivity->GetNumberOfConnectivityEntries() - aNbCells; + aPoints->Allocate( aNbPoints ); + + typedef typename VISU::TL::TEnum2VTKArrayType< points_type >::TResult TPointsDataArray; + TPointsDataArray* anInputPointsArray = TPointsDataArray::SafeDownCast( anInputPoints->GetData() ); + TPointsDataArray* anOutputPointsArray = TPointsDataArray::SafeDownCast( aPoints->GetData() ); + vtkPointData *anInputPointData = theInput->GetPointData(); + vtkPointData *anOutputPointData = theOutput->GetPointData(); + anOutputPointData->Allocate( aNbPoints ); + + typedef typename VISU::TL::TEnum2VTKArrayType< elno_type >::TResult TElnoDataArray; + TElnoDataArray* anElnoPointDataArray = TElnoDataArray::New(); + anElnoPointDataArray->SetNumberOfComponents( aGetElnoNodeData.getNbComp() ); + anElnoPointDataArray->SetNumberOfTuples( aNbPoints ); + anElnoPointDataArray->SetName( "ELNO_POINT_DATA" ); + + TPointsDataArray* anElnoPointCoords = TPointsDataArray::New(); + anElnoPointCoords->SetNumberOfComponents( 3 ); + anElnoPointCoords->SetNumberOfTuples( aNbPoints ); + anElnoPointCoords->SetName( "ELNO_POINT_COORDS" ); + + if ( theShrinkFactor > 0.0 ) { + ShrunkExecute< points_type, TPointsDataArray, elno_type, TElnoDataArray > + ( aConnectivity, + anInputPointData, + anOutputPointData, + anInputPointsArray, + anOutputPointsArray, + aGetElnoNodeData, + anElnoPointDataArray, + anElnoPointCoords, + theShrinkFactor ); + } else { + SimpleExecute< points_type, TPointsDataArray, elno_type, TElnoDataArray > + ( aConnectivity, + anInputPointData, + anOutputPointData, + anInputPointsArray, + anOutputPointsArray, + aGetElnoNodeData, + anElnoPointDataArray, + anElnoPointCoords ); + } + theOutput->SetPoints( aPoints ); theOutput->SetCells( theInput->GetCellTypesArray(), @@ -158,11 +251,11 @@ namespace aCellData->RemoveArray( "ELNO_COMPONENT_MAPPER" ); aCellData->RemoveArray( "ELNO_FIELD" ); - aPointData->AddArray( anElnoPointsDataArray ); - anElnoPointsDataArray->Delete(); + anOutputPointData->AddArray( anElnoPointDataArray ); + anElnoPointDataArray->Delete(); - aPointData->AddArray( aPointCoords ); - aPointCoords->Delete(); + anOutputPointData->AddArray( anElnoPointCoords ); + anElnoPointCoords->Delete(); return 1; } @@ -173,17 +266,22 @@ namespace int Execute( vtkUnstructuredGrid *theInput, vtkUnstructuredGrid *theOutput, vtkDataArray *theElnoDataArray, - vtkDataArray *theElnoDataMapper ) + vtkDataArray *theElnoDataMapper, + vtkFloatingPointType theShrinkFactor ) { switch( theElnoDataArray->GetDataType() ){ case VTK_DOUBLE: - return Execute2< points_type, VTK_DOUBLE >( theInput, theOutput, theElnoDataArray, theElnoDataMapper ); + return Execute2< points_type, VTK_DOUBLE > + ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor ); case VTK_FLOAT: - return Execute2< points_type, VTK_FLOAT >( theInput, theOutput, theElnoDataArray, theElnoDataMapper ); + return Execute2< points_type, VTK_FLOAT > + ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor ); case VTK_INT: - return Execute2< points_type, VTK_INT >( theInput, theOutput, theElnoDataArray, theElnoDataMapper ); + return Execute2< points_type, VTK_INT > + ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor ); case VTK_LONG: - return Execute2< points_type, VTK_LONG >( theInput, theOutput, theElnoDataArray, theElnoDataMapper ); + return Execute2< points_type, VTK_LONG > + ( theInput, theOutput, theElnoDataArray, theElnoDataMapper, theShrinkFactor ); default: break; } @@ -217,13 +315,13 @@ int VISU_ElnoDisassembleFilter::RequestData( vtkInformation *vtkNotUsed(request) vtkPoints *aPoints = anInput->GetPoints(); switch( aPoints->GetDataType() ){ case VTK_DOUBLE: - return ::Execute< VTK_DOUBLE >( anInput, anOutput, anElnoDataArray, anElnoDataMapper ); + return ::Execute< VTK_DOUBLE >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor ); case VTK_FLOAT: - return ::Execute< VTK_FLOAT >( anInput, anOutput, anElnoDataArray, anElnoDataMapper ); + return ::Execute< VTK_FLOAT >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor ); case VTK_INT: - return ::Execute< VTK_INT >( anInput, anOutput, anElnoDataArray, anElnoDataMapper ); + return ::Execute< VTK_INT >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor ); case VTK_LONG: - return ::Execute< VTK_LONG >( anInput, anOutput, anElnoDataArray, anElnoDataMapper ); + return ::Execute< VTK_LONG >( anInput, anOutput, anElnoDataArray, anElnoDataMapper, myShrinkFactor ); default: break; } diff --git a/src/PIPELINE/VISU_ElnoDisassembleFilter.hxx b/src/PIPELINE/VISU_ElnoDisassembleFilter.hxx index 730f3c55..0053c13e 100644 --- a/src/PIPELINE/VISU_ElnoDisassembleFilter.hxx +++ b/src/PIPELINE/VISU_ElnoDisassembleFilter.hxx @@ -29,12 +29,17 @@ public: static VISU_ElnoDisassembleFilter *New(); + void SetShrinkFactor( vtkFloatingPointType theValue ); + vtkFloatingPointType GetShrinkFactor(); + protected: VISU_ElnoDisassembleFilter(); ~VISU_ElnoDisassembleFilter(); int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *); + vtkFloatingPointType myShrinkFactor; + private: VISU_ElnoDisassembleFilter(const VISU_ElnoDisassembleFilter&); // Not implemented. void operator=(const VISU_ElnoDisassembleFilter&); // Not implemented. diff --git a/src/PIPELINE/VISU_ElnoWarpVector.hxx b/src/PIPELINE/VISU_ElnoWarpVector.hxx index 3047a726..29c1662a 100644 --- a/src/PIPELINE/VISU_ElnoWarpVector.hxx +++ b/src/PIPELINE/VISU_ElnoWarpVector.hxx @@ -38,6 +38,7 @@ protected: ~VISU_ElnoWarpVector(); int RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *); + vtkFloatingPointType myScaleFactor; private: