]> SALOME platform Git repositories - modules/visu.git/blob - src/VISU_I/VISU_StreamLines_i.cc
Salome HOME
Merge from PortingMED3 07Apr11
[modules/visu.git] / src / VISU_I / VISU_StreamLines_i.cc
1 //  Copyright (C) 2007-2010  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
23 //  VISU OBJECT : interactive object for VISU entities implementation
24 //  File   : VISU_PrsObject_i.cxx
25 //  Author : Alexey PETROV
26 //  Module : VISU
27 //
28 #include "VISU_StreamLines_i.hh"
29 #include "VISU_Prs3dUtils.hh"
30 #include "VISU_Result_i.hh"
31
32 #include "VISU_Actor.h"
33 #include "VISU_StreamLinesPL.hxx"
34 #include "VISU_Convertor.hxx"
35
36 #include "SUIT_ResourceMgr.h"
37 #include "SALOME_Event.h"
38
39 #include <vtkDataSetMapper.h>
40 #include <vtkAppendFilter.h>
41 #include <vtkUnstructuredGrid.h>
42
43 #ifdef _DEBUG_
44 static int MYDEBUG = 0;
45 #else
46 static int MYDEBUG = 0;
47 #endif
48
49 using namespace std;
50
51 //---------------------------------------------------------------
52 size_t
53 VISU::StreamLines_i
54 ::IsPossible(Result_i* theResult, 
55              const std::string& theMeshName, 
56              VISU::Entity theEntity,
57              const std::string& theFieldName, 
58              CORBA::Long theTimeStampNumber,
59              bool theIsMemoryCheck)
60 {
61   try{
62     if(!TSuperClass::IsPossible(theResult, theMeshName, theEntity, theFieldName, theTimeStampNumber, false)) 
63       return 0;
64
65     VISU::Result_i::PInput anInput = theResult->GetInput(theMeshName,
66                                                          theEntity,
67                                                          theFieldName,
68                                                          theTimeStampNumber);
69     VISU::PUnstructuredGridIDMapper anIDMapper = 
70       anInput->GetTimeStampOnMesh(theMeshName,
71                                   VISU::TEntity(theEntity),
72                                   theFieldName,
73                                   theTimeStampNumber);
74
75     vtkUnstructuredGrid* aDataSet = anIDMapper->GetUnstructuredGridOutput();
76     size_t aResult = VISU_StreamLinesPL::IsPossible(aDataSet);
77     MESSAGE("StreamLines_i::IsPossible - aResult = "<<aResult);
78     return aResult;
79   }catch(std::exception& exc){
80     INFOS("Follow exception was occured :\n"<<exc.what());
81   }catch(...){
82     INFOS("Unknown exception was occured!");
83   }
84   return 0;
85 }
86
87 //---------------------------------------------------------------
88 int VISU::StreamLines_i::myNbPresent = 0;
89
90 //---------------------------------------------------------------
91 QString
92 VISU::StreamLines_i
93 ::GenerateName() 
94 {
95   return VISU::GenerateName("StreamLines",myNbPresent++);
96 }
97
98 //---------------------------------------------------------------
99 const string VISU::StreamLines_i::myComment = "STREAMLINES";
100
101 //---------------------------------------------------------------
102 const char* 
103 VISU::StreamLines_i
104 ::GetComment() const 
105
106   return myComment.c_str();
107 }
108
109 //---------------------------------------------------------------
110 const char*
111 VISU::StreamLines_i
112 ::GetIconName()
113 {
114   if (!IsGroupsUsed())
115     return "ICON_TREE_STREAM_LINES";
116   else
117     return "ICON_TREE_STREAM_LINES_GROUPS";
118 }
119
120 //---------------------------------------------------------------
121 VISU::StreamLines_i
122 ::StreamLines_i(EPublishInStudyMode thePublishInStudyMode) :
123   ColoredPrs3d_i(thePublishInStudyMode),
124   ScalarMap_i(thePublishInStudyMode),
125   MonoColorPrs_i(thePublishInStudyMode),
126   myStreamLinesPL(NULL),
127   myAppendFilter(vtkAppendFilter::New())
128 {}
129
130
131 //---------------------------------------------------------------
132 void
133 VISU::StreamLines_i
134 ::SameAs(const Prs3d_i* theOrigin)
135 {
136   TSuperClass::SameAs(theOrigin);
137
138   if(const StreamLines_i* aPrs3d = dynamic_cast<const StreamLines_i*>(theOrigin)) {
139     StreamLines_i* anOrigin = const_cast<StreamLines_i*>(aPrs3d);
140     SetSource(anOrigin->GetSource());
141   }
142 }
143
144
145 //---------------------------------------------------------------
146 VISU::Storable* 
147 VISU::StreamLines_i
148 ::Create(const std::string& theMeshName, 
149          VISU::Entity theEntity,
150          const std::string& theFieldName, 
151          CORBA::Long theTimeStampNumber)
152 {
153   return TSuperClass::Create(theMeshName,theEntity,theFieldName,theTimeStampNumber);
154 }
155
156
157 //---------------------------------------------------------------
158 VISU::Storable* 
159 VISU::StreamLines_i
160 ::Restore(SALOMEDS::SObject_ptr theSObject,
161           const Storable::TRestoringMap& theMap)
162 {
163   if(!TSuperClass::Restore(theSObject, theMap))
164     return NULL;
165
166   double anIntegrationStep = VISU::Storable::FindValue(theMap,"myIntegrationStep").toDouble();
167   double aPropagationTime = VISU::Storable::FindValue(theMap,"myPropagationTime").toDouble();
168   double aStepLength = VISU::Storable::FindValue(theMap,"myStepLength").toDouble();
169   int aDirection = VISU::StreamLines::Direction(VISU::Storable::FindValue(theMap,"myDirection").toInt());
170   double aPercents = VISU::Storable::FindValue(theMap,"myPercents").toDouble();
171   SetParams(anIntegrationStep,
172             aPropagationTime,
173             aStepLength,
174             VISU::Prs3d::_nil(),
175             aPercents,
176             VISU::StreamLines::Direction(aDirection));
177   mySourceEntry = (const char*)VISU::Storable::FindValue(theMap,"mySourceEntry").toLatin1();
178
179   return this;
180 }
181
182
183 //---------------------------------------------------------------
184 void
185 VISU::StreamLines_i
186 ::ToStream(std::ostringstream& theStr)
187 {
188   TSuperClass::ToStream(theStr);
189
190   Storable::DataToStream( theStr, "myIntegrationStep", GetIntegrationStep());
191   Storable::DataToStream( theStr, "myPropagationTime", GetPropagationTime());
192   Storable::DataToStream( theStr, "myStepLength", GetStepLength());
193
194   Storable::DataToStream( theStr, "myDirection", int(GetDirection()));
195   Storable::DataToStream( theStr, "myPercents", GetUsedPoints());
196
197   Storable::DataToStream( theStr, "mySourceEntry", mySourceEntry.c_str());
198 }
199
200
201 //---------------------------------------------------------------
202 VISU::StreamLines_i
203 ::~StreamLines_i()
204 {
205   if(MYDEBUG) MESSAGE("StreamLines_i::~StreamLines_i()");
206   myAppendFilter->Delete();
207 }
208
209
210 //---------------------------------------------------------------
211 CORBA::Boolean
212 VISU::StreamLines_i
213 ::SetParams(CORBA::Double theIntStep,
214             CORBA::Double thePropogationTime,
215             CORBA::Double theStepLength,
216             VISU::Prs3d_ptr thePrs3d,
217             CORBA::Double thePercents,
218             VISU::StreamLines::Direction theDirection)
219 {
220   struct TEvent: public SALOME_Event 
221   {
222     CORBA::Double myIntStep;
223     CORBA::Double myPropogationTime;
224     CORBA::Double myStepLength;
225     VISU::Prs3d_ptr myPrs3d;
226     CORBA::Double myPercents;
227     VISU::StreamLines::Direction myDirection;
228     
229     VISU_StreamLinesPL* myPipeLine;
230     VISU::Prs3d_i*& myPrs3dServant;
231     vtkAppendFilter* myAppendFilter;
232     size_t& myIsAccepted;
233     
234     TEvent(CORBA::Double theIntStep,
235            CORBA::Double thePropogationTime,
236            CORBA::Double theStepLength,
237            VISU::Prs3d_ptr thePrs3d,
238            CORBA::Double thePercents,
239            VISU::StreamLines::Direction theDirection,
240            VISU_StreamLinesPL* thePipeLine,
241            VISU::Prs3d_i*& thePrs3dServant,
242            vtkAppendFilter* theAppendFilter,
243            size_t& theIsAccepted):
244       myIntStep(theIntStep),
245       myPropogationTime(thePropogationTime),
246       myStepLength(theStepLength),
247       myPrs3d(thePrs3d),
248       myPercents(thePercents),
249       myDirection(theDirection),
250       myPipeLine(thePipeLine),
251       myPrs3dServant(thePrs3dServant),
252       myAppendFilter(theAppendFilter),
253       myIsAccepted(theIsAccepted)
254     {}
255
256     virtual
257     void
258     Execute()
259     {
260       myPrs3dServant = NULL;
261       vtkPointSet* aSource = NULL;
262       if (!CORBA::is_nil(myPrs3d)){
263         myPrs3dServant = dynamic_cast<VISU::Prs3d_i*>(VISU::GetServant(myPrs3d).in());
264         if(myPrs3dServant){
265           myAppendFilter->RemoveAllInputs();
266           myAppendFilter->AddInput(myPrs3dServant->GetPipeLine()->GetMapper()->GetInput());
267           aSource = myAppendFilter->GetOutput();
268         }
269       }
270       myIsAccepted = myPipeLine->SetParams(myIntStep,
271                                            myPropogationTime,
272                                            myStepLength,
273                                            aSource,
274                                            myPercents,
275                                            myDirection);
276     }
277   };
278
279   size_t anIsAccepted = false;
280   VISU::Prs3d_i* aPrs3dServant = NULL;
281   ProcessVoidEvent(new TEvent(theIntStep,
282                               thePropogationTime,
283                               theStepLength,
284                               thePrs3d,
285                               thePercents,
286                               theDirection,
287                               GetSpecificPL(),
288                               aPrs3dServant,
289                               myAppendFilter,
290                               anIsAccepted));
291   if(anIsAccepted)
292     SetSource(aPrs3dServant);
293
294   return anIsAccepted;
295 }
296
297
298 //---------------------------------------------------------------
299 void
300 VISU::StreamLines_i
301 ::SetSource(VISU::Prs3d_ptr thePrs3d)
302 {
303   if(!thePrs3d->_is_nil()){
304     VISU::Prs3d_i* aPrs3di = dynamic_cast<VISU::Prs3d_i*>(VISU::GetServant(thePrs3d).in());
305     SetSource(aPrs3di);
306   }
307 }
308
309 //---------------------------------------------------------------
310 void
311 VISU::StreamLines_i
312 ::SetSource(VISU::Prs3d_i* thePrs3d)
313 {
314   mySourceEntry = "";
315   if(!thePrs3d)
316     return;
317
318   SALOMEDS::SObject_var aSObject = thePrs3d->GetSObject();
319   CORBA::String_var aString = aSObject->GetID();
320   if(mySourceEntry == aString.in())
321     return;
322   
323   VISU::TSetModified aModified(this);
324   
325   mySourceEntry = aString.in();
326   myParamsTime.Modified();
327 }
328
329 //---------------------------------------------------------------
330 void
331 VISU::StreamLines_i
332 ::SetSource()
333 {
334   if(!myStreamLinesPL->GetSource() && mySourceEntry == "") 
335     return;
336
337   if(myStreamLinesPL->GetSource() == myAppendFilter->GetOutput()) 
338     return;
339
340   VISU::Prs3d_var aPrs3d = GetSource();
341   SetParams(GetIntegrationStep(),
342             GetPropagationTime(),
343             GetStepLength(),
344             aPrs3d,
345             GetUsedPoints(),
346             GetDirection());
347 }
348
349
350 //---------------------------------------------------------------
351 CORBA::Double 
352 VISU::StreamLines_i
353 ::GetIntegrationStep() 
354
355   return myStreamLinesPL->GetIntegrationStep();
356 }
357
358 //---------------------------------------------------------------
359 CORBA::Double
360 VISU::StreamLines_i
361 ::GetPropagationTime() 
362
363   return myStreamLinesPL->GetPropagationTime();
364 }
365
366 //---------------------------------------------------------------
367 CORBA::Double
368 VISU::StreamLines_i
369 ::GetStepLength() 
370
371   return myStreamLinesPL->GetStepLength();
372 }
373
374 //---------------------------------------------------------------
375 VISU::StreamLines::Direction
376 VISU::StreamLines_i
377 ::GetDirection() 
378
379   return VISU::StreamLines::Direction(myStreamLinesPL->GetDirection());
380 }
381
382
383 //---------------------------------------------------------------
384 VISU::Prs3d_ptr
385 VISU::StreamLines_i
386 ::GetSource()
387 {
388   VISU::Prs3d_var aPrs3d;
389   if(MYDEBUG) MESSAGE("StreamLines_i::GetSource() mySourceEntry = '"<<mySourceEntry<<"'");
390   if(mySourceEntry != ""){
391     SALOMEDS::SObject_var aSObject = GetStudyDocument()->FindObjectID(mySourceEntry.c_str());
392     CORBA::Object_var anObj = SObjectToObject(aSObject);
393     if(!CORBA::is_nil(anObj)) aPrs3d = VISU::Prs3d::_narrow(anObj);
394   }
395   return aPrs3d._retn();
396 }
397
398 //---------------------------------------------------------------
399 CORBA::Double 
400 VISU::StreamLines_i
401 ::GetUsedPoints() 
402
403   return myStreamLinesPL->GetUsedPoints();
404 }
405
406
407 //---------------------------------------------------------------
408 void 
409 VISU::StreamLines_i
410 ::CreatePipeLine(VISU_PipeLine* thePipeLine)
411 {
412   if(!thePipeLine){
413     myStreamLinesPL = VISU_StreamLinesPL::New();
414   }else
415     myStreamLinesPL = dynamic_cast<VISU_StreamLinesPL*>(thePipeLine);
416
417   TSuperClass::CreatePipeLine(myStreamLinesPL);
418 }
419
420
421 //----------------------------------------------------------------------------
422 bool
423 VISU::StreamLines_i
424 ::CheckIsPossible() 
425 {
426   return IsPossible(GetCResult(),GetCMeshName(),GetEntity(),GetCFieldName(),GetTimeStampNumber(),true);
427 }
428
429 //---------------------------------------------------------------
430 void
431 VISU::StreamLines_i
432 ::Update() 
433 {
434   SetSource();
435   TSuperClass::Update();
436 }
437
438
439 //---------------------------------------------------------------
440 VISU_Actor* 
441 VISU::StreamLines_i
442 ::CreateActor() 
443 {
444   if(VISU_Actor* anActor = TSuperClass::CreateActor(true)){
445     anActor->SetVTKMapping(true);
446     SUIT_ResourceMgr* aResourceMgr = VISU::GetResourceMgr();
447     int  aDispMode = aResourceMgr->integerValue("VISU", "stream_lines_represent", 1);
448     anActor->SetRepresentation(aDispMode);
449     return anActor;
450   }
451   return NULL;
452 }
453
454
455 //---------------------------------------------------------------
456 void
457 VISU::StreamLines_i
458 ::UpdateActor(VISU_Actor* theActor) 
459 {
460   TSuperClass::UpdateActor(theActor);
461 }