]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
To improve memory size calculation
authorapo <apo@opencascade.com>
Wed, 29 Nov 2006 14:00:05 +0000 (14:00 +0000)
committerapo <apo@opencascade.com>
Wed, 29 Nov 2006 14:00:05 +0000 (14:00 +0000)
src/OBJECT/VISU_Actor.cxx
src/PIPELINE/VISU_DeformedShapePL.cxx
src/PIPELINE/VISU_DeformedShapePL.hxx
src/PIPELINE/VISU_PipeLine.cxx
src/PIPELINE/VISU_PipeLine.hxx
src/PIPELINE/VISU_ScalarMapPL.cxx
src/PIPELINE/VISU_ScalarMapPL.hxx
src/VISU_I/VISU_Prs3d_i.cc

index 6d945d86d526573d4c36b7886b6dab6433f753dd..be9dbf2932c9a625ef12221f52e2a22b179af19d 100644 (file)
@@ -419,8 +419,17 @@ VISU_Actor
 ::GetMemorySize()
 {
   vtkDataSet* aDataSet = GetMapper()->GetInput();
-  unsigned long int aSize = aDataSet->GetActualMemorySize();
-  return aSize * 1024;
+  unsigned long int aSize = aDataSet->GetActualMemorySize() * 1024;
+
+  aDataSet = myGeomFilter->GetOutput();
+  aSize += aDataSet->GetActualMemorySize() * 1024;
+
+  if(IsShrunk()){
+    aDataSet = myShrinkFilter->GetOutput();
+    aSize += aDataSet->GetActualMemorySize() * 1024;
+  }
+
+  return aSize;
 }
 
 //----------------------------------------------------------------------------
index 43314db9b910c8c28d50a718016feece030331d1..eef5dbe2b05983abcf3504c9980b2e16c35f6879 100644 (file)
@@ -99,6 +99,7 @@ VISU_DeformedShapePL
   return myScaleFactor;
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_DeformedShapePL
 ::Init()
@@ -117,6 +118,7 @@ VISU_DeformedShapePL
     SetScale(0.0);
 }
 
+//----------------------------------------------------------------------------
 VISU_ScalarMapPL::THook* 
 VISU_DeformedShapePL
 ::DoHook()
@@ -125,6 +127,7 @@ VISU_DeformedShapePL
   return myWarpVector->GetOutput();
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_DeformedShapePL
 ::Update()
@@ -132,6 +135,24 @@ VISU_DeformedShapePL
   VISU_ScalarMapPL::Update();
 }
 
+//----------------------------------------------------------------------------
+unsigned long int
+VISU_DeformedShapePL
+::GetMemorySize()
+{
+  vtkDataSet* aDataSet = myCellDataToPointData->GetOutput();
+  unsigned long int aSize = aDataSet->GetActualMemorySize() * 1024;
+  
+  aDataSet = myWarpVector->GetOutput();
+  aSize += aDataSet->GetActualMemorySize() * 1024;
+
+  aSize += Superclass::GetMemorySize();
+
+  return aSize;
+}
+
+
+//----------------------------------------------------------------------------
 void
 VISU_DeformedShapePL
 ::SetMapScale(vtkFloatingPointType theMapScale)
index 96c3d0ebe7ea62e18f21561abd8b88404930a2d5..4049c668e0613a0e19403c5eb79c5e43798d8de2 100644 (file)
@@ -70,6 +70,11 @@ public:
   void
   Update();
 
+  //! Gets memory size used by the instance (bytes).
+  virtual
+  unsigned long int
+  GetMemorySize();
+
   virtual
   void
   SetMapScale(vtkFloatingPointType theMapScale = 1.0);
index 3905fdb4a541e6578668560751b48af0d9f108ce..d244c466e26fb7b69a72c79eafc9a7138393e477 100644 (file)
@@ -53,6 +53,7 @@ static int MYDEBUG = 0;
 static int MYDEBUG = 0;
 #endif
 
+//----------------------------------------------------------------------------
 VISU_PipeLine
 ::VISU_PipeLine():
   myMapper(vtkDataSetMapper::New()),
@@ -75,12 +76,14 @@ VISU_PipeLine
   myIsShrinkable = false;
 }
 
+//----------------------------------------------------------------------------
 VISU_PipeLine
 ::~VISU_PipeLine()
 {
   if(MYDEBUG) MESSAGE("VISU_PipeLine::~VISU_PipeLine - "<<this);
 }
 
+//----------------------------------------------------------------------------
 unsigned long int 
 VISU_PipeLine
 ::GetMTime()
@@ -89,6 +92,7 @@ VISU_PipeLine
   return aTime;
 }
 
+//----------------------------------------------------------------------------
 // Turn debugging output on.
 void
 VISU_PipeLine
@@ -98,6 +102,7 @@ VISU_PipeLine
   Superclass::DebugOn();
 }
 
