Salome HOME
Fix for Bug IPAL8945
[modules/visu.git] / src / CONVERTOR / VISUConvertor.cxx
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:    VISUConvertor.cxx
24 //  Author:  Alexey PETROV
25 //  Module : VISU
26
27 #include "VISU_Convertor.hxx"
28 #include "VISU_ConvertorUtils.hxx"
29
30 #include <fstream>      
31 #include <strstream>
32 #include <vtkCellType.h>
33 #include <qdir.h>
34 #include <qfileinfo.h>
35 #include <qstringlist.h>
36 #include <memory>       
37 #include "VISU_Convertor_impl.hxx"
38
39 #include <vtkUnstructuredGrid.h>
40
41 using namespace std;
42
43 #ifdef DEBUG
44 static int MYDEBUG = 1;
45 #else
46 static int MYDEBUG = 0;
47 #endif
48
49 //#define _DEXCEPT_
50
51 typedef vtkUnstructuredGrid TOutput;
52
53 void parseFile(const char* theFileName) 
54 {
55 #ifndef _DEXCEPT_
56   try{
57 #endif
58     MSG(MYDEBUG,"'"<<theFileName<<"'...");
59     auto_ptr<VISU_Convertor> aCon(CreateConvertor(theFileName));
60     //aCon->GetSize();
61     //return;
62     aCon->BuildEntities();
63     aCon->BuildFields();
64     aCon->BuildMinMax();
65     const VISU::TMeshMap& aMeshMap = aCon->GetMeshMap();
66     //return;
67     VISU::TMeshMap::const_iterator aMeshMapIter = aMeshMap.begin();
68     for(; aMeshMapIter != aMeshMap.end(); aMeshMapIter++){
69       //continue;
70
71       const string& aMeshName = aMeshMapIter->first;
72       const VISU::PMesh& aMesh = aMeshMapIter->second;
73       const VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh->myMeshOnEntityMap;
74       VISU::TMeshOnEntityMap::const_iterator aMeshOnEntityMapIter;
75
76       //Import fields
77       aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
78       for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
79         const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
80         const VISU::PMeshOnEntity& aMeshOnEntity = aMeshOnEntityMapIter->second;
81         const VISU::TFieldMap& aFieldMap = aMeshOnEntity->myFieldMap;
82         VISU::TFieldMap::const_reverse_iterator aFieldMapIter = aFieldMap.rbegin();
83         for(; aFieldMapIter != aFieldMap.rend(); aFieldMapIter++){
84           const string& aFieldName = aFieldMapIter->first;
85           const VISU::PField& aField = aFieldMapIter->second;
86           const VISU::TValField& aValField = aField->myValField;
87           VISU::TValField::const_iterator aValFieldIter = aValField.begin();
88           for(; aValFieldIter != aValField.end(); aValFieldIter++){
89             int aTimeStamp = aValFieldIter->first;
90
91             if(anEntity != VISU::NODE_ENTITY){
92               VISU::PGaussPtsIDMapper aGaussMesh = 
93                 aCon->GetTimeStampOnGaussPts(aMeshName,anEntity,aFieldName,aTimeStamp);
94               VISU::TVTKOutput* aDataSet = aGaussMesh->GetVTKOutput();
95               /*
96               int aNbCells = aDataSet->GetNumberOfCells();
97               for(int anCellId = 0; anCellId < aNbCells; anCellId++){
98                 VISU::TGaussPointID anObjID = aGaussMesh->GetObjID(anCellId);
99                 cout<<anObjID.first<<"; "<<anObjID.second<<endl;
100               }
101               */
102             }else{
103               VISU::PIDMapper anIDMapper = 
104                 aCon->GetTimeStampOnMesh(aMeshName,anEntity,aFieldName,aTimeStamp);
105               VISU::TVTKOutput* aDataSet = anIDMapper->GetVTKOutput();
106               /*
107               int aNbCells = aDataSet->GetNumberOfCells();
108               for(int anCellId = 0; anCellId < aNbCells; anCellId++){
109                 int anObjID = anIDMapper->GetElemObjID(anCellId);
110                 int aVTKID  = anIDMapper->GetElemVTKID(anObjID);
111                 cout<<anObjID<<"; "<<aVTKID<<endl;
112               }
113               */
114             }
115             //goto OK;
116           }
117         }
118       }
119
120       //Importing groups
121       const VISU::TGroupMap& aGroupMap = aMesh->myGroupMap;
122       VISU::TGroupMap::const_iterator aGroupMapIter = aGroupMap.begin();
123       for(; aGroupMapIter != aGroupMap.end(); aGroupMapIter++){
124         const string& aGroupName = aGroupMapIter->first;
125         aCon->GetMeshOnGroup(aMeshName,aGroupName);
126       }
127
128       //continue;
129
130       //Import mesh on entity
131       aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
132       for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
133         const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
134         VISU::PIDMapper anIDMapper = aCon->GetMeshOnEntity(aMeshName,anEntity);
135         VISU::TVTKOutput* aDataSet = anIDMapper->GetVTKOutput();
136         {
137           /*
138           int aNbCells, anCellId, anObjID, aVTKID;
139           aNbCells = aDataSet->GetNumberOfCells();
140           for(anCellId = 0; anCellId < aNbCells; anCellId++){
141             anObjID = anIDMapper->GetElemObjID(anCellId);
142             aVTKID  = anIDMapper->GetElemVTKID(anObjID);
143             cout<<anObjID<<"; "<<aVTKID<<endl;
144           }
145           */
146         }
147       }
148
149       //Import families
150       aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
151       for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
152         const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
153         const VISU::PMeshOnEntity& aMeshOnEntity = aMeshOnEntityMapIter->second;
154         //aCon->GetMeshOnEntity(aMeshName,anEntity);
155         const VISU::TFamilyMap& aFamilyMap = aMeshOnEntity->myFamilyMap;
156         VISU::TFamilyMap::const_iterator aFamilyMapIter = aFamilyMap.begin();
157         for(; aFamilyMapIter != aFamilyMap.end(); aFamilyMapIter++){
158           const string& aFamilyName = aFamilyMapIter->first;
159           aCon->GetFamilyOnEntity(aMeshName,anEntity,aFamilyName);
160         }
161       }
162     }
163   OK:
164     MSG(MYDEBUG,"OK");
165 #ifndef _DEXCEPT_
166   }catch(std::exception& exc){
167     MSG(MYDEBUG,"Follow exception was occured in file:"<<theFileName<<"\n"<<exc.what());
168   }catch(...){
169     MSG(MYDEBUG,"Unknown exception was occured in VISU_Convertor_impl in file:"<<theFileName);
170   } 
171 #endif
172 }
173
174 int main(int argc, char** argv){ 
175   if(argc > 1){
176     QFileInfo fi(argv[1]);
177     for(int i = 0; i < 1; i++){
178       if(fi.exists()){
179         if(fi.isDir()){
180           QDir aDir(fi.absFilePath());
181           QStringList aStringList = aDir.entryList("*.med",QDir::Files);
182           int jEnd = aStringList.count();
183           for(int j = 0; j < jEnd; j++){
184             parseFile(aDir.filePath(aStringList[j]).latin1());
185           }
186         }else{
187           parseFile(argv[1]);
188         }
189       }
190     }
191     return 0;
192   }
193   return 1;
194 }