Salome HOME
Crash after pre-viewing mesh translation
[modules/smesh.git] / src / SMESHGUI / SMESHGUI_MeshEditPreview.cxx
1 // Copyright (C) 2007-2016  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 // SMESH SMESHGUI : GUI for SMESH component
24 // File   : SMESHGUI_MeshEditPreview.cxx
25 // Author : Open CASCADE S.A.S.
26 // SMESH includes
27 //
28 #include "SMESHGUI_MeshEditPreview.h"
29
30 #include "SMESHGUI_VTKUtils.h"
31 #include "SMESH_Actor.h"
32 #include "SMESH_ActorUtils.h"
33
34 // SALOME GUI includes
35 #include <SVTK_Renderer.h>
36 #include <SVTK_ViewWindow.h>
37 #include <VTKViewer_CellLocationsArray.h>
38
39 // VTK includes
40 #include <vtkCellArray.h>
41 #include <vtkCoordinate.h>
42 #include <vtkDataSetMapper.h>
43 #include <vtkIdList.h>
44 #include <vtkPoints.h>
45 #include <vtkProperty.h>
46 #include <vtkRenderer.h>
47 #include <vtkTextActor.h>
48 #include <vtkTextMapper.h>
49 #include <vtkUnsignedCharArray.h>
50 #include <vtkUnstructuredGrid.h>
51
52 // Qt includes
53 #include <QColor>
54
55 // IDL includes
56 #include <SALOMEconfig.h>
57 #include CORBA_SERVER_HEADER(SMESH_MeshEditor)
58
59 #include <gp_Ax3.hxx>
60
61 //================================================================================
62 /*!
63  * \brief Constructor
64  */
65 //================================================================================
66
67 SMESHGUI_MeshEditPreview::SMESHGUI_MeshEditPreview(SVTK_ViewWindow* theViewWindow):
68   myViewWindow(theViewWindow)
69 {
70   myGrid = vtkUnstructuredGrid::New();
71
72   // Create and display actor
73   vtkDataSetMapper* aMapper = vtkDataSetMapper::New();
74   aMapper->SetInputData( myGrid );
75
76   myPreviewActor = SALOME_Actor::New();
77   myPreviewActor->SetInfinitive(true);
78   myPreviewActor->VisibilityOn();
79   myPreviewActor->PickableOff();
80
81   double aFactor,aUnits;
82   myPreviewActor->SetResolveCoincidentTopology(true);
83   myPreviewActor->GetPolygonOffsetParameters(aFactor,aUnits);
84   myPreviewActor->SetPolygonOffsetParameters(aFactor,0.2*aUnits);
85
86   double anRGB[3];
87   SMESH::GetColor( "SMESH", "selection_element_color", anRGB[0], anRGB[1], anRGB[2], QColor( 0, 170, 255 ) );
88   SetColor( anRGB[0], anRGB[1], anRGB[2] );
89
90   myPreviewActor->SetMapper( aMapper );
91   aMapper->Delete();
92
93   myViewWindow->AddActor(myPreviewActor);
94 }
95
96 //================================================================================
97 /*!
98  * \brief Destroy
99  */
100 //================================================================================
101
102 SMESHGUI_MeshEditPreview::~SMESHGUI_MeshEditPreview()
103 {
104   myGrid->Delete();
105
106   for ( size_t iA = 0; iA < myLabelActors.size(); ++iA )
107     if ( myLabelActors[iA] )
108     {
109       myPreviewActor->GetRenderer()->RemoveActor( myLabelActors[iA] );
110       myLabelActors[iA]->Delete();
111     }
112
113   myViewWindow->RemoveActor(myPreviewActor);
114   myPreviewActor->Delete();
115 }
116
117 //================================================================================
118 /*!
119  * \brief Returns vtk cell type
120  */
121 //================================================================================
122
123 vtkIdType getCellType( const SMDSAbs_ElementType theType,
124                        const bool                thePoly,
125                        const int                 theNbNodes )
126 {
127   switch( theType ) 
128   {
129   case SMDSAbs_Ball:              return VTK_VERTEX;
130   case SMDSAbs_Node:              return VTK_VERTEX;
131   case SMDSAbs_Edge:
132     if( theNbNodes == 2 )         return VTK_LINE;
133     else if ( theNbNodes == 3 )   return VTK_QUADRATIC_EDGE;
134     else return VTK_EMPTY_CELL;
135
136   case SMDSAbs_Face  :
137     if (thePoly && theNbNodes>2 ) return VTK_POLYGON;
138     else if ( theNbNodes == 3 )   return VTK_TRIANGLE;
139     else if ( theNbNodes == 4 )   return VTK_QUAD;
140     else if ( theNbNodes == 6 )   return VTK_QUADRATIC_TRIANGLE;
141     else if ( theNbNodes == 7 )   return VTK_BIQUADRATIC_TRIANGLE;
142     else if ( theNbNodes == 8 )   return VTK_QUADRATIC_QUAD;
143     else if ( theNbNodes == 9 )   return VTK_BIQUADRATIC_QUAD;
144     else return VTK_EMPTY_CELL;
145
146   case SMDSAbs_Volume:
147     if (thePoly && theNbNodes>3 ) return VTK_CONVEX_POINT_SET;
148     else if ( theNbNodes == 4 )   return VTK_TETRA;
149     else if ( theNbNodes == 5 )   return VTK_PYRAMID;
150     else if ( theNbNodes == 6 )   return VTK_WEDGE;
151     else if ( theNbNodes == 8 )   return VTK_HEXAHEDRON;
152     else if ( theNbNodes == 10 )  return VTK_QUADRATIC_TETRA;
153     else if ( theNbNodes == 20 )  return VTK_QUADRATIC_HEXAHEDRON;
154     else if ( theNbNodes == 27 )  return VTK_TRIQUADRATIC_HEXAHEDRON;
155     else if ( theNbNodes == 15  ) return VTK_QUADRATIC_WEDGE;
156     else if ( theNbNodes == 18  ) return VTK_BIQUADRATIC_QUADRATIC_WEDGE;
157     else if ( theNbNodes == 13  ) return VTK_QUADRATIC_PYRAMID;//VTK_CONVEX_POINT_SET;
158     else return VTK_EMPTY_CELL;
159
160   default: return VTK_EMPTY_CELL;
161   }
162 }
163
164 //================================================================================
165 /*!
166  * \brief Set preview data
167  */
168 //================================================================================
169
170 void SMESHGUI_MeshEditPreview::SetData (const SMESH::MeshPreviewStruct& previewData)
171 {
172   // Create points
173   const SMESH::nodes_array& aNodesXYZ = previewData.nodesXYZ;
174   vtkPoints* aPoints = vtkPoints::New();
175   aPoints->SetNumberOfPoints(aNodesXYZ.length());
176
177   for ( size_t i = 0; i < aNodesXYZ.length(); i++ ) {
178     aPoints->SetPoint( i, aNodesXYZ[i].x, aNodesXYZ[i].y, aNodesXYZ[i].z );
179   }
180   myGrid->SetPoints(aPoints);
181
182   aPoints->Delete();
183
184   // Create cells
185   const SMESH::long_array&  anElemConnectivity = previewData.elementConnectivities;
186   const SMESH::types_array& anElemTypes = previewData.elementTypes;
187
188   vtkIdType aCellsSize = anElemConnectivity.length() + anElemTypes.length();
189   vtkIdType aNbCells = anElemTypes.length();
190
191   vtkCellArray* aConnectivity = vtkCellArray::New();
192   aConnectivity->Allocate( aCellsSize, 0 );
193
194   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
195   aCellTypesArray->SetNumberOfComponents( 1 );
196   aCellTypesArray->Allocate( aNbCells * aCellTypesArray->GetNumberOfComponents() );
197
198   vtkIdList *anIdList = vtkIdList::New();
199   int aNodePos = 0;
200
201   for ( size_t i = 0; i < anElemTypes.length(); i++ ) {
202     const SMESH::ElementSubType& anElementSubType = anElemTypes[i];
203     SMDSAbs_ElementType aType = SMDSAbs_ElementType(anElementSubType.SMDS_ElementType);
204     vtkIdType aNbNodes = anElementSubType.nbNodesInElement;
205     anIdList->SetNumberOfIds( aNbNodes );
206
207     for ( vtkIdType aNodeId = 0; aNodeId < aNbNodes; aNodeId++ ){
208       anIdList->SetId( aNodeId, anElemConnectivity[aNodePos] );
209       aNodePos++;
210     }
211
212     aConnectivity->InsertNextCell( anIdList );
213     aCellTypesArray->InsertNextValue( getCellType( aType,
214                                                    anElemTypes[i].isPoly, 
215                                                    aNbNodes ) );
216   }
217   anIdList->Delete();
218
219   // Insert cells in grid
220   VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
221   aCellLocationsArray->SetNumberOfComponents( 1 );
222   aCellLocationsArray->SetNumberOfTuples( aNbCells );
223
224   aConnectivity->InitTraversal();
225   for( vtkIdType idType = 0, *pts, npts; aConnectivity->GetNextCell( npts, pts ); idType++ )
226     aCellLocationsArray->SetValue( idType, aConnectivity->GetTraversalLocation( npts ) );
227
228   myGrid->SetCells( aCellTypesArray, aCellLocationsArray, aConnectivity );
229
230   myPreviewActor->GetMapper()->Update();
231
232   aCellTypesArray->Delete();
233   aCellLocationsArray->Delete();
234   aConnectivity->Delete();
235
236   SetVisibility(true);
237 }
238
239 //================================================================================
240 /*!
241  * \brief Set shape of an arrow of a unit length and nb of arrows
242  */
243 //================================================================================
244
245 void SMESHGUI_MeshEditPreview::SetArrowShapeAndNb( int         nbArrows,
246                                                    double      headLength,
247                                                    double      headRadius,
248                                                    double      start,
249                                                    const char* labels)
250 {
251   const int theNbPoints = 10; // in one arrow
252   myUnitArrowPnts.reserve( theNbPoints );
253   myUnitArrowPnts.clear();
254
255   // unit arrow || OZ
256
257   for ( int i = 0; i < theNbPoints - 2; ++i )
258   {
259     double angle = i * 2 * M_PI / ( theNbPoints - 2 );
260     myUnitArrowPnts.push_back( gp_Pnt( headRadius * Cos( angle ),
261                                        headRadius * Sin( angle ),
262                                        1. - headLength ));
263   }
264   myUnitArrowPnts.push_back( gp_Pnt( 0, 0, start ));
265   myUnitArrowPnts.push_back( gp_Pnt( 0, 0, 1 ));
266
267
268   // nodes of all arrows
269
270   vtkPoints* aPoints = vtkPoints::New();
271   aPoints->SetNumberOfPoints( theNbPoints * nbArrows );
272   for ( int iP = 0, iA = 0; iA < nbArrows; ++iA )
273     for ( int i = 0; i < theNbPoints; ++i, ++iP )
274       aPoints->SetPoint( iP,
275                          myUnitArrowPnts[i].X(),
276                          myUnitArrowPnts[i].Y(),
277                          myUnitArrowPnts[i].Z()  );
278   myGrid->SetPoints(aPoints);
279   aPoints->Delete();
280
281   // connectivity of all arrows
282
283   const int theNbCells = ( theNbPoints - 1 ); // in one arrow
284   myGrid->Allocate( theNbCells  * nbArrows );
285   for ( int nP = 0, iA = 0; iA < nbArrows; ++iA, nP += theNbPoints )
286   {
287     vtkIdType conn[3] = { theNbPoints - 1 + nP, // arrow end
288                           theNbPoints - 3 + nP, // point on a circle
289                           nP };                 // point on a circle
290     for ( int i = 0; i < theNbCells-1; ++i )
291     {
292       myGrid->InsertNextCell( VTK_TRIANGLE, 3, conn );
293       conn[1] = conn[2];
294       conn[2] = conn[2] + 1;
295     }
296     conn[1] = theNbPoints - 2 + nP;
297     myGrid->InsertNextCell( VTK_LINE, 2, conn );
298   }
299
300   myLabelActors.resize( nbArrows, ( vtkTextActor*) NULL );
301   char label[] = "X";
302   if ( labels )
303     for ( int iA = 0; iA < nbArrows; ++iA )
304     {
305       label[0] = labels[iA];
306       vtkTextMapper* text = vtkTextMapper::New();
307       text->SetInput( label );
308       vtkCoordinate* coord = vtkCoordinate::New();
309
310       myLabelActors[iA] = vtkTextActor::New();
311       //myLabelActors[iA]->SetMapper( text );
312       myLabelActors[iA]->SetInput( label );
313       myLabelActors[iA]->SetTextScaleModeToNone();
314       myLabelActors[iA]->PickableOff();
315       myLabelActors[iA]->GetPositionCoordinate()->SetReferenceCoordinate( coord );
316
317       text->Delete();
318       coord->Delete();
319
320       myPreviewActor->GetRenderer()->AddActor(myLabelActors[iA]);
321     }
322 }
323
324 //================================================================================
325 /*!
326  * \brief Set data to show moved/rotated/scaled arrows
327  *  \param [in] axes - location and direction of the arrows
328  *  \param [in] length - length of arrows
329  */
330 //================================================================================
331
332 void SMESHGUI_MeshEditPreview::SetArrows( const gp_Ax1* axes,
333                                           double        length )
334 {
335   vtkPoints* aPoints = myGrid->GetPoints();
336
337   for ( int iP = 0, iA = 0; iA < (int) myLabelActors.size(); ++iA )
338   {
339     gp_Trsf trsf;
340     trsf.SetTransformation( gp_Ax3( axes[iA].Location(), axes[iA].Direction() ), gp::XOY() );
341
342     for ( size_t i = 0; i < myUnitArrowPnts.size(); ++i, ++iP )
343     {
344       gp_Pnt p = myUnitArrowPnts[i].Scaled( gp::Origin(), length );
345       p.Transform( trsf );
346       aPoints->SetPoint( iP, p.X(), p.Y(), p.Z() );
347     }
348     if ( myLabelActors[iA] )
349       if ( vtkCoordinate* aCoord =
350            myLabelActors[iA]->GetPositionCoordinate()->GetReferenceCoordinate() )
351       {
352         double p[3];
353         aPoints->GetPoint( iP-1, p );
354         aCoord->SetValue( p );
355       }
356   }
357
358   myGrid->Modified();
359 }
360
361 //================================================================================
362 /*!
363  * \brief Set visibility
364  */
365 //================================================================================
366
367 void SMESHGUI_MeshEditPreview::SetVisibility (bool theVisibility)
368 {
369   myPreviewActor->SetVisibility(theVisibility);
370   for ( size_t iA = 0; iA < myLabelActors.size(); ++iA )
371     if ( myLabelActors[iA] )
372       myLabelActors[iA]->SetVisibility(theVisibility);
373   SMESH::RepaintCurrentView();
374 }
375
376 //================================================================================
377 /*!
378  * \brief Set preview color
379  */
380 //================================================================================
381
382 void SMESHGUI_MeshEditPreview::SetColor(double R, double G, double B)
383 {
384   myPreviewActor->SetColor( R, G, B );
385 }
386
387 //================================================================================
388 /*!
389  * \brief Get preview actor
390  */
391 //================================================================================
392
393 SALOME_Actor* SMESHGUI_MeshEditPreview::GetActor() const
394
395   return myPreviewActor;
396 }
397
398 //================================================================================
399 /*!
400  * \brief Returns the priewed vtkUnstructuredGrid
401  */
402 //================================================================================
403
404 vtkUnstructuredGrid* SMESHGUI_MeshEditPreview::GetGrid() const
405 {
406   return myGrid;
407 }
408
409 //================================================================================
410 /*!
411  * \brief Returns myViewWindow
412  */
413 //================================================================================
414
415 SVTK_ViewWindow* SMESHGUI_MeshEditPreview::GetViewWindow() const
416 {
417   return myViewWindow;
418 }