]> SALOME platform Git repositories - modules/visu.git/blob - src/PIPELINE/VISU_StreamLinesPL.cxx
Salome HOME
This commit was generated by cvs2git to track changes on a CVS vendor
[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 "SALOME_GeometryFilter.h"
32
33 #include <vtkPointSet.h>
34 #include <vtkStreamLine.h>
35
36 #ifdef _DEBUG_
37 static int MYDEBUG = 0;
38 #else
39 static int MYDEBUG = 0;
40 #endif
41
42 static float EPS = 1.0e-20;
43 int VISU_StreamLinesPL::myMaxIncrementMemorySize = 64;
44
45 void SetStreamerSource(vtkStreamer *theStreamer, vtkPointSet* theDataSet, float thePercents){
46   theDataSet->Update();
47   VISU_UsedPointsFilter* anUsedPointsFilter = VISU_UsedPointsFilter::New();
48   anUsedPointsFilter->SetInput(theDataSet);
49   anUsedPointsFilter->SetPercentsOfUsedPoints(thePercents);
50   theStreamer->SetSource(anUsedPointsFilter->GetOutput());
51   anUsedPointsFilter->Register(theStreamer);
52   anUsedPointsFilter->Delete();
53 }
54
55
56 vtkStandardNewMacro(VISU_StreamLinesPL);
57
58 VISU_StreamLinesPL::VISU_StreamLinesPL(){
59   myStream = vtkStreamLine::New();
60   myGeomFilter = SALOME_GeometryFilter::New();
61   mySource = NULL;
62 }
63
64 VISU_StreamLinesPL::~VISU_StreamLinesPL(){
65   myStream->Delete();
66   myGeomFilter->Delete();
67   SetSource(NULL);
68 }
69
70 void VISU_StreamLinesPL::ShallowCopy(VISU_PipeLine *thePipeLine){
71   VISU_DeformedShapePL::ShallowCopy(thePipeLine);
72   if(VISU_StreamLinesPL *aPipeLine = dynamic_cast<VISU_StreamLinesPL*>(thePipeLine)){
73
74     //Order of setting of the values are important 
75     SetUsedPoints(aPipeLine->GetUsedPoints());
76     SetSource(aPipeLine->GetSource());
77
78     SetDirection(aPipeLine->GetDirection());
79
80     //Order of setting of the values are important 
81     SetIntegrationStep(aPipeLine->GetIntegrationStep());
82     SetPropagationTime(aPipeLine->GetPropagationTime());
83     SetStepLength(aPipeLine->GetStepLength());
84   }
85 }
86
87 void VISU_StreamLinesPL::SetSource(vtkPointSet* theSource){
88   if (mySource != NULL) mySource->UnRegister(this);
89   mySource = theSource;
90   if (mySource != NULL) mySource->Register(this);
91   if(myInput && myInput->GetCellData()->GetNumberOfArrays()){
92     vtkCellDataToPointData *aFilter = vtkCellDataToPointData::New();
93     aFilter->SetInput(myFieldTransform->GetUnstructuredGridOutput());
94     aFilter->PassCellDataOn();
95     myStream->SetInput(aFilter->GetOutput());
96     aFilter->Register(myStream);
97     aFilter->Delete();
98     vtkCellCenters *centers = vtkCellCenters::New(); // for vectors on cells
99     centers->SetInput(myFieldTransform->GetUnstructuredGridOutput());
100     centers->VertexCellsOn();
101     mySource = mySource? mySource: centers->GetOutput();
102     SetStreamerSource(myStream,mySource,myPercents);
103     centers->Register(myStream);
104     centers->Delete();
105   }else{
106     myStream->SetInput(myFieldTransform->GetUnstructuredGridOutput());
107     mySource = mySource? mySource: myFieldTransform->GetUnstructuredGridOutput();
108     SetStreamerSource(myStream,mySource,myPercents);
109   }
110   myStream->Modified();
111   Modified();
112 }
113
114
115 float VISU_StreamLinesPL::GetMaxIntegrationStep(vtkPointSet* theDataSet) {
116   if(!theDataSet) return 0.0;
117   theDataSet->Update();
118   float* aBounds = theDataSet->GetBounds();
119   float aMaxSizeY = aBounds[3] - aBounds[2];
120   float aMaxSizeZ = aBounds[5] - aBounds[4];
121   float aMinMax = aBounds[1] - aBounds[0];
122   if (aMinMax < EPS || (aMaxSizeY < aMinMax && aMaxSizeY > EPS)) aMinMax = aMaxSizeY;
123   if (aMinMax < EPS || (aMaxSizeZ < aMinMax && aMaxSizeZ > EPS)) aMinMax = aMaxSizeZ;
124   return aMinMax / 2.0;
125
126
127 float VISU_StreamLinesPL::GetMinIntegrationStep(vtkPointSet* theDataSet) {
128   if(!theDataSet) return 0.0;
129   float aVolume = 1.0, tmp;
130   int degree = 0;
131   theDataSet->Update();
132   float* aBounds = theDataSet->GetBounds();
133   for(int i = 0, j = 0; i < 3; ++i, j = 2*i){
134     float tmp = aBounds[j+1] - aBounds[j];
135     if (tmp > EPS ) {
136       aVolume *= tmp;
137       degree += 1;
138     }
139   }
140   if (degree < 1) return 0.0; // absolutely empty object
141   unsigned long aSize = 1024*theDataSet->GetActualMemorySize();
142   // we have to use no more, than myMaxIncrementMemorySize*aSize amount of memory
143   unsigned long aRealSize = GetAvailableMemory(myMaxIncrementMemorySize*aSize);
144   if(MYDEBUG) cout<<"GetMinIntegrationStep - GetActualMemorysize() = "<<aRealSize<<"; "<<aSize<<endl;
145   float anAverageVolume = aVolume / aRealSize; 
146   float aStep = pow(double(anAverageVolume), double(1.0/double(degree)));
147   return aStep;
148 }
149
150 float VISU_StreamLinesPL::GetBasePropagationTime(vtkPointSet* theDataSet) {
151   if(!theDataSet) return 0.0;
152   theDataSet->Update();
153   float* aScalarRange = theDataSet->GetScalarRange();
154   return theDataSet->GetLength()/fabs(aScalarRange[1] + aScalarRange[0])*2.0;
155 }
156
157
158 int VISU_StreamLinesPL::GetDirection(){
159   return myStream->GetIntegrationDirection();
160 }
161 void VISU_StreamLinesPL::SetDirection(int theDirection){
162   myStream->SetIntegrationDirection(theDirection);
163   Modified();
164 }
165
166
167 float VISU_StreamLinesPL::GetPropagationTime() { 
168   return myStream->GetMaximumPropagationTime();
169 }
170 void VISU_StreamLinesPL::SetPropagationTime(float theTime) { 
171   if(myStream->GetMaximumPropagationTime() == theTime) return;
172   myStream->SetMaximumPropagationTime(theTime);
173   Modified();
174 }
175
176
177 float VISU_StreamLinesPL::GetIntegrationStep(){
178   return myStream->GetIntegrationStepLength();
179 }
180 float VISU_StreamLinesPL::SetIntegrationStep(float theStep){ 
181   if(myStream->GetIntegrationStepLength() == theStep) return theStep;
182   myExtractor->Update();
183   vtkPointSet* aDataSet = GetSource();
184   aDataSet->Update();
185   float anIntegrationStep = theStep;
186   float aMinIntegrationStep = GetMinIntegrationStep(aDataSet);
187   float aMaxIntegrationStep = GetMaxIntegrationStep(aDataSet);
188   if(aMinIntegrationStep < anIntegrationStep && anIntegrationStep < aMaxIntegrationStep){
189     myStream->SetIntegrationStepLength(theStep);
190     Modified();
191   }else{
192     anIntegrationStep = aMaxIntegrationStep / 10.0;
193     float aMinMax = aDataSet->GetLength()/myInput->GetNumberOfPoints();
194     if(aMinMax > anIntegrationStep) 
195       anIntegrationStep = (anIntegrationStep*9.0+aMinMax)/10.0;
196     if(aMinIntegrationStep > anIntegrationStep) 
197       anIntegrationStep = aMinIntegrationStep;
198   }
199   return anIntegrationStep;
200 }
201
202
203 float VISU_StreamLinesPL::GetStepLength() { 
204   return myStream->GetStepLength();
205 }
206 float VISU_StreamLinesPL::GetMinStepLength(float theIntegrationStep){
207   return theIntegrationStep*2.0;
208 }
209 float VISU_StreamLinesPL::GetMaxStepLength(float thePropagationTime){
210   return thePropagationTime/100.0;
211 }
212 float VISU_StreamLinesPL::CorrectStepLength(float theStep, float theIntegrationStep, float thePropagationTime){
213   float aMinStep = GetMinStepLength(theIntegrationStep);
214   if(theStep < aMinStep) theStep = aMinStep;
215   float aMaxStep = GetMaxStepLength(thePropagationTime);
216   if(theStep > aMaxStep) theStep = aMaxStep;
217   return theStep;
218 }
219 float VISU_StreamLinesPL::SetStepLength(float theStep){
220   if(myStream->GetStepLength() == theStep) return theStep;
221   float aStepLength = CorrectStepLength(theStep,GetIntegrationStep(),GetPropagationTime());
222   if(aStepLength == theStep){
223     myStream->SetStepLength(theStep);
224     Modified();
225   }
226   return aStepLength;
227 }
228
229
230 void VISU_StreamLinesPL::Init(){
231   VISU_ScalarMapPL::Init();
232
233   SetDirection(VTK_INTEGRATE_BOTH_DIRECTIONS);
234   vtkPointSet* aDataSet = myExtractor->GetOutput();                                                     
235   //Order of setting of the values are important 
236   myPercents = 0.3;
237   SetSource(NULL);
238
239   SetIntegrationStep(SetIntegrationStep(GetMaxIntegrationStep(aDataSet)/10.0));
240   SetPropagationTime(GetBasePropagationTime(aDataSet));
241   SetStepLength(SetStepLength(GetPropagationTime()/20.));
242 }
243
244 VISU_ScalarMapPL::THook* VISU_StreamLinesPL::DoHook(){
245   SetSource(NULL);
246   myInput->Update();
247   float *aBounds = myInput->GetBounds();
248   myGeomFilter->SetExtent(aBounds);
249   myGeomFilter->ExtentClippingOn();
250   myGeomFilter->SetInput(myStream->GetOutput());
251   return myGeomFilter->GetOutput();
252 }
253
254 void VISU_StreamLinesPL::Update(){
255   VISU_ScalarMapPL::Update();
256 }
257
258 void VISU_StreamLinesPL::SetMapScale(float theMapScale){
259   VISU_ScalarMapPL::SetMapScale(theMapScale);
260 }