Salome HOME
Merge from V6_main_20120808 08Aug12
[tools/medcoupling.git] / src / INTERP_KERNEL / VTKNormalizedUnstructuredMesh.txx
1 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19 #ifndef __VTKNORMALIZEDUNSTRUCTUREDMESH_TXX__
20 #define __VTKNORMALIZEDUNSTRUCTUREDMESH_TXX__
21
22 #include "VTKNormalizedUnstructuredMesh.hxx"
23
24 #include "vtkUnstructuredGrid.h"
25 #include "vtkCellArray.h"
26 #include "vtkPoints.h"
27
28 template<int MESHDIM>
29 VTKNormalizedUnstructuredMesh<MESHDIM>::VTKNormalizedUnstructuredMesh(vtkUnstructuredGrid *mesh):_mesh_in_vtk_mode(mesh),
30                                                                                                  _tmp_index_array(0)
31 {
32 }
33
34 template<int MESHDIM>
35 VTKNormalizedUnstructuredMesh<MESHDIM>::~VTKNormalizedUnstructuredMesh()
36 {
37   _mesh_in_vtk_mode->Delete();
38   releaseTempArrays();
39 }
40
41 template<int MESHDIM>
42 void VTKNormalizedUnstructuredMesh<MESHDIM>::getBoundingBox(double *boundingBox) const
43 {
44   double tmp[6];
45   _mesh_in_vtk_mode->GetBounds(tmp);
46   for(unsigned i=0;i<3;i++)
47     {
48       boundingBox[i]=tmp[2*i];
49       boundingBox[3+i]=tmp[2*i+1];
50     }
51 }
52
53 template<int MESHDIM>
54 NormalizedCellType VTKNormalizedUnstructuredMesh<MESHDIM>::getTypeOfElement(vtkIdType eltId) const
55 {
56   int cellType=_mesh_in_vtk_mode->GetCellType(eltId);
57   int convTab[30]={0,0,0,0,0,(int)NORM_TRI3,0,(int)NORM_POLYGON,0,(int)NORM_QUAD4,(int)NORM_TETRA4,0,(int)NORM_HEXA8
58                    0,(int)NORM_PYRA5,0,0,0,(int)NORM_TRI6,(int)NORM_QUAD8,};
59 }
60
61 template<int MESHDIM>
62 unsigned long VTKNormalizedUnstructuredMesh<MESHDIM>::getNumberOfElements() const
63 {
64   return _mesh_in_vtk_mode->GetNumberOfCells();
65 }
66
67 template<int MESHDIM>
68 unsigned long VTKNormalizedUnstructuredMesh<MESHDIM>::getNumberOfNodes() const
69 {
70   return _mesh_in_vtk_mode->GetNumberOfPoints();
71 }
72
73 template<int MESHDIM>
74 const vtkIdType *VTKNormalizedUnstructuredMesh<MESHDIM>::getConnectivityPtr() const
75 {
76   vtkIdType *ret=_mesh_in_vtk_mode->GetCells()->GetPointer();
77   if(_tmp_index_array)
78     return ret;
79   else
80     {
81       putinMEDFormat();
82       return ret;
83     }
84 }
85
86 template<int MESHDIM>
87 const double *VTKNormalizedUnstructuredMesh<MESHDIM>::getCoordinatesPtr() const
88 {
89   return (const double *)_mesh_in_vtk_mode->GetPoints()->GetVoidPointer(0);
90 }
91
92 template<int MESHDIM>
93 const vtkIdType *VTKNormalizedUnstructuredMesh<MESHDIM>::getConnectivityIndexPtr() const
94 {
95   if(_tmp_index_array)
96     return _tmp_index_array;
97   else
98     {
99       putinMEDFormat();
100       return _tmp_index_array;
101     }
102 }
103
104 template<int MESHDIM>
105 void VTKNormalizedUnstructuredMesh<MESHDIM>::putinMEDFormat() const
106 {
107   long nbOfElem=getNumberOfElements();
108   _tmp_index_array=new vtkIdType[nbOfElem+1];
109   _tmp_index_array[0]=0;
110   vtkIdType *coarseConn=_mesh_in_vtk_mode->GetCells()->GetPointer();
111   long ptInCC=0;
112   vtkIdType *finalConn=coarseConn;
113   for(long i=0;i<nbOfElem;i++)
114     {
115       vtkIdType cellLgth=coarseConn[ptInCC];
116       for(vtkIdType j=0;j<cellLgth;j++)
117         *finalConn++=coarseConn[ptInCC+j+1];
118       _tmp_index_array[i+1]=_tmp_index_array[i]+cellLgth;
119       ptInCC+=cellLgth+1;
120     }
121   int gh=0;
122   gh++;
123 }
124
125 template<int MESHDIM>
126 void VTKNormalizedUnstructuredMesh<MESHDIM>::releaseTempArrays()
127 {
128   delete [] _tmp_index_array;
129   _tmp_index_array=0;
130 }
131
132 #endif