]> SALOME platform Git repositories - modules/visu.git/blob - src/PIPELINE/VISU_DeformedShapeAndScalarMapPL.cxx
Salome HOME
Merge from V6_3_BR 06/06/2011
[modules/visu.git] / src / PIPELINE / VISU_DeformedShapeAndScalarMapPL.cxx
1 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 //  VISU DeformedShapeAndScalarMapPL
21 // File:    VISU_DeformedShapeAndScalarMapPL.cxx
22 // Author:  Eugeny Nikolaev
23 // Module : VISU
24 //
25 #include "VISU_DeformedShapeAndScalarMapPL.hxx"
26 #include "VISU_FieldTransform.hxx"
27 #include "VISU_Extractor.hxx"
28 #include "VISU_LookupTable.hxx"
29 #include "VISU_DeformedShapePL.hxx"
30 #include "VTKViewer_TransformFilter.h"
31 #include "VTKViewer_Transform.h"
32 #include "VISU_MergeFilter.hxx"
33 #include "VISU_ElnoDisassembleFilter.hxx"
34 #include "VISU_PipeLineUtils.hxx"
35 #include "SALOME_ExtractGeometry.h"
36
37 #include <vtkPlane.h>
38 #include <vtkWarpVector.h>
39 #include <vtkImplicitBoolean.h>
40 #include <vtkImplicitFunction.h>
41 #include <vtkUnstructuredGrid.h>
42 #include <vtkPointDataToCellData.h>
43 #include <vtkImplicitFunctionCollection.h>
44
45 #ifdef WNT
46 #include <float.h>
47 #define isnan _isnan
48 #endif
49
50 //----------------------------------------------------------------------------
51 vtkStandardNewMacro(VISU_DeformedShapeAndScalarMapPL)
52
53 //----------------------------------------------------------------------------
54 /*!
55  * Constructor. Creating new instances of vtkWarpVector,vtkMergeFilter,vtkUnstructuredGrid
56  * Where:
57  * \li myDeformVectors is vtkWarpVector  - deformation vectors data
58  * \li myScalarsMergeFilter   is vtkMergeFilter - merge filter.
59  * Merge filter which unify the deformation and scalars
60  * \li myScalars is vtk shared pointer to vtkUnstructuredGrid - scalars data
61 */
62 VISU_DeformedShapeAndScalarMapPL
63 ::VISU_DeformedShapeAndScalarMapPL():
64   myScaleFactor(1.0),
65   myMapScaleFactor(1.0)
66 {
67   myWarpVector = vtkWarpVector::New();
68
69   myScalarsMergeFilter = VISU_MergeFilter::New();
70   myScalarsMergeFilter->SetMergingInputs(true);
71
72   myScalarsExtractor = VISU_Extractor::New();
73
74   myScalarsFieldTransform = VISU_FieldTransform::New();
75
76   myCellDataToPointData = VISU_CellDataToPointData::New();
77   myScalarsElnoDisassembleFilter = VISU_ElnoDisassembleFilter::New();
78
79   vtkImplicitBoolean* anImplicitBoolean = vtkImplicitBoolean::New();
80   anImplicitBoolean->SetOperationTypeToIntersection();
81
82   myExtractGeometry = SALOME_ExtractGeometry::New();
83   myExtractGeometry->SetImplicitFunction(anImplicitBoolean);
84 }
85
86 //----------------------------------------------------------------------------
87 /*!
88  * Destructor.
89  * Delete all fields.
90 */
91 VISU_DeformedShapeAndScalarMapPL
92 ::~VISU_DeformedShapeAndScalarMapPL()
93 {
94   myWarpVector->Delete();
95
96   myScalarsMergeFilter->Delete();
97   
98   myScalarsExtractor->Delete();
99
100   myScalarsFieldTransform->Delete();
101
102   myCellDataToPointData->Delete();
103 }
104
105 //----------------------------------------------------------------------------
106 /*!
107  * Initial method
108  */
109 void
110 VISU_DeformedShapeAndScalarMapPL
111 ::Init()
112 {
113   Superclass::Init();
114   
115   SetScale(VISU_DeformedShapePL::GetDefaultScale(this));
116 }
117
118 //----------------------------------------------------------------------------
119 /*!
120  * Build method
121  * Building of deformation and puts result to merge filter.
122  */
123 void
124 VISU_DeformedShapeAndScalarMapPL
125 ::Build()
126 {
127   Superclass::Build();
128 }
129
130
131 //----------------------------------------------------------------------------
132 vtkDataSet* 
133 VISU_DeformedShapeAndScalarMapPL
134 ::InsertCustomPL()
135 {
136   GetMapper()->SetColorModeToMapScalars();
137   GetMapper()->ScalarVisibilityOn();
138
139   VISU::CellDataToPoint(myWarpVector,
140                         myCellDataToPointData,
141                         GetMergedInput());
142   
143   myScalars = vtkUnstructuredGrid::SafeDownCast(GetMergedInput());
144
145   UpdateScalars();
146
147   myScalarsFieldTransform->SetInput(myScalarsExtractor->GetOutput());
148
149   // Sets geometry for merge filter
150   myScalarsMergeFilter->SetGeometry(myWarpVector->GetUnstructuredGridOutput());
151
152   vtkDataSet* aScalarsDataSet = myScalarsFieldTransform->GetOutput();
153   myScalarsMergeFilter->SetScalars(aScalarsDataSet);
154   myScalarsMergeFilter->AddField("VISU_CELLS_MAPPER", aScalarsDataSet);
155   myScalarsMergeFilter->AddField("VISU_POINTS_MAPPER", aScalarsDataSet);
156
157   return myScalarsMergeFilter->GetOutput();
158 }
159
160
161 //----------------------------------------------------------------------------
162 /*!
163  *  Update method
164  */
165 void
166 VISU_DeformedShapeAndScalarMapPL
167 ::Update()
168 {
169   Superclass::Update();
170   //{
171   //  std::string aFileName = std::string(getenv("HOME"))+"/"+getenv("USER")+"-myScalarsExtractor.vtk";
172   //  VISU::WriteToFile(myScalarsExtractor->GetUnstructuredGridOutput(), aFileName);
173   //}
174   //{
175   //  std::string aFileName = std::string(getenv("HOME"))+"/"+getenv("USER")+"-myWarpVector.vtk";
176   //  VISU::WriteToFile(myWarpVector->GetUnstructuredGridOutput(), aFileName);
177   //}
178   //{
179   //  std::string aFileName = std::string(getenv("HOME"))+"/"+getenv("USER")+"-myScalarsMergeFilter.vtk";
180   //  VISU::WriteToFile(myScalarsMergeFilter->GetUnstructuredGridOutput(), aFileName);
181   //}
182 }
183
184 //----------------------------------------------------------------------------
185 unsigned long int
186 VISU_DeformedShapeAndScalarMapPL
187 ::GetMemorySize()
188 {
189   unsigned long int aSize = Superclass::GetMemorySize();
190
191   if(vtkDataSet* aDataSet = myWarpVector->GetOutput())
192     aSize += aDataSet->GetActualMemorySize() * 1024;
193   
194   if(vtkDataSet* aDataSet = myScalarsExtractor->GetOutput())
195     aSize += aDataSet->GetActualMemorySize() * 1024;
196
197   if(vtkDataSet* aDataSet = myScalarsMergeFilter->GetOutput())
198     aSize += aDataSet->GetActualMemorySize() * 1024;
199
200   if(myCellDataToPointData->GetInput())
201     if(vtkDataSet* aDataSet = myCellDataToPointData->GetOutput())
202       aSize += aDataSet->GetActualMemorySize() * 1024;
203
204   return aSize;
205 }
206
207 //----------------------------------------------------------------------------
208 /*!
209  * Update scalars method.
210  * Put scalars to merge filter.
211  */
212 void
213 VISU_DeformedShapeAndScalarMapPL
214 ::UpdateScalars()
215 {
216   vtkDataSet* aScalars = GetScalars();
217   myScalarsElnoDisassembleFilter->SetInput(aScalars);
218   myExtractGeometry->SetInput(myScalarsElnoDisassembleFilter->GetOutput());
219   myScalarsExtractor->SetInput(myExtractGeometry->GetOutput());
220
221   if(VISU::IsDataOnCells(myScalarsElnoDisassembleFilter->GetOutput()))
222     GetMapper()->SetScalarModeToUseCellData();
223   else
224     GetMapper()->SetScalarModeToUsePointData();
225 }
226
227 //----------------------------------------------------------------------------
228 /*!
229  * Copy information about pipline.
230  * Copy scale and scalars.
231  */
232 void
233 VISU_DeformedShapeAndScalarMapPL
234 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
235                 bool theIsCopyInput)
236 {
237   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
238
239   if(VISU_DeformedShapeAndScalarMapPL *aPipeLine = dynamic_cast<VISU_DeformedShapeAndScalarMapPL*>(thePipeLine)){
240      SetImplicitFunction(aPipeLine->GetImplicitFunction());
241      SetScale(aPipeLine->GetScale());
242      SetScalars(aPipeLine->GetScalars());
243   }
244 }
245
246 //----------------------------------------------------------------------------
247 /*!
248  * Set scalars.
249  * Sets vtkDataSet with scalars values to VISU_Extractor filter for scalars extraction.
250  */
251 void
252 VISU_DeformedShapeAndScalarMapPL
253 ::SetScalars(vtkDataSet *theScalars)
254 {
255   if(GetScalars() == theScalars)
256     return;
257   
258   myScalars = vtkUnstructuredGrid::SafeDownCast(theScalars);
259   UpdateScalars();
260 }
261
262 //----------------------------------------------------------------------------
263 /*!
264  * Get pointer to input scalars.
265  */
266 vtkDataSet* 
267 VISU_DeformedShapeAndScalarMapPL
268 ::GetScalars()
269 {
270   return myScalars.GetPointer();
271 }
272
273 //----------------------------------------------------------------------------
274 /*!
275  * Removes all clipping planes (for myScalars)
276  */
277 void
278 VISU_DeformedShapeAndScalarMapPL
279 ::RemoveAllClippingPlanes()
280 {
281   if(vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()){
282     vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
283     aFunction->RemoveAllItems();
284     aBoolean->Modified();
285   }
286   Superclass::RemoveAllClippingPlanes();
287 }
288
289 //----------------------------------------------------------------------------
290 /*!
291  * Removes a clipping plane (for myScalars)
292  */
293 void
294 VISU_DeformedShapeAndScalarMapPL
295 ::RemoveClippingPlane(vtkIdType theID)
296 {
297   if(vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()){
298     vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
299     if(theID >= 0 && theID < aFunction->GetNumberOfItems())
300       aFunction->RemoveItem(theID);
301     aBoolean->Modified();
302   }
303   Superclass::RemoveClippingPlane(theID);
304 }
305
306 //----------------------------------------------------------------------------
307 /*!
308  * Adds a clipping plane (for myScalars)
309  */
310 bool 
311 VISU_DeformedShapeAndScalarMapPL
312 ::AddClippingPlane(vtkPlane* thePlane)
313 {
314   if (thePlane) {
315     if (vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()) {
316       vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
317       aFunction->AddItem(thePlane);
318
319       // Check, that at least one cell present after clipping.
320       // This check was introduced because of bug IPAL8849.
321       vtkDataSet* aClippedDataSet = GetClippedInput();
322       if(aClippedDataSet->GetNumberOfCells() < 1)
323         return false;
324     }
325   }
326   return Superclass::AddClippingPlane(thePlane);
327 }
328
329 //----------------------------------------------------------------------------
330 /*!
331  * Sets implicit function of clipping
332  */
333 void
334 VISU_DeformedShapeAndScalarMapPL
335 ::SetImplicitFunction(vtkImplicitFunction *theFunction)
336 {
337   myExtractGeometry->SetImplicitFunction(theFunction);
338
339
340 //----------------------------------------------------------------------------
341 /*!
342  * Gets implicit function of clipping
343  */
344 vtkImplicitFunction * 
345 VISU_DeformedShapeAndScalarMapPL
346 ::GetImplicitFunction()
347 {
348   return myExtractGeometry->GetImplicitFunction();
349 }
350
351 //----------------------------------------------------------------------------
352 /*!
353  * Sets scale for deformed shape
354  */
355 void
356 VISU_DeformedShapeAndScalarMapPL
357 ::SetScale(vtkFloatingPointType theScale) 
358 {
359   if(VISU::CheckIsSameValue(myScaleFactor, theScale))
360     return;
361
362   myScaleFactor = theScale;
363   myWarpVector->SetScaleFactor(theScale*myMapScaleFactor);
364 }
365
366 //----------------------------------------------------------------------------
367 /*!
368  * Gets scale of deformed shape.
369  */
370 vtkFloatingPointType
371 VISU_DeformedShapeAndScalarMapPL
372 ::GetScale() 
373 {
374   return myScaleFactor;
375 }
376
377 //----------------------------------------------------------------------------
378 /*!
379  * Set scale factor of deformation.
380  */
381 void
382 VISU_DeformedShapeAndScalarMapPL
383 ::SetMapScale(vtkFloatingPointType theMapScale)
384 {
385   myMapScaleFactor = theMapScale;
386   Superclass::SetMapScale(theMapScale);
387   myWarpVector->SetScaleFactor(myScaleFactor*theMapScale);
388 }
389
390 //----------------------------------------------------------------------------
391 /*!
392  * Gets scalar mode.
393  */
394 int
395 VISU_DeformedShapeAndScalarMapPL
396 ::GetScalarMode()
397 {
398   return myScalarsExtractor->GetScalarMode();
399 }
400
401 //----------------------------------------------------------------------------
402 /*!
403  * Sets scalar mode.
404  */
405 void
406 VISU_DeformedShapeAndScalarMapPL
407 ::SetScalarMode(int theScalarMode)
408 {
409   VISU_ScalarMapPL::SetScalarMode(theScalarMode, GetScalars(), myScalarsExtractor);
410 }
411
412 //----------------------------------------------------------------------------
413 void
414 VISU_DeformedShapeAndScalarMapPL
415 ::SetScaling(int theScaling) 
416 {
417   if(GetScaling() == theScaling)
418     return;
419
420   GetBarTable()->SetScale(theScaling);
421
422   if(theScaling == VTK_SCALE_LOG10)
423     myScalarsFieldTransform->SetScalarTransform(&(VISU_FieldTransform::Log10));
424   else
425     myScalarsFieldTransform->SetScalarTransform(&(VISU_FieldTransform::Ident));
426 }
427
428
429 //----------------------------------------------------------------------------
430 void
431 VISU_DeformedShapeAndScalarMapPL
432 ::SetScalarRange(vtkFloatingPointType theRange[2])
433 {
434   if (isnan(theRange[0]) || isnan(theRange[1]))
435     throw std::runtime_error("NAN values in the presentation");
436
437   if(VISU::CheckIsSameRange(theRange, GetScalarRange()))
438     return;
439
440   myScalarsFieldTransform->SetScalarRange(theRange);
441   GetBarTable()->SetRange(theRange);
442 }
443
444
445 //----------------------------------------------------------------------------
446 vtkFloatingPointType* 
447 VISU_DeformedShapeAndScalarMapPL
448 ::GetScalarRange() 
449 {
450   return myScalarsFieldTransform->GetScalarRange();
451 }
452
453
454 //----------------------------------------------------------------------------
455 /*!
456  * Gets ranges of extracted scalars
457  * \param theRange[2] - output values
458  * \li theRange[0] - minimum value
459  * \li theRange[1] - maximum value
460  */
461 void 
462 VISU_DeformedShapeAndScalarMapPL
463 ::GetSourceRange(vtkFloatingPointType theRange[2])
464 {
465   myScalarsExtractor->Update();
466   myScalarsExtractor->GetUnstructuredGridOutput()->GetScalarRange(theRange);
467
468   if (isnan(theRange[0]) || isnan(theRange[1]))
469     throw std::runtime_error("NAN values in the presentation");
470 }
471
472
473 //----------------------------------------------------------------------------
474 void
475 VISU_DeformedShapeAndScalarMapPL
476 ::SetGaussMetric(VISU::TGaussMetric theGaussMetric) 
477 {
478   if(GetGaussMetric() == theGaussMetric)
479     return;
480
481   myScalarsExtractor->SetGaussMetric(theGaussMetric);
482 }
483
484
485 //----------------------------------------------------------------------------
486 VISU::TGaussMetric
487 VISU_DeformedShapeAndScalarMapPL
488 ::GetGaussMetric() 
489 {
490   return myScalarsExtractor->GetGaussMetric();
491 }