Salome HOME
Fix selection of cells with delegateToVtk
authorChristophe Bourcier <christophe.bourcier@cea.fr>
Mon, 14 Feb 2022 15:08:23 +0000 (16:08 +0100)
committerrnv <rnv@opencascade.com>
Tue, 22 Feb 2022 14:13:33 +0000 (17:13 +0300)
src/VTKViewer/VTKViewer_GeometryFilter.cxx
src/VTKViewer/VTKViewer_GeometryFilter.h

index 577206a2d0f93ed717ac2f26305c1b701fb93a09..e05ec36f281b97b6094bfb9fd570fc5c5facc8e8 100644 (file)
@@ -108,8 +108,8 @@ VTKViewer_GeometryFilter
   static int forceDelegateToVtk = -1;
   if ( forceDelegateToVtk < 0 )
   {
-    QString env = Qtx::getenv( "SALOME_ACTOR_DELEGATE_TO_VTK" );
-    forceDelegateToVtk = (int)(env == "1");
+    QString env = Qtx::getenv( "SALOME_ACTOR_DO_NOT_DELEGATE_TO_VTK" );
+    forceDelegateToVtk = (int)(env != "1");
   }
   delegateToVtk = forceDelegateToVtk > 0;
 }
@@ -166,6 +166,25 @@ struct vtkExcludedFaces
   ~vtkExcludedFaces() { delete this->Links; }
 };
 
+// fill myVTK2ObjIds to get the correspondence between vtk ids (displayed edges and faces
+// computed in vtkGeometryFilter::UnstructuredGridExecute) and original cell ids (mesh cells)
+void VTKViewer_GeometryFilter
+::FillVTK2ObjIds(vtkPolyData *output) {
+
+  vtkDataArray* vtkOriginalCellIds = output->GetCellData()->GetArray("vtkOriginalCellIds");
+
+  if (vtkOriginalCellIds == nullptr)
+    throw std::runtime_error("vtkOriginalCellIds is null. Something is wrong.");
+
+  const vtkIdType numTuples = vtkOriginalCellIds->GetNumberOfTuples();
+  myVTK2ObjIds.resize(numTuples);
+
+  // copy ids of vtkOriginalCellIds into myVTK2ObjIds
+  vtkIdTypeArray *vtkOriginalCellIdsInt(vtkIdTypeArray::SafeDownCast(vtkOriginalCellIds));
+  vtkIdType* origIds(vtkOriginalCellIdsInt->GetPointer(0));
+  myVTK2ObjIds.assign(origIds, origIds+numTuples);
+}
+
 int
 VTKViewer_GeometryFilter
 ::RequestData(
@@ -217,9 +236,16 @@ VTKViewer_GeometryFilter
           return this->vtkGeometryFilter::PolyDataExecute(input, output, &exc);
        case VTK_UNSTRUCTURED_GRID:
          {
-          vtkUnstructuredGrid* inputUnstrctured = static_cast<vtkUnstructuredGrid*>(input);
+          vtkUnstructuredGrid* inputUnstructured = static_cast<vtkUnstructuredGrid*>(input);
+
+          // The "info", is passed to provide information about the unstructured grid. This
+          // is done to avoid repeated evaluations of "info" in the vtkGeometryFilter and
+          // the vtkDataSetSurfaceFilter (which vtkGeometryFilter  might delagate to in case
+          // of nonlinear data).
+          vtkGeometryFilterHelper* info = vtkGeometryFilterHelper::CharacterizeUnstructuredGrid(inputUnstructured);
+
           bool NotFitForDelegation = false;
-          if ( vtkUnsignedCharArray* types = inputUnstrctured->GetCellTypesArray() )
+          if ( vtkUnsignedCharArray* types = inputUnstructured->GetCellTypesArray() )
             {
              std::set<vtkIdType> ElementsNotFitToDelegate;
 
@@ -251,7 +277,21 @@ VTKViewer_GeometryFilter
          if ( NotFitForDelegation )
             return this->UnstructuredGridExecute(input, output, outInfo);
          else
-            return this->vtkGeometryFilter::UnstructuredGridExecute(input, output, nullptr, &exc);
+           {
+           int ret;
+            if ( myStoreMapping ) {
+              // pass through cell ids to get original cell ids
+              this->PassThroughCellIds = true;
+              ret = this->vtkGeometryFilter::UnstructuredGridExecute(input, output, info, &exc);
+              FillVTK2ObjIds(output);
+            }
+            else {
+              // no need to get original cell ids
+              this->PassThroughCellIds = false;
+              ret = this->vtkGeometryFilter::UnstructuredGridExecute(input, output, info, &exc);
+            }
+            return ret;
+           }
          }
      }
 
index af65b0267d040111d3addd2072a6a0f403cb9d8b..119427545d1ee734735199cbe1d9053e6ab124fd 100644 (file)
@@ -130,7 +130,7 @@ protected:
   //special cases for performance
   
   /*! \fn void UnstructuredGridExecute();
-   * \brief Filter culculation method for data object type is VTK_UNSTRUCTURED_GRID.
+   * \brief Filter calculation method for data object type is VTK_UNSTRUCTURED_GRID.
    */
   int UnstructuredGridExecute (vtkDataSet *, vtkPolyData *, vtkInformation *);
 
@@ -141,6 +141,12 @@ protected:
                          TMapOfVectorId& theDimension2VTK2ObjIds,
                          bool triangulate = false);
 
+  /*! \fn void FillVTK2ObjIds(vtkPolyData *output);
+   * \brief fill myVTK2ObjIds to get the correspondence between vtk ids (displayed edges and faces
+   * computed in vtkGeometryFilter::UnstructuredGridExecute) and original cell ids (mesh cells)
+   */
+  void FillVTK2ObjIds(vtkPolyData *output);
+
   // Delegate VTKViewer_GeometryFilter::UnstructuredGridExecute to vtkGeometryFilter::UnstructuredGridExecute
   bool delegateToVtk = false;