]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
VISU impruvements T1.5 - Implement Clipping Planes functionality
authorapo <apo@opencascade.com>
Wed, 13 Apr 2005 12:03:23 +0000 (12:03 +0000)
committerapo <apo@opencascade.com>
Wed, 13 Apr 2005 12:03:23 +0000 (12:03 +0000)
src/PIPELINE/Makefile.in
src/PIPELINE/SALOME_ExtractGeometry.cxx [new file with mode: 0755]
src/PIPELINE/SALOME_ExtractGeometry.h [new file with mode: 0755]
src/PIPELINE/VISU_PipeLine.cxx
src/PIPELINE/VISU_PipeLine.hxx
src/VISU_I/VISU_Prs3d_i.cc
src/VISU_I/VISU_Prs3d_i.hh

index c0bf576f5143f92b4516b60df1d92a0a4b111e7a..f4ba169472c82dda3e7bed67b2a80069f6881431 100644 (file)
@@ -51,7 +51,8 @@ LIB_SRC = VISU_PipeLine.cxx VISU_PipeLineUtils.cxx \
        VISU_VectorsPL.cxx VISU_StreamLinesPL.cxx \
        VISU_LookupTable.cxx VISU_ScalarBarActor.cxx \
        VISU_Extractor.cxx VISU_FieldTransform.cxx \
-       VISU_UsedPointsFilter.cxx
+       VISU_UsedPointsFilter.cxx \
+       SALOME_ExtractGeometry.cxx
 
 # Executables targets
 
