From: sln Date: Tue, 9 Dec 2008 12:25:04 +0000 (+0000) Subject: 0013557: field values display X-Git-Tag: Phase8_Part1_16122008~11 X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=c5caea948257c461d6f5e2997410551833d8914e;p=modules%2Fvisu.git 0013557: field values display Now you can display values applied to the cells or nodes of 3D presentation intended for visualization of calculation data. New methods are provided in actor. virtual void SetValuesLabeled( const bool theIsValLabeled ); virtual bool GetValuesLabeled() const; virtual vtkDataSet* GetValLabelsInput() const; vtkTextProperty* GetsValLabelsProps() const; --- diff --git a/src/OBJECT/VISU_Actor.cxx b/src/OBJECT/VISU_Actor.cxx index f001de9a..cb8245af 100644 --- a/src/OBJECT/VISU_Actor.cxx +++ b/src/OBJECT/VISU_Actor.cxx @@ -60,6 +60,12 @@ #include #include #include +#include +#include +#include +#include +#include +#include #include #include @@ -139,6 +145,44 @@ VISU_Actor aPickingSettings->AddObserver(VISU::UpdatePickingSettingsEvent, myEventCallbackCommand.GetPointer(), myPriority); + + //Definition of values labeling pipeline + + myValLblDataSet = vtkUnstructuredGrid::New(); + + myValCellCenters = vtkCellCenters::New(); + myValCellCenters->SetInput(myValLblDataSet); + + myValMaskPoints = vtkMaskPoints::New(); + myValMaskPoints->SetInput(myValCellCenters->GetOutput()); + myValMaskPoints->SetOnRatio(1); + + myValSelectVisiblePoints = vtkSelectVisiblePoints::New(); + myValSelectVisiblePoints->SetInput(myValMaskPoints->GetOutput()); + myValSelectVisiblePoints->SelectInvisibleOff(); + myValSelectVisiblePoints->SetTolerance(0.1); + + myValLabeledDataMapper = vtkLabeledDataMapper::New(); + myValLabeledDataMapper->SetInput(myValSelectVisiblePoints->GetOutput()); + myValLabeledDataMapper->SetLabelFormat("%g"); + myValLabeledDataMapper->SetLabelModeToLabelScalars(); + + vtkTextProperty* aClsTextProp = vtkTextProperty::New(); + aClsTextProp->SetFontFamilyToTimes(); + static int aCellsFontSize = 12; + aClsTextProp->SetFontSize(aCellsFontSize); + aClsTextProp->SetBold(1); + aClsTextProp->SetItalic(0); + aClsTextProp->SetShadow(0); + myValLabeledDataMapper->SetLabelTextProperty(aClsTextProp); + aClsTextProp->Delete(); + + myIsValLabeled = false; + + myValLabels = vtkActor2D::New(); + myValLabels->SetMapper(myValLabeledDataMapper); + myValLabels->GetProperty()->SetColor(0,1,0); + myValLabels->SetVisibility( myIsValLabeled ); } //---------------------------------------------------------------------------- @@ -182,6 +226,15 @@ VISU_Actor VISU_Actor ::~VISU_Actor() { + // Deleting of values labeling pipeline + myValLblDataSet->Delete(); + myValLabeledDataMapper->RemoveAllInputs(); + myValLabeledDataMapper->Delete(); + myValSelectVisiblePoints->Delete(); + myValMaskPoints->Delete(); + myValCellCenters->Delete(); + myValLabels->Delete(); + if(MYDEBUG) MESSAGE("~VISU_Actor() - this = "<AddActor(myAnnotationActor.GetPointer()); theRenderer->AddActor(myTextActor.GetPointer()); + + myValSelectVisiblePoints->SetRenderer( theRenderer ); + theRenderer->AddActor2D( myValLabels ); + } //================================================================== @@ -487,10 +544,20 @@ VISU_Actor { theRenderer->RemoveActor(myAnnotationActor.GetPointer()); theRenderer->RemoveActor(myTextActor.GetPointer()); + theRenderer->RemoveActor(myValLabels); Superclass::RemoveFromRender(theRenderer); myDestroySignal(this); } +//---------------------------------------------------------------------------- +void +VISU_Actor +::SetVisibility(int theMode) +{ + Superclass::SetVisibility( theMode ); + myValLabels->SetVisibility( myIsValLabeled && theMode ); +} + //---------------------------------------------------------------------------- void VISU_Actor @@ -1101,3 +1168,92 @@ VISU_Actor Update(); } + +// --------------------------------------------------------------- + +void VISU_Actor::SetValuesLabeled( const bool theIsValLabeled ) +{ + vtkDataSet* aGrid = GetValLabelsInput(); + if ( !aGrid ) + return; + + bool isOnPnt = VISU::IsDataOnPoints( aGrid ); + bool isOnCell = VISU::IsDataOnCells( aGrid ); + if ( !isOnPnt && !isOnCell ) + { + // try to specify location of scalars "manually" + vtkCellData* aCData = aGrid->GetCellData(); + if ( aCData ) + { + vtkDataArray* anArr = aCData->GetScalars(); + if ( anArr && anArr->GetNumberOfTuples() ) + isOnCell = true; + } + + if ( !isOnCell ) + { + vtkPointData* aPData = aGrid->GetPointData(); + if ( aPData ) + { + vtkDataArray* anArr = aPData->GetScalars(); + if ( anArr && anArr->GetNumberOfTuples() ) + isOnPnt = true; + } + } + + if ( !isOnPnt && !isOnCell ) + { + myValLabels->SetVisibility( false ); + return; + } + } + + myIsValLabeled = theIsValLabeled; + + if ( myIsValLabeled ) + { + vtkDataSet* aDataSet = aGrid; + + if ( isOnCell ) + { + myValCellCenters->SetInput( aDataSet ); + myValMaskPoints->SetInput( myValCellCenters->GetOutput() ); + } + else if ( isOnPnt ) + myValMaskPoints->SetInput( aDataSet ); + + myValLabels->SetVisibility( GetVisibility() ); + } + else + myValLabels->SetVisibility( false ); + + Modified(); +} + +//---------------------------------------------------------------------------- + +bool VISU_Actor::GetValuesLabeled() const +{ + return myIsValLabeled; +} + +//---------------------------------------------------------------------------- + +vtkTextProperty* VISU_Actor::GetsValLabelsProps() const +{ + return myValLabeledDataMapper->GetLabelTextProperty(); +} + +//---------------------------------------------------------------------------- + +vtkDataSet* VISU_Actor::GetValLabelsInput() +{ + vtkDataSet* aDataSet = 0; + VISU_PipeLine* aPL = GetPipeLine(); + if ( aPL ) + aDataSet = aPL->GetOutput(); + if ( !aDataSet ) + aDataSet = GetInput(); + return aDataSet; +} + diff --git a/src/OBJECT/VISU_Actor.h b/src/OBJECT/VISU_Actor.h index 34c1dcc8..69fcf141 100644 --- a/src/OBJECT/VISU_Actor.h +++ b/src/OBJECT/VISU_Actor.h @@ -49,6 +49,12 @@ class VISU_PipeLine; class vtkPlane; class vtkImplicitFunctionCollection; class vtkFeatureEdges; +class vtkTextProperty; +class vtkCellCenters; +class vtkSelectVisiblePoints; +class vtkLabeledDataMapper; +class vtkMaskPoints; +class vtkActor2D; class VISU_FramedTextActor; @@ -107,6 +113,10 @@ class VTKOCC_EXPORT VISU_Actor : public VISU_ActorBase virtual void RemoveFromRender(); + + virtual + void + SetVisibility(int theMode); //---------------------------------------------------------------------------- virtual @@ -308,6 +318,22 @@ class VTKOCC_EXPORT VISU_Actor : public VISU_ActorBase UpdatePickingSettings(); //---------------------------------------------------------------------------- + //! Methods for values labeling + virtual + void + SetValuesLabeled( const bool theIsValLabeled ); + + virtual + bool + GetValuesLabeled() const; + + virtual + vtkDataSet* + GetValLabelsInput(); + + vtkTextProperty* + GetsValLabelsProps() const; + protected: VISU_Actor(); @@ -364,6 +390,15 @@ class VTKOCC_EXPORT VISU_Actor : public VISU_ActorBase Selection_Mode myLastSelectionMode; bool myIsSubElementsHighlighted; + + // fields for values labeling + bool myIsValLabeled; + vtkDataSet* myValLblDataSet; + vtkActor2D* myValLabels; + vtkMaskPoints* myValMaskPoints; + vtkCellCenters* myValCellCenters; + vtkLabeledDataMapper* myValLabeledDataMapper; + vtkSelectVisiblePoints* myValSelectVisiblePoints; }; #endif //VISU_ACTOR_H