]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
IPAL8849,8850: problems with empty data set; MEN_FILE replaced with MEN_DESK_FILE...
authorjfa <jfa@opencascade.com>
Tue, 19 Jul 2005 13:12:52 +0000 (13:12 +0000)
committerjfa <jfa@opencascade.com>
Tue, 19 Jul 2005 13:12:52 +0000 (13:12 +0000)
src/PIPELINE/VISU_PipeLine.cxx
src/PIPELINE/VISU_PipeLine.hxx
src/VISUGUI/VISU_msg_en.po
src/VISUGUI/VisuGUI.cxx
src/VISUGUI/VisuGUI_ClippingDlg.cxx
src/VISU_I/VISU_Prs3d_i.cc
src/VISU_I/VISU_Prs3d_i.hh

index 21137aabe7d7e572ea18b781081473f598bb20ea..d670a681c9628bf6bfccf38c7e55871f2591c04e 100644 (file)
@@ -164,14 +164,22 @@ float VISU_PipeLine::GetAvailableMemory(float theSize, float theMinSize)
 
 //------------------------ Clipping planes -----------------------------------
 
-void VISU_PipeLine::AddClippingPlane(vtkPlane* thePlane)
+bool VISU_PipeLine::AddClippingPlane(vtkPlane* thePlane)
 {
-  if(thePlane){
-    if(vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()){
+  if (thePlane) {
+    if (vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()) {
       vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
       aFunction->AddItem(thePlane);
+
+      // Check, that at least one cell present after clipping.
+      // This check was introduced because of bug IPAL8849.
+      vtkUnstructuredGrid* aClippedGrid = GetInput2();
+      if (aClippedGrid->GetNumberOfCells() < 1) {
+        return false;
+      }
     }
   }
+  return true;
 }
 
 vtkPlane* VISU_PipeLine::GetClippingPlane(vtkIdType theID) const
index 5bdd0f2ec612e854cfb02986746ef01e0e959ee8..4dd6db6ea26d3770e796b46295a50c63775051a8 100644 (file)
@@ -93,7 +93,7 @@ public:
   // Clipping planes
   void RemoveAllClippingPlanes();
   vtkIdType GetNumberOfClippingPlanes() const;
-  void AddClippingPlane(vtkPlane* thePlane);
+  bool AddClippingPlane(vtkPlane* thePlane);
   vtkPlane* GetClippingPlane(vtkIdType theID) const;
 
   void SetPlaneParam(float theDir[3], float theDist, vtkPlane* thePlane);
index a2c4f76e76d8a4d02be430fea5bb7284c19f2550..19fe079eee6e220436ef2878194d5fb9956effc9 100644 (file)
@@ -974,6 +974,9 @@ msgstr "Plane# %1"
 msgid "VisuGUI_ClippingDlg::PLANES_COMBO_ITEM_no"
 msgstr "No planes"
 
+msgid "VisuGUI_ClippingDlg::WRN_EMPTY_RESULTING_PRS"
+msgstr "Impossible to use given clipping planes because of VTK restrictions. \n Please, provide non-empty resulting presentation."
+
 # -------------- Plot 3D --------------
 
 msgid "VisuGUI_Plot3DDlg::TITLE"
index 13b4a38a4557f114017b2c07f7b889620ef4e3e8..0d31866cc5e5156056dd4aa3af6cc6c77e4373f0 100644 (file)
@@ -2305,7 +2305,7 @@ createMenus()
 {
   // Add actions to menus
   int aMenuId;
-  aMenuId = createMenu( tr( "MEN_FILE" ), -1 );
+  aMenuId = createMenu( tr( "MEN_DESK_FILE" ), -1 );
   createMenu( separator(), aMenuId, -1, 10 );
   createMenu( VISU_IMPORT_FROM_FILE, aMenuId, 10 ); // import from file
   createMenu( VISU_EXPLORE_MED, aMenuId, 10 ); // explore MED file
index b62289281f3ed0f1bff5bed3350b498292a22bea..68f4b0c218557da672fc638da5c46f953fc055ac 100644 (file)
@@ -7,12 +7,15 @@
 #include "VISU_Prs3d_i.hh"
 #include "VISU_Result_i.hh"
 
+#include "VISU_PipeLine.hxx"
+
 #include "SalomeApp_SelectionMgr.h"
 
 #include "SVTK_ViewWindow.h"
 
 #include "SUIT_Session.h"
 #include "SUIT_Desktop.h"
+#include "SUIT_MessageBox.h"
 #include "SUIT_ResourceMgr.h"
 #include "SUIT_OverrideCursor.h"
 
@@ -331,27 +334,26 @@ VisuGUI_ClippingDlg::VisuGUI_ClippingDlg (VisuGUI* theModule,
   onSelectionChanged();
 
   // signals and slots connections :
-  connect(ComboBoxPlanes, SIGNAL(activated(int)), this, SLOT(onSelectPlane(int)));
-  connect(buttonNew, SIGNAL(clicked()), this, SLOT(ClickOnNew()));
-  connect(buttonDelete, SIGNAL(clicked()), this, SLOT(ClickOnDelete()));
-  connect(ComboBoxOrientation, SIGNAL(activated(int)), this, SLOT(onSelectOrientation(int)));
-  connect(SpinBoxDistance, SIGNAL(valueChanged(double)), this, SLOT(SetCurrentPlaneParam()));
-  connect(SpinBoxRot1, SIGNAL(valueChanged(double)), this, SLOT(SetCurrentPlaneParam()));
-  connect(SpinBoxRot2, SIGNAL(valueChanged(double)), this, SLOT(SetCurrentPlaneParam()));
-  connect(ButtonGroupIJKAxis, SIGNAL(clicked(int)), this, SLOT(onIJKAxisChanged(int)));
-  connect(SpinBoxIJKIndex, SIGNAL(valueChanged(int)), this, SLOT(SetCurrentPlaneIJKParam()));
-  connect(CheckBoxIJKPlaneReverse, SIGNAL(toggled(bool)), this, SLOT(SetCurrentPlaneIJKParam()));
-  connect(TabPane, SIGNAL(currentChanged (QWidget*)), this, SLOT(onTabChanged(QWidget*)));
-
-  connect(PreviewCheckBox, SIGNAL(toggled(bool)), this, SLOT(OnPreviewToggle(bool)));
+  connect(ComboBoxPlanes         , SIGNAL(activated(int))           , this, SLOT(onSelectPlane(int)));
+  connect(buttonNew              , SIGNAL(clicked())                , this, SLOT(ClickOnNew()));
+  connect(buttonDelete           , SIGNAL(clicked())                , this, SLOT(ClickOnDelete()));
+  connect(ComboBoxOrientation    , SIGNAL(activated(int))           , this, SLOT(onSelectOrientation(int)));
+  connect(SpinBoxDistance        , SIGNAL(valueChanged(double))     , this, SLOT(SetCurrentPlaneParam()));
+  connect(SpinBoxRot1            , SIGNAL(valueChanged(double))     , this, SLOT(SetCurrentPlaneParam()));
+  connect(SpinBoxRot2            , SIGNAL(valueChanged(double))     , this, SLOT(SetCurrentPlaneParam()));
+  connect(ButtonGroupIJKAxis     , SIGNAL(clicked(int))             , this, SLOT(onIJKAxisChanged(int)));
+  connect(SpinBoxIJKIndex        , SIGNAL(valueChanged(int))        , this, SLOT(SetCurrentPlaneIJKParam()));
+  connect(CheckBoxIJKPlaneReverse, SIGNAL(toggled(bool))            , this, SLOT(SetCurrentPlaneIJKParam()));
+  connect(TabPane                , SIGNAL(currentChanged (QWidget*)), this, SLOT(onTabChanged(QWidget*)));
+
+  connect(PreviewCheckBox  , SIGNAL(toggled(bool)), this, SLOT(OnPreviewToggle(bool)));
   connect(AutoApplyCheckBox, SIGNAL(toggled(bool)), this, SLOT(ClickOnApply()));
-  connect(buttonOk, SIGNAL(clicked()), this, SLOT(ClickOnOk()));
+
+  connect(buttonOk    , SIGNAL(clicked()), this, SLOT(ClickOnOk()));
+  connect(buttonApply , SIGNAL(clicked()), this, SLOT(ClickOnApply()));
   connect(buttonCancel, SIGNAL(clicked()), this, SLOT(ClickOnCancel()));
-  connect(buttonApply, SIGNAL(clicked()), this, SLOT(ClickOnApply()));
-  connect(myVisuGUI, SIGNAL(SignalCloseAllDialogs()), this, SLOT(ClickOnOk()));
+
   connect(mySelectionMgr, SIGNAL(currentSelectionChanged()), this, SLOT(onSelectionChanged()));
-  /* to close dialog if study frame change */
-  connect(myVisuGUI, SIGNAL(SignalStudyFrameChanged()), this, SLOT(ClickOnCancel()));
 
   this->show();
 }
@@ -459,16 +461,57 @@ void VisuGUI_ClippingDlg::ClickOnApply()
   if (SVTK_ViewWindow* aViewWindow = VISU::GetViewWindow(myVisuGUI)) {
     SUIT_OverrideCursor wc;
 
+    // Save clipping planes, currently applied to the presentation
+    // to enable restoring this state in case of failure.
+    // Refer to bugs IPAL8849, IPAL8850 for more information.
+    typedef vtkSmartPointer<vtkPlane> TPln;
+    typedef std::vector<TPln> TPlns;
+    bool isFailed = false;
+    TPlns anOldPlanes;
+    int iopl = 0, nbOldPlanes = myPrs3d->GetNumberOfClippingPlanes();
+    for (; iopl < nbOldPlanes; iopl++) {
+      anOldPlanes.push_back(myPrs3d->GetClippingPlane(iopl));
+    }
+
+    // Try to apply new clipping
     myPrs3d->RemoveAllClippingPlanes();
 
     VISU::TPlanes::iterator anIter = myPlanes.begin();
-    for (;anIter != myPlanes.end();anIter++) {
+    for (; anIter != myPlanes.end(); anIter++) {
       OrientedPlane* anOrientedPlane = OrientedPlane::New(aViewWindow);
       anOrientedPlane->ShallowCopy(anIter->GetPointer());
-      myPrs3d->AddClippingPlane(anOrientedPlane);
+      if (!myPrs3d->AddClippingPlane(anOrientedPlane)) {
+        isFailed = true;
+      }
       anOrientedPlane->Delete();
     }
 
+    // Check contents of the resulting (clipped) presentation data
+    if (!isFailed) {
+      VISU_PipeLine* aPL = myPrs3d->GetPL();
+      VISU_PipeLine::TMapper* aM = aPL->GetMapper();
+      vtkDataSet* aPrsData = aM->GetInput();
+      aPrsData->Update();
+      if (aPrsData->GetNumberOfCells() < 1) {
+        isFailed = true;
+      }
+    }
+
+    if (isFailed) {
+      // Restore previous clipping state because of failure.
+      myPrs3d->RemoveAllClippingPlanes();
+
+      TPlns::iterator anOldIter = anOldPlanes.begin();
+      for (; anOldIter != anOldPlanes.end(); anOldIter++) {
+        myPrs3d->AddClippingPlane(anOldIter->GetPointer());
+      }
+
+      SUIT_MessageBox::warn1(VISU::GetDesktop(myVisuGUI),
+                             tr("WRN_VISU"),
+                             tr("WRN_EMPTY_RESULTING_PRS"),
+                             tr("BUT_OK") );
+    }
+
     VISU::RenderViewWindow(aViewWindow);
   }
 }
index 699e683f50ab2aac97288a0979bc54a55bfb81f0..5e39d80e3799988ba81fef1fdf9bbb2c69e7e4de 100644 (file)
@@ -174,8 +174,8 @@ vtkIdType VISU::Prs3d_i::GetNumberOfClippingPlanes() const{
   return myPipeLine->GetNumberOfClippingPlanes();
 }
 
-void VISU::Prs3d_i::AddClippingPlane(vtkPlane* thePlane){
-  myPipeLine->AddClippingPlane(thePlane);
+bool VISU::Prs3d_i::AddClippingPlane(vtkPlane* thePlane){
+  return myPipeLine->AddClippingPlane(thePlane);
 }
 
 vtkPlane* VISU::Prs3d_i::GetClippingPlane(vtkIdType theID) const{
index cc0d0a37b8e0e7d83481e1746f7d5f3b684d2c69..124a97ff75ca206d0f9ecea8bb68695a355542f1 100644 (file)
@@ -105,7 +105,7 @@ namespace VISU{
     // Clipping planes
     void RemoveAllClippingPlanes();
     vtkIdType GetNumberOfClippingPlanes() const;
-    void AddClippingPlane(vtkPlane* thePlane);
+    bool AddClippingPlane(vtkPlane* thePlane);
     vtkPlane* GetClippingPlane(vtkIdType theID) const;
 
     void SetPlaneParam(float theDir[3], float theDist, vtkPlane* thePlane);