diff --git a/src/PIPELINE/SALOME_ExtractGeometry.cxx b/src/PIPELINE/SALOME_ExtractGeometry.cxx
new file mode 100755 (executable)
index 0000000..d4b03f4
--- /dev/null
@@ -0,0 +1,267 @@
+//  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+//  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
+// 
+//  This library is free software; you can redistribute it and/or 
+//  modify it under the terms of the GNU Lesser General Public 
+//  License as published by the Free Software Foundation; either 
+//  version 2.1 of the License. 
+// 
+//  This library is distributed in the hope that it will be useful, 
+//  but WITHOUT ANY WARRANTY; without even the implied warranty of 
+//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
+//  Lesser General Public License for more details. 
+// 
+//  You should have received a copy of the GNU Lesser General Public 
+//  License along with this library; if not, write to the Free Software 
+//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
+// 
+//  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
+
+
+#include "SALOME_ExtractGeometry.h"
+
+#include <vtkCell.h>
+#include <vtkCellData.h>
+#include <vtkFloatArray.h>
+#include <vtkIdList.h>
+#include <vtkImplicitFunction.h>
+#include <vtkObjectFactory.h>
+#include <vtkPointData.h>
+#include <vtkUnstructuredGrid.h>
+
+#include <vtkImplicitBoolean.h>
+#include <vtkImplicitFunctionCollection.h>
+
+using namespace std;
+
+
+vtkStandardNewMacro(SALOME_ExtractGeometry);
+
+
+SALOME_ExtractGeometry::SALOME_ExtractGeometry()
+{}
+
+
+SALOME_ExtractGeometry::~SALOME_ExtractGeometry(){}
+
+
+vtkIdType SALOME_ExtractGeometry::GetElemObjId(int theID){
+  if(myElemVTK2ObjIds.empty() || theID > myElemVTK2ObjIds.size()) 
+    return theID;
+  return myElemVTK2ObjIds[theID];
+}
+
+
+vtkIdType SALOME_ExtractGeometry::GetNodeObjId(int theID){
+  if(myNodeVTK2ObjIds.empty() || theID > myNodeVTK2ObjIds.size()) 
+    return theID;
+  return myNodeVTK2ObjIds[theID];
+}
+
+
+void SALOME_ExtractGeometry::SetImplicitBoolean(vtkImplicitBoolean* theImplicitBoolean)
+{
+  myImplicitBoolean = theImplicitBoolean;
+  SetImplicitFunction(theImplicitBoolean);
+}
+
+
+void SALOME_ExtractGeometry::SetStoreMapping(bool theStoreMapping)
+{
+  myStoreMapping = theStoreMapping;
+  Modified();
+}
+
+
+void SALOME_ExtractGeometry::Execute()
+{
+  if(myImplicitBoolean.GetPointer()){
+    if(vtkImplicitFunctionCollection* aFunction = myImplicitBoolean->GetFunction()){
+      if(aFunction->GetNumberOfItems() == 0){
+       vtkDebugMacro(<< "Extracting geometry - ShallowCopy");
+       GetOutput()->ShallowCopy(GetInput());
+       return;
+      }
+    }
+  }
+  Execute2();
+}
+
+void SALOME_ExtractGeometry::Execute2()
+{
+  vtkIdType ptId, numPts, numCells, i, cellId, newCellId, newId, *pointMap;
+  vtkIdList *cellPts;
+  vtkCell *cell;
+  int numCellPts;
+  float *x;
+  float multiplier;
+  vtkPoints *newPts;
+  vtkIdList *newCellPts;
+  vtkDataSet *input = this->GetInput();
+  vtkPointData *pd = input->GetPointData();
+  vtkCellData *cd = input->GetCellData();
+  vtkUnstructuredGrid *output = this->GetOutput();
+  vtkPointData *outputPD = output->GetPointData();
+  vtkCellData *outputCD = output->GetCellData();
+  int npts;
+  numCells = input->GetNumberOfCells();
+  numPts = input->GetNumberOfPoints();
+  
+  vtkDebugMacro(<< "Extracting geometry");
+
+  if ( ! this->ImplicitFunction )
+    {
+    vtkErrorMacro(<<"No implicit function specified");
+    return;
+    }
+
+  newCellPts = vtkIdList::New();
+  newCellPts->Allocate(VTK_CELL_SIZE);
+
+  if ( this->ExtractInside )
+    {
+    multiplier = 1.0;
+    }
+  else 
+    {
+    multiplier = -1.0;
+    }
+
+  // Loop over all points determining whether they are inside the
+  // implicit function. Copy the points and point data if they are.
+  //
+  pointMap = new vtkIdType[numPts]; // maps old point ids into new
+  for (i=0; i < numPts; i++)
+    {
+    pointMap[i] = -1;
+    }
+
+  output->Allocate(numCells/4); //allocate storage for geometry/topology
+  newPts = vtkPoints::New();
+  newPts->Allocate(numPts/4,numPts);
+  outputPD->CopyAllocate(pd);
+  outputCD->CopyAllocate(cd);
+  vtkFloatArray *newScalars = NULL;
+  
+  if(myStoreMapping){
+    myElemVTK2ObjIds.clear();
+    myElemVTK2ObjIds.reserve(numCells);
+    myNodeVTK2ObjIds.clear();
+    myNodeVTK2ObjIds.reserve(numPts);
+  }
+
+  if ( ! this->ExtractBoundaryCells )
+    {
+    for ( ptId=0; ptId < numPts; ptId++ )
+      {
+      x = input->GetPoint(ptId);
+      if ( (this->ImplicitFunction->FunctionValue(x)*multiplier) < 0.0 )
+        {
+        newId = newPts->InsertNextPoint(x);
+        pointMap[ptId] = newId;
+       myNodeVTK2ObjIds.push_back(ptId);
+        outputPD->CopyData(pd,ptId,newId);
+        }
+      }
+    }
+  else
+    {
+    // To extract boundary cells, we have to create supplemental information
+    if ( this->ExtractBoundaryCells )
+      {
+      float val;
+      newScalars = vtkFloatArray::New();
+      newScalars->SetNumberOfValues(numPts);
+
+      for (ptId=0; ptId < numPts; ptId++ )
+        {
+        x = input->GetPoint(ptId);
+        val = this->ImplicitFunction->FunctionValue(x) * multiplier;
+        newScalars->SetValue(ptId, val);
+        if ( val < 0.0 )
+          {
+          newId = newPts->InsertNextPoint(x);
+          pointMap[ptId] = newId;
+         myNodeVTK2ObjIds.push_back(ptId);
+          outputPD->CopyData(pd,ptId,newId);
+          }
+        }
+      }
+    }
+
+  // Now loop over all cells to see whether they are inside implicit
+  // function (or on boundary if ExtractBoundaryCells is on).
+  //
+  for (cellId=0; cellId < numCells; cellId++)
+    {
+    cell = input->GetCell(cellId);
+    cellPts = cell->GetPointIds();
+    numCellPts = cell->GetNumberOfPoints();
+
+    newCellPts->Reset();
+    if ( ! this->ExtractBoundaryCells ) //requires less work
+      {
+      for ( npts=0, i=0; i < numCellPts; i++, npts++)
+        {
+        ptId = cellPts->GetId(i);
+        if ( pointMap[ptId] < 0 )
+          {
+          break; //this cell won't be inserted
+          }
+        else
+          {
+          newCellPts->InsertId(i,pointMap[ptId]);
+          }
+        }
+      } //if don't want to extract boundary cells
+    
+    else //want boundary cells
+      {
+      for ( npts=0, i=0; i < numCellPts; i++ )
+        {
+        ptId = cellPts->GetId(i);
+        if ( newScalars->GetValue(ptId) <= 0.0 )
+          {
+          npts++;
+          }
+        }
+      if ( npts > 0 )
+        {
+        for ( i=0; i < numCellPts; i++ )
+          {
+          ptId = cellPts->GetId(i);
+          if ( pointMap[ptId] < 0 )
+            {
+            x = input->GetPoint(ptId);
+            newId = newPts->InsertNextPoint(x);
+            pointMap[ptId] = newId;
+           myNodeVTK2ObjIds.push_back(ptId);
+            outputPD->CopyData(pd,ptId,newId);
+            }
+          newCellPts->InsertId(i,pointMap[ptId]);
+          }
+        }//a boundary or interior cell
+      }//if mapping boundary cells
+      
+    if ( npts >= numCellPts || (this->ExtractBoundaryCells && npts > 0) )
+      {
+      newCellId = output->InsertNextCell(cell->GetCellType(),newCellPts);
+      myElemVTK2ObjIds.push_back(cellId);
+      outputCD->CopyData(cd,cellId,newCellId);
+      }
+    }//for all cells
+
+  // Update ourselves and release memory
+  //
+  delete [] pointMap;
+  newCellPts->Delete();
+  output->SetPoints(newPts);
+  newPts->Delete();
+  
+  if ( this->ExtractBoundaryCells )
+    {
+    newScalars->Delete();
+    }
+
+  output->Squeeze();
+}
diff --git a/src/PIPELINE/SALOME_ExtractGeometry.h b/src/PIPELINE/SALOME_ExtractGeometry.h
new file mode 100755 (executable)
index 0000000..ebbaec5
--- /dev/null
@@ -0,0 +1,70 @@
+//  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
+//  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
+// 
+//  This library is free software; you can redistribute it and/or 
+//  modify it under the terms of the GNU Lesser General Public 
+//  License as published by the Free Software Foundation; either 
+//  version 2.1 of the License. 
+// 
+//  This library is distributed in the hope that it will be useful, 
+//  but WITHOUT ANY WARRANTY; without even the implied warranty of 
+//  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
+//  Lesser General Public License for more details. 
+// 
+//  You should have received a copy of the GNU Lesser General Public 
+//  License along with this library; if not, write to the Free Software 
+//  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
+// 
+//  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
+
+#ifndef SALOME_ExtractGeometry_H
+#define SALOME_ExtractGeometry_H
+
+#include <vtkExtractGeometry.h>
+#include <vtkSmartPointer.h>
+#include <vector>
+
+class vtkImplicitBoolean;
+
+class SALOME_ExtractGeometry : public vtkExtractGeometry{
+public:
+  vtkTypeMacro(SALOME_ExtractGeometry,vtkExtractGeometry);
+
+  static SALOME_ExtractGeometry *New();
+
+  void SetImplicitBoolean(vtkImplicitBoolean* theImplicitBoolean);
+  vtkImplicitBoolean* GetImplicitBoolean() const { 
+    return myImplicitBoolean.GetPointer();
+  }
+
+  void SetStoreMapping(bool theStoreMapping);
+  bool GetStoreMapping() const { 
+    return myStoreMapping;
+  }
+
+  virtual vtkIdType GetNodeObjId(int theID);
+  virtual vtkIdType GetElemObjId(int theID);
+
+protected:
+  SALOME_ExtractGeometry();
+  ~SALOME_ExtractGeometry();
+
+  virtual void Execute();
+  void Execute2();
+
+private:
+  bool myStoreMapping;   
+  typedef std::vector<vtkIdType> TVectorId;
+  TVectorId myElemVTK2ObjIds;
+  TVectorId myNodeVTK2ObjIds;
+
+  vtkSmartPointer<vtkImplicitBoolean> myImplicitBoolean;
+
+  SALOME_ExtractGeometry(const SALOME_ExtractGeometry&);  // Not implemented.
+  void operator=(const SALOME_ExtractGeometry&);  // Not implemented.
+};
+
+
+#endif
+
+
index 2524787c457e5aeabdf6eb98cb5a38ff8c47c630..9cb57f20cf9fa3b7d84b47a106fbda1ae2fb870b 100644 (file)
@@ -29,7 +29,7 @@
 #include "VISU_PipeLine.hxx"
 #include "VISU_PipeLineUtils.hxx"
 
