]> SALOME platform Git repositories - modules/gui.git/commitdiff
Salome HOME
update of the package
authorsrn <srn@opencascade.com>
Mon, 6 Feb 2006 10:13:35 +0000 (10:13 +0000)
committersrn <srn@opencascade.com>
Mon, 6 Feb 2006 10:13:35 +0000 (10:13 +0000)
src/VTKViewer/Makefile.in
src/VTKViewer/VTKViewer_AppendFilter.cxx [new file with mode: 0644]
src/VTKViewer/VTKViewer_AppendFilter.h [new file with mode: 0644]
src/VTKViewer/VTKViewer_ConvexTool.cxx
src/VTKViewer/VTKViewer_GeometryFilter.cxx
src/VTKViewer/VTKViewer_GeometryFilter.h
src/VTKViewer/resources/VTKViewerM_images.po [new file with mode: 0644]
src/VTKViewer/resources/VTKViewerM_msg_en.po [new file with mode: 0644]
src/VTKViewer/resources/view_graduated_axes.png [new file with mode: 0755]
src/VTKViewer/resources/view_scaling.png [new file with mode: 0644]

index 4c72c0d5c1f1c3e37e003fbb754a70418b995ecd..f36e7362052accdec68b871154d859570ed5663c 100755 (executable)
@@ -18,6 +18,7 @@ EXPORT_HEADERS= VTKViewer_Actor.h \
                VTKViewer_ConvexTool.h \
                VTKViewer_Filter.h \
                VTKViewer_GeometryFilter.h \
+               VTKViewer_AppendFilter.h \
                VTKViewer_Algorithm.h \
                VTKViewer.h \
                VTKViewer_InteractorStyle.h \
@@ -36,9 +37,12 @@ EXPORT_HEADERS= VTKViewer_Actor.h \
                VTKViewer_ViewWindow.h \
                VTKViewer_Functor.h
 
-PO_FILES = VTKViewer_images.po \
-          VTKViewer_msg_en.po 
-                    
+PO_FILES = \
+       VTKViewer_images.po \
+       VTKViewer_msg_en.po \
+       VTKViewerM_images.po \
+       VTKViewerM_msg_en.po 
+
 # Libraries targets
 LIB = libVTKViewer.la
 
@@ -47,6 +51,7 @@ LIB_SRC= VTKViewer_Actor.cxx \
         VTKViewer_ExtractUnstructuredGrid.cxx \
         VTKViewer_Filter.cxx \
         VTKViewer_GeometryFilter.cxx \
+        VTKViewer_AppendFilter.cxx \
         VTKViewer_InteractorStyle.cxx \
         VTKViewer_PassThroughFilter.cxx \
         VTKViewer_RectPicker.cxx \
@@ -70,23 +75,6 @@ LIB_MOC = \
          VTKViewer_ViewManager.h \
          VTKViewer_ViewModel.h \
          VTKViewer_ViewWindow.h 
-RESOURCES_FILES = \
-       view_back.png \
-       view_bottom.png \
-       view_camera_dump.png \
-       view_fitall.png \
-       view_fitarea.png \
-       view_front.png \
-       view_glpan.png \
-       view_left.png \
-       view_pan.png \
-       view_reset.png \
-       view_right.png \
-       view_rotate.png \
-       view_top.png \
-       view_triedre.png \
-       view_zoom.png
 
 CPPFLAGS+=$(QT_INCLUDES) $(VTK_INCLUDES) $(OCC_INCLUDES)
 
