Salome HOME
Merge from V6_main (04/10/2012)
[modules/med.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 // Author : Anthony Geay (CEA/DEN)
20 #ifndef __VTKNORMALIZEDUNSTRUCTUREDMESH_TXX__
21 #define __VTKNORMALIZEDUNSTRUCTUREDMESH_TXX__
22
23 #include "VTKNormalizedUnstructuredMesh.hxx"
24
25 #include "vtkUnstructuredGrid.h"
26 #include "vtkCellArray.h"
27 #include "vtkPoints.h"
28
29 template<int MESHDIM>
30 VTKNormalizedUnstructuredMesh<MESHDIM>::VTKNormalizedUnstructuredMesh(vtkUnstructuredGrid *mesh):_mesh_in_vtk_mode(mesh),
31                                                                                                  _tmp_index_array(0)
32 {
33 }
34
35 template<int MESHDIM>
36 VTKNormalizedUnstructuredMesh<MESHDIM>::~VTKNormalizedUnstructuredMesh()
37 {
38   _mesh_in_vtk_mode->Delete();
39   releaseTempArrays();
40 }
41
42 template<int MESHDIM>
43 void VTKNormalizedUnstructuredMesh<MESHDIM>::getBoundingBox(double *boundingBox) const
44 {
45   double tmp[6];
46   _mesh_in_vtk_mode->GetBounds(tmp);
47   for(unsigned i=0;i<3;i++)
48     {
49       boundingBox[i]=tmp[2*i];
50       boundingBox[3+i]=tmp[2*i+1];
51     }
52 }
53
54 template<int MESHDIM>
55 NormalizedCellType VTKNormalizedUnstructuredMesh<MESHDIM>::getTypeOfElement(vtkIdType eltId) const
56 {
57   int cellType=_mesh_in_vtk_mode->GetCellType(eltId);
58   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
59                    0,(int)NORM_PYRA5,0,0,0,(int)NORM_TRI6,(int)NORM_QUAD8,};
60 }
61
62 template<int MESHDIM>
63 unsigned long VTKNormalizedUnstructuredMesh<MESHDIM>::getNumberOfElements() const
64 {
65   return _mesh_in_vtk_mode->GetNumberOfCells();
66 }
67
68 template<int MESHDIM>
69 unsigned long VTKNormalizedUnstructuredMesh<MESHDIM>::getNumberOfNodes() const
70 {
71   return _mesh_in_vtk_mode->GetNumberOfPoints();
72 }
73
74 template<int MESHDIM>
75 const vtkIdType *VTKNormalizedUnstructuredMesh<MESHDIM>::getConnectivityPtr() const
76 {
77   vtkIdType *ret=_mesh_in_vtk_mode->GetCells()->GetPointer();
78   if(_tmp_index_array)
79     return ret;
80   else
81     {
82       putinMEDFormat();
83       return ret;
84     }
85 }
86
87 template<int MESHDIM>
88 const double *VTKNormalizedUnstructuredMesh<MESHDIM>::getCoordinatesPtr() const
89 {
90   return (const double *)_mesh_in_vtk_mode->GetPoints()->GetVoidPointer(0);
91 }
92
93 template<int MESHDIM>
94 const vtkIdType *VTKNormalizedUnstructuredMesh<MESHDIM>::getConnectivityIndexPtr() const
95 {
96   if(_tmp_index_array)
97     return _tmp_index_array;
98   else
99     {
100       putinMEDFormat();
101       return _tmp_index_array;
102     }
103 }
104
105 template<int MESHDIM>
106 void VTKNormalizedUnstructuredMesh<MESHDIM>::putinMEDFormat() const
107 {
108   long nbOfElem=getNumberOfElements();
109   _tmp_index_array=new vtkIdType[nbOfElem+1];
110   _tmp_index_array[0]=0;
111   vtkIdType *coarseConn=_mesh_in_vtk_mode->GetCells()->GetPointer();
112   long ptInCC=0;
113   vtkIdType *finalConn=coarseConn;
114   for(long i=0;i<nbOfElem;i++)
115     {
116       vtkIdType cellLgth=coarseConn[ptInCC];
117       for(vtkIdType j=0;j<cellLgth;j++)
118         *finalConn++=coarseConn[ptInCC+j+1];
119       _tmp_index_array[i+1]=_tmp_index_array[i]+cellLgth;
120       ptInCC+=cellLgth+1;
121     }
122   int gh=0;
123   gh++;
124 }
125
126 template<int MESHDIM>
127 void VTKNormalizedUnstructuredMesh<MESHDIM>::releaseTempArrays()
128 {
129   delete [] _tmp_index_array;
130   _tmp_index_array=0;
131 }
132
133 #endif