]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_Convertor.hxx
Salome HOME
DCQ : Merge with Ecole Ete a6.
[modules/visu.git] / src / CONVERTOR / VISU_Convertor.hxx
1 //  VISU CONVERTOR :
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 <numeric>
37 #include <set>
38 #include <utility>
39 #include <vector>
40 #include <string>
41 #include <stdexcept>
42
43 #include <vtkPoints.h>
44 #include <vtkSystemIncludes.h>
45 #include <vtkUnstructuredGrid.h>
46 #include <vtkMergeDataObjectFilter.h>
47 #include <vtkFieldDataToAttributeDataFilter.h>
48
49 #include "VISU_ExtractUnstructuredGrid.hxx"
50
51 #include <vtkSmartPointer.h>
52
53 namespace VISU{
54   enum  TEntity {NODE_ENTITY, EDGE_ENTITY, FACE_ENTITY, CELL_ENTITY};
55
56   typedef vtkSmartPointer<vtkPoints> TVTKPoints;
57   typedef vtkSmartPointer<vtkUnstructuredGrid> TVTKSource;
58   typedef vtkSmartPointer<vtkMergeDataObjectFilter> TVTKMergetFilter;
59   typedef vtkSmartPointer<VISU_ExtractUnstructuredGrid> TVTKExtractFilter;
60   typedef vtkSmartPointer<vtkFieldDataToAttributeDataFilter> TVTKAttribyteFilter;
61   typedef std::set<std::string> TBindGroups;
62
63   struct TFamily{
64     TVTKSource myStorage;
65     vtkIdType myId;
66     std::string myName;
67     TEntity myEntity;
68     TBindGroups myGroups;
69     vtkIdType myNbCells, myCellsSize;
70     typedef std::set<vtkIdType> TSubMeshOnCellType;
71     typedef std::map<vtkIdType,TSubMeshOnCellType> TSubMesh;
72     TSubMesh mySubMesh;
73     TFamily() : myNbCells(0), myCellsSize(0) {}
74   };
75   typedef std::map<std::string,TFamily> TFamilyMap;
76
77   struct TField{
78     TVTKExtractFilter myExtractFilter;
79     vtkIdType myId;
80     std::string myName;
81     TEntity myEntity;
82     std::string myMeshName;
83     vtkIdType myNbComp, myNbValField, myDataSize, myIsTrimmed;
84     typedef std::vector<float> TValForCellsWithType;
85     typedef std::map<vtkIdType,TValForCellsWithType> TValForCells;
86     typedef std::pair<double,std::string> TTime;
87     typedef std::vector<std::string> TCompNames;
88     typedef std::vector<std::string> TUnitNames;
89     struct TValForTime{
90       TVTKAttribyteFilter myAttribyteFilter;
91       TVTKMergetFilter myMergeFilter;
92       TVTKSource myStorage;
93       vtkIdType myId;
94       std::string myMeshName;
95       TEntity myEntity;
96       std::string myFieldName;
97       vtkIdType myNbComp;
98       TTime myTime;
99       TValForCells myValForCells;
100       TValForTime() : myNbComp(0) {}
101       ~TValForTime() { 
102         if(myMergeFilter.GetPointer())
103           myMergeFilter->UnRegisterAllOutputs();
104         if(myAttribyteFilter.GetPointer())
105           myAttribyteFilter->UnRegisterAllOutputs();
106       }
107     };
108     typedef std::map<vtkIdType,TValForTime> TValField;
109     TValField myValField;
110     TCompNames myCompNames;
111     TUnitNames myUnitNames;
112     TField() : myNbComp(0), myNbValField(0), myDataSize(0), myIsTrimmed(0) {}
113     void ShallowCopy(const TField& aField);
114     ~TField() { 
115       if(myExtractFilter.GetPointer())
116         myExtractFilter->UnRegisterAllOutputs();
117     }
118   };
119   typedef std::map<std::string,TField> TFieldMap;
120  
121   struct TMeshOnEntity{
122     TVTKSource myStorage;
123     std::string myMeshName;
124     TEntity myEntity;
125     vtkIdType myNbCells, myCellsSize;
126     typedef std::vector<vtkIdType> TConnect;
127     typedef std::vector<TConnect> TConnForCellType;
128     typedef std::map<vtkIdType,TConnForCellType> TCellsConn;
129     TCellsConn myCellsConn;
130     TFamilyMap myFamilyMap;
131     TFieldMap myFieldMap;
132     TMeshOnEntity() : myNbCells(0), myCellsSize(0) {}
133     std::pair<vtkIdType,vtkIdType> GetCellsDims(const std::string& theFamilyName = "") const;
134       
135   };
136   typedef std::map<TEntity,TMeshOnEntity> TMeshOnEntityMap;
137   const TFamily* GetFamily(const VISU::TMeshOnEntity& theMeshOnEntity, 
138                      const std::string& theFamilyName);
139
140   TFamily* GetFamily(VISU::TMeshOnEntity& theMeshOnEntity, 
141                      const std::string& theFamilyName);
142
143   typedef std::pair<std::string,TEntity> TFamilyAndEntity;
144   typedef std::set<TFamilyAndEntity> TFamilyAndEntitySet;
145   struct TGroup{
146     TVTKSource myStorage;
147     std::string myName;
148     std::string myMeshName;
149     vtkIdType myNbCells, myCellsSize;
150     TGroup() : myNbCells(0), myCellsSize(0) {}
151     TFamilyAndEntitySet myFamilyAndEntitySet;
152   };
153   typedef std::map<std::string,TGroup> TGroupMap;
154
155   struct TMesh{
156     TVTKPoints myPoints;
157     vtkIdType myDim, myNbPoints;
158     std::string myName;
159     typedef float TCoord;
160     typedef std::vector<TCoord> TPointsCoord;
161     TPointsCoord myPointsCoord;
162     TMeshOnEntityMap myMeshOnEntityMap;
163     TGroupMap myGroupMap;
164     TMesh() : myDim(0), myNbPoints(0) {}
165     const TField* GetField(const std::string& theFieldName) const;
166   };
167   typedef std::map<std::string,TMesh> TMeshMap;
168   void WriteToFile(vtkUnstructuredGrid* theDataSet, const std::string& theFileName);
169 };
170
171 class VISU_Convertor{
172 protected:
173   std::string myName;
174   VISU::TMeshMap myMeshMap;
175   int myIsDone;
176 public:
177   virtual ~VISU_Convertor(){};
178   virtual const std::string& GetName() { return myName;}
179   virtual int IsDone() const { return myIsDone;}
180   typedef vtkUnstructuredGrid TOutput;
181
182   virtual VISU_Convertor* Build() = 0;
183   virtual const VISU::TMeshMap& GetMeshMap() ;
184   virtual float GetSize() = 0;
185
186   virtual TOutput* GetMeshOnEntity(const std::string& theMeshName, 
187                                    const VISU::TEntity& theEntity,
188                                    const std::string& theFamilyName = "") = 0;
189     
190   virtual float GetMeshOnEntitySize(const std::string& theMeshName, 
191                                     const VISU::TEntity& theEntity,
192                                     const std::string& theFamilyName = "") = 0;
193   
194   virtual TOutput* GetMeshOnGroup(const std::string& theMeshName, 
195                                   const std::string& theGroupName) = 0;
196
197   virtual float GetMeshOnGroupSize(const std::string& theMeshName, 
198                                    const std::string& theGroupName) = 0;
199
200   virtual TOutput* GetTimeStampOnMesh(const std::string& theMeshName, 
201                                       const VISU::TEntity& theEntity,
202                                       const std::string& theFieldName,
203                                       int theStampsNum) = 0;
204    
205   virtual float GetTimeStampSize(const std::string& theMeshName, 
206                                  const VISU::TEntity& theEntity,
207                                  const std::string& theFieldName,
208                                  int theStampsNum) = 0;
209     
210   virtual float GetFieldOnMeshSize(const std::string& theMeshName, 
211                                    const VISU::TEntity& theEntity,
212                                    const std::string& theFieldName) = 0;
213  
214   virtual const VISU::TField& GetField(const std::string& theMeshName, 
215                                        VISU::TEntity theEntity, 
216                                        const std::string& theFieldName) = 0;
217
218   virtual const VISU::TField::TValForTime& GetTimeStamp(const std::string& theMeshName, 
219                                                         const VISU::TEntity& theEntity,
220                                                         const std::string& theFieldName,
221                                                         int theStampsNum) = 0;
222     
223   static std::string GenerateName(const VISU::TField::TTime& aTime);
224   static std::string GenerateName(const std::string& theName, unsigned int theTimeId);
225   static void WriteToFile(vtkUnstructuredGrid* theDataSet, const std::string& theFileName);
226 };
227
228 extern "C"{
229   VISU_Convertor* CreateMedConvertor(const std::string& theFileName) ;
230   VISU_Convertor* CreateDatConvertor(const std::string& theFileName) ;
231   VISU_Convertor* CreateConvertor(const std::string& theFileName) ;
232 };
233
234 #endif