Version: 6.3.1

src/SMDS/SMDS_UnstructuredGrid.hxx

Go to the documentation of this file.
00001 // Copyright (C) 2010-2011  CEA/DEN, EDF R&D, OPEN CASCADE
00002 //
00003 // This library is free software; you can redistribute it and/or
00004 // modify it under the terms of the GNU Lesser General Public
00005 // License as published by the Free Software Foundation; either
00006 // version 2.1 of the License.
00007 //
00008 // This library is distributed in the hope that it will be useful,
00009 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00010 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00011 // Lesser General Public License for more details.
00012 //
00013 // You should have received a copy of the GNU Lesser General Public
00014 // License along with this library; if not, write to the Free Software
00015 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00016 //
00017 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00018 //
00019 
00020 // File:    SMDS_UnstructuredGrid.hxx
00021 // Author:  prascle
00022 // Created: September 16, 2009, 10:28 PM
00023 
00024 #ifndef _SMDS_UNSTRUCTUREDGRID_HXX
00025 #define _SMDS_UNSTRUCTUREDGRID_HXX
00026 
00027 #include <vtkUnstructuredGrid.h>
00028 #include <vtkCellLinks.h>
00029 #include "chrono.hxx"
00030 
00031 #include <vector>
00032 #include <set>
00033 #include <map>
00034 
00035 //#define VTK_HAVE_POLYHEDRON
00036 //#ifdef VTK_HAVE_POLYHEDRON
00037 #define VTK_MAXTYPE VTK_POLYHEDRON
00038 //#else
00039 //  #define VTK_MAXTYPE VTK_QUADRATIC_PYRAMID
00040 //#endif
00041 
00042 #define NBMAXNEIGHBORS 100
00043 
00044 // allow very huge polyhedrons in tests
00045 #define NBMAXNODESINCELL 5000
00046 
00047 class SMDS_Downward;
00048 class SMDS_Mesh;
00049 class SMDS_MeshVolume;
00050 
00051 class SMDS_CellLinks: public vtkCellLinks
00052 {
00053 public:
00054   vtkCellLinks::Link* ResizeL(vtkIdType sz);
00055   vtkIdType GetLinksSize();
00056   static SMDS_CellLinks* New();
00057 protected:
00058   SMDS_CellLinks();
00059   ~SMDS_CellLinks();
00060 };
00061 
00062 class SMDS_UnstructuredGrid: public vtkUnstructuredGrid
00063 {
00064 public:
00065   void setSMDS_mesh(SMDS_Mesh *mesh);
00066   void compactGrid(std::vector<int>& idNodesOldToNew, int newNodeSize, std::vector<int>& idCellsOldToNew,
00067                    int newCellSize);
00068 
00069   virtual unsigned long GetMTime();
00070   virtual void Update();
00071   virtual void UpdateInformation();
00072   virtual vtkPoints *GetPoints();
00073 
00074   //#ifdef VTK_HAVE_POLYHEDRON
00075   int InsertNextLinkedCell(int type, int npts, vtkIdType *pts);
00076   //#endif
00077 
00078   int CellIdToDownId(int vtkCellId);
00079   void setCellIdToDownId(int vtkCellId, int downId);
00080   void CleanDownwardConnectivity();
00081   void BuildDownwardConnectivity(bool withEdges);
00082   int GetNeighbors(int* neighborsVtkIds, int* downIds, unsigned char* downTypes, int vtkId);
00083   int GetParentVolumes(int* volVtkIds, int vtkId);
00084   int GetParentVolumes(int* volVtkIds, int downId, unsigned char downType);
00085   void GetNodeIds(std::set<int>& nodeSet, int downId, unsigned char downType);
00086   void ModifyCellNodes(int vtkVolId, std::map<int, int> localClonedNodeIds);
00087   int getOrderedNodesOfFace(int vtkVolId, std::vector<vtkIdType>& orderedNodes);
00088   void BuildLinks();
00089   SMDS_MeshVolume* extrudeVolumeFromFace(int vtkVolId, int domain1, int domain2, std::set<int>& originalNodes,
00090                                          std::map<int, std::map<int, int> >& nodeDomains,
00091                                          std::map<int, std::map<long,int> >& nodeQuadDomains);
00092   vtkCellLinks* GetLinks()
00093   {
00094     return Links;
00095   }
00096   SMDS_Downward* getDownArray(unsigned char vtkType)
00097   {
00098     return _downArray[vtkType];
00099   }
00100   static SMDS_UnstructuredGrid* New();
00101   SMDS_Mesh *_mesh;
00102 protected:
00103   SMDS_UnstructuredGrid();
00104   ~SMDS_UnstructuredGrid();
00105   void copyNodes(vtkPoints *newPoints, std::vector<int>& idNodesOldToNew, int& alreadyCopied, int start, int end);
00106   void copyBloc(vtkUnsignedCharArray *newTypes, std::vector<int>& idCellsOldToNew, std::vector<int>& idNodesOldToNew,
00107                 vtkCellArray* newConnectivity, vtkIdTypeArray* newLocations, vtkIdType* pointsCell, int& alreadyCopied,
00108                 int start, int end);
00109 
00110   std::vector<int> _cellIdToDownId; 
00111   std::vector<unsigned char> _downTypes;
00112   std::vector<SMDS_Downward*> _downArray;
00113 };
00114 
00115 #endif  /* _SMDS_UNSTRUCTUREDGRID_HXX */
00116 
Copyright © 2007-2011 CEA/DEN, EDF R&D, OPEN CASCADE
Copyright © 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN, CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS