]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISUConvertor.cxx
Salome HOME
Join modifications from BR_Dev_For_4_0 tag V4_1_1.
[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.salome-platform.org/ or email : webmaster.salome@opencascade.com
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 _DEBUG_ID_MAPPING_
50 #define _DEXCEPT_
51
52 typedef vtkUnstructuredGrid TOutput;
53
54 void parseFile(const char* theFileName) 
55 {
56 #ifndef _DEXCEPT_
57   try{
58 #endif
59     MSG(MYDEBUG,"'"<<theFileName<<"'...");
60     //theFileName = "Apointe.med";
61     auto_ptr<VISU_Convertor> aCon(CreateConvertor(theFileName));
62     //aCon->GetSize();
63     //return;
64     aCon->BuildEntities();
65     aCon->BuildGroups();
66     aCon->BuildFields();
67     aCon->BuildMinMax();
68     const VISU::TMeshMap& aMeshMap = aCon->GetMeshMap();
69     //return;
70     VISU::TMeshMap::const_iterator aMeshMapIter = aMeshMap.begin();
71     for(; aMeshMapIter != aMeshMap.end(); aMeshMapIter++){
72       //continue;
73
74       const string& aMeshName = aMeshMapIter->first;
75       const VISU::PMesh& aMesh = aMeshMapIter->second;
76       const VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh->myMeshOnEntityMap;
77       VISU::TMeshOnEntityMap::const_iterator aMeshOnEntityMapIter;
78
79       //Import fields
80       aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
81       for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
82         const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
83         const VISU::PMeshOnEntity& aMeshOnEntity = aMeshOnEntityMapIter->second;
84         const VISU::TFieldMap& aFieldMap = aMeshOnEntity->myFieldMap;
85         VISU::TFieldMap::const_reverse_iterator aFieldMapIter = aFieldMap.rbegin();
86         for(; aFieldMapIter != aFieldMap.rend(); aFieldMapIter++){
87           const string& aFieldName = aFieldMapIter->first;
88           const VISU::PField& aField = aFieldMapIter->second;
89           const VISU::TValField& aValField = aField->myValField;
90           VISU::TValField::const_iterator aValFieldIter = aValField.begin();
91           for(; aValFieldIter != aValField.end(); aValFieldIter++){
92             int aTimeStamp = aValFieldIter->first;
93
94             if(anEntity != VISU::NODE_ENTITY){
95               VISU::PGaussPtsIDMapper aGaussMesh = 
96                 aCon->GetTimeStampOnGaussPts(aMeshName,anEntity,aFieldName,aTimeStamp);
97 #ifdef _DEBUG_ID_MAPPING_
98               vtkDataSet* aDataSet = aGaussMesh->GetOutput();
99               aDataSet->Update();
100               int aNbCells = aDataSet->GetNumberOfCells();
101               cout<<"aNbCells = "<<aNbCells<<endl;
102               for(int anCellId = 0; anCellId < aNbCells; anCellId++){
103                 VISU::TGaussPointID anObjID = aGaussMesh->GetObjID(anCellId);
104                 cout<<anObjID.first<<"; "<<anObjID.second<<"; "<<aGaussMesh->GetNodeVTKID(anObjID.first)<<endl;
105                 vtkFloatingPointType* aCoord = aGaussMesh->GetNodeCoord(anCellId);
106                 cout<<aCoord[0]<<"; "<<aCoord[1]<<"; "<<aCoord[2]<<endl;
107               }
108 #endif
109             }else{
110               //continue;
111               VISU::PIDMapper anIDMapper = 
112                 aCon->GetTimeStampOnMesh(aMeshName,anEntity,aFieldName,aTimeStamp);
113 #ifdef _DEBUG_ID_MAPPING_
114               vtkDataSet* aDataSet = anIDMapper->GetOutput();
115               aDataSet->Update();
116               int aNbCells = aDataSet->GetNumberOfCells();
117               for(int anCellId = 0; anCellId < aNbCells; anCellId++){
118                 int anObjID = anIDMapper->GetElemObjID(anCellId);
119                 int aVTKID  = anIDMapper->GetElemVTKID(anObjID);
120                 cout<<anObjID<<"; "<<aVTKID<<endl;
121               }
122 #endif
123             }
124             //goto OK;
125           }
126         }
127       }
128
129       //Importing groups
130       const VISU::TGroupMap& aGroupMap = aMesh->myGroupMap;
131       VISU::TGroupMap::const_iterator aGroupMapIter = aGroupMap.begin();
132       for(; aGroupMapIter != aGroupMap.end(); aGroupMapIter++){
133         const string& aGroupName = aGroupMapIter->first;
134         aCon->GetMeshOnGroup(aMeshName,aGroupName);
135       }
136
137       //continue;
138
139       //Import mesh on entity
140       aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
141       for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
142         const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
143         VISU::PIDMapper anIDMapper = aCon->GetMeshOnEntity(aMeshName,anEntity);
144 #ifdef _DEBUG_ID_MAPPING_
145         vtkDataSet* aDataSet = anIDMapper->GetOutput();
146         int aNbCells, anCellId, anObjID, aVTKID;
147         aNbCells = aDataSet->GetNumberOfCells();
148         for(anCellId = 0; anCellId < aNbCells; anCellId++){
149           anObjID = anIDMapper->GetElemObjID(anCellId);
150           aVTKID  = anIDMapper->GetElemVTKID(anObjID);
151           cout<<anObjID<<"; "<<aVTKID<<endl;
152         }
153 #endif
154       }
155
156       //Import families
157       aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
158       for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
159         const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
160         const VISU::PMeshOnEntity& aMeshOnEntity = aMeshOnEntityMapIter->second;
161         //aCon->GetMeshOnEntity(aMeshName,anEntity);
162         const VISU::TFamilyMap& aFamilyMap = aMeshOnEntity->myFamilyMap;
163         VISU::TFamilyMap::const_iterator aFamilyMapIter = aFamilyMap.begin();
164         for(; aFamilyMapIter != aFamilyMap.end(); aFamilyMapIter++){
165           const string& aFamilyName = aFamilyMapIter->first;
166           aCon->GetFamilyOnEntity(aMeshName,anEntity,aFamilyName);
167         }
168       }
169     }
170   OK:
171     MSG(MYDEBUG,"OK");
172 #ifndef _DEXCEPT_
173   }catch(std::exception& exc){
174     MSG(MYDEBUG,"Follow exception was occured in file:"<<theFileName<<"\n"<<exc.what());
175   }catch(...){
176     MSG(MYDEBUG,"Unknown exception was occured in VISU_Convertor_impl in file:"<<theFileName);
177   } 
178 #endif
179 }
180
181 int main(int argc, char** argv){ 
182   if(argc > 1){
183     QFileInfo fi(argv[1]);
184     for(int i = 0; i < 1; i++){
185       if(fi.exists()){
186         if(fi.isDir()){
187           QDir aDir(fi.absFilePath());
188           QStringList aStringList = aDir.entryList("*.med",QDir::Files);
189           int jEnd = aStringList.count();
190           for(int j = 0; j < jEnd; j++){
191             parseFile(aDir.filePath(aStringList[j]).latin1());
192           }
193         }else{
194           parseFile(argv[1]);
195         }
196       }
197     }
198     return 0;
199   }
200   return 1;
201 }