]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_Convertor.hxx
Salome HOME
sources v1.2c
[modules/visu.git] / src / CONVERTOR / VISU_Convertor.hxx
1 //  VISU OBJECT : interactive object for VISU entities implementation
2 //
3 //  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
5 // 
6 //  This library is free software; you can redistribute it and/or 
7 //  modify it under the terms of the GNU Lesser General Public 
8 //  License as published by the Free Software Foundation; either 
9 //  version 2.1 of the License. 
10 // 
11 //  This library is distributed in the hope that it will be useful, 
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of 
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
14 //  Lesser General Public License for more details. 
15 // 
16 //  You should have received a copy of the GNU Lesser General Public 
17 //  License along with this library; if not, write to the Free Software 
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
19 // 
20 //  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
21 //
22 //
23 //  File   : VISU_Convertor.hxx
24 //  Author : Alexey PETROV
25 //  Module : VISU
26
27 #ifndef VISU_Convertor_HeaderFile
28 #define VISU_Convertor_HeaderFile
29
30 #include <stdio.h>
31
32 #include <algorithm>
33 #include <iterator>
34 #include <list>
35 #include <map>
36 #include <memory>
37 #include <numeric>
38 #include <set>
39 #include <utility>
40 #include <vector>
41 #include <string>
42 #include <stdexcept>
43
44 #ifndef MESSAGE
45 #define MESSAGE(msg) std::cout<<msg<<endl
46 #endif
47
48 #include <vtkUnstructuredGrid.h>
49 #include <vtkSystemIncludes.h>
50 #include <vtkPoints.h>
51
52 namespace VISU{
53   class TVtkCellInfo{
54   public:
55     TVtkCellInfo() {}
56     TVtkCellInfo(const char* theName, vtkIdType theSize) : 
57       myName(theName), mySize(theSize) {}
58     const char* myName;
59     vtkIdType mySize;
60   };
61   typedef std::map<vtkIdType,TVtkCellInfo> TVtkCellInfoMap;
62   extern TVtkCellInfoMap aVtkCellInfoMap;
63
64   enum  TEntity {NODE_ENTITY, EDGE_ENTITY, FACE_ENTITY, CELL_ENTITY};
65   template <class _Tp> class vtk_ptr {
66   private:
67     _Tp* _M_ptr;
68   public:
69     typedef _Tp element_type;
70     explicit vtk_ptr(_Tp* __p = 0) : _M_ptr(__p) {}
71     vtk_ptr(const vtk_ptr& __a) : _M_ptr(0) {}
72     vtk_ptr(vtk_ptr& __a) : _M_ptr(__a.release()) {}
73     vtk_ptr& operator=(const vtk_ptr& __a) { return *this;}
74     vtk_ptr& operator=(vtk_ptr& __a) {
75       if (&__a != this) {
76         if(_M_ptr) _M_ptr->Delete(); 
77         _M_ptr = __a.release();
78       }
79       return *this;
80     }
81     ~vtk_ptr() { 
82       if(_M_ptr) _M_ptr->Delete(); 
83       _M_ptr = 0;
84     }
85     _Tp& operator*() const { return *_M_ptr;}
86     _Tp* operator->() const { return _M_ptr;}
87     _Tp* get() const { return _M_ptr;}
88     _Tp* release() {
89       _Tp* __tmp = _M_ptr;
90       _M_ptr = 0;
91       return __tmp;
92     }
93     void reset(_Tp* __p = 0) {
94       if(_M_ptr) _M_ptr->Delete();
95       _M_ptr = __p;
96     }
97   };
98   typedef vtk_ptr<vtkUnstructuredGrid> TVTKSource;
99   typedef vtk_ptr<vtkPoints> TVTKPoints;
100
101   typedef std::set<std::string> TBindGroups;
102
103   struct TFamily{
104     TVTKSource myStorage;
105     vtkIdType myId;
106     std::string myName;
107     TEntity myEntity;
108     TBindGroups myGroups;
109     vtkIdType myNbCells, myCellsSize;
110     typedef std::set<vtkIdType> TSubMeshOnCellType;
111     typedef std::map<vtkIdType,TSubMeshOnCellType> TSubMesh;
112     TSubMesh mySubMesh;
113     TFamily() : myNbCells(0), myCellsSize(0) {}
114   };
115   typedef std::map<std::string,TFamily> TFamilyMap;
116
117   struct TField{
118     vtkIdType myId;
119     std::string myName;
120     TEntity myEntity;
121     std::string myMeshName;
122     vtkIdType myNbComp, myNbValField, myDataSize;
123     typedef std::vector<float> TValForCellsWithType;
124     typedef std::map<vtkIdType,TValForCellsWithType> TValForCells;
125     typedef std::pair<double,std::string> TTime;
126     typedef std::vector<std::string> TCompNames;
127     typedef std::vector<std::string> TUnitNames;
128     struct TValForTime{
129       TVTKSource myStorage;
130       vtkIdType myId;
131       std::string myMeshName;
132       TEntity myEntity;
133       std::string myFieldName;
134       vtkIdType myNbComp;
135       TTime myTime;
136       TValForCells myValForCells;
137       TValForTime() : myNbComp(0) {}
138     };
139     typedef std::map<vtkIdType,TValForTime> TValField;
140     TValField myValField;
141     TCompNames myCompNames;
142     TUnitNames myUnitNames;
143     TField() : myNbComp(0), myNbValField(0), myDataSize(0) {}
144     void ShallowCopy(const TField& aField);
145   };
146   typedef std::map<std::string,TField> TFieldMap;
147
148   struct TMeshOnEntity{
149     TVTKSource myStorage;
150     std::string myMeshName;
151     TEntity myEntity;
152     vtkIdType myNbCells, myCellsSize;
153     typedef std::vector<vtkIdType> TConnect;
154     typedef std::vector<TConnect> TConnForCellType;
155     typedef std::map<vtkIdType,TConnForCellType> TCellsConn;
156     TCellsConn myCellsConn;
157     TFamilyMap myFamilyMap;
158     TFieldMap myFieldMap;
159     TMeshOnEntity() : myNbCells(0), myCellsSize(0) {}
160     std::pair<vtkIdType,vtkIdType> GetCellsDims(const std::string& theFamilyName = "") const
161       throw(std::runtime_error&); 
162   };
163   typedef std::map<TEntity,TMeshOnEntity> TMeshOnEntityMap;
164   const TFamily* GetFamily(const VISU::TMeshOnEntity& theMeshOnEntity, 
165                      const std::string& theFamilyName)
166     throw(std::runtime_error&); 
167   TFamily* GetFamily(VISU::TMeshOnEntity& theMeshOnEntity, 
168                      const std::string& theFamilyName)
169     throw(std::runtime_error&); 
170
171   typedef std::pair<std::string,TEntity> TFamilyAndEntity;
172   typedef std::set<TFamilyAndEntity> TFamilyAndEntitySet;
173   struct TGroup{
174     TVTKSource myStorage;
175     std::string myName;
176     std::string myMeshName;
177     vtkIdType myNbCells, myCellsSize;
178     TGroup() : myNbCells(0), myCellsSize(0) {}
179     TFamilyAndEntitySet myFamilyAndEntitySet;
180   };
181   typedef std::map<std::string,TGroup> TGroupMap;
182
183   struct TMesh{
184     TVTKPoints myPoints;
185     vtkIdType myDim, myNbPoints;
186     std::string myName;
187     typedef float TCoord;
188     typedef std::vector<TCoord> TPointsCoord;
189     TPointsCoord myPointsCoord;
190     TMeshOnEntityMap myMeshOnEntityMap;
191     TGroupMap myGroupMap;
192     TMesh() : myDim(0), myNbPoints(0) {}
193     const TField* GetField(const std::string& theFieldName) const;
194   };
195   typedef std::map<std::string,TMesh> TMeshMap;
196   void WriteToFile(vtkUnstructuredGrid* theDataSet, const std::string& theFileName);
197 };
198
199 class VISU_Convertor{
200 protected:
201   std::string myName;
202   VISU::TMeshMap myMeshMap;
203   int myIsDone;
204 public:
205   virtual ~VISU_Convertor(){};
206   virtual const std::string& GetName() { return myName;}
207   virtual int IsDone() const { return myIsDone;}
208   typedef vtkUnstructuredGrid TOutput;
209
210   virtual VISU_Convertor* Build() throw (std::runtime_error&) = 0;
211   virtual const VISU::TMeshMap& GetMeshMap() throw(std::runtime_error&);
212
213
214   virtual TOutput* GetMeshOnEntity(const std::string& theMeshName, 
215                                    const VISU::TEntity& theEntity,
216                                    const std::string& theFamilyName = "")
217     throw(std::runtime_error&) = 0;
218   virtual vtkIdType GetMeshOnEntitySize(const std::string& theMeshName, 
219                                         const VISU::TEntity& theEntity,
220                                         const std::string& theFamilyName = "")
221     throw (std::runtime_error&) = 0;
222   
223   
224   virtual TOutput* GetMeshOnGroup(const std::string& theMeshName, 
225                                   const std::string& theGroupName)
226     throw(std::runtime_error&) = 0;
227   virtual vtkIdType GetMeshOnGroupSize(const std::string& theMeshName, 
228                                        const std::string& theGroupName)
229     throw(std::runtime_error&) = 0;
230   
231
232   virtual TOutput* GetTimeStampOnMesh(const std::string& theMeshName, 
233                                       const VISU::TEntity& theEntity,
234                                       const std::string& theFieldName,
235                                       int theStampsNum)
236     throw(std::runtime_error&) = 0;
237   virtual vtkIdType GetTimeStampSize(const std::string& theMeshName, 
238                                      const VISU::TEntity& theEntity,
239                                      const std::string& theFieldName,
240                                      int theStampsNum)
241     throw(std::runtime_error&) = 0;
242   virtual vtkIdType GetFieldOnMeshSize(const std::string& theMeshName, 
243                                        const VISU::TEntity& theEntity,
244                                        const std::string& theFieldName)
245     throw(std::runtime_error&) = 0;
246   virtual const VISU::TField& GetField(const std::string& theMeshName, 
247                                        VISU::TEntity theEntity, 
248                                        const std::string& theFieldName)
249     throw(std::runtime_error&) = 0;
250   virtual const VISU::TField::TValForTime& GetTimeStamp(const std::string& theMeshName, 
251                                                         const VISU::TEntity& theEntity,
252                                                         const std::string& theFieldName,
253                                                         int theStampsNum)
254     throw(std::runtime_error&) = 0;
255
256   static std::string GenerateName(const VISU::TField::TTime& aTime);
257   static std::string GenerateName(const std::string& theName, unsigned int theTimeId);
258   static void WriteToFile(vtkUnstructuredGrid* theDataSet, const std::string& theFileName);
259 };
260
261 extern "C"{
262   VISU_Convertor* CreateMedConvertor(const std::string& theFileName) throw(std::runtime_error&);
263   VISU_Convertor* CreateDatConvertor(const std::string& theFileName) throw(std::runtime_error&);
264   VISU_Convertor* CreateConvertor(const std::string& theFileName) throw(std::runtime_error&);
265 };
266
267 #endif