diff --git a/src/VTKViewer/VTKViewer_AppendFilter.cxx b/src/VTKViewer/VTKViewer_AppendFilter.cxx
new file mode 100644 (file)
index 0000000..dce05c4
--- /dev/null
@@ -0,0 +1,363 @@
+//  SALOME OBJECT : kernel of SALOME component
+//
+//  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 
+//
+//
+//
+//  File   : VTKViewer_GeometryFilter.cxx
+//  Author : 
+//  Module : SALOME
+//  $Header$
+
+#include "VTKViewer_AppendFilter.h"
+
+#include "VTKViewer_ConvexTool.h"
+
+#include <vtkSmartPointer.h>
+#include <vtkCellArray.h>
+#include <vtkCellData.h>
+#include <vtkGenericCell.h>
+#include <vtkHexahedron.h>
+#include <vtkMergePoints.h>
+#include <vtkObjectFactory.h>
+#include <vtkPointData.h>
+#include <vtkPolyData.h>
+#include <vtkPyramid.h>
+#include <vtkStructuredGrid.h>
+#include <vtkTetra.h>
+#include <vtkUnsignedCharArray.h>
+#include <vtkUnstructuredGrid.h>
+#include <vtkVoxel.h>
+#include <vtkWedge.h>
+#include <vtkDataSetCollection.h>
+
+#include <vector>
+#include <map>
+using namespace std;
+
+
+#ifdef _DEBUG_
+//static int MYDEBUG = 0;
+//static int MYDEBUGWITHFILES = 0;
+#else
+//static int MYDEBUG = 0;
+//static int MYDEBUGWITHFILES = 0;
+#endif
+
+#if defined __GNUC__
+  #if __GNUC__ == 2
+    #define __GNUC_2__
+  #endif
+#endif
+
+vtkCxxRevisionMacro(VTKViewer_AppendFilter, "$Revision$");
+vtkStandardNewMacro(VTKViewer_AppendFilter);
+
+VTKViewer_AppendFilter
+::VTKViewer_AppendFilter() 
+{
+  myDoMappingFlag=false;
+}
+
+VTKViewer_AppendFilter
+::~VTKViewer_AppendFilter()
+{}
+
+void
+VTKViewer_AppendFilter
+::SetDoMappingFlag(const bool theFlag)
+{
+  myDoMappingFlag=theFlag;
+}
+
+bool 
+VTKViewer_AppendFilter
+::DoMappingFlag() const
+{
+  return myDoMappingFlag;
+}
+
+void
+VTKViewer_AppendFilter
+::SetPoints(vtkPoints* thePoints)
+{
+  myPoints = thePoints;
+}
+
+vtkPoints*
+VTKViewer_AppendFilter
+::GetPoints()
+{
+  return myPoints.GetPointer();
+}
+
+void
+VTKViewer_AppendFilter
+::Execute()
+{
+  if (myPoints.GetPointer()) {
+    MakeOutput();
+  }
+  else {
+    vtkAppendFilter::Execute();
+  }
+  if (myDoMappingFlag){
+    DoMapping();
+  }
+}
+
+void 
+VTKViewer_AppendFilter
+::Reset()
+{
+  myNodeIds.clear();
+  myCellIds.clear();
+  myNodeRanges.clear();
+  myCellRanges.clear();
+  myNodeMapObjIDVtkID.clear();
+  myCellMapObjIDVtkID.clear();
+}
+//==================================================================
+// function: DoMapping
+// purpose :
+//==================================================================
+void 
+VTKViewer_AppendFilter
+::DoMapping()
+{
+  int i, j, i1, i2, iNodeCnt, iCellCnt; 
+  IteratorOfDataMapOfIntegerInteger aMapIt;
+  vtkIdType aNbPnts, aNbCells, aId;
+  vtkDataSet *pDS;
+  //
+  Reset();
+  //
+  iNodeCnt=0;
+  iCellCnt=0;
+  for (i=0; i<NumberOfInputs; ++i) {
+    pDS=(vtkDataSet *)Inputs[i];
+    //
+    // Nodes
+    if (!myPoints.GetPointer()) {
+      aNbPnts=pDS->GetNumberOfPoints();
+      i1=myNodeIds.size();
+      i2=i1+aNbPnts-1;
+      myNodeRanges.push_back(i1);
+      myNodeRanges.push_back(i2);
+      //
+      for(j=0; j<aNbPnts; ++j) {
+       aId=(vtkIdType)j;
+       myNodeIds.push_back(aId);
+       //
+       aMapIt=myNodeMapObjIDVtkID.find(aId);
+       if (aMapIt==myNodeMapObjIDVtkID.end()) {
+         // if not found
+         myNodeMapObjIDVtkID[aId]=iNodeCnt;
+       }
+       ++iNodeCnt;
+      }
+    }
+    //
+    // Cells
+    aNbCells=pDS->GetNumberOfCells();
+    i1=myCellIds.size();
+    i2=i1+aNbCells-1;
+    myCellRanges.push_back(i1);
+    myCellRanges.push_back(i2);
+    for(j=0; j<aNbCells; ++j) {
+      aId=(vtkIdType)j;
+      myCellIds.push_back(aId);
+      //
+      aMapIt=myCellMapObjIDVtkID.find(aId);
+      if (aMapIt==myCellMapObjIDVtkID.end()) {
+       // if not found
+       myCellMapObjIDVtkID[aId]=iCellCnt;
+      }
+      ++iCellCnt;
+    }
+  }
+}
+
+//---------------------------------------------------------------
+vtkIdType
+VTKViewer_AppendFilter
+::GetPointOutputID(vtkIdType theInputID)
+{
+  if (myPoints.GetPointer()) {
+    return theInputID;
+  }
+  //
+  int aVtkID=-1;
+  IteratorOfDataMapOfIntegerInteger aMapIt;
+  //
+  aMapIt=myNodeMapObjIDVtkID.find(theInputID);
+  if (aMapIt!=myNodeMapObjIDVtkID.end()) {
+    // found
+    PairOfDataMapOfIntegerInteger& aPair=(*aMapIt);
+    aVtkID=aPair.second;
+  }
+  return aVtkID;
+}
+
+
+//---------------------------------------------------------------
+vtkIdType 
+VTKViewer_AppendFilter
+::GetCellOutputID(vtkIdType theInputID)
+{
+  int aVtkID=-1;
+  IteratorOfDataMapOfIntegerInteger aMapIt;
+  //
+  aMapIt=myCellMapObjIDVtkID.find(theInputID);
+  if (aMapIt!=myCellMapObjIDVtkID.end()) {
+    // found
+    PairOfDataMapOfIntegerInteger& aPair=(*aMapIt);
+    aVtkID=aPair.second;
+  }
+  return aVtkID;
+}
+
+
+//---------------------------------------------------------------
+vtkIdType 
+VTKViewer_AppendFilter
+::GetPointInputID(vtkIdType theOutputID, 
+                 vtkIdType& theInputDataSetID)
+{
+  if (myPoints.GetPointer()) {
+    theInputDataSetID=0;
+    return theOutputID;
+  }
+  //
+  int aNb, aNbRanges, aRetID, i, i1, i2, j;
+  //
+  aRetID=-1;
+  theInputDataSetID=-1;
+  //
+  aNb=myNodeIds.size();
+  if (theOutputID<0 ||  theOutputID>=aNb) {
+    return aRetID;
+  }
+  //
+  aRetID=(int)myNodeIds[theOutputID];
+  //
+  aNbRanges=myNodeRanges.size()/2;
+  for (i=0; i<aNbRanges; ++i) {
+    j=2*i;
+    i1=myNodeRanges[j];
+    i2=myNodeRanges[j+1];
+    if (theOutputID>=i1 && theOutputID<=i2) {
+      theInputDataSetID=i;
+    }
+  }
+  //
+  return aRetID;
+}
+
+
+//---------------------------------------------------------------
+vtkIdType 
+VTKViewer_AppendFilter
+::GetCellInputID(vtkIdType theOutputID, 
+                vtkIdType& theInputDataSetID)
+{
+  int aNb, aNbRanges, aRetID, i, i1, i2, j;
+  //
+  aRetID=-1;
+  theInputDataSetID=-1;
+  //
+  aNb=myCellIds.size();
+  if (theOutputID<0 ||  theOutputID>=aNb) {
+    return aRetID;
+  }
+  //
+  aRetID=(int)myCellIds[theOutputID];
+  //
+  aNbRanges=myCellRanges.size()/2;
+  for (i=0; i<aNbRanges; ++i) {
+    j=2*i;
+    i1=myCellRanges[j];
+    i2=myCellRanges[j+1];
+    if (theOutputID>=i1 && theOutputID<=i2) {
+      theInputDataSetID=i;
+    }
+  }
+  //
+  return aRetID;
+}
+
+
+//---------------------------------------------------------------
+void 
+VTKViewer_AppendFilter
+::MakeOutput()
+{
+  int idx;
+  vtkIdType numPts, numCells, newCellId, cellId;
+  vtkCellData *cd;
+  vtkIdList *ptIds;
+  vtkDataSet *ds;
+  vtkUnstructuredGrid *output = this->GetOutput();
+  //
+  numPts = myPoints->GetNumberOfPoints();
+  if (numPts < 1) {
+    return;
+  }
+  //
+  numCells = 0;
+  for (idx = 0; idx < this->NumberOfInputs; ++idx) {
+    ds = (vtkDataSet *)(this->Inputs[idx]);
+    if (ds != NULL)  {
+      if ( ds->GetNumberOfPoints() <= 0 && ds->GetNumberOfCells() <= 0 )  {
+        continue; //no input, just skip
+      }
+      numCells += ds->GetNumberOfCells();
+    }//if non-empty dataset
+  }//for all inputs
+  if (numCells < 1) {
+    return;
+  }
+  //
+  // Now can allocate memory
+  output->Allocate(numCells); 
+  ptIds = vtkIdList::New(); 
+  ptIds->Allocate(VTK_CELL_SIZE);
+  //
+  // Append each input dataset together
+  //
+  // 1.points
+  output->SetPoints(myPoints.GetPointer());
+  // 2.cells
+  for (idx = 0; idx < this->NumberOfInputs; ++idx) {
+    ds = (vtkDataSet *)(this->Inputs[idx]);
+    if (ds != NULL) {
+      numCells = ds->GetNumberOfCells(); 
+      cd = ds->GetCellData();
+      // copy cell and cell data
+      for (cellId=0; cellId<numCells; cellId++)  {
+        ds->GetCellPoints(cellId, ptIds);
+        newCellId = output->InsertNextCell(ds->GetCellType(cellId), ptIds);
+      }
+    }
+  }
+  //
+  ptIds->Delete();
+}
+
diff --git a/src/VTKViewer/VTKViewer_AppendFilter.h b/src/VTKViewer/VTKViewer_AppendFilter.h
new file mode 100644 (file)
index 0000000..f9c6b24
--- /dev/null
@@ -0,0 +1,91 @@
+#ifndef VTKVIEWER_APPENDFILTER_H
+#define VTKVIEWER_APPENDFILTER_H
+
+#include "VTKViewer.h"
+
+#include <vtkAppendFilter.h>
+#include <vtkSmartPointer.h>
+
+#include <vector>
+#include <map>
+
+class vtkPoints;
+
+/*! \brief This class used same as vtkAppendFilter. See documentation on VTK for more information.
+ */
+class VTKVIEWER_EXPORT VTKViewer_AppendFilter : public vtkAppendFilter 
+{
+public:
+  /*! \fn static VTKViewer_AppendFilter *New()
+   */
+  static VTKViewer_AppendFilter *New();
+  
+  /*! \fn vtkTypeRevisionMacro(VTKViewer_AppendFilter, vtkAppendFilter)
+   *  \brief VTK type revision macros.
+   */
+  vtkTypeRevisionMacro(VTKViewer_AppendFilter, vtkAppendFilter);
+
+  void SetDoMappingFlag(const bool theFlag);
+
+  bool DoMappingFlag() const;
+
+  void
+  SetPoints(vtkPoints* thePoints);
+
+  vtkPoints*
+  GetPoints();
+
+  vtkIdType
+  GetPointOutputID(vtkIdType theInputID);
+
+  vtkIdType
+  GetCellOutputID(vtkIdType theInputID);
+
+  vtkIdType 
+  GetPointInputID(vtkIdType theOutputID, 
+                 vtkIdType& theInputDataSetID);
+
+  vtkIdType
+  GetCellInputID(vtkIdType theOutputID, 
+                vtkIdType& theInputDataSetID);
+
+protected:
+  /*! \fn VTKViewer_AppendFilter();
+   * \brief Constructor
+   */
+  VTKViewer_AppendFilter();
+  /*! \fn ~VTKViewer_AppendFilter();
+   * \brief Destructor.
+   */
+  ~VTKViewer_AppendFilter();
+  /*! \fn void Execute();
+   * \brief Filter culculation method.
+   */
+  virtual void Execute();
+  //
+  void DoMapping();
+
+  void Reset();
+
+  void MakeOutput();
+
+  //
+  vtkSmartPointer<vtkPoints> myPoints;
+
+private:
+  typedef std::vector<vtkIdType> TVectorId;
+  typedef std::vector<int> VectorInt;
+  typedef std::map <int,int>                  DataMapOfIntegerInteger;
+  typedef DataMapOfIntegerInteger::iterator   IteratorOfDataMapOfIntegerInteger;
+  typedef DataMapOfIntegerInteger::value_type PairOfDataMapOfIntegerInteger;
+private:
+  bool      myDoMappingFlag;
+  TVectorId myNodeIds;
+  TVectorId myCellIds;
+  VectorInt myNodeRanges;
+  VectorInt myCellRanges;
+  DataMapOfIntegerInteger myNodeMapObjIDVtkID;
+  DataMapOfIntegerInteger myCellMapObjIDVtkID;
+};
+
+#endif
index a4dded2cc20d6fa4fde630c48a97e3783a6b5b4f..679e484289b2c5d38d31ad1f015b39ceeeeb15ca 100644 (file)
@@ -48,6 +48,10 @@ typedef std::map<vtkIdType,TUIDS> TPTOIDS; // id points -> unique ids
 
 namespace CONVEX_TOOL
 {
+  // all pairs
+  typedef std::pair<vtkIdType,vtkIdType> TPair;
+  typedef std::set<TPair> TSet;
+  
   void
   WriteToFile(vtkUnstructuredGrid* theDataSet, const std::string& theFileName)
   {
@@ -384,7 +388,7 @@ static void GetAndRemoveIdsOnOneLine(vtkPoints* points,
   float P[3][3];
   vtkIdType current_points_ids[2];
   if(MYDEBUG_REMOVE) 
-    if(input_two_points_ids.size()!=2) cout << "Error. Must be two ids in variable input_two_points_ids"<<endl;
+    if(input_two_points_ids.size()!=2) cout << "Error. Must be two ids in variable input_two_points_ids="<<input_two_points_ids.size()<<endl;
   TUIDS::const_iterator aInPointsIter = input_two_points_ids.begin();
   for(int i=0;i<2 && aInPointsIter!=input_two_points_ids.end();aInPointsIter++,i++){
     current_points_ids[i] = *aInPointsIter;
@@ -518,15 +522,13 @@ static void GetAndRemoveIdsOnOneLine(vtkPoints* points,
   out_two_points_ids.insert(current_points_ids[1]);
 }
 
-static vtkConvexPointSet* RemoveAllUnneededPoints(vtkConvexPointSet* convex){
+static vtkSmartPointer<vtkConvexPointSet> RemoveAllUnneededPoints(vtkConvexPointSet* convex){
   vtkSmartPointer<vtkConvexPointSet> out = vtkConvexPointSet::New();
   
   TUIDS two_points,input_points,out_two_points_ids,removed_points_ids,loc_removed_points_ids;
   vtkIdList* aPointIds = convex->GetPointIds();
   int numIds = aPointIds->GetNumberOfIds();
-  // all pairs
-  typedef std::pair<vtkIdType,vtkIdType> TPair;
-  typedef std::set<TPair> TSet;
+  if (numIds<2) return out;
   TSet good_point_ids;
   TSet aLists[numIds-2];
   for(int i=0;i<numIds-2;i++){
@@ -624,7 +626,7 @@ static vtkConvexPointSet* RemoveAllUnneededPoints(vtkConvexPointSet* convex){
   out->Modified();
   out->Initialize();
   
-  return out.GetPointer();
+  return out;
 }
 
 void GetPolygonalFaces(vtkUnstructuredGrid* theGrid,int cellId,TCellArray &outputCellArray)
@@ -632,7 +634,7 @@ void GetPolygonalFaces(vtkUnstructuredGrid* theGrid,int cellId,TCellArray &outpu
   if (theGrid->GetCellType(cellId) == VTK_CONVEX_POINT_SET){
     // get vtkCell type = VTK_CONVEX_POINT_SET
     if(vtkConvexPointSet* convex_in = static_cast<vtkConvexPointSet*>(theGrid->GetCell(cellId))){
-      vtkConvexPointSet* convex = RemoveAllUnneededPoints(convex_in);
+      vtkSmartPointer<vtkConvexPointSet> convex = RemoveAllUnneededPoints(convex_in);
       TCellArray f2points;
       float convex_center[3]; // convex center point coorinat
       int aNbFaces = convex->GetNumberOfFaces();
index 10dec386c4228c636cdbca7446a41d14b3356cc2..d1baa5b115af77e360ce622ebad33b0a1a1a811e 100755 (executable)
 using namespace std;
 
 
-#ifdef _DEBUG_
-static int MYDEBUG = 0;
-#else
-static int MYDEBUG = 0;
-#endif
-
 #if defined __GNUC__
   #if __GNUC__ == 2
     #define __GNUC_2__
   #endif
 #endif
 
+//----------------------------------------------------------------------------
 vtkCxxRevisionMacro(VTKViewer_GeometryFilter, "$Revision$");
 vtkStandardNewMacro(VTKViewer_GeometryFilter);
 
 VTKViewer_GeometryFilter::VTKViewer_GeometryFilter(): 
   myShowInside(0),
-  myStoreMapping(0)
+  myStoreMapping(0),
+  myIsWireframeMode(0)
 {}
 
 
@@ -76,6 +72,7 @@ VTKViewer_GeometryFilter::~VTKViewer_GeometryFilter()
 {}
 
 
+//----------------------------------------------------------------------------
 void VTKViewer_GeometryFilter::Execute()
 {
   vtkDataSet *input= this->GetInput();
@@ -94,22 +91,7 @@ void VTKViewer_GeometryFilter::Execute()
 }
 
 
-void VTKViewer_GeometryFilter::SetStoreMapping(int theStoreMapping){
-  myStoreMapping = theStoreMapping;
-  this->Modified();
-}
-
-
-vtkIdType VTKViewer_GeometryFilter::GetElemObjId(int theVtkID){
-  if(myVTK2ObjIds.empty() || theVtkID > myVTK2ObjIds.size()) return -1;
-#if defined __GNUC_2__
-  return myVTK2ObjIds[theVtkID];
-#else
-  return myVTK2ObjIds.at(theVtkID);
-#endif
-}
-
-
+//----------------------------------------------------------------------------
 void VTKViewer_GeometryFilter::UnstructuredGridExecute()
 {
   vtkUnstructuredGrid *input= (vtkUnstructuredGrid *)this->GetInput();
@@ -128,7 +110,6 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
   vtkPointData *outputPD = output->GetPointData();
   
   vtkCellData *outputCD = output->GetCellData();
-  //vtkCellArray *Verts, *Lines, *Polys, *Strips;
   vtkIdList *cellIds, *faceIds;
   char *cellVis;
   vtkIdType newCellId;
@@ -190,14 +171,6 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
   outputCD->CopyAllocate(cd,numCells,numCells/2);
 
   output->Allocate(numCells/4+1,numCells);
-  //Verts = vtkCellArray::New();
-  //Verts->Allocate(numCells/4+1,numCells);
-  //Lines = vtkCellArray::New();
-  //Lines->Allocate(numCells/4+1,numCells);
-  //Polys = vtkCellArray::New();
-  //Polys->Allocate(numCells/4+1,numCells);
-  //Strips = vtkCellArray::New();
-  //Strips->Allocate(numCells/4+1,numCells);
   
   // Loop over the cells determining what's visible
   if (!allVisible)
@@ -278,9 +251,8 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
         case VTK_LINE: 
         case VTK_POLY_LINE:
           newCellId = output->InsertNextCell(aCellType,npts,pts);
-         if(myStoreMapping){
-           myVTK2ObjIds.push_back(cellId); //apo
-         }
+         if(myStoreMapping)
+           myVTK2ObjIds.push_back(cellId);
           outputCD->CopyData(cd,cellId,newCellId);
           break;
 
@@ -288,25 +260,22 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
         case VTK_QUAD:
         case VTK_POLYGON:
           newCellId = output->InsertNextCell(aCellType,npts,pts);
-         if(myStoreMapping){
-           myVTK2ObjIds.push_back(cellId); //apo
-         }
+         if(myStoreMapping)
+           myVTK2ObjIds.push_back(cellId);
           outputCD->CopyData(cd,cellId,newCellId);
           break;
 
         case VTK_TRIANGLE_STRIP:
           newCellId = output->InsertNextCell(aCellType,npts,pts);
-         if(myStoreMapping){
-           myVTK2ObjIds.push_back(cellId); //apo
-         }
+         if(myStoreMapping)
+           myVTK2ObjIds.push_back(cellId);
           outputCD->CopyData(cd,cellId,newCellId);
           break;
 
         case VTK_PIXEL:
           newCellId = output->InsertNextCell(aCellType,npts,pts);
-         if(myStoreMapping){
-           myVTK2ObjIds.push_back(cellId); //apo
-         }
+         if(myStoreMapping)
+           myVTK2ObjIds.push_back(cellId);
          outputCD->CopyData(cd,cellId,newCellId);
           break;
          
@@ -350,20 +319,17 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
             faceIds->InsertNextId(pts[faceVerts[0]]);
             faceIds->InsertNextId(pts[faceVerts[1]]);
             faceIds->InsertNextId(pts[faceVerts[2]]);
-            numFacePts = 3;
            aCellType = VTK_TRIANGLE;
+            numFacePts = 3;
             input->GetCellNeighbors(cellId, faceIds, cellIds);
-            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside == 1 ||
+            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside ||
                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
               {
               for ( i=0; i < numFacePts; i++)
-                {
                 aNewPts[i] = pts[faceVerts[i]];
-                }
               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-             if(myStoreMapping){
-               myVTK2ObjIds.push_back(cellId); //apo
-             }
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
               outputCD->CopyData(cd,cellId,newCellId);
               }
             }
@@ -378,20 +344,17 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
             faceIds->InsertNextId(pts[faceVerts[1]]);
             faceIds->InsertNextId(pts[faceVerts[2]]);
             faceIds->InsertNextId(pts[faceVerts[3]]);
-            numFacePts = 4;
            aCellType = VTK_QUAD;
+            numFacePts = 4;
             input->GetCellNeighbors(cellId, faceIds, cellIds);
-            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside == 1 || 
+            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside || 
                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
               {
               for ( i=0; i < numFacePts; i++)
-                {
                 aNewPts[i] = pts[faceVerts[PixelConvert[i]]];
-                }
               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-             if(myStoreMapping){
-               myVTK2ObjIds.push_back(cellId); //apo
-             }
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
               outputCD->CopyData(cd,cellId,newCellId);
               }
             }
@@ -406,20 +369,17 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
             faceIds->InsertNextId(pts[faceVerts[1]]);
             faceIds->InsertNextId(pts[faceVerts[2]]);
             faceIds->InsertNextId(pts[faceVerts[3]]);
-            numFacePts = 4;
            aCellType = VTK_QUAD;
+            numFacePts = 4;
             input->GetCellNeighbors(cellId, faceIds, cellIds);
-            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside == 1 ||
+            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside ||
                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
               {
               for ( i=0; i < numFacePts; i++)
-                {
                 aNewPts[i] = pts[faceVerts[i]];
-                }
               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-             if(myStoreMapping){
-               myVTK2ObjIds.push_back(cellId); //apo
-             }
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
               outputCD->CopyData(cd,cellId,newCellId);
               }
             }
@@ -433,26 +393,23 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
             faceIds->InsertNextId(pts[faceVerts[0]]);
             faceIds->InsertNextId(pts[faceVerts[1]]);
             faceIds->InsertNextId(pts[faceVerts[2]]);
-            numFacePts = 3;
            aCellType = VTK_TRIANGLE;
+            numFacePts = 3;
             if (faceVerts[3] >= 0)
               {
               faceIds->InsertNextId(pts[faceVerts[3]]);
-              numFacePts = 4;
              aCellType = VTK_QUAD;
+              numFacePts = 4;
               }
             input->GetCellNeighbors(cellId, faceIds, cellIds);
-            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside == 1 || 
+            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside || 
                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
               {
               for ( i=0; i < numFacePts; i++)
-                {
                 aNewPts[i] = pts[faceVerts[i]];
-                }
               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-             if(myStoreMapping){
-               myVTK2ObjIds.push_back(cellId); //apo
-             }
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
               outputCD->CopyData(cd,cellId,newCellId);
               }
             }
@@ -466,223 +423,319 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
             faceIds->InsertNextId(pts[faceVerts[0]]);
             faceIds->InsertNextId(pts[faceVerts[1]]);
             faceIds->InsertNextId(pts[faceVerts[2]]);
-            numFacePts = 3;
            aCellType = VTK_TRIANGLE;
+            numFacePts = 3;
             if (faceVerts[3] >= 0)
               {
               faceIds->InsertNextId(pts[faceVerts[3]]);
-              numFacePts = 4;
              aCellType = VTK_QUAD;
+              numFacePts = 4;
               }
             input->GetCellNeighbors(cellId, faceIds, cellIds);
-            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside == 1 || 
+            if ( cellIds->GetNumberOfIds() <= 0 || myShowInside || 
                  (!allVisible && !cellVis[cellIds->GetId(0)]) )
               {
               for ( i=0; i < numFacePts; i++)
-                {
                 aNewPts[i] = pts[faceVerts[i]];
-                }
               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-             if(myStoreMapping){
-               myVTK2ObjIds.push_back(cellId); //apo
-             }
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
               outputCD->CopyData(cd,cellId,newCellId);
               }
             }
           break;
        }
         //Quadratic cells
-        case VTK_QUADRATIC_EDGE: {
-         newCellId = output->InsertNextCell(VTK_POLY_LINE,npts,pts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-
-          break;
-       }
-        case VTK_QUADRATIC_TRIANGLE: {
-         numFacePts = 6;
-         aCellType = VTK_POLYGON;
-         
-         aNewPts[0] = pts[0];
-         aNewPts[1] = pts[3];
-         aNewPts[2] = pts[1];
-         aNewPts[3] = pts[4];
-         aNewPts[4] = pts[2];
-         aNewPts[5] = pts[5];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-         break;
-       }
-        case VTK_QUADRATIC_QUAD: {
-         numFacePts = 8;
-         aCellType = VTK_POLYGON;
-         
-         aNewPts[0] = pts[0];
-         aNewPts[1] = pts[4];
-         aNewPts[2] = pts[1];
-         aNewPts[3] = pts[5];
-         aNewPts[4] = pts[2];
-         aNewPts[5] = pts[6];
-         aNewPts[6] = pts[3];
-         aNewPts[7] = pts[7];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-         break;
-       }
-        case VTK_QUADRATIC_TETRA: {
-         numFacePts = 8;
-         aCellType = VTK_POLYGON;
-         
-         aNewPts[0] = pts[0];
-         aNewPts[1] = pts[4];
-         aNewPts[2] = pts[1];
-         aNewPts[3] = pts[5];
-         aNewPts[4] = pts[2];
-         aNewPts[5] = pts[6];
-         aNewPts[6] = pts[3];
-         aNewPts[7] = pts[7];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-         break;
-       }
-        case VTK_QUADRATIC_HEXAHEDRON: {
-         numFacePts = 8;
-         aCellType = VTK_POLYGON;
-         
-         //---------------------------------------------------------------
-         aNewPts[0] = pts[0];
-         aNewPts[1] = pts[8];
-         aNewPts[2] = pts[1];
-         aNewPts[3] = pts[17];
-         aNewPts[4] = pts[5];
-         aNewPts[5] = pts[12];
-         aNewPts[6] = pts[4];
-         aNewPts[7] = pts[16];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-
-         //---------------------------------------------------------------
-         aNewPts[0] = pts[1];
-         aNewPts[1] = pts[9];
-         aNewPts[2] = pts[2];
-         aNewPts[3] = pts[18];
-         aNewPts[4] = pts[6];
-         aNewPts[5] = pts[13];
-         aNewPts[6] = pts[5];
-         aNewPts[7] = pts[17];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-
-         //---------------------------------------------------------------
-         aNewPts[0] = pts[2];
-         aNewPts[1] = pts[10];
-         aNewPts[2] = pts[3];
-         aNewPts[3] = pts[19];
-         aNewPts[4] = pts[7];
-         aNewPts[5] = pts[14];
-         aNewPts[6] = pts[6];
-         aNewPts[7] = pts[18];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-
-         //---------------------------------------------------------------
-         aNewPts[0] = pts[3];
-         aNewPts[1] = pts[11];
-         aNewPts[2] = pts[0];
-         aNewPts[3] = pts[16];
-         aNewPts[4] = pts[4];
-         aNewPts[5] = pts[15];
-         aNewPts[6] = pts[7];
-         aNewPts[7] = pts[19];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-
-         //---------------------------------------------------------------
-         aNewPts[0] = pts[0];
-         aNewPts[1] = pts[8];
-         aNewPts[2] = pts[1];
-         aNewPts[3] = pts[9];
-         aNewPts[4] = pts[2];
-         aNewPts[5] = pts[10];
-         aNewPts[6] = pts[3];
-         aNewPts[7] = pts[11];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-
-         //---------------------------------------------------------------
-         aNewPts[0] = pts[4];
-         aNewPts[1] = pts[12];
-         aNewPts[2] = pts[5];
-         aNewPts[3] = pts[13];
-         aNewPts[4] = pts[6];
-         aNewPts[5] = pts[14];
-         aNewPts[6] = pts[7];
-         aNewPts[7] = pts[15];
-
-         newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
-         if(myStoreMapping)
-           myVTK2ObjIds.push_back(cellId);
-
-          outputCD->CopyData(cd,cellId,newCellId);
-
-         break;
-       }
-        } //switch
+        case VTK_QUADRATIC_EDGE:
+        case VTK_QUADRATIC_TRIANGLE:
+        case VTK_QUADRATIC_QUAD:
+        case VTK_QUADRATIC_TETRA:
+        case VTK_QUADRATIC_HEXAHEDRON:
+         if(!myIsWireframeMode){
+           vtkGenericCell *cell = vtkGenericCell::New();
+           input->GetCell(cellId,cell);
+           vtkIdList *pts = vtkIdList::New();  
+           vtkPoints *coords = vtkPoints::New();
+           vtkIdList *cellIds = vtkIdList::New();
+           vtkIdType newCellId;
+           
+           if ( cell->GetCellDimension() == 1 ) {
+             aCellType = VTK_LINE;
+             numFacePts = 2;
+             cell->Triangulate(0,pts,coords);
+             for (i=0; i < pts->GetNumberOfIds(); i+=2) {
+               aNewPts[0] = pts->GetId(i);
+               aNewPts[1] = pts->GetId(i+1);
+               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+               if(myStoreMapping)
+                 myVTK2ObjIds.push_back(cellId);
+               outputCD->CopyData(cd,cellId,newCellId);
+              }
+            }
+           else if ( cell->GetCellDimension() == 2 ) {
+             aCellType = VTK_TRIANGLE;
+             numFacePts = 3;
+             cell->Triangulate(0,pts,coords);
+             for (i=0; i < pts->GetNumberOfIds(); i+=3) {
+               aNewPts[0] = pts->GetId(i);
+               aNewPts[1] = pts->GetId(i+1);
+               aNewPts[2] = pts->GetId(i+2);
+               newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+               if(myStoreMapping)
+                 myVTK2ObjIds.push_back(cellId);
+               outputCD->CopyData(cd,cellId,newCellId);
+              }
+            } 
+           else //3D nonlinear cell
+            {
+             aCellType = VTK_TRIANGLE;
+             numFacePts = 3;
+             for (int j=0; j < cell->GetNumberOfFaces(); j++){
+               vtkCell *face = cell->GetFace(j);
+               input->GetCellNeighbors(cellId, face->PointIds, cellIds);
+               if ( cellIds->GetNumberOfIds() <= 0 || myShowInside ) {
+                 face->Triangulate(0,pts,coords);
+                 for (i=0; i < pts->GetNumberOfIds(); i+=3) {
+                   aNewPts[0] = pts->GetId(i);
+                   aNewPts[1] = pts->GetId(i+1);
+                   aNewPts[2] = pts->GetId(i+2);
+                   newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+                   if(myStoreMapping)
+                     myVTK2ObjIds.push_back(cellId);
+                   outputCD->CopyData(cd,cellId,newCellId);
+                  }
+                }
+              }
+            } //3d cell
+           cellIds->Delete();
+           coords->Delete();
+           pts->Delete();
+           cell->Delete();
+           break;
+          }else{
+           switch(aCellType){
+           case VTK_QUADRATIC_EDGE: {
+             aCellType = VTK_POLY_LINE;
+             numFacePts = 3;
+             
+             aNewPts[0] = pts[0];
+             aNewPts[2] = pts[1];
+             aNewPts[1] = pts[2];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             break;
+           }
+           case VTK_QUADRATIC_TRIANGLE: {
+             aCellType = VTK_POLYGON;
+             numFacePts = 6;
+             
+             aNewPts[0] = pts[0];
+             aNewPts[1] = pts[3];
+             aNewPts[2] = pts[1];
+             aNewPts[3] = pts[4];
+             aNewPts[4] = pts[2];
+             aNewPts[5] = pts[5];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             break;
+           }
+           case VTK_QUADRATIC_QUAD: {
+             aCellType = VTK_POLYGON;
+             numFacePts = 8;
+             
+             aNewPts[0] = pts[0];
+             aNewPts[1] = pts[4];
+             aNewPts[2] = pts[1];
+             aNewPts[3] = pts[5];
+             aNewPts[4] = pts[2];
+             aNewPts[5] = pts[6];
+             aNewPts[6] = pts[3];
+             aNewPts[7] = pts[7];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             break;
+           }
+           case VTK_QUADRATIC_TETRA: {
+             aCellType = VTK_POLYGON;
+             numFacePts = 6;
+             
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[0];
+             aNewPts[1] = pts[4];
+             aNewPts[2] = pts[1];
+             aNewPts[3] = pts[5];
+             aNewPts[4] = pts[2];
+             aNewPts[5] = pts[6];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[0];
+             aNewPts[1] = pts[7];
+             aNewPts[2] = pts[3];
+             aNewPts[3] = pts[8];
+             aNewPts[4] = pts[1];
+             aNewPts[5] = pts[4];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[1];
+             aNewPts[1] = pts[8];
+             aNewPts[2] = pts[3];
+             aNewPts[3] = pts[9];
+             aNewPts[4] = pts[2];
+             aNewPts[5] = pts[5];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[2];
+             aNewPts[1] = pts[9];
+             aNewPts[2] = pts[3];
+             aNewPts[3] = pts[7];
+             aNewPts[4] = pts[0];
+             aNewPts[5] = pts[6];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+
+             break;
+           }
+           case VTK_QUADRATIC_HEXAHEDRON: {
+             aCellType = VTK_POLYGON;
+             numFacePts = 8;
+             
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[0];
+             aNewPts[1] = pts[8];
+             aNewPts[2] = pts[1];
+             aNewPts[3] = pts[17];
+             aNewPts[4] = pts[5];
+             aNewPts[5] = pts[12];
+             aNewPts[6] = pts[4];
+             aNewPts[7] = pts[16];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[1];
+             aNewPts[1] = pts[9];
+             aNewPts[2] = pts[2];
+             aNewPts[3] = pts[18];
+             aNewPts[4] = pts[6];
+             aNewPts[5] = pts[13];
+             aNewPts[6] = pts[5];
+             aNewPts[7] = pts[17];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[2];
+             aNewPts[1] = pts[10];
+             aNewPts[2] = pts[3];
+             aNewPts[3] = pts[19];
+             aNewPts[4] = pts[7];
+             aNewPts[5] = pts[14];
+             aNewPts[6] = pts[6];
+             aNewPts[7] = pts[18];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[3];
+             aNewPts[1] = pts[11];
+             aNewPts[2] = pts[0];
+             aNewPts[3] = pts[16];
+             aNewPts[4] = pts[4];
+             aNewPts[5] = pts[15];
+             aNewPts[6] = pts[7];
+             aNewPts[7] = pts[19];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[0];
+             aNewPts[1] = pts[8];
+             aNewPts[2] = pts[1];
+             aNewPts[3] = pts[9];
+             aNewPts[4] = pts[2];
+             aNewPts[5] = pts[10];
+             aNewPts[6] = pts[3];
+             aNewPts[7] = pts[11];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             
+             //---------------------------------------------------------------
+             aNewPts[0] = pts[4];
+             aNewPts[1] = pts[12];
+             aNewPts[2] = pts[5];
+             aNewPts[3] = pts[13];
+             aNewPts[4] = pts[6];
+             aNewPts[5] = pts[14];
+             aNewPts[6] = pts[7];
+             aNewPts[7] = pts[15];
+             
+             newCellId = output->InsertNextCell(aCellType,numFacePts,aNewPts);
+             if(myStoreMapping)
+               myVTK2ObjIds.push_back(cellId);
+             
+             outputCD->CopyData(cd,cellId,newCellId);
+             
+             break;
+           }}
+         }
+       } //switch
       } //if visible
     } //for all cells
   
-  if(MYDEBUG && myStoreMapping){
-    for(int i = 0, iEnd = myVTK2ObjIds.size(); i < iEnd; i++){
-      cout<<myVTK2ObjIds[i]<<", ";
-    }
-    cout<<"\n";
-  }
-
-  // Update ourselves and release memory
-  //
-  //output->SetVerts(Verts);
-  //Verts->Delete();
-  //output->SetLines(Lines);
-  //Lines->Delete();
-  //output->SetPolys(Polys);
-  //Polys->Delete();
-  //output->SetStrips(Strips);
-  //Strips->Delete();
-  
   output->Squeeze();
 
   vtkDebugMacro(<<"Extracted " << input->GetNumberOfPoints() << " points,"
@@ -697,11 +750,72 @@ void VTKViewer_GeometryFilter::UnstructuredGridExecute()
 }
 
 
-void VTKViewer_GeometryFilter::SetInside(int theShowInside){
-  if(myShowInside == theShowInside) return;
+//----------------------------------------------------------------------------
+void
+VTKViewer_GeometryFilter
+::SetInside(int theShowInside)
+{
+  if(myShowInside == theShowInside) 
+    return;
+
   myShowInside = theShowInside;
   this->Modified();
 }
-int VTKViewer_GeometryFilter::GetInside(){
+
+int
+VTKViewer_GeometryFilter
+::GetInside()
+{
   return myShowInside;
 }
+
+
+//----------------------------------------------------------------------------
+void 
+VTKViewer_GeometryFilter
+::SetWireframeMode(int theIsWireframeMode)
+{
+  if(myIsWireframeMode == theIsWireframeMode)
+    return;
+
+  myIsWireframeMode = theIsWireframeMode;
+  this->Modified();
+}
+
+int
+VTKViewer_GeometryFilter
+::GetWireframeMode()
+{
+  return myIsWireframeMode;
+}
+
+
+//----------------------------------------------------------------------------
+void
+VTKViewer_GeometryFilter
+::SetStoreMapping(int theStoreMapping)
+{
+  if(myStoreMapping == theStoreMapping) 
+    return;
+
+  myStoreMapping = theStoreMapping;
+  this->Modified();
+}
+
+int
+VTKViewer_GeometryFilter
+::GetStoreMapping()
+{
+  return myStoreMapping;
+}
+
+
+//----------------------------------------------------------------------------
+vtkIdType VTKViewer_GeometryFilter::GetElemObjId(int theVtkID){
+  if(myVTK2ObjIds.empty() || theVtkID > myVTK2ObjIds.size()) return -1;
+#if defined __GNUC_2__
+  return myVTK2ObjIds[theVtkID];
+#else
+  return myVTK2ObjIds.at(theVtkID);
+#endif
+}
index b302dd04a69b91be98068d7dad45a9a2f8dd8adc..3290640361cbd6e7b476f5986e1878ddac85f5dd 100755 (executable)
@@ -47,6 +47,16 @@ public:
    * \retval myShowInside
    */
   int GetInside();
+  /*! \fn void SetWireframeMode(int theIsWireframeMode)
+   * \brief Sets \a myIsWireframeMode flag. \a myIsWireframeMode is changed, call this->Modified().
+   * \param theIsWireframeMode - used for changing value of \a myIsWireframeMode variable.
+   */
+  void SetWireframeMode(int theIsWireframeMode);
+  /*! \fn int GetWireframeMode()
+   * \brief Return value of \a myIsWireframeMode
+   * \retval myIsWireframeMode
+   */
+  int GetWireframeMode();
   /*! \fn void SetStoreMapping(int theStoreMapping);
    * \brief Sets \a myStoreMapping flag and call this->Modified()
    * \param theStoreMapping - used for changing value of \a myStoreMapping variable.
@@ -56,7 +66,7 @@ public:
    * \brief Return value of \a myStoreMapping
    * \retval myStoreMapping
    */
-  int GetStoreMapping(){ return myStoreMapping;}
+  int GetStoreMapping();
   /*! \fn virtual vtkIdType GetNodeObjId(int theVtkID)
    * \brief Return input value theVtkID
    * \retval theVtkID
@@ -93,6 +103,7 @@ private:
   TVectorId myVTK2ObjIds;
   int       myShowInside;
   int       myStoreMapping;
+  int       myIsWireframeMode;
 };
 
 #endif
diff --git a/src/VTKViewer/resources/VTKViewerM_images.po b/src/VTKViewer/resources/VTKViewerM_images.po
new file mode 100644 (file)
index 0000000..ce2c6c7
--- /dev/null
@@ -0,0 +1,40 @@
+#  VISU VISUGUI : GUI of VISU component
+#
+#  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 
+#
+#
+#
+#  File   : 
+#  Module : 
+
+msgid ""
+msgstr ""
+"Project-Id-Version: PROJECT VERSION\n"
+"POT-Creation-Date: 2002-05-28 10:57:43 AM CEST\n"
+"PO-Revision-Date: 2005-05-10 15:20+0400\n"
+"Last-Translator: FULLNAME <EMAIL@ADDRESS>\n"
+"Content-Type: text/plain; charset=iso-8859-1\n"
+
+
+msgid "ICON_SVTK_SCALING"
+msgstr "view_scaling.png"
+
+msgid "ICON_GRADUATED_AXES"
+msgstr "view_graduated_axes.png"
diff --git a/src/VTKViewer/resources/VTKViewerM_msg_en.po b/src/VTKViewer/resources/VTKViewerM_msg_en.po
new file mode 100644 (file)
index 0000000..7f2c72c
--- /dev/null
@@ -0,0 +1,109 @@
+msgid ""
+msgstr ""
+"Project-Id-Version: PROJECT VERSION\n"
+"POT-Creation-Date: 2002-02-22 16:56:46 CET\n"
+"PO-Revision-Date: 2005-06-27 12:38+0400\n"
+"Last-Translator: FULLNAME <EMAIL@ADDRESS>\n"
+"Content-Type: text/plain; charset=iso-8859-1\n"
+
+
+#: SVTK_NonIsometricDlg.cxx
+
+msgid "SVTK_NonIsometricDlg::MEN_SCALING"
+msgstr "Scaling"
+
+msgid "SVTK_NonIsometricDlg::O&K"
+msgstr ""
+
+msgid "SVTK_NonIsometricDlg::&Apply"
+msgstr ""
+
+msgid "SVTK_NonIsometricDlg::&Cancel"
+msgstr ""
+
+msgid "SVTK_NonIsometricDlg::&Reset"
+msgstr ""
+
+msgid "SVTK_NonIsometricDlg::DLG_TITLE"
+msgstr "Scaling"
+
+msgid "SVTK_NonIsometricDlg::LBL_X"
+msgstr "X :"
+
+msgid "SVTK_NonIsometricDlg::LBL_Y"
+msgstr "Y :"
+
+msgid "SVTK_NonIsometricDlg::LBL_Z"
+msgstr "Z :"
+
+msgid "SVTK_MainWindow::MNU_SVTK_SCALING"
+msgstr "Scaling"
+
+msgid "SVTK_MainWindow::DSC_SVTK_SCALING"
+msgstr "Scaling"
+
+msgid "SVTK_FontWidget::ARIAL"
+msgstr "Arial"
+
+msgid "SVTK_FontWidget::COURIER"
+msgstr "Courier"
+
+msgid "SVTK_FontWidget::TIMES"
+msgstr "Times"
+
+msgid "SVTK_FontWidget::BOLD"
+msgstr "Bold"
+
+msgid "SVTK_FontWidget::ITALIC"
+msgstr "Italic"
+
+msgid "SVTK_FontWidget::SHADOW"
+msgstr "Shadow"
+
+msgid "SVTK_CubeAxesDlg::CAPTION"
+msgstr "Graduated axes"
+
+msgid "SVTK_CubeAxesDlg::X_AXIS"
+msgstr "X axis"
+
+msgid "SVTK_CubeAxesDlg::Y_AXIS"
+msgstr "Y axis"
+
+msgid "SVTK_CubeAxesDlg::Z_AXIS"
+msgstr "Z axis"
+
+msgid "SVTK_CubeAxesDlg::IS_VISIBLE"
+msgstr "Is visible"
+
+msgid "SVTK_AxisWidget::AXIS_NAME"
+msgstr "Axis name"
+
+msgid "SVTK_AxisWidget::IS_VISIBLE"
+msgstr "Is visible"
+
+msgid "SVTK_AxisWidget::NAME"
+msgstr "Name"
+
+msgid "SVTK_AxisWidget::FONT"
+msgstr "Font"
+
+msgid "SVTK_AxisWidget::LABELS"
+msgstr "Labels"
+
+msgid "SVTK_AxisWidget::NUMBER"
+msgstr "Number"
+
+msgid "SVTK_AxisWidget::OFFSET"
+msgstr "Offset"
+
+msgid "SVTK_AxisWidget::TICK_MARKS"
+msgstr "Tick marks"
+
+msgid "SVTK_AxisWidget::LENGTH"
+msgstr "Length"
+
+msgid "SVTK_MainWindow::MNU_SVTK_GRADUATED_AXES"
+msgstr "Graduated axes"
+
+msgid "SVTK_MainWindow::DSC_SVTK_GRADUATED_AXES"
+msgstr "Graduated axes"
\ No newline at end of file
diff --git a/src/VTKViewer/resources/view_graduated_axes.png b/src/VTKViewer/resources/view_graduated_axes.png
new file mode 100755 (executable)
index 0000000..9c9f3d2
Binary files /dev/null and b/src/VTKViewer/resources/view_graduated_axes.png differ
diff --git a/src/VTKViewer/resources/view_scaling.png b/src/VTKViewer/resources/view_scaling.png
new file mode 100644 (file)
index 0000000..5d34436
Binary files /dev/null and b/src/VTKViewer/resources/view_scaling.png differ