]> SALOME platform Git repositories - modules/visu.git/blob - src/VISU_I/VISU_Convertor.hxx
Salome HOME
NRI : Remove dependence with KERNEL.
[modules/visu.git] / src / VISU_I / VISU_Convertor.hxx
1 // File:        VISU_Convertor.hxx
2 // Created:     Fri Jan 10 11:37:50 2003
3 // Author:      Alexey PETROV
4 //              <apo@ivanox.nnov.matra-dtv.fr>
5
6 #ifndef VISU_Convertor_HeaderFile
7 #define VISU_Convertor_HeaderFile
8
9 #include <stdio.h>
10
11 #include <algorithm>
12 #include <iterator>
13 #include <list>
14 #include <map>
15 #include <memory>
16 #include <numeric>
17 #include <set>
18 #include <utility>
19 #include <vector>
20 #include <string>
21 #include <stdexcept>
22
23 #include <qstring.h>
24
25 #include "utilities.h"
26 #ifndef MESSAGE
27 #define MESSAGE(msg) cout<<msg<<endl
28 #endif
29
30
31 //class vtkUnstructuredGridReader;
32 #include "vtkUnstructuredGridReader.h"
33
34 namespace VISU{
35   class TVtkCellInfo{
36   public:
37     TVtkCellInfo() {}
38     TVtkCellInfo(const char* theName, int theSize) : 
39       myName(theName), mySize(theSize) {}
40     const char* myName;
41     int mySize;
42   };
43   typedef std::map<int,TVtkCellInfo> TVtkCellInfoMap;
44   extern TVtkCellInfoMap aVtkCellInfoMap;
45
46   enum  TEntity {NODE_ENTITY, EDGE_ENTITY, FACE_ENTITY, CELL_ENTITY};
47   template <class _Tp> class vtk_ptr {
48   private:
49     _Tp* _M_ptr;
50   public:
51     typedef _Tp element_type;
52     explicit vtk_ptr(_Tp* __p = 0) : _M_ptr(__p) {}
53     vtk_ptr(const vtk_ptr& __a) : _M_ptr(0) {}
54     vtk_ptr(vtk_ptr& __a) : _M_ptr(__a.release()) {}
55     vtk_ptr& operator=(const vtk_ptr& __a) { return *this;}
56     vtk_ptr& operator=(vtk_ptr& __a) {
57       if (&__a != this) {
58         if(_M_ptr) _M_ptr->Delete(); 
59         _M_ptr = __a.release();
60       }
61       return *this;
62     }
63     ~vtk_ptr() { 
64       if(_M_ptr) _M_ptr->Delete(); 
65       _M_ptr = 0;
66     }
67     _Tp& operator*() const { return *_M_ptr;}
68     _Tp* operator->() const { return _M_ptr;}
69     _Tp* get() const { return _M_ptr;}
70     _Tp* release() {
71       _Tp* __tmp = _M_ptr;
72       _M_ptr = 0;
73       return __tmp;
74     }
75     void reset(_Tp* __p = 0) {
76       if(_M_ptr) _M_ptr->Delete();
77       _M_ptr = __p;
78     }
79   };
80   typedef vtk_ptr<vtkUnstructuredGridReader> TVTKReader;
81
82   typedef std::set<std::string> TBindGroups;
83
84   struct TFamily{
85     TVTKReader myStorage;
86     int myId;
87     string myName;
88     TEntity myEntity;
89     TBindGroups myGroups;
90     typedef std::set<int> TSubMeshOnCellType;
91     typedef std::map<int,TSubMeshOnCellType> TSubMesh;
92     TSubMesh mySubMesh;
93   };
94   typedef std::map<std::string,TFamily> TFamilyMap;
95
96   struct TField{
97     int myId;
98     string myName;
99     TEntity myEntity;
100     string myMeshName;
101     int myNbComp;
102     typedef std::vector<float> TValForCellsWithType;
103     typedef std::map<int,TValForCellsWithType> TValForCells;
104     typedef std::pair<double,std::string> TTime;
105     struct TValForTime{
106       TVTKReader myStorage;
107       int myId;
108       TTime myTime;
109       TValForCells myValForCells;
110     };
111     typedef std::map<int,TValForTime> TValField;
112     TValField myValField;
113     void ShallowCopy(const TField& aField);
114   };
115   typedef map<string,TField> TFieldMap;
116
117   struct TMeshOnEntity{
118     TVTKReader myStorage;
119     string myMeshName;
120     TEntity myEntity;
121     typedef vector<int> TConnect;
122     typedef vector<TConnect> TConnForCellType;
123     typedef map<int,TConnForCellType> TCellsConn;
124     TCellsConn myCellsConn;
125     TFamilyMap myFamilyMap;
126     TFieldMap myFieldMap;
127     pair<int,int> GetCellsDims(const string& theFamilyName = "") const
128       throw(std::runtime_error&); 
129   };
130   typedef std::map<TEntity,TMeshOnEntity> TMeshOnEntityMap;
131   const TFamily* GetFamily(const VISU::TMeshOnEntity& theMeshOnEntity, 
132                      const string& theFamilyName)
133     throw(std::runtime_error&); 
134   TFamily* GetFamily(VISU::TMeshOnEntity& theMeshOnEntity, 
135                      const string& theFamilyName)
136     throw(std::runtime_error&); 
137
138   typedef std::pair<std::string,TEntity> TFamilyAndEntity;
139   typedef std::set<TFamilyAndEntity> TFamilyAndEntitySet;
140   struct TGroup{
141     TVTKReader myStorage;
142     string myName;
143     string myMeshName;
144     TFamilyAndEntitySet myFamilyAndEntitySet;
145   };
146   typedef std::map<std::string,TGroup> TGroupMap;
147
148   struct TMesh{
149     int myDim;
150     string myName;
151     typedef vector<float> TPointsCoord;
152     TPointsCoord myPointsCoord;
153     TMeshOnEntityMap myMeshOnEntityMap;
154     TGroupMap myGroupMap;
155     void CreateMeshOnNodes();
156     const TField* GetField(const string& theFieldName) const;
157   };
158   typedef std::map<std::string,TMesh> TMeshMap;
159 };
160
161 class VISU_Convertor{
162 protected:
163   std::string myName;
164   VISU::TMeshMap myMeshMap;
165 public:
166   virtual ~VISU_Convertor(){};
167   typedef vtkUnstructuredGridReader OutputType;
168   const char* GetName() { return myName.c_str();}
169
170   virtual OutputType* GetMeshOnEntity(const string& theMeshName, 
171                                       const VISU::TEntity& theEntity,
172                                       const string& theFamilyName = "")
173     throw(std::runtime_error&) = 0;
174
175   virtual OutputType* GetMeshOnGroup(const string& theMeshName, 
176                                      const string& theGroupName)
177     throw(std::runtime_error&) = 0;
178
179   virtual OutputType* GetFieldOnMesh(const string& theMeshName, 
180                                      const VISU::TEntity& theEntity,
181                                      const string& theFieldName,
182                                      int theStampsNum)
183     throw(std::runtime_error&) = 0;
184
185   virtual int ToMedFile(const string& theFileName) throw(std::runtime_error&) = 0;
186   virtual int ToDatFile(const string& theFileName) throw(std::runtime_error&) = 0;
187   virtual int ToFile(const string& theFileName) throw(std::runtime_error&) = 0;
188
189   const VISU::TMeshMap* GetMeshMap() const { return &myMeshMap;}
190   const VISU::TMesh* GetMesh(const string& theMeshName = "") const;
191   const VISU::TField* GetField(const string& theMeshName, VISU::TEntity theEntity, 
192                                const string& theFieldName) const;
193   static string GenerateName(const VISU::TField::TTime& aTime);
194   static string GenerateName(const string& theName, unsigned int theTimeId);
195 };
196
197 template<class T> std::string dtos(const string& fmt, T val){
198   static QString aString;
199   aString.sprintf(fmt.c_str(),val);
200   return aString.latin1();
201 }
202 extern "C"{
203   VISU_Convertor* CreateMedConvertor(const std::string& theFileName) throw(std::runtime_error&);
204   VISU_Convertor* CreateDatConvertor(const std::string& theFileName) throw(std::runtime_error&);
205   VISU_Convertor* CreateConvertor(const std::string& theFileName) throw(std::runtime_error&);
206 };
207
208 #endif