-#include "SALOME_PassThroughFilter.h"
+#include "SALOME_ExtractGeometry.h"
 
 #include <limits.h>
 
@@ -47,89 +47,75 @@ static int MYVTKDEBUG = 0;
 
 #ifdef _DEBUG_
 static int MYDEBUG = 0;
-static int MYDEBUGWITHFILES = 0;
 #else
 static int MYDEBUG = 0;
-static int MYDEBUGWITHFILES = 0;
 #endif
 
 VISU_PipeLine::VISU_PipeLine()
 {
+  if(MYDEBUG) MESSAGE("VISU_PipeLine - "<<this);
   // Clipping planes
-  myPassFilter0 = SALOME_PassThroughFilter::New();
-
-  myExtractGeometry = vtkExtractGeometry::New();
-  myExtractGeometry->SetReleaseDataFlag(true);
-  myIsImplicitFunctionUsed = false;
+  myExtractGeometry = SALOME_ExtractGeometry::New();
+  //myExtractGeometry->SetReleaseDataFlag(true);
+  myExtractGeometry->Delete();
+  //myExtractGeometry->DebugOn();
 
-  myImplicitBoolean = vtkImplicitBoolean::New();
-  myImplicitBoolean->SetOperationTypeToIntersection();
-  myExtractGeometry->SetImplicitFunction(myImplicitBoolean);
+  vtkImplicitBoolean* anImplicitBoolean = vtkImplicitBoolean::New();
+  myExtractGeometry->SetImplicitBoolean(anImplicitBoolean);
+  anImplicitBoolean->SetOperationTypeToIntersection();
+  anImplicitBoolean->Delete();
 
   // Mapper
   myMapper = TMapper::New();
-
   myInput = NULL;
+
   SetDebug(MYVTKDEBUG);
 }
 
 VISU_PipeLine::~VISU_PipeLine()
 {
-  if (MYDEBUG) MESSAGE("~VISU_PipeLine - myInput = "<<myInput->GetReferenceCount());
-  SetInput(NULL);
-  myMapper->RemoveAllInputs();
+  if(MYDEBUG) MESSAGE("~VISU_PipeLine - "<<this);
   myMapper->Delete();
-
-  //myPassFilter0->UnRegisterAllOutputs(); 
-  myPassFilter0->Delete();
-  //myExtractGeometry->UnRegisterAllOutputs();
-  myExtractGeometry->Delete();
-  myImplicitBoolean->Delete();
 }
 
-void VISU_PipeLine::ShallowCopy (VISU_PipeLine *thePipeLine) {
-  myIsImplicitFunctionUsed = thePipeLine->myIsImplicitFunctionUsed;
-  myCippingPlaneCont = thePipeLine->myCippingPlaneCont;
-
-  myImplicitBoolean->GetFunction()->RemoveAllItems();
-  for (int ipln = 0; ipln < myCippingPlaneCont.size(); ipln++) {
-    myImplicitBoolean->GetFunction()->AddItem(myCippingPlaneCont[ipln].Get());
-  }
-  myImplicitBoolean->GetFunction()->Modified(); // VTK bug
-
+void VISU_PipeLine::ShallowCopy(VISU_PipeLine *thePipeLine){
   SetInput(thePipeLine->GetInput());
-  SetImplicitFunctionUsed(myIsImplicitFunctionUsed);
   myMapper->ShallowCopy(thePipeLine->GetMapper());
+  myExtractGeometry->SetImplicitBoolean(thePipeLine->myExtractGeometry->GetImplicitBoolean());
   Build();
 }
 
-TInput* VISU_PipeLine::GetInput2() {
-  return myPassFilter0->GetUnstructuredGridOutput();
+void VISU_PipeLine::SameAs(VISU_PipeLine *thePipeLine){
+  ShallowCopy(thePipeLine);
+  myExtractGeometry->SetImplicitBoolean(vtkImplicitBoolean::New());
+  myExtractGeometry->GetImplicitBoolean()->Delete();
 }
 
-void VISU_PipeLine::SetInput (TInput* theInput)
+TInput* VISU_PipeLine::GetInput() const
 {
-  if (myInput != theInput) {
-    if (myInput != NULL) myInput->UnRegister(this);
-    myInput = theInput;
-    if (myInput != NULL) {
-      myInput->Register(this);
-
-      // Clipping planes
-      SetImplicitFunctionUsed(myIsImplicitFunctionUsed);
-
-      myInput->Update();
-    } else {
-      myMapper->SetInput(NULL);
-    }
-    Modified();
-  }
+  return myInput;
+}
+
+TInput* VISU_PipeLine::GetInput2() const
+{
+  vtkUnstructuredGrid* aDataSet = myExtractGeometry->GetOutput();
+  aDataSet->Update();
+  return aDataSet;
+}
+
+void VISU_PipeLine::SetInput(TInput* theInput)
+{
+  myExtractGeometry->SetInput(theInput);
+  if((myInput = theInput))
+    myInput->Update();
+
+  Modified();
 }
 
 VISU_PipeLine::TMapper* VISU_PipeLine::GetMapper() { 
-  if (myInput) {
-    if (!myMapper->GetInput()) {
-      myInput->Update();
+  if(GetInput()){
+    if(!myMapper->GetInput()){
+      GetInput2()->Update();
       Build();
     }
     myMapper->Update();
@@ -141,7 +127,8 @@ void VISU_PipeLine::Update(){
   myMapper->Update();
 }
 
-int VISU_PipeLine::CheckAvailableMemory(const float& theSize){
+int VISU_PipeLine::CheckAvailableMemory(const float& theSize)
+{
   try{
     if(theSize > ULONG_MAX) return 0;
     size_t aSize = size_t(theSize);
@@ -173,56 +160,48 @@ float VISU_PipeLine::GetAvailableMemory(float theSize, float theMinSize){
 
 //------------------------ Clipping planes -----------------------------------
 
-bool VISU_PipeLine::IsImplicitFunctionUsed() const
+void VISU_PipeLine::AddClippingPlane(vtkPlane* thePlane)
 {
-  return myIsImplicitFunctionUsed;
-}
-
-void VISU_PipeLine::SetImplicitFunctionUsed (bool theIsImplicitFunctionUsed)
-{
-  if (myInput != NULL) {
-    if (theIsImplicitFunctionUsed) {
-      myExtractGeometry->SetInput(myInput);
-      myPassFilter0->SetInput(myExtractGeometry->GetOutput());
-    } else {
-      myPassFilter0->SetInput(myInput);
+  if(thePlane){
+    if(vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()){
+      vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
+      aFunction->AddItem(thePlane);
     }
   }
-
-  Modified();
-  myIsImplicitFunctionUsed = theIsImplicitFunctionUsed;
-//  SetStoreClippingMapping(myStoreClippingMapping);
 }
 
-vtkIdType VISU_PipeLine::AddClippingPlane (vtkPlane* thePlane)
+vtkPlane* VISU_PipeLine::GetClippingPlane(vtkIdType theID) const
 {
-  if (thePlane) {
-    myImplicitBoolean->GetFunction()->AddItem(thePlane);
-    myCippingPlaneCont.push_back(thePlane);
-    if (!IsImplicitFunctionUsed())
-      SetImplicitFunctionUsed(true);
+  vtkPlane* aPlane = NULL;
+  if(theID >= 0 && theID < GetNumberOfClippingPlanes()){
+    if(vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()){
+      vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
+      vtkImplicitFunction* aFun = NULL;
+      aFunction->InitTraversal();
+      for(vtkIdType anID = 0; anID <= theID; anID++)
+       aFun = aFunction->GetNextItem();
+      aPlane = dynamic_cast<vtkPlane*>(aFun);
+    }
   }
-  return myCippingPlaneCont.size();
+  return aPlane;
 }
 
 void VISU_PipeLine::RemoveAllClippingPlanes()
 {
-  myImplicitBoolean->GetFunction()->RemoveAllItems();
-  myImplicitBoolean->GetFunction()->Modified(); // VTK bug
-  myCippingPlaneCont.clear();
-  SetImplicitFunctionUsed(false);
-}
-
-vtkIdType VISU_PipeLine::GetNumberOfClippingPlanes()
-{
-  return myCippingPlaneCont.size();
+  if(vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()){
+    vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
+    aFunction->RemoveAllItems();
+    aBoolean->Modified(); // VTK bug
+  }
 }
 
-vtkPlane* VISU_PipeLine::GetClippingPlane (vtkIdType theID)
+vtkIdType VISU_PipeLine::GetNumberOfClippingPlanes() const
 {
-  if (theID >= myCippingPlaneCont.size())
-    return NULL;
-  return myCippingPlaneCont[theID].Get();
+  if(vtkImplicitBoolean* aBoolean = myExtractGeometry->GetImplicitBoolean()){
+    vtkImplicitFunctionCollection* aFunction = aBoolean->GetFunction();
+    return aFunction->GetNumberOfItems();
+  }
+  return 0;
 }
 
 static void ComputeBoundsParam (vtkDataSet* theDataSet,
@@ -293,7 +272,7 @@ void VISU_PipeLine::SetPlaneParam (float theDir[3], float theDist, vtkPlane* the
 {
   thePlane->SetNormal(theDir);
   float anOrigin[3];
-  ::DistanceToPosition(myInput,theDir,theDist,anOrigin);
+  ::DistanceToPosition(GetInput(),theDir,theDist,anOrigin);
   thePlane->SetOrigin(anOrigin);
 }
 
@@ -303,5 +282,5 @@ void VISU_PipeLine::GetPlaneParam (float theDir[3], float& theDist, vtkPlane* th
 
   float anOrigin[3];
   thePlane->GetOrigin(anOrigin);
-  ::PositionToDistance(myInput,theDir,anOrigin,theDist);
+  ::PositionToDistance(GetInput(),theDir,anOrigin,theDist);
 }
index ebbe5d54729dba2d34ce36e7b0445fcf684f45e7..0f1670ced333d1f24f4f81e50f239778ff723288 100644 (file)
@@ -56,22 +56,24 @@ class vtkUnstructuredGrid;
 class vtkExtractGeometry;
 class vtkImplicitBoolean;
 class vtkPlane;
-class SALOME_PassThroughFilter;
+
+class SALOME_ExtractGeometry;
 
 typedef vtkUnstructuredGrid TInput;
 
 class VISU_PipeLine : public vtkObject{
 protected:
-  VISU_PipeLine();
-  VISU_PipeLine(const VISU_PipeLine&);
+
 public:
   vtkTypeMacro(VISU_PipeLine,vtkObject);
   virtual ~VISU_PipeLine();
+
   virtual void ShallowCopy(VISU_PipeLine *thePipeLine);
+  virtual void SameAs(VISU_PipeLine *thePipeLine);
 
 public:  
   virtual void SetInput(TInput* theInput);
-  virtual TInput* GetInput() { return myInput; }
+  virtual TInput* GetInput() const;
 
   typedef vtkDataSetMapper TMapper;
   virtual TMapper* GetMapper();
@@ -84,36 +86,26 @@ public:
                                  float theMinSize = 1024*1024.0);
 
   // Clipping planes
-  virtual void SetPlaneParam (float theDir[3], float theDist, vtkPlane* thePlane);
-  virtual void GetPlaneParam (float theDir[3], float& theDist, vtkPlane* thePlane);
-
-  virtual void RemoveAllClippingPlanes();
-  virtual vtkIdType GetNumberOfClippingPlanes();
-  virtual vtkPlane* GetClippingPlane (vtkIdType theID);
-  virtual vtkIdType AddClippingPlane (vtkPlane* thePlane); 
+  void RemoveAllClippingPlanes();
+  vtkIdType GetNumberOfClippingPlanes() const;
+  void AddClippingPlane(vtkPlane* thePlane); 
+  vtkPlane* GetClippingPlane(vtkIdType theID) const;
 
-  void SetImplicitFunctionUsed (bool theIsImplicitFunctionUsed);
-  bool IsImplicitFunctionUsed() const;
+  void SetPlaneParam(float theDir[3], float theDist, vtkPlane* thePlane);
+  void GetPlaneParam(float theDir[3], float& theDist, vtkPlane* thePlane);
 
 protected:
+  VISU_PipeLine();
+  VISU_PipeLine(const VISU_PipeLine&);
+
+  virtual TInput* GetInput2() const;
   virtual void Build() = 0;
-  virtual TInput* GetInput2();
   
   TMapper *myMapper;
+  TInput *myInput;
 
   // Clipping planes
-  bool myIsImplicitFunctionUsed;
-
-  vtkImplicitBoolean* myImplicitBoolean;
-  typedef TVTKSmartPtr<vtkPlane> TPlanePtr;
-  typedef std::vector<TPlanePtr> TCippingPlaneCont;
-  TCippingPlaneCont myCippingPlaneCont;
-
-  SALOME_PassThroughFilter* myPassFilter0;
-  vtkExtractGeometry* myExtractGeometry;
-
-private:
-  TInput *myInput;
+  TVTKSmartPtr<SALOME_ExtractGeometry> myExtractGeometry;
 };
 
 #endif
index e0b42eec92841068787284aafa2f08548a1b73ec..c34170a33098f4615ad912f2d87333640c848758 100644 (file)
@@ -55,8 +55,8 @@ VISU::Prs3d_i::Prs3d_i(Result_i* theResult, SALOMEDS::SObject_ptr theSObject) :
 
 void VISU::Prs3d_i::SameAs(const Prs3d_i* theOrigin)
 {
-  Prs3d_i* aOrigin = const_cast<Prs3d_i*>(theOrigin);
-  myPipeLine->ShallowCopy(aOrigin->GetPL());
+  if(Prs3d_i* aOrigin = const_cast<Prs3d_i*>(theOrigin))
+    myPipeLine->SameAs(aOrigin->GetPL());
 }
 
 VISU::Prs3d_i::~Prs3d_i() {
@@ -136,24 +136,24 @@ void VISU::Prs3d_i::GetBounds(float aBounds[6]){
 
 // Clipping planes
 
-void VISU::Prs3d_i::SetPlaneParam (float theDir[3], float theDist, vtkPlane* thePlane) {
-  myPipeLine->SetPlaneParam(theDir, theDist, thePlane);
-}
-
-void VISU::Prs3d_i::RemoveAllClippingPlanes() {
+void VISU::Prs3d_i::RemoveAllClippingPlanes(){
   myPipeLine->RemoveAllClippingPlanes();
 }
 
-vtkIdType VISU::Prs3d_i::GetNumberOfClippingPlanes() {
+vtkIdType VISU::Prs3d_i::GetNumberOfClippingPlanes() const{
   return myPipeLine->GetNumberOfClippingPlanes();
 }
 
-vtkPlane* VISU::Prs3d_i::GetClippingPlane (vtkIdType theID) {
+void VISU::Prs3d_i::AddClippingPlane(vtkPlane* thePlane){
+  myPipeLine->AddClippingPlane(thePlane);
+}
+
+vtkPlane* VISU::Prs3d_i::GetClippingPlane(vtkIdType theID) const{
   return myPipeLine->GetClippingPlane(theID);
 }
 
-vtkIdType VISU::Prs3d_i::AddClippingPlane (vtkPlane* thePlane) {
-  return myPipeLine->AddClippingPlane(thePlane);
+void VISU::Prs3d_i::SetPlaneParam (float theDir[3], float theDist, vtkPlane* thePlane) {
+  myPipeLine->SetPlaneParam(theDir, theDist, thePlane);
 }
 
 VISU::Result_i* VISU::GetResult(SALOMEDS::SObject_ptr theSObject){
index feb1a8abe659005347c1e23ce9f901d03d34e1d2..d2a01dc3396910121d54abb1004a30c94a8b0963 100644 (file)
@@ -42,6 +42,7 @@ class vtkUnstructuredGrid;
 
 namespace VISU{
   class Result_i;
+
   class Prs3d_i : 
     public virtual POA_VISU::Prs3d,
     public virtual SALOME::GenericObj_i,
@@ -89,14 +90,14 @@ namespace VISU{
     virtual SALOMEDS::SObject_var GetSObject();
 
     // Clipping planes
-    virtual void SetPlaneParam (float theDir[3], float theDist, vtkPlane* thePlane);
-    //virtual void GetPlaneParam (float theDir[3], float& theDist, vtkPlane* thePlane);
+    void RemoveAllClippingPlanes();
+    vtkIdType GetNumberOfClippingPlanes() const;
+    void AddClippingPlane(vtkPlane* thePlane); 
+    vtkPlane* GetClippingPlane(vtkIdType theID) const;
 
-    virtual void RemoveAllClippingPlanes();
-    virtual vtkIdType GetNumberOfClippingPlanes();
-    virtual vtkPlane* GetClippingPlane (vtkIdType theID);
-    virtual vtkIdType AddClippingPlane (vtkPlane* thePlane); 
+    void SetPlaneParam(float theDir[3], float theDist, vtkPlane* thePlane);
   };
+
   Result_i* GetResult(SALOMEDS::SObject_ptr theSObject);
   template<class TPrs3d>
   Storable* Restore (SALOMEDS::SObject_ptr theSObject,