+//----------------------------------------------------------------------------
 // Turn debugging output off.
 void
 VISU_PipeLine
@@ -107,6 +112,7 @@ VISU_PipeLine
   Superclass::DebugOff();
 }
 
+//----------------------------------------------------------------------------
 void 
 VISU_PipeLine
 ::ShallowCopy(VISU_PipeLine *thePipeLine)
@@ -121,6 +127,7 @@ VISU_PipeLine
   Build();
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_PipeLine
 ::SameAs(VISU_PipeLine *thePipeLine)
@@ -130,6 +137,7 @@ VISU_PipeLine
   GetImplicitFunction()->Delete();
 }
 
+//----------------------------------------------------------------------------
 TInput* 
 VISU_PipeLine
 ::GetInput() const
@@ -137,6 +145,7 @@ VISU_PipeLine
   return myInput.GetPointer();
 }
 
+//----------------------------------------------------------------------------
 vtkDataSet* 
 VISU_PipeLine
 ::GetOutput()
@@ -144,6 +153,7 @@ VISU_PipeLine
   return GetMapper()->GetInput();
 }
 
+//----------------------------------------------------------------------------
 TInput* 
 VISU_PipeLine
 ::GetInput2() const
@@ -153,6 +163,7 @@ VISU_PipeLine
   return aDataSet;
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_PipeLine
 ::SetInput(TInput* theInput)
@@ -167,6 +178,7 @@ VISU_PipeLine
   myInput = theInput;
 }
 
+//----------------------------------------------------------------------------
 VISU_PipeLine::TMapper* 
 VISU_PipeLine
 ::GetMapper()
@@ -181,6 +193,7 @@ VISU_PipeLine
   return myMapper.GetPointer();
 }
 
+//----------------------------------------------------------------------------
 void 
 VISU_PipeLine
 ::Update()
@@ -188,6 +201,7 @@ VISU_PipeLine
   myMapper->Update();
 }
 
+//----------------------------------------------------------------------------
 size_t
 VISU_PipeLine
 ::CheckAvailableMemory(size_t theSize)
@@ -205,6 +219,7 @@ VISU_PipeLine
   return 0;
 }
 
+//----------------------------------------------------------------------------
 size_t
 VISU_PipeLine
 ::GetAvailableMemory(size_t theSize, 
@@ -226,6 +241,23 @@ VISU_PipeLine
   return aMax;
 }
 
+//----------------------------------------------------------------------------
+unsigned long int
+VISU_PipeLine
+::GetMemorySize()
+{
+  vtkDataSet* aDataSet = myExtractGeometry->GetInput();
+  unsigned long int aSize = aDataSet->GetActualMemorySize() * 1024;
+  
+  aDataSet = myExtractGeometry->GetOutput();
+  aSize += aDataSet->GetActualMemorySize() * 1024;
+
+  aDataSet = GetMapper()->GetInput();
+  aSize += aDataSet->GetActualMemorySize() * 1024;
+
+  return aSize;
+}
+
 //------------------------ Clipping planes -----------------------------------
 bool 
 VISU_PipeLine
@@ -247,6 +279,7 @@ VISU_PipeLine
   return true;
 }
 
+//----------------------------------------------------------------------------
 vtkPlane* 
 VISU_PipeLine
 ::GetClippingPlane(vtkIdType theID) const
@@ -265,6 +298,7 @@ VISU_PipeLine
   return aPlane;
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_PipeLine
 ::RemoveAllClippingPlanes()
@@ -276,6 +310,7 @@ VISU_PipeLine
   }
 }
 
+//----------------------------------------------------------------------------
 vtkIdType
 VISU_PipeLine
 ::GetNumberOfClippingPlanes() const
@@ -287,6 +322,7 @@ VISU_PipeLine
   return 0;
 }
 
+//----------------------------------------------------------------------------
 static
 void
 ComputeBoundsParam (vtkDataSet* theDataSet,
@@ -335,6 +371,7 @@ ComputeBoundsParam (vtkDataSet* theDataSet,
   theMinPnt[2] = aMinPnt[2];
 }
 
+//----------------------------------------------------------------------------
 static
 void
 DistanceToPosition(vtkDataSet* theDataSet,
@@ -350,6 +387,7 @@ DistanceToPosition(vtkDataSet* theDataSet,
   thePos[2] = aMinPnt[2]-theDirection[2]*aLength;
 }
 
+//----------------------------------------------------------------------------
 static
 void
 PositionToDistance (vtkDataSet* theDataSet,
@@ -363,6 +401,7 @@ PositionToDistance (vtkDataSet* theDataSet,
   theDist = (aPrj-aMinBoundPrj)/(aMaxBoundPrj-aMinBoundPrj);
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_PipeLine
 ::SetPlaneParam(vtkFloatingPointType theDir[3], 
@@ -375,6 +414,7 @@ VISU_PipeLine
   thePlane->SetOrigin(anOrigin);
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_PipeLine
 ::GetPlaneParam(vtkFloatingPointType theDir[3], 
@@ -388,10 +428,7 @@ VISU_PipeLine
   ::PositionToDistance(GetInput(),theDir,anOrigin,theDist);
 }
 
-//=======================================================================
-//function : IsPlanarInput
-//purpose  :
-//=======================================================================
+//----------------------------------------------------------------------------
 bool
 VISU_PipeLine
 ::IsPlanarInput() const
index e10275260d67fea3e7b9c75450c14c3e63bc7556..95bf380edb705383b3f37c5ab67647ed2024307d 100644 (file)
@@ -147,9 +147,14 @@ public:
   
   static
   size_t
-  GetAvailableMemory(size_t theSize = 16*1024*1024,
+  GetAvailableMemory(size_t theSize,
                     size_t theMinSize = 1024*1024);
 
+  //! Gets memory size used by the instance (bytes).
+  virtual
+  unsigned long int
+  GetMemorySize();
+
   // Clipping planes
   void 
   RemoveAllClippingPlanes();
index 0c6b990328b57d3bff5dabf721b35ae914757bd5..755e45753552f2da6e0dc9847ccbf48299734f4b 100644 (file)
@@ -204,6 +204,7 @@ VISU_ScalarMapPL
 }
 
 
+//----------------------------------------------------------------------------
 void
 VISU_ScalarMapPL
 ::Init()
@@ -212,6 +213,7 @@ VISU_ScalarMapPL
   SetSourceRange();
 }
 
+//----------------------------------------------------------------------------
 void
 VISU_ScalarMapPL
 ::Build() 
@@ -220,6 +222,24 @@ VISU_ScalarMapPL
 }
 
 
+//----------------------------------------------------------------------------
+unsigned long int
+VISU_ScalarMapPL
+::GetMemorySize()
+{
+  vtkDataSet* aDataSet = myExtractor->GetOutput();
+  unsigned long int aSize = aDataSet->GetActualMemorySize() * 1024;
+  
+  aDataSet = myFieldTransform->GetOutput();
+  aSize += aDataSet->GetActualMemorySize() * 1024;
+
+  aSize += Superclass::GetMemorySize();
+
+  return aSize;
+}
+
+
+//----------------------------------------------------------------------------
 void
 VISU_ScalarMapPL
 ::Update() 
index de07ce2a7953a61ebd02466fe26b6ce3bca324de..fb7be40a889f2e60102b705c6fa858dc833a279b 100644 (file)
@@ -114,6 +114,11 @@ public:
   void
   Update();
   
+  //! Gets memory size used by the instance (bytes).
+  virtual
+  unsigned long int
+  GetMemorySize();
+
   virtual
   VISU_LookupTable*
   GetMapperTable();
index 1ee0f33de3dd903c436ea49c3b7618504652c020..d319ad259d659aa62f53103e272c955987285ec6 100644 (file)
@@ -601,20 +601,23 @@ CORBA::Float
 VISU::Prs3d_i
 ::GetMemorySize()
 {
-  static int INCMEMORY = 4;
-  vtkDataSet* aDataSet = GetPipeLine()->GetMapper()->GetInput();
-  CORBA::Float aSize = aDataSet->GetActualMemorySize() * INCMEMORY * 1024.0;
-  //cout<<"Prs3d_i::GetMemorySize - "<<this<<"; aDataSet = "<<aSize / (1024.0 * 1024.0)<<endl;
+  // To calculate memory used by VISU Converter
+  static int INCMEMORY = 2;
+  vtkDataSet* aDataSet = GetPipeLine()->GetInput();
+  CORBA::Float aSize = aDataSet->GetActualMemorySize() * INCMEMORY / 1024.0;
 
+  // To calculate memory used by VISU PipeLine
+  aSize += GetPipeLine()->GetMemorySize();
+
+  // To calculate memory used by VISU Actos
   int anEnd = myActorCollection->GetNumberOfItems();
   for(int anId = 0; anId < anEnd; anId++)
     if(vtkObject* anObject = myActorCollection->GetItemAsObject(anId))
       if(VISU_Actor* anActor = dynamic_cast<VISU_Actor*>(anObject)){
        aSize += anActor->GetMemorySize();
-       //cout<<"Prs3d_i::GetMemorySize - "<<this<<"; anActor = "<<anActor->GetMemorySize() / (1024.0 * 1024.0)<<endl;
       }
 
-  //cout<<"Prs3d_i::GetMemorySize - "<<this<<"; aSize = "<<aSize / (1024.0 * 1024.0)<<endl;
+  // Convert to mega bytes
   return aSize / (1024.0 * 1024.0);
 }