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