From: Afeef Date: Wed, 14 Apr 2021 11:20:11 +0000 (+0200) Subject: Addition to non-delegation list: X-Git-Url: http://git.salome-platform.org/gitweb/?a=commitdiff_plain;h=fa5fc5a9213f3a5337dd14f3823e7b3dd5144b94;p=modules%2Fgui.git Addition to non-delegation list: Add all 3d quadratic, biquadratic and triquadratic elements to the list of elements which have to be handled natively. --- diff --git a/src/VTKViewer/VTKViewer_GeometryFilter.cxx b/src/VTKViewer/VTKViewer_GeometryFilter.cxx index 169d7d803..f9984da96 100644 --- a/src/VTKViewer/VTKViewer_GeometryFilter.cxx +++ b/src/VTKViewer/VTKViewer_GeometryFilter.cxx @@ -211,19 +211,35 @@ VTKViewer_GeometryFilter if ( vtkUnsignedCharArray* types = inputUnstrctured->GetCellTypesArray() ) { std::set ElementsNotFitToDelegate; + + // All quadratic, biquadratic, and triquadratic elements not fit for delegation + // as SMESH has special display with curves mode for meshes containing these + // elements hence such meshes are not handled by "vtkGeometryFilter" instead + // the native VTKViewer_GeometryFilter::UnstructuredGridExecute is used. ElementsNotFitToDelegate.insert( VTK_QUADRATIC_EDGE ); ElementsNotFitToDelegate.insert( VTK_QUADRATIC_TRIANGLE ); ElementsNotFitToDelegate.insert( VTK_BIQUADRATIC_TRIANGLE ); ElementsNotFitToDelegate.insert( VTK_QUADRATIC_QUAD ); ElementsNotFitToDelegate.insert( VTK_BIQUADRATIC_QUAD ); ElementsNotFitToDelegate.insert( VTK_QUADRATIC_POLYGON ); + ElementsNotFitToDelegate.insert( VTK_QUADRATIC_TETRA ); + ElementsNotFitToDelegate.insert( VTK_QUADRATIC_HEXAHEDRON ); + ElementsNotFitToDelegate.insert( VTK_TRIQUADRATIC_HEXAHEDRON ); + ElementsNotFitToDelegate.insert( VTK_QUADRATIC_WEDGE ); + ElementsNotFitToDelegate.insert( VTK_BIQUADRATIC_QUADRATIC_WEDGE ); + ElementsNotFitToDelegate.insert( VTK_QUADRATIC_PYRAMID ); + + // Some openMP tests reveal that meshes with polyhedrons can sometimes cause + // problems as such we avoide delegation = ElementsNotFitToDelegate. It would be + // nice to investigate and resolve the problem with multi-therding in future. ElementsNotFitToDelegate.insert( VTK_POLYHEDRON ); + for ( int i = 0; i < types->GetNumberOfTuples() && !NotFitForDelegation; ++i ) NotFitForDelegation = ElementsNotFitToDelegate.count( types->GetValue(i) ); } - if ( NotFitForDelegation ) + if ( NotFitForDelegation ) return this->UnstructuredGridExecute(input, output, outInfo); - else + else return this->vtkGeometryFilter::UnstructuredGridExecute(input, output, nullptr, &exc); } }