Version: 6.3.1

src/SMDS/SMDS_MeshElementIDFactory.cxx

Go to the documentation of this file.
00001 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
00002 //
00003 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
00004 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
00005 //
00006 // This library is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 2.1 of the License.
00010 //
00011 // This library is distributed in the hope that it will be useful,
00012 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00013 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00014 // Lesser General Public License for more details.
00015 //
00016 // You should have received a copy of the GNU Lesser General Public
00017 // License along with this library; if not, write to the Free Software
00018 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
00019 //
00020 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
00021 //
00022 
00023 //  SMESH SMDS : implementaion of Salome mesh data structure
00024 //  File   : SMDS_MeshElementIDFactory.cxx
00025 //  Author : Jean-Michel BOULCOURT
00026 //  Module : SMESH
00027 //
00028 #ifdef _MSC_VER
00029 #pragma warning(disable:4786)
00030 #endif
00031 
00032 #include "SMDS_MeshElementIDFactory.hxx"
00033 #include "SMDS_MeshElement.hxx"
00034 #include "SMDS_Mesh.hxx"
00035 
00036 #include "utilities.h"
00037 
00038 #include "SMDS_UnstructuredGrid.hxx"
00039 #include <vtkCellType.h>
00040 
00041 #include <climits>
00042 
00043 using namespace std;
00044 
00045 //=======================================================================
00046 //function : SMDS_MeshElementIDFactory
00047 //purpose  : 
00048 //=======================================================================
00049 SMDS_MeshElementIDFactory::SMDS_MeshElementIDFactory():
00050   SMDS_MeshNodeIDFactory()
00051 {
00052 //    myIDElements.clear();
00053 //    myVtkIndex.clear();
00054     myVtkCellTypes.clear();
00055     myVtkCellTypes.reserve(SMDSEntity_Last);
00056     myVtkCellTypes[SMDSEntity_Node]            = VTK_VERTEX;
00057     myVtkCellTypes[SMDSEntity_0D]              = VTK_VERTEX;
00058     myVtkCellTypes[SMDSEntity_Edge]            = VTK_LINE;
00059     myVtkCellTypes[SMDSEntity_Quad_Edge]       = VTK_QUADRATIC_EDGE;
00060     myVtkCellTypes[SMDSEntity_Triangle]        = VTK_TRIANGLE;
00061     myVtkCellTypes[SMDSEntity_Quad_Triangle]   = VTK_QUADRATIC_TRIANGLE;
00062     myVtkCellTypes[SMDSEntity_Quadrangle]      = VTK_QUAD;
00063     myVtkCellTypes[SMDSEntity_Quad_Quadrangle] = VTK_QUADRATIC_TRIANGLE;
00064     myVtkCellTypes[SMDSEntity_Polygon]         = VTK_POLYGON;
00065     myVtkCellTypes[SMDSEntity_Quad_Polygon]    = VTK_POLYGON; // -PR- verifer
00066     myVtkCellTypes[SMDSEntity_Tetra]           = VTK_TETRA;
00067     myVtkCellTypes[SMDSEntity_Quad_Tetra]      = VTK_QUADRATIC_TETRA;
00068     myVtkCellTypes[SMDSEntity_Pyramid]         = VTK_PYRAMID;
00069     myVtkCellTypes[SMDSEntity_Quad_Pyramid]    = VTK_CONVEX_POINT_SET;
00070     myVtkCellTypes[SMDSEntity_Hexa]            = VTK_HEXAHEDRON;
00071     myVtkCellTypes[SMDSEntity_Quad_Hexa]       = VTK_QUADRATIC_HEXAHEDRON;
00072     myVtkCellTypes[SMDSEntity_Penta]           = VTK_WEDGE;
00073     myVtkCellTypes[SMDSEntity_Quad_Penta]      = VTK_QUADRATIC_WEDGE;
00074 //#ifdef VTK_HAVE_POLYHEDRON
00075     myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_POLYHEDRON;
00076 //#else
00077 //    myVtkCellTypes[SMDSEntity_Polyhedra]       = VTK_CONVEX_POINT_SET;
00078 //#endif
00079     myVtkCellTypes[SMDSEntity_Quad_Polyhedra]  = VTK_CONVEX_POINT_SET;
00080 }
00081 
00082 int SMDS_MeshElementIDFactory::SetInVtkGrid(SMDS_MeshElement * elem)
00083 {
00084    // --- retrieve nodes ID
00085 
00086   SMDS_MeshCell *cell = dynamic_cast<SMDS_MeshCell*>(elem);
00087   assert(cell);
00088   vector<vtkIdType> nodeIds;
00089   SMDS_ElemIteratorPtr it = elem->nodesIterator();
00090   while(it->more())
00091   {
00092       int nodeId = (static_cast<const SMDS_MeshNode*>(it->next()))->getVtkId();
00093       MESSAGE("   node in cell " << cell->getVtkId() << " : " << nodeId)
00094       nodeIds.push_back(nodeId);
00095   }
00096 
00097   // --- insert cell in vtkUnstructuredGrid
00098 
00099   vtkUnstructuredGrid * grid = myMesh->getGrid();
00100   //int locType = elem->GetType();
00101   int typ = VTK_VERTEX;//GetVtkCellType(locType);
00102   int cellId = grid->InsertNextLinkedCell(typ, nodeIds.size(), &nodeIds[0]);
00103   cell->setVtkId(cellId); 
00104   //MESSAGE("SMDS_MeshElementIDFactory::SetInVtkGrid " << cellId);
00105   return cellId;
00106 }
00107 
00108 //=======================================================================
00109 //function : BindID
00110 //purpose  : 
00111 //=======================================================================
00112 
00113 bool SMDS_MeshElementIDFactory::BindID(int ID, SMDS_MeshElement * elem)
00114 {
00115   MESSAGE("SMDS_MeshElementIDFactory::BindID " << ID);
00116   SetInVtkGrid(elem);
00117   return myMesh->registerElement(ID, elem);
00118 }
00119 
00120 //=======================================================================
00121 //function : MeshElement
00122 //purpose  : 
00123 //=======================================================================
00124 SMDS_MeshElement* SMDS_MeshElementIDFactory::MeshElement(int ID)
00125 {
00126   if ((ID<1) || (ID>=myMesh->myCells.size()))
00127     return NULL;
00128   const SMDS_MeshElement* elem = GetMesh()->FindElement(ID);
00129   return (SMDS_MeshElement*)(elem);
00130 }
00131 
00132 //=======================================================================
00133 //function : GetFreeID
00134 //purpose  : 
00135 //=======================================================================
00136 
00137 int SMDS_MeshElementIDFactory::GetFreeID()
00138 {
00139   int ID;
00140   do {
00141     ID = SMDS_MeshIDFactory::GetFreeID();
00142   } while ( MeshElement( ID ));
00143   return ID;
00144 }
00145 
00146 //=======================================================================
00147 //function : ReleaseID
00148 //purpose  : 
00149 //=======================================================================
00150 void SMDS_MeshElementIDFactory::ReleaseID(int ID, int vtkId)
00151 {
00152   if (ID < 1) // TODO check case ID == O
00153     {
00154       MESSAGE("~~~~~~~~~~~~~~ SMDS_MeshElementIDFactory::ReleaseID ID = " << ID);
00155       return;
00156     }
00157   //MESSAGE("~~~~~~~~~~~~~~ SMDS_MeshElementIDFactory::ReleaseID smdsId vtkId " << ID << " " << vtkId);
00158   if (vtkId >= 0)
00159     {
00160       assert(vtkId < myMesh->myCellIdVtkToSmds.size());
00161       myMesh->myCellIdVtkToSmds[vtkId] = -1;
00162       myMesh->setMyModified();
00163     }
00164   SMDS_MeshIDFactory::ReleaseID(ID);
00165   if (ID == myMax)
00166     myMax = 0;
00167   if (ID == myMin)
00168     myMax = 0;
00169 }
00170 
00171 //=======================================================================
00172 //function : updateMinMax
00173 //purpose  : 
00174 //=======================================================================
00175 
00176 void SMDS_MeshElementIDFactory::updateMinMax() const
00177 {
00178   myMin = INT_MAX;
00179   myMax = 0;
00180   for (int i = 0; i < myMesh->myCells.size(); i++)
00181     {
00182       if (myMesh->myCells[i])
00183         {
00184           int id = myMesh->myCells[i]->GetID();
00185           if (id > myMax)
00186             myMax = id;
00187           if (id < myMin)
00188             myMin = id;
00189         }
00190     }
00191   if (myMin == INT_MAX)
00192     myMin = 0;
00193 }
00194 
00195 //=======================================================================
00196 //function : elementsIterator
00197 //purpose  : Return an iterator on elements of the factory
00198 //=======================================================================
00199 
00200 SMDS_ElemIteratorPtr SMDS_MeshElementIDFactory::elementsIterator() const
00201 {
00202     return myMesh->elementsIterator(SMDSAbs_All);
00203 }
00204 
00205 void SMDS_MeshElementIDFactory::Clear()
00206 {
00207   //myMesh->myCellIdSmdsToVtk.clear();
00208   myMesh->myCellIdVtkToSmds.clear();
00209   myMin = myMax = 0;
00210   SMDS_MeshIDFactory::Clear();
00211 }
00212 
00213 int SMDS_MeshElementIDFactory::GetVtkCellType(int SMDSType)
00214 {
00215     assert((SMDSType >=0) && (SMDSType< SMDSEntity_Last));
00216     return myVtkCellTypes[SMDSType];
00217 }
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