]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_ConvertorUtils.cxx
Salome HOME
Fix for the "0051899: curves are not shown in opened study" issue.
[modules/visu.git] / src / CONVERTOR / VISU_ConvertorUtils.cxx
1 // Copyright (C) 2007-2013  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.
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 //  VISU OBJECT : interactive object for VISU entities implementation
24 //  File   : VISU_Convertor_impl.cxx
25 //  Author : Alexey PETROV
26 //  Module : VISU
27 //
28 #include "VISU_ConvertorUtils.hxx"
29
30 #include <vtkCellType.h>
31
32 #include <vtkUnstructuredGridWriter.h>
33 #include <vtkUnstructuredGrid.h>
34
35 #include <vtkPolyDataWriter.h>
36 #include <vtkPolyData.h>
37
38 #include <vtkInformationVector.h>
39 #include <vtkInformation.h>
40 #include <vtkExecutive.h>
41
42 #include <vtkTimerLog.h>
43 #include <vtkPointData.h>
44 #include <vtkCellData.h>
45 #include <vtkDataSet.h>
46 #include <vtkIdList.h>
47
48 #include <vtkIntArray.h>
49 #include <algorithm>
50
51 #ifdef _DEBUG_
52 static int MYDEBUG = 0;
53 #else
54 static int MYDEBUG = 0;
55 #endif
56
57 namespace VISU
58 {
59   //---------------------------------------------------------------
60   vtkIdType
61   VISUGeom2NbNodes(EGeometry theGeom)
62   { 
63     switch(theGeom){
64 #ifndef VISU_ENABLE_QUADRATIC
65     case VISU::eSEG3: 
66       return 2;
67     case VISU::eTRIA6: 
68       return 3;
69     case VISU::eQUAD8: 
70       return 4;
71     case VISU::eQUAD9:
72       return 4;
73     case VISU::eTETRA10: 
74       return 4;
75     case VISU::eHEXA20: 
76       return 8;
77     case VISU::eHEXA27: 
78       return 8;
79     case VISU::ePENTA15: 
80       return 6;
81     case VISU::ePYRA13: 
82       return 5;
83 #endif
84     case VISU::ePOLYGONE: 
85     case VISU::ePOLYEDRE: 
86       return -1;
87     default:
88       return theGeom % 100;
89     }
90   }
91
92
93   //---------------------------------------------------------------
94   vtkIdType
95   VISUGeom2VTK(EGeometry theGeom)
96   { 
97     switch(theGeom){
98     case VISU::ePOINT1: 
99       return VTK_VERTEX;
100     case VISU::eSEG2: 
101       return VTK_LINE;
102     case VISU::eTRIA3: 
103       return VTK_TRIANGLE;
104     case VISU::eQUAD4: 
105       return VTK_QUAD;
106     case VISU::eTETRA4: 
107       return VTK_TETRA;
108     case VISU::eHEXA8: 
109       return VTK_HEXAHEDRON;
110     case VISU::ePENTA6: 
111       return VTK_WEDGE;
112     case VISU::ePYRA5: 
113       return VTK_PYRAMID;
114     case VISU::eOCTA12: 
115       return VTK_HEXAGONAL_PRISM;
116     case VISU::ePOLYGONE: 
117       return VTK_POLYGON;
118     case VISU::ePOLYEDRE: 
119       return VTK_CONVEX_POINT_SET;
120
121 #ifndef VISU_ENABLE_QUADRATIC
122     case VISU::eSEG3: 
123       return VTK_LINE;
124     case VISU::eTRIA6: 
125       return VTK_TRIANGLE;
126     case VISU::eQUAD8: 
127       return VTK_QUAD;
128     case VISU::eQUAD9: 
129       return VTK_QUAD;
130     case VISU::eTETRA10: 
131       return VTK_TETRA;
132     case VISU::eHEXA20: 
133       return VTK_HEXAHEDRON;
134     case VISU::eHEXA27: 
135       return VTK_TRIQUADRATIC_HEXAHEDRON;
136     case VISU::ePENTA15: 
137       return VTK_WEDGE;
138     case VISU::ePYRA13: 
139       return VTK_PYRAMID;
140
141 #else
142
143     case VISU::eSEG3: 
144 #if defined(VISU_USE_VTK_QUADRATIC)
145       return VTK_QUADRATIC_EDGE;
146 #else
147       return VTK_POLY_LINE;
148 #endif
149
150     case VISU::eTRIA6: 
151 #if defined(VISU_USE_VTK_QUADRATIC)
152       return VTK_QUADRATIC_TRIANGLE;
153 #else
154       return VTK_POLYGON;
155 #endif
156
157     case VISU::eQUAD8: 
158 #if defined(VISU_USE_VTK_QUADRATIC)
159       return VTK_QUADRATIC_QUAD;
160 #else
161       return VTK_POLYGON;
162 #endif
163
164     case VISU::eQUAD9:
165 #if defined(VISU_USE_VTK_QUADRATIC)
166       return VTK_BIQUADRATIC_QUAD;
167 #else
168       return VTK_POLYGON;
169 #endif
170
171     case VISU::eTETRA10: 
172 #if defined(VISU_USE_VTK_QUADRATIC)
173       return VTK_QUADRATIC_TETRA;
174 #else
175       return VTK_CONVEX_POINT_SET;
176 #endif
177
178     case VISU::eHEXA20: 
179 #if defined(VISU_USE_VTK_QUADRATIC)
180       return VTK_QUADRATIC_HEXAHEDRON;
181 #else
182       return VTK_CONVEX_POINT_SET;
183 #endif
184
185     case VISU::eHEXA27: 
186 #if defined(VISU_USE_VTK_QUADRATIC)
187       return VTK_TRIQUADRATIC_HEXAHEDRON;
188 #else
189       return VTK_CONVEX_POINT_SET;
190 #endif
191
192     case VISU::ePENTA15: 
193 #if defined(VISU_USE_VTK_QUADRATIC)
194       return VTK_QUADRATIC_WEDGE;
195 #else
196       return VTK_CONVEX_POINT_SET;
197 #endif
198
199     case VISU::ePYRA13: 
200 #if defined(VISU_USE_VTK_QUADRATIC)
201       return VTK_QUADRATIC_PYRAMID;
202 #else
203       return VTK_CONVEX_POINT_SET;
204 #endif
205
206 #endif //VISU_ENABLE_QUADRATIC
207
208     default:
209       return -1;
210     }
211   }
212
213
214   //---------------------------------------------------------------
215   void 
216   WriteToFile(vtkUnstructuredGrid* theDataSet, 
217               const std::string& theFileName)
218   {
219     vtkUnstructuredGridWriter* aWriter = vtkUnstructuredGridWriter::New();
220     //aWriter->SetFileType(VTK_BINARY);
221     aWriter->SetFileName(theFileName.c_str());
222     aWriter->SetInputData(theDataSet);
223     aWriter->Write();
224     aWriter->Delete();
225   }
226
227
228   //---------------------------------------------------------------
229   void 
230   WriteToFile(vtkPolyData* theDataSet, 
231               const std::string& theFileName)
232   {
233     vtkPolyDataWriter* aWriter = vtkPolyDataWriter::New();
234     //aWriter->SetFileType(VTK_BINARY);
235     aWriter->SetFileName(theFileName.c_str());
236     aWriter->SetInputData(theDataSet);
237     aWriter->Write();
238     aWriter->Delete();
239   }
240
241
242   //---------------------------------------------------------------
243   bool 
244   IsDataOnPoints(vtkDataSet* theDataSet)
245   {
246     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData();
247     return aDataSetAttributes->GetArray("VISU_FIELD") != NULL;
248   }
249
250
251   //---------------------------------------------------------------
252   bool 
253   IsDataOnCells(vtkDataSet* theDataSet)
254   {
255     //theDataSet->Update();
256     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
257     return aDataSetAttributes->GetArray("VISU_FIELD") != NULL;
258   }
259
260   //---------------------------------------------------------------
261   bool 
262   IsElnoData(vtkDataSet* theDataSet)
263   {
264     //theDataSet->Update();
265
266     if ( vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData() )
267       if ( aDataSetAttributes->GetArray( "ELNO_FIELD" )  != NULL )
268         return true;
269
270     if ( vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData() )
271       if ( aDataSetAttributes->GetArray( "ELNO_POINT_COORDS" )  != NULL )
272         return true;
273
274     return false;
275   }
276
277
278   //---------------------------------------------------------------
279   vtkIdType
280   GetVTKID(vtkDataArray *theIDDataArray, vtkIdType theID, int theEntity)
281   {
282     if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(theIDDataArray)){
283       int aNbTuples = anIntArray->GetNumberOfTuples();
284       int* aPointer = anIntArray->GetPointer(0);
285       for(int aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
286         if(*aPointer == theID && *(aPointer + 1) == theEntity){
287           return aTupleId;
288         }
289         aPointer += 2;
290       }
291     }
292     return -1;
293   }
294
295
296   //---------------------------------------------------------------
297   vtkIdType
298   GetObjectID(vtkDataArray *theIDDataArray, vtkIdType theID)
299   {
300     if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(theIDDataArray)){
301       int aNbComp = anIntArray->GetNumberOfComponents();
302       int* aPointer = anIntArray->GetPointer(theID*aNbComp);
303       return *aPointer;
304     }
305     return -1;
306   }
307
308
309   //---------------------------------------------------------------
310   vtkIdType
311   GetElemVTKID(vtkDataSet *theDataSet, vtkIdType theID, int theEntity)
312   {
313     //theDataSet->Update();
314     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
315     if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_CELLS_MAPPER")){
316       if(theEntity < 0){
317         {
318           vtkIdType anID = GetVTKID(aDataArray, theID, VISU::CELL_ENTITY);
319           if(anID != -1)
320             return anID;
321         }
322         {
323           vtkIdType anID = GetVTKID(aDataArray, theID, VISU::FACE_ENTITY);
324           if(anID != -1)
325             return anID;
326         }
327         {
328           vtkIdType anID = GetVTKID(aDataArray, theID, VISU::EDGE_ENTITY);
329           if(anID != -1)
330             return anID;
331         }
332         {
333           vtkIdType anID = GetVTKID(aDataArray, theID, VISU::NODE_ENTITY);
334           if(anID != -1)
335             return anID;
336         }
337       }else
338         return GetVTKID(aDataArray, theID, theEntity);
339     }
340     return -1;
341   }
342
343
344   //---------------------------------------------------------------
345   vtkIdType
346   GetElemObjID(vtkDataSet *theDataSet, vtkIdType theID)
347   {
348     //theDataSet->Update();
349     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
350     if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_CELLS_MAPPER"))
351       return GetObjectID(aDataArray, theID);
352
353     return -1;
354   }
355
356
357   //---------------------------------------------------------------
358   vtkCell* 
359   GetElemCell(vtkDataSet *theDataSet, vtkIdType  theObjID)
360   {
361     vtkIdType aVTKID = GetElemVTKID(theDataSet, theObjID);
362     return theDataSet->GetCell(aVTKID);
363   }
364
365
366   //---------------------------------------------------------------
367   vtkIdType
368   GetNodeVTKID(vtkDataSet *theDataSet, vtkIdType theID)
369   {
370     //theDataSet->Update();
371     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData();
372     if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_POINTS_MAPPER"))
373       return GetVTKID(aDataArray, theID, VISU::NODE_ENTITY);
374
375     return -1;
376   }
377
378
379   //---------------------------------------------------------------
380   vtkIdType
381   GetNodeObjID(vtkDataSet *theDataSet, vtkIdType theID)
382   {
383     //theDataSet->Update();
384     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData();
385     if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_POINTS_MAPPER"))
386       return GetObjectID(aDataArray, theID);
387
388     return -1;
389   }
390
391
392   //---------------------------------------------------------------
393   double* 
394   GetNodeCoord(vtkDataSet *theDataSet, vtkIdType theObjID)
395   {
396     vtkIdType aVTKID = GetNodeVTKID(theDataSet, theObjID);
397     if ( aVTKID >= 0 ) return theDataSet->GetPoint(aVTKID);
398     return 0;
399   }
400
401
402   //---------------------------------------------------------------
403   TGaussPointID
404   GetObjID(vtkDataSet *theDataSet, vtkIdType theID)
405   {
406     //theDataSet->Update();
407     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
408     if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_CELLS_MAPPER")){
409       if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
410         vtkIdType anID = 2 * theID;
411         TCellID aCellID = anIntArray->GetValue(anID);
412         TLocalPntID aLocalPntID = anIntArray->GetValue(anID + 1);
413         return TGaussPointID(aCellID, aLocalPntID);
414       }
415     }
416     return TGaussPointID();
417   }
418
419
420   //---------------------------------------------------------------
421   TInputCellID
422   GetInputCellID(vtkDataSet *theDataSet, vtkIdType theObjID)
423   {
424     //theDataSet->Update();
425     vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetCellData();
426     if(vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_INPUTS_MAPPER")){
427       if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
428         vtkIdType aVTKID = GetElemVTKID(theDataSet, theObjID);
429         vtkIdType aTupleID = 2 * aVTKID;
430         TCellID aCellID = anIntArray->GetValue(aTupleID);
431         TInputID anInputID = anIntArray->GetValue(aTupleID + 1);
432         return TInputCellID(anInputID, aCellID);
433       }
434     }
435     return TInputCellID();
436   }
437
438
439   //---------------------------------------------------------------
440   vtkDataSet*
441   GetInput(vtkInformationVector **theInputVector, 
442            vtkIdType theInputId)
443   {
444     if(vtkInformation* anInformation = theInputVector[0]->GetInformationObject(theInputId))
445       return vtkDataSet::SafeDownCast(anInformation->Get(vtkDataObject::DATA_OBJECT()));
446     return NULL;
447   }
448
449
450   //---------------------------------------------------------------
451   vtkDataSet*
452   GetOutput(vtkInformationVector *theOutputVector)
453   {
454     if(vtkInformation* anInformation = theOutputVector->GetInformationObject(0))
455       return vtkDataSet::SafeDownCast(anInformation->Get(vtkDataObject::DATA_OBJECT()));
456     return NULL;
457   }
458
459   //---------------------------------------------------------------
460   TElnoPoints
461   GetElnoPoints(vtkDataSet *theDataSet, vtkIdType theNodeObjID)
462   {
463     TElnoPoints aResult;
464     if(theDataSet && IsElnoData(theDataSet)) {
465       vtkDataSetAttributes *aDataSetAttributes = theDataSet->GetPointData();
466       vtkDataArray *aDataArray = aDataSetAttributes->GetArray("VISU_POINTS_MAPPER");
467       if(aDataArray){
468         if(vtkIntArray *anIntArray = dynamic_cast<vtkIntArray*>(aDataArray)){
469           int aNbTuples = anIntArray->GetNumberOfTuples();
470           int* aPointer = anIntArray->GetPointer(0);
471           for(int aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
472             if( *aPointer == theNodeObjID ) {
473               vtkIdList *aCellIds = vtkIdList::New();
474               theDataSet->GetPointCells(aTupleId,aCellIds);
475               if(aCellIds->GetNumberOfIds() == 1){
476                 aResult.push_back(TElnoPointID(aTupleId,aCellIds->GetId(0)));
477               }
478             }
479             aPointer += 2;
480           }
481         }
482       }
483     }
484     return aResult;
485   }
486
487   //---------------------------------------------------------------
488   TTimerLog
489   ::TTimerLog(int theIsDebug,
490               const std::string& theName):
491     myIsDebug(MYDEBUG + theIsDebug),
492     myTimerLog(vtkTimerLog::New()),
493     myPrefixPrinter(myIsDebug == 1),
494     myName(theName)
495   {
496     myCPUTime = myTimerLog->GetCPUTime();
497     BEGMSG(myIsDebug > 1,"{\n");
498   }
499
500   //---------------------------------------------------------------
501   TTimerLog
502   ::~TTimerLog()
503   {
504     myCPUTime = myTimerLog->GetCPUTime() - myCPUTime;
505
506     if(myIsDebug > 1){
507       BEGMSG(myIsDebug,"} = "<<myCPUTime<<" secs ("<<myName<<")\n");
508     }else{
509       BEGMSG(myIsDebug,myName<<" takes "<<myCPUTime<<" secs\n");
510     }
511     
512     myTimerLog->Delete();
513     myTimerLog = NULL;
514   }
515
516
517   //---------------------------------------------------------------
518 }