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