Salome HOME
Update Help for VISU module.
[modules/visu.git] / src / PIPELINE / VISU_StreamLinesPL.cxx
1 //  VISU OBJECT : interactive object for VISU entities implementation
2 //
3 //  Copyright (C) 2003  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.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org
21 //
22 //
23 // File:    VISU_PipeLine.cxx
24 // Author:  Alexey PETROV
25 // Module : VISU
26
27
28 #include "VISU_StreamLinesPL.hxx"
29 #include "VISU_PipeLineUtils.hxx"
30 #include "VISU_UsedPointsFilter.hxx"
31 #include "VTKViewer_GeometryFilter.h"
32
33 #include <algo.h>
34
35 #include <vtkCell.h>
36 #include <vtkPointSet.h>
37 #include <vtkStreamLine.h>
38
39 #ifdef _DEBUG_
40 static int MYDEBUG = 1;
41 #else
42 static int MYDEBUG = 0;
43 #endif
44
45 static float EPS = 1.0e-7;
46 static float aMinNbOfSteps = 1.0E+2;
47 //static float aMaxNbOfSteps = 1.0E+3;
48 static float aCoeffOfIntStep = 1.0E+1;
49
50
51 vtkStandardNewMacro(VISU_StreamLinesPL);
52
53 VISU_StreamLinesPL::VISU_StreamLinesPL(){
54   myStream = vtkStreamLine::New();
55   myCenters = vtkCellCenters::New();
56   myGeomFilter = VTKViewer_GeometryFilter::New();
57   myPointsFilter = VISU_UsedPointsFilter::New();
58   myPercents = 0.3;
59   mySource = NULL;
60 }
61
62 VISU_StreamLinesPL::~VISU_StreamLinesPL(){
63   myPointsFilter->UnRegisterAllOutputs();
64   myPointsFilter->Delete();
65
66   myCenters->UnRegisterAllOutputs();
67   myCenters->Delete();
68
69   myGeomFilter->UnRegisterAllOutputs();
70   myGeomFilter->Delete();
71
72   myStream->UnRegisterAllOutputs();
73   myStream->Delete();
74 }
75
76 void VISU_StreamLinesPL::ShallowCopy(VISU_PipeLine *thePipeLine){
77   VISU_DeformedShapePL::ShallowCopy(thePipeLine);
78   if(VISU_StreamLinesPL *aPipeLine = dynamic_cast<VISU_StreamLinesPL*>(thePipeLine))
79     SetParams(aPipeLine->GetIntegrationStep(),
80               aPipeLine->GetPropagationTime(),
81               aPipeLine->GetStepLength(),
82               aPipeLine->GetSource(),
83               aPipeLine->GetUsedPoints(),
84               aPipeLine->GetDirection());
85 }
86
87
88 float VISU_StreamLinesPL::GetNecasseryMemorySize(vtkIdType theNbOfPoints, float theStepLength,
89                                                  float thePropogationTime, float thePercents)
90 {
91   static float aStreamPointSize = sizeof(float)*15 + sizeof(vtkIdType)*2;
92   static float aStreamArraySize = aStreamPointSize*1024; // == 69632
93
94   float aNbCells = thePercents*theNbOfPoints*2.0;
95   float aNbPointsPerCell = thePropogationTime/theStepLength;
96   float aCellsSize = aNbCells*(1+aNbPointsPerCell);
97   float aPointsSize = aCellsSize*3.0*sizeof(float);
98
99   float aConnectivitySize = aCellsSize*sizeof(vtkIdType);
100   float aTypesSize = aNbCells*sizeof(char);
101   float aLocationsSize = aNbCells*sizeof(int);
102   //float aNbCellsPerPoint = aCellsSize / aNbCells - 1;
103   float aMeshSize = aPointsSize + aConnectivitySize + aTypesSize + aLocationsSize;
104
105   float anAssignedDataSize = aCellsSize*4.0*sizeof(float);
106   float anOutputDataSetSize = aMeshSize + anAssignedDataSize;
107
108   float aResult = aStreamArraySize*aNbCells + anOutputDataSetSize;
109   return aResult;
110 }
111
112 int VISU_StreamLinesPL::FindPossibleParams(vtkPointSet* theDataSet, float& theStepLength,
113                                            float& thePropogationTime, float& thePercents)
114 {
115   static float aPercentsDecrease = 3.0, aStepLengthIncrease = 9.0;
116   vtkIdType aNbOfPoints = theDataSet->GetNumberOfPoints();
117   float aSize = GetNecasseryMemorySize(aNbOfPoints,theStepLength,thePropogationTime,thePercents);
118   int isPoss = CheckAvailableMemory(aSize);
119   if(!isPoss){
120     float aMaxStepLength = max(GetMaxStepLength(theDataSet),thePropogationTime);
121     float aMinStepLength = GetMinStepLength(theDataSet);
122     float aDeltaStepLength = (aMaxStepLength - aMinStepLength)/aStepLengthIncrease;
123     for(int i = 2, aStepChanged = 1, aPerecentsChanged = 1; aStepChanged || aPerecentsChanged; i++){
124       float aStepLength = theStepLength + aDeltaStepLength;
125       if(aStepLength < aMaxStepLength) theStepLength = aStepLength;
126       else if(aStepChanged){
127         aStepLength = aMaxStepLength;
128         aStepChanged = 0;
129       }
130       float aPercents = thePercents /= aPercentsDecrease;
131       if(aPercents*aNbOfPoints > 1) thePercents = aPercents;
132       else if(aPerecentsChanged) {
133         thePercents = 1.1 / aNbOfPoints;
134         aPerecentsChanged = 0;
135       }
136       aSize = GetNecasseryMemorySize(aNbOfPoints,theStepLength,thePropogationTime,thePercents);
137       if(CheckAvailableMemory(aSize)){
138         isPoss = i;
139         break;
140       }
141     }
142   }
143   if(MYDEBUG) MESSAGE("FindPossibleParams - aSize = "<<aSize<<"; isPoss = "<<isPoss);
144   return isPoss;
145 }
146
147
148 int VISU_StreamLinesPL::SetParams(float theIntStep,
149                                   float thePropogationTime,
150                                   float theStepLength,
151                                   vtkPointSet* theSource,
152                                   float thePercents,
153                                   int theDirection,
154                                   int isOnlyTry)
155 {
156   vtkPointSet* aDataSet = theSource? theSource: myFieldTransform->GetUnstructuredGridOutput();
157   aDataSet->Update();
158   vtkIdType aNbOfPoints = aDataSet->GetNumberOfPoints();
159   vtkPointSet* aPointSet = myExtractor->GetOutput();
160   if(thePercents*aNbOfPoints < 1) thePercents = 2.0/aNbOfPoints;
161   theIntStep = CorrectIntegrationStep(theIntStep,aPointSet,thePercents);
162   thePropogationTime = CorrectPropagationTime(thePropogationTime,aPointSet);
163   theStepLength = CorrectStepLength(theStepLength,aPointSet);
164   int isAccepted = FindPossibleParams(aPointSet,theStepLength,thePropogationTime,thePercents);
165   if((!isOnlyTry && isAccepted) || (isOnlyTry && isAccepted == 1)){
166     mySource = theSource;
167     myPercents = thePercents;
168     if(GetInput2()->GetCellData()->GetNumberOfArrays()){
169       myCenters->SetInput(aDataSet);
170       myCenters->VertexCellsOn();
171       aDataSet = myCenters->GetOutput();
172     }
173     myPointsFilter->SetInput(aDataSet);
174     myPointsFilter->SetPercentsOfUsedPoints(thePercents);
175     aDataSet = myPointsFilter->GetOutput();
176     myStream->SetSource(aDataSet);
177     myStream->SetIntegrationStepLength(theIntStep);
178     myStream->SetMaximumPropagationTime(thePropogationTime);
179     myStream->SetStepLength(theStepLength);
180     myStream->SetSavePointInterval(theIntStep*aMinNbOfSteps);
181     myStream->SetIntegrationDirection(theDirection);
182     myStream->Modified();
183     Modified();
184   }
185   return isAccepted;
186 }
187
188
189 vtkPointSet* VISU_StreamLinesPL::GetSource() {
190   return mySource;
191 }
192
193 float VISU_StreamLinesPL::GetUsedPoints() {
194   return myPercents;
195 }
196
197 vtkDataSet* VISU_StreamLinesPL::GetStreamerSource(){
198   return myStream->GetSource();
199 }
200
201 float VISU_StreamLinesPL::GetVelocityCoeff(){
202   return GetVelocityCoeff(myExtractor->GetOutput());
203 }
204
205 float VISU_StreamLinesPL::GetVelocityCoeff(vtkPointSet* theDataSet){
206   float* aScalarRange = theDataSet->GetScalarRange();
207   return (fabs(aScalarRange[1]) + fabs(aScalarRange[0]))/2.0;
208 }
209
210
211 int VISU_StreamLinesPL::IsPossible(vtkPointSet* theDataSet, float thePercents){
212   float aStepLength = GetBaseStepLength(theDataSet);
213   float aBasePropTime = GetBasePropagationTime(theDataSet);
214   VISU_UsedPointsFilter *aPointsFilter = VISU_UsedPointsFilter::New();
215   aPointsFilter->SetInput(theDataSet);
216   vtkPointSet* aDataSet = aPointsFilter->GetOutput();
217   aDataSet->Update();
218   int aRes = FindPossibleParams(aDataSet,aStepLength,aBasePropTime,thePercents);
219   aPointsFilter->UnRegisterAllOutputs();
220   aPointsFilter->Delete();
221   return aRes;
222 }
223
224
225 float VISU_StreamLinesPL::GetIntegrationStep(){
226   return myStream->GetIntegrationStepLength();
227 }
228 float VISU_StreamLinesPL::GetStepLength() {
229   return myStream->GetStepLength();
230 }
231 float VISU_StreamLinesPL::GetPropagationTime() {
232   return myStream->GetMaximumPropagationTime();
233 }
234 int VISU_StreamLinesPL::GetDirection(){
235   return myStream->GetIntegrationDirection();
236 }
237
238
239 float VISU_StreamLinesPL::GetMinIntegrationStep(vtkPointSet* theDataSet, float thePercents) {
240   if(!theDataSet) return -1.0;
241   float aVolume = 1.0;
242   int degree = 0;
243   theDataSet->Update();
244   float* aBounds = theDataSet->GetBounds();
245   for(int i = 0, j = 0; i < 3; ++i, j = 2*i){
246     float tmp = aBounds[j+1] - aBounds[j];
247     if (tmp > EPS ) {
248       aVolume *= tmp;
249       degree += 1;
250     }
251   }
252   if (degree < 1) return 0.0; // absolutely empty object
253   float anStepLength = GetMaxIntegrationStep(theDataSet)/aCoeffOfIntStep;
254   float aBasePropTime = GetBasePropagationTime(theDataSet)/GetVelocityCoeff(theDataSet);
255   thePercents = 1.0;
256   vtkIdType aNbOfPoints = theDataSet->GetNumberOfPoints();
257   float aSize = GetNecasseryMemorySize(aNbOfPoints,anStepLength,aBasePropTime,thePercents);
258   float aRealSize = GetAvailableMemory(aSize);
259   float anAverageVolume = aVolume / aRealSize;
260   float aStep = pow(double(anAverageVolume), double(1.0/double(degree)));
261   return aStep;
262 }
263 float VISU_StreamLinesPL::GetMinIntegrationStep(){
264   return GetMinIntegrationStep(myExtractor->GetOutput(),GetUsedPoints());
265 }
266
267
268 float VISU_StreamLinesPL::GetMaxIntegrationStep(vtkPointSet* theDataSet) {
269   if(!theDataSet) return -1.0;
270   theDataSet->Update();
271   float aLength = theDataSet->GetLength();
272   float* aBounds = theDataSet->GetBounds();
273   float aMaxSizeY = (aBounds[3]-aBounds[2])/aLength;
274   float aMaxSizeZ = (aBounds[5]-aBounds[4])/aLength;
275   float aMinMax = (aBounds[1] - aBounds[0])/aLength;
276   if (aMinMax < EPS || (aMaxSizeY < aMinMax && aMaxSizeY > EPS)) aMinMax = aMaxSizeY;
277   if (aMinMax < EPS || (aMaxSizeZ < aMinMax && aMaxSizeZ > EPS)) aMinMax = aMaxSizeZ;
278   return aMinMax*aLength/2.0;
279 }
280 float VISU_StreamLinesPL::GetMaxIntegrationStep(){
281   return GetMaxIntegrationStep(myExtractor->GetOutput());
282 }
283
284 float VISU_StreamLinesPL::GetBaseIntegrationStep(vtkPointSet* theDataSet, float thePercents){
285   theDataSet->Update();
286   float aMinIntegrationStep = GetMinIntegrationStep(theDataSet,thePercents);
287   float aMaxIntegrationStep = GetMaxIntegrationStep(theDataSet);
288   float anIntegrationStep = aMaxIntegrationStep / aCoeffOfIntStep;
289   float aMinMax = theDataSet->GetLength()/theDataSet->GetNumberOfPoints();
290   if(aMinMax > anIntegrationStep)
291     anIntegrationStep = (anIntegrationStep*aCoeffOfIntStep*0.9+aMinMax)/aCoeffOfIntStep;
292   if(aMinIntegrationStep > anIntegrationStep)
293     anIntegrationStep = aMinIntegrationStep;
294   return anIntegrationStep;
295 }
296
297 float VISU_StreamLinesPL::CorrectIntegrationStep(float theStep, vtkPointSet* theDataSet, float thePercents){
298   theDataSet->Update();
299   float aMinIntegrationStep = GetMinIntegrationStep(theDataSet,thePercents);
300   float aMaxIntegrationStep = GetMaxIntegrationStep(theDataSet);
301   if(aMinIntegrationStep > theStep)
302     theStep = aMinIntegrationStep;
303   if(aMaxIntegrationStep < theStep)
304     theStep = aMaxIntegrationStep;
305   return theStep;
306 }
307
308
309 float VISU_StreamLinesPL::GetMinPropagationTime(vtkPointSet* theDataSet){
310   if(!theDataSet) return -1.0;
311   return GetMinStepLength(theDataSet);
312 }
313 float VISU_StreamLinesPL::GetMinPropagationTime(){
314   return GetMinPropagationTime(myExtractor->GetOutput());
315 }
316
317 float VISU_StreamLinesPL::GetMaxPropagationTime(vtkPointSet* theDataSet){
318   if(!theDataSet) return -1.0;
319   return GetBasePropagationTime(theDataSet)*aMinNbOfSteps;
320 }
321 float VISU_StreamLinesPL::GetMaxPropagationTime(){
322   return GetMaxPropagationTime(myExtractor->GetOutput());
323 }
324
325 float VISU_StreamLinesPL::CorrectPropagationTime(float thePropagationTime, vtkPointSet* theDataSet){
326   float aMinPropagationTime = GetMinPropagationTime(theDataSet);
327   float aMaxPropagationTime = GetMaxPropagationTime(theDataSet);
328   if(aMinPropagationTime > thePropagationTime)
329     thePropagationTime = aMinPropagationTime;
330   if(aMaxPropagationTime < thePropagationTime)
331     thePropagationTime = aMaxPropagationTime;
332   return thePropagationTime;
333 }
334 float VISU_StreamLinesPL::GetBasePropagationTime(vtkPointSet* theDataSet){
335   if(!theDataSet) return -1.0;
336   theDataSet->Update();
337   float aPropagationTime = theDataSet->GetLength()/GetVelocityCoeff(theDataSet);
338   return aPropagationTime;
339 }
340 float VISU_StreamLinesPL::GetBasePropagationTime(){
341   return GetBasePropagationTime(myExtractor->GetOutput());
342 }
343
344
345 float VISU_StreamLinesPL::GetMinStepLength(vtkPointSet* theDataSet){
346   static float aNbOfStepsOfIntStep = 1.0E+1;
347   float anIntStep = GetMinIntegrationStep(theDataSet);
348   float aStepLength = anIntStep*aNbOfStepsOfIntStep/GetVelocityCoeff(theDataSet);
349   return aStepLength;
350 }
351 float VISU_StreamLinesPL::GetMinStepLength(){
352   return GetMinStepLength(myExtractor->GetOutput());
353 }
354
355 float VISU_StreamLinesPL::GetMaxStepLength(vtkPointSet* theDataSet){
356   float aStepLength = GetBasePropagationTime(theDataSet);
357   return aStepLength;
358 }
359 float VISU_StreamLinesPL::GetMaxStepLength(){
360   return GetMaxStepLength(myExtractor->GetOutput());
361 }
362
363 float VISU_StreamLinesPL::CorrectStepLength(float theStep, vtkPointSet* theDataSet){
364   float aMinStep = GetMinStepLength(theDataSet);
365   if(theStep < aMinStep) theStep = aMinStep;
366   float aMaxStep = GetMaxStepLength(theDataSet);
367   if(theStep > aMaxStep) theStep = aMaxStep;
368   return theStep;
369 }
370 float VISU_StreamLinesPL::GetBaseStepLength(vtkPointSet* theDataSet){
371   static float anAvgNbOfSteps = 1.0E+2;
372   float aPropagationTime = GetBasePropagationTime(theDataSet);
373   float aStepLength = aPropagationTime/anAvgNbOfSteps;
374   aStepLength = CorrectStepLength(aStepLength,theDataSet);
375   return aStepLength;
376 }
377
378
379 void VISU_StreamLinesPL::Init(){
380   VISU_ScalarMapPL::Init();
381   vtkPointSet* aDataSet = myExtractor->GetOutput();
382   float anIntStep = GetBaseIntegrationStep(aDataSet);
383   float aPropagationTime = GetBasePropagationTime(aDataSet);
384   float aStepLength = GetBaseStepLength(aDataSet);
385   SetParams(anIntStep,aPropagationTime,aStepLength);
386 }
387
388 VISU_ScalarMapPL::THook* VISU_StreamLinesPL::DoHook(){
389   GetInput2()->Update();
390   VISU::CellDataToPoint(myStream,myCellDataToPointData,GetInput2(),myFieldTransform);
391   float *aBounds = GetInput2()->GetBounds();
392   myGeomFilter->SetExtent(aBounds);
393   myGeomFilter->ExtentClippingOn();
394   myGeomFilter->SetInput(myStream->GetOutput());
395   return myGeomFilter->GetOutput();
396 }
397
398 void VISU_StreamLinesPL::Update(){
399   VISU_ScalarMapPL::Update();
400 }
401
402 void VISU_StreamLinesPL::SetMapScale(float theMapScale){
403   VISU_ScalarMapPL::SetMapScale(theMapScale);
404 }