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