]> SALOME platform Git repositories - modules/visu.git/blob - src/PIPELINE/VISU_DeformedShapePL.cxx
Salome HOME
b3c93c87467bdae7534b9e5246205bad8636d0eb
[modules/visu.git] / src / PIPELINE / VISU_DeformedShapePL.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 OBJECT : interactive object for VISU entities implementation
23 // File:    VISU_PipeLine.cxx
24 // Author:  Alexey PETROV
25 // Module : VISU
26 //
27 #include "VISU_DeformedShapePL.hxx"
28 #include "VISU_PipeLineUtils.hxx"
29 #include "VTKViewer_Transform.h"
30
31 #include <vtkWarpVector.h>
32
33
34 //----------------------------------------------------------------------------
35 vtkStandardNewMacro(VISU_DeformedShapePL);
36
37
38 //----------------------------------------------------------------------------
39 VISU_DeformedShapePL
40 ::VISU_DeformedShapePL():
41   myScaleFactor(0.0),
42   myMapScaleFactor(1.0)
43 {
44   SetIsShrinkable(true);
45   SetIsFeatureEdgesAllowed(true);
46
47   myWarpVector = vtkWarpVector::New();
48   myCellDataToPointData = vtkCellDataToPointData::New();
49 }
50
51
52 //----------------------------------------------------------------------------
53 VISU_DeformedShapePL
54 ::~VISU_DeformedShapePL()
55 {
56   myWarpVector->Delete();
57
58   myCellDataToPointData->Delete();
59 }
60
61
62 //----------------------------------------------------------------------------
63 unsigned long int 
64 VISU_DeformedShapePL
65 ::GetMTime()
66 {
67   unsigned long int aTime = Superclass::GetMTime();
68
69   aTime = std::max(aTime, myWarpVector->GetMTime());
70   aTime = std::max(aTime, myCellDataToPointData->GetMTime());
71
72   return aTime;
73 }
74
75
76 //----------------------------------------------------------------------------
77 void
78 VISU_DeformedShapePL
79 ::DoShallowCopy(VISU_PipeLine *thePipeLine,
80                 bool theIsCopyInput)
81 {
82   Superclass::DoShallowCopy(thePipeLine, theIsCopyInput);
83
84   if(VISU_DeformedShapePL *aPipeLine = dynamic_cast<VISU_DeformedShapePL*>(thePipeLine)){
85     SetScale(aPipeLine->GetScale());
86   }
87 }
88
89
90 //----------------------------------------------------------------------------
91 vtkFloatingPointType
92 VISU_DeformedShapePL
93 ::GetScaleFactor(vtkDataSet* theDataSet)
94 {
95   if(!theDataSet)
96     return 0.0;
97
98   theDataSet->Update();
99
100   int aNbCells = theDataSet->GetNumberOfCells();
101   int aNbPoints = theDataSet->GetNumberOfPoints();
102   int aNbElem = aNbCells? aNbCells: aNbPoints;
103
104   vtkFloatingPointType* aBounds = theDataSet->GetBounds();
105   vtkFloatingPointType aVolume = 1, aVol, idim = 0;
106   for(int i = 0; i < 6; i += 2){
107     aVol = fabs(aBounds[i+1] - aBounds[i]);
108     if(aVol > 0) {
109       idim++;
110       aVolume *= aVol;
111     }
112   }
113   if( aNbElem == 0 || fabs(idim) < 1.0 / VTK_LARGE_FLOAT )
114     return 0.0; // to avoid division by zero
115   aVolume /= aNbElem;
116   return pow(aVolume, vtkFloatingPointType(1.0/idim));
117 }
118
119
120 //----------------------------------------------------------------------------
121 vtkFloatingPointType
122 VISU_DeformedShapePL
123 ::GetDefaultScale(VISU_ScalarMapPL* theScalarMapPL)
124 {
125   vtkFloatingPointType aSourceRange[2];
126   theScalarMapPL->GetSourceRange(aSourceRange);
127   
128   static vtkFloatingPointType EPS = 1.0 / VTK_LARGE_FLOAT;
129   if(fabs(aSourceRange[1]) > EPS){
130     vtkDataSet* aDataSet = theScalarMapPL->GetMergedInput();
131     vtkFloatingPointType aScaleFactor = VISU_DeformedShapePL::GetScaleFactor(aDataSet);
132     return aScaleFactor / aSourceRange[1];
133   }
134   return 0.0;
135 }
136
137
138 //----------------------------------------------------------------------------
139 void
140 VISU_DeformedShapePL
141 ::SetScale(vtkFloatingPointType theScale) 
142 {
143   if(VISU::CheckIsSameValue(myWarpVector->GetScaleFactor(), theScale))
144     return;
145   
146   myWarpVector->SetScaleFactor(theScale*myMapScaleFactor);
147   myScaleFactor = theScale;
148 }
149
150
151 //----------------------------------------------------------------------------
152 vtkFloatingPointType
153 VISU_DeformedShapePL
154 ::GetScale() 
155 {
156   return myScaleFactor;
157 }
158
159
160 //----------------------------------------------------------------------------
161 void
162 VISU_DeformedShapePL
163 ::Init()
164 {
165   Superclass::Init();
166
167   SetScale(VISU_DeformedShapePL::GetDefaultScale(this));
168 }
169
170
171 //----------------------------------------------------------------------------
172 void
173 VISU_DeformedShapePL
174 ::Update()
175 {
176   Superclass::Update();
177   //{
178   //  std::string aFileName = std::string(getenv("HOME"))+"/"+getenv("USER")+"-myWarpVector.vtk";
179   //  VISU::WriteToFile(myWarpVector->GetUnstructuredGridOutput(), aFileName);
180   //}
181 }
182
183
184 //----------------------------------------------------------------------------
185 vtkDataSet* 
186 VISU_DeformedShapePL
187 ::InsertCustomPL()
188 {
189   VISU::CellDataToPoint(myWarpVector,
190                         myCellDataToPointData,
191                         GetMergedInput());
192
193   return myWarpVector->GetOutput();
194 }
195
196
197 //----------------------------------------------------------------------------
198 unsigned long int
199 VISU_DeformedShapePL
200 ::GetMemorySize()
201 {
202   unsigned long int aSize = Superclass::GetMemorySize();
203
204   if(myWarpVector->GetInput())
205     if(vtkDataSet* aDataSet = myWarpVector->GetOutput())
206       aSize += aDataSet->GetActualMemorySize() * 1024;
207   
208   if(myCellDataToPointData->GetInput())
209     if(vtkDataSet* aDataSet = myCellDataToPointData->GetOutput())
210       aSize += aDataSet->GetActualMemorySize() * 1024;
211
212   return aSize;
213 }
214
215
216 //----------------------------------------------------------------------------
217 void
218 VISU_DeformedShapePL
219 ::SetMapScale(vtkFloatingPointType theMapScale)
220 {
221   myMapScaleFactor = theMapScale;
222   Superclass::SetMapScale(theMapScale);
223
224   vtkFloatingPointType aMapScale = myScaleFactor * theMapScale;
225   if(VISU::CheckIsSameValue(myWarpVector->GetScaleFactor(), aMapScale))
226     return;
227
228   myWarpVector->SetScaleFactor( aMapScale );
229 }
230
231
232 //----------------------------------------------------------------------------