Salome HOME
8f6ed718c8b525d776c91fd5b71ac67a4aa1d0d5
[modules/visu.git] / src / CONVERTOR / VISU_Convertor_impl.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   : VISU_Convertor_impl.cxx
24 //  Author : Alexey PETROV
25 //  Module : VISU
26
27 #include "VISU_Convertor_impl.hxx"
28
29 #include <vtkIdList.h>
30 #include <vtkCellType.h>
31 #include <vtkIntArray.h>
32 #include <vtkCellArray.h>
33 #include <vtkFloatArray.h>
34 #include <vtkUnsignedCharArray.h>
35 #include <vtkPointData.h>
36 #include <vtkCellData.h>
37
38 #include <qstring.h>
39 #include <qfileinfo.h>
40
41 #include <valarray>     
42 #include <memory>
43
44 using namespace std;
45
46 #ifdef _DEBUG_
47 static int MYDEBUG = 0;
48 static int MYDEBUGWITHFILES = 0;
49 #else
50 static int MYDEBUG = 0;
51 static int MYDEBUGWITHFILES = 0;
52 #endif
53
54
55 VISU_Convertor_impl::VISU_Convertor_impl() {
56   myIsDone = false;
57 }
58
59 VISU_Convertor_impl::~VISU_Convertor_impl() {}
60
61 VISU_Convertor::TOutput* 
62 VISU_Convertor_impl::GetMeshOnEntity(const string& theMeshName, 
63                                      const VISU::TEntity& theEntity,
64                                      const string& theFamilyName)
65   throw (std::runtime_error&)
66 {
67   if(MYDEBUG) 
68     MESSAGE("GetMeshOnEntity - theMeshName = '"<<theMeshName<<
69             "'; theEntity = "<<theEntity<<"; theFamilyName = '"<<theFamilyName<<"'");
70   //Cheching possibility do the query
71   VISU::TMesh* pMesh = NULL;
72   VISU::TFamily* pFamily = NULL;
73   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
74   FindMeshOnEntity(theMeshName,pMesh,theEntity,pMeshOnEntity,theFamilyName,pFamily);
75   VISU::TMesh& aMesh = *pMesh;
76   VISU::TMeshOnEntity& aMeshOnEntity = *pMeshOnEntity;
77   VISU::TVTKSource* pSource;
78   if(pFamily != NULL)
79     pSource = &(pFamily->myStorage);
80   else
81     pSource = &(aMeshOnEntity.myStorage);
82   VISU::TVTKSource& aSource = *pSource;
83   //Main part of code
84   if(aSource.get() == NULL){
85     aSource.reset(TOutput::New());
86     LoadMeshOnEntity(aMeshOnEntity,theFamilyName);
87     GetPoints(aSource,aMesh);
88     GetCellsOnEntity(aSource,aMeshOnEntity,theFamilyName);
89     if(MYDEBUGWITHFILES){
90       string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
91       string aFamilyName = QString(theFamilyName.c_str()).simplifyWhiteSpace().latin1();
92       string aFileName = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
93       aFileName += aMeshName + dtos("-%d-",int(theEntity)) + aFamilyName + "-Conv.vtk";
94       VISU::WriteToFile(aSource.get(),aFileName);
95     }
96   }
97   return aSource.get();
98 }
99
100 VISU_Convertor::TOutput* 
101 VISU_Convertor_impl::GetMeshOnGroup(const string& theMeshName, 
102                                     const string& theGroupName)
103      throw(std::runtime_error&)
104 {
105   if(MYDEBUG) MESSAGE("GetMeshOnGroup - theMeshName = '"<<theMeshName<<
106                       "'; theGroupName = '"<<theGroupName<<"'");
107   //Cheching possibility do the query
108   VISU::TMesh* pMesh = NULL;
109   VISU::TGroup* pGroup = NULL;
110   FindMeshOnGroup(theMeshName,pMesh,theGroupName,pGroup);
111   VISU::TMesh& aMesh = *pMesh;
112   VISU::TGroup& aGroup = *pGroup;
113   const VISU::TFamilyAndEntitySet& aFamilyAndEntitySet = aGroup.myFamilyAndEntitySet;
114   VISU::TVTKSource& aSource = aGroup.myStorage;
115   //Main part of code
116   if(aSource.get() == NULL){
117     aSource.reset(TOutput::New());
118     LoadMeshOnGroup(aMesh,aFamilyAndEntitySet);
119     GetPoints(aSource,aMesh);
120     GetCellsOnGroup(aSource,aMesh,aFamilyAndEntitySet);
121     if(MYDEBUGWITHFILES){
122       string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
123       string aGroupName = QString(theGroupName.c_str()).simplifyWhiteSpace().latin1();
124       string aFileName = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
125       aFileName += aMeshName + "-" + aGroupName + "-Conv.vtk";
126       VISU::WriteToFile(aSource.get(),aFileName);
127     }
128   }
129   return aSource.get();
130 }
131
132 VISU_Convertor::TOutput* 
133 VISU_Convertor_impl::GetTimeStampOnMesh(const string& theMeshName, 
134                                         const VISU::TEntity& theEntity,
135                                         const string& theFieldName,
136                                         int theStampsNum)
137      throw(std::runtime_error&)
138 {
139   if(MYDEBUG){
140     MESSAGE("GetTimeStampOnMesh - theMeshName = '"<<theMeshName<<"; theEntity = "<<theEntity);
141     MESSAGE("GetTimeStampOnMesh - theFieldName = '"<<theFieldName<<"'; theStampsNum = "<<theStampsNum);
142   }
143   //Cheching possibility do the query
144   VISU::TMesh* pMesh = NULL;
145   VISU::TField* pField = NULL;
146   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
147   VISU::TField::TValForTime* pValForTime = NULL;
148   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
149                 theFieldName,pField,theStampsNum,pValForTime);
150   VISU::TMesh& aMesh = *pMesh;
151   VISU::TMeshOnEntity& aMeshOnEntity = *pMeshOnEntity;
152   VISU::TMeshOnEntity& aVTKMeshOnEntity = *pVTKMeshOnEntity;
153   VISU::TField& aField = *pField;
154   VISU::TField::TValForTime& aValForTime = *pValForTime;
155   VISU::TVTKSource& aSource = aValForTime.myStorage;
156   //Main part of code
157   if(aSource.get() == NULL){
158     aSource.reset(TOutput::New());
159     LoadFieldOnMesh(aMesh,aMeshOnEntity,aField,aValForTime);
160     try{
161       LoadMeshOnEntity(aVTKMeshOnEntity);
162     }catch(std::runtime_error& exc){
163       aVTKMeshOnEntity = aMeshOnEntity;
164       MESSAGE("Follow exception was accured :\n"<<exc.what());
165     }catch(...){
166       aVTKMeshOnEntity = aMeshOnEntity;
167       MESSAGE("Unknown exception was accured!");
168     }
169     GetMeshOnEntity(aVTKMeshOnEntity.myMeshName,aVTKMeshOnEntity.myEntity);
170     aSource->ShallowCopy(aVTKMeshOnEntity.myStorage.get());
171     GetField(aSource,aMesh,aVTKMeshOnEntity,aField,aValForTime);
172     if(MYDEBUGWITHFILES){
173       string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
174       string aFieldName = QString(theFieldName.c_str()).simplifyWhiteSpace().latin1();
175       string aFileName = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
176       aFileName += aMeshName + dtos("-%d-",int(theEntity)) + aFieldName + dtos("-%d",theStampsNum) + "-Conv.vtk";
177       VISU::WriteToFile(aSource.get(),aFileName);
178     }
179   }
180   return aSource.get();
181 }
182
183 inline void PrintCells(int& theStartId,
184                        vtkCellArray* theConnectivity, 
185                        const VISU::TMeshOnEntity::TConnect& theVector)
186 {
187   vtkIdList *anIdList = vtkIdList::New();
188   int kEnd = theVector.size();
189   anIdList->SetNumberOfIds(kEnd);
190   for(int k = 0; k < kEnd; k++){
191     anIdList->SetId(k,theVector[k]);
192     //anIdList->InsertNextId(theVector[k]);
193   }
194   theConnectivity->InsertNextCell(anIdList);
195   anIdList->Delete();
196 }
197
198 void VISU_Convertor_impl::GetCellsOnEntity(VISU::TVTKSource& theStorage,
199                                            const VISU::TMeshOnEntity& theMeshOnEntity, 
200                                            const string& theFamilyName) 
201   const throw (std::runtime_error&)
202 {
203   //Check on existing family
204   const VISU::TFamily* pFamily = VISU::GetFamily(theMeshOnEntity,theFamilyName);
205   bool isFamilyPresent = (pFamily != NULL);
206   const VISU::TFamily& aFamily = *pFamily;
207   //Main part of code
208   pair<int,int> aCellsDim = theMeshOnEntity.GetCellsDims(theFamilyName);
209   int aNbCells = aCellsDim.first, aCellsSize = aCellsDim.second;
210   vtkCellArray* aConnectivity = vtkCellArray::New();
211   //vtkIdType *anIdArray = aConnectivity->WritePointer(0,aCellsSize);
212   aConnectivity->Allocate(aCellsSize,0);
213   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
214   aCellTypesArray->SetNumberOfComponents(1);
215   aCellTypesArray->SetNumberOfTuples(aNbCells);
216   if(MYDEBUG) MESSAGE("GetCellsOnEntity - isFamilyPresent = "<<isFamilyPresent);
217   const VISU::TMeshOnEntity::TCellsConn &aCellsConn = theMeshOnEntity.myCellsConn;
218   VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
219   for(int i = 0, j = 0; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
220     const VISU::TMeshOnEntity::TConnForCellType& anArray = aCellsConnIter->second;
221     int aVtkType = aCellsConnIter->first;
222     if(MYDEBUG) MESSAGE("GetCellsOnEntity - aVtkType = "<<aVtkType<<"; anArray.size() = "<<anArray.size());
223     if(!isFamilyPresent)
224       for(int k = 0, kEnd = anArray.size(); k < kEnd; k++, i++){
225         PrintCells(i,aConnectivity,anArray[k]);
226         aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
227         //aCellTypesArray->InsertNextValue((unsigned char)aVtkType);
228       }
229     else{
230       const VISU::TFamily::TSubMesh& aSubMesh = aFamily.mySubMesh;
231       if(aSubMesh.empty()) throw std::runtime_error("GetCells >> There is no elements on the family !!!");
232       VISU::TFamily::TSubMesh::const_iterator aSubMeshIter = aSubMesh.find(aVtkType);
233       if(aSubMeshIter == aSubMesh.end()) continue;
234       const VISU::TFamily::TSubMeshOnCellType& aSubMeshOnCellType = aSubMeshIter->second;
235       if(MYDEBUG) MESSAGE("GetCellsOnEntity - aSubMeshOnCellType.size() = "<<aSubMeshOnCellType.size());
236       VISU::TFamily::TSubMeshOnCellType::const_iterator aSubMeshOnCellTypeIter = aSubMeshOnCellType.begin();
237       for(; aSubMeshOnCellTypeIter != aSubMeshOnCellType.end(); aSubMeshOnCellTypeIter++, i++){
238         PrintCells(i,aConnectivity,anArray[*aSubMeshOnCellTypeIter]);
239         aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
240         //aCellTypesArray->InsertNextValue((unsigned char)aVtkType);
241       }
242     }
243   }
244   vtkIdType *pts = 0, npts = 0;
245   vtkIntArray* aCellLocationsArray = vtkIntArray::New();
246   aCellLocationsArray->SetNumberOfComponents(1);
247   aCellLocationsArray->SetNumberOfTuples(aNbCells);
248   aConnectivity->InitTraversal();
249   for(int i=0; aConnectivity->GetNextCell(npts,pts); i++){
250     aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
251     //aCellLocationsArray->InsertNextValue(aConnectivity->GetTraversalLocation(npts));
252   }
253   theStorage->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
254
255
256
257 void VISU_Convertor_impl::GetCellsOnGroup(VISU::TVTKSource& theStorage,
258                                           const VISU::TMesh& theMesh,
259                                           const VISU::TFamilyAndEntitySet& theFamilyAndEntitySet) 
260   const throw (std::runtime_error&)
261 {
262   //Calculate dimentions of the group
263   int aNbCells = 0, aCellsSize = 0;
264   VISU::TFamilyAndEntitySet::const_iterator aFamilyAndEntitySetIter = theFamilyAndEntitySet.begin();
265   for(; aFamilyAndEntitySetIter != theFamilyAndEntitySet.end(); aFamilyAndEntitySetIter++){
266     const VISU::TFamilyAndEntity& aFamilyAndEntity = *aFamilyAndEntitySetIter;
267     const string& aFamilyName = aFamilyAndEntity.first;
268     const VISU::TEntity& anEntity = aFamilyAndEntity.second;
269     const VISU::TMeshOnEntity& aMeshOnEntity = theMesh.myMeshOnEntityMap.find(anEntity)->second;
270     pair<int,int> aCellsDim = aMeshOnEntity.GetCellsDims(aFamilyName);
271     aNbCells += aCellsDim.first;
272     aCellsSize += aCellsDim.second;
273   }
274   vtkCellArray* aConnectivity = vtkCellArray::New();
275   //vtkIdType *anIdArray = aConnectivity->WritePointer(0,aCellsSize);
276   aConnectivity->Allocate(aCellsSize,0);
277   vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
278   aCellTypesArray->SetNumberOfComponents(1);
279   aCellTypesArray->SetNumberOfTuples(aNbCells);
280   aFamilyAndEntitySetIter = theFamilyAndEntitySet.begin();
281   for(int i = 0, j = 0; aFamilyAndEntitySetIter != theFamilyAndEntitySet.end(); aFamilyAndEntitySetIter++){
282     const VISU::TFamilyAndEntity& aFamilyAndEntity = *aFamilyAndEntitySetIter;
283     const string& aFamilyName = aFamilyAndEntity.first;
284     const VISU::TEntity& anEntity = aFamilyAndEntity.second;
285     const VISU::TMeshOnEntity& aMeshOnEntity = theMesh.myMeshOnEntityMap.find(anEntity)->second;
286     const VISU::TFamily& aFamily = *(VISU::GetFamily(aMeshOnEntity,aFamilyName));
287     const VISU::TMeshOnEntity::TCellsConn &aCellsConn = aMeshOnEntity.myCellsConn;
288     VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
289     for(; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
290       const VISU::TMeshOnEntity::TConnForCellType& anArray = aCellsConnIter->second;
291       int aVtkType = aCellsConnIter->first;
292       if(MYDEBUG) MESSAGE("GetCellsOnGroup - aVtkType = "<<aVtkType<<"; anArray.size() = "<<anArray.size());
293       const VISU::TFamily::TSubMesh& aSubMesh = aFamily.mySubMesh;
294       if(aSubMesh.empty()) throw std::runtime_error("GetCells >> There is no elements on the family !!!");
295       VISU::TFamily::TSubMesh::const_iterator aSubMeshIter = aSubMesh.find(aVtkType);
296       if(aSubMeshIter == aSubMesh.end()) continue;
297       const VISU::TFamily::TSubMeshOnCellType& aSubMeshOnCellType = aSubMeshIter->second;
298       if(MYDEBUG) MESSAGE("GetCellsOnGroup - aSubMeshOnCellType.size() = "<<aSubMeshOnCellType.size());
299       VISU::TFamily::TSubMeshOnCellType::const_iterator aSubMeshOnCellTypeIter = aSubMeshOnCellType.begin();
300       for(; aSubMeshOnCellTypeIter != aSubMeshOnCellType.end(); aSubMeshOnCellTypeIter++, i++){
301         PrintCells(i,aConnectivity,anArray[*aSubMeshOnCellTypeIter]);
302         aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
303         //aCellTypesArray->InsertNextValue((unsigned char)aVtkType);
304       }
305     }
306   }
307   vtkIdType *pts = 0, npts = 0;
308   vtkIntArray* aCellLocationsArray = vtkIntArray::New();
309   aCellLocationsArray->SetNumberOfComponents(1);
310   aCellLocationsArray->SetNumberOfTuples(aNbCells);
311   aConnectivity->InitTraversal();
312   for(int i = 0; aConnectivity->GetNextCell(npts,pts); i++){
313     aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
314     //aCellLocationsArray->InsertNextValue(aConnectivity->GetTraversalLocation(npts));
315   }
316   theStorage->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
317
318
319
320 void VISU_Convertor_impl::GetPoints(VISU::TVTKSource& theStorage, const VISU::TMesh& theMesh) 
321   const throw (std::runtime_error&)
322 {
323   vtkPoints* aPoints = theMesh.myPoints.get();
324   if(!aPoints){
325     aPoints = vtkPoints::New();
326     const VISU::TMesh::TPointsCoord& anArray = theMesh.myPointsCoord;
327     vtkIdType iEnd = theMesh.myPointsCoord.size();
328     vtkIdType aNbPoints = iEnd / theMesh.myDim;
329     aPoints->SetNumberOfPoints(aNbPoints);
330     if(MYDEBUG) 
331       MESSAGE("GetPoints - aNbPoints = "<<aNbPoints<<"; myDim = "<<theMesh.myDim);
332     switch(theMesh.myDim) {
333     case 1:
334       for (vtkIdType i = 0, j = 0; i < iEnd; i += theMesh.myDim, j++) 
335         aPoints->SetPoint(j,anArray[i],0.0,0.0);
336       break;
337     case 2:
338       for (vtkIdType i = 0, j = 0; i < iEnd; i += theMesh.myDim, j++) 
339         aPoints->SetPoint(j,anArray[i],anArray[i+1],0.0);
340       break;
341     case 3: 
342       for (vtkIdType i = 0, j = 0; i < iEnd; i += theMesh.myDim, j++) 
343         aPoints->SetPoint(j,anArray[i],anArray[i+1],anArray[i+2]);
344       break;
345     }
346   }
347   theStorage->SetPoints(aPoints);
348 }
349
350 void VISU_Convertor_impl::GetField(VISU::TVTKSource& theStorage,
351                                    const VISU::TMesh& theMesh,
352                                    const VISU::TMeshOnEntity& theMeshOnEntity, 
353                                    const VISU::TField& theField, 
354                                    const VISU::TField::TValForTime& theValForTime)
355   const throw (std::runtime_error&)
356 {
357   if(MYDEBUG) MESSAGE("GetField - aTime = "<<theValForTime.myId<<"; theField.myName = "<<theField.myName);
358   int aNumberOfTuples;
359   vtkDataSetAttributes* aDataSetAttributes;
360   switch(theField.myEntity) {
361   case VISU::NODE_ENTITY : 
362     {
363       int aNbPoints = theMesh.myPointsCoord.size()/theMesh.myDim;
364       aNumberOfTuples = aNbPoints;
365       aDataSetAttributes = theStorage->GetPointData();
366       break;
367     }
368   default: 
369     {
370       pair<int,int> aCellsDim = theMeshOnEntity.GetCellsDims();
371       int aNbCells = aCellsDim.first, aCellsSize = aCellsDim.second;
372       aNumberOfTuples = aNbCells;
373       aDataSetAttributes = theStorage->GetCellData();
374     }
375   }
376   vtkFloatArray *aFloatArray = vtkFloatArray::New();
377   switch(theField.myNbComp) {
378   case 1:
379     aFloatArray->SetNumberOfComponents(1);
380     aDataSetAttributes->SetScalars(aFloatArray);
381     break;
382   default:
383     aFloatArray->SetNumberOfComponents(3);
384     aDataSetAttributes->SetVectors(aFloatArray);
385   }
386   aFloatArray->SetNumberOfTuples(aNumberOfTuples);
387   //const VISU::TField::TTime& aTime = theValForTime.myTime;
388   //string aFieldName = theField.myMeshName + ", " + theField.myName + ": " + GenerateName(aTime);
389   //aFloatArray->SetName(aFieldName.c_str());
390   if(MYDEBUG) MESSAGE("GetField - aNumberOfTuples = "<<aNumberOfTuples);
391   const VISU::TField::TValForCells& aValForCells = theValForTime.myValForCells;
392   VISU::TField::TValForCells::const_iterator aValForCellsIter = aValForCells.begin();
393   for(int k = 0; aValForCellsIter != aValForCells.end(); aValForCellsIter++) {
394     const VISU::TField::TValForCellsWithType& anArray = aValForCellsIter->second;
395     int iEnd = anArray.size()/theField.myNbComp;
396     int aVtkType = aValForCellsIter->first;
397     if(MYDEBUG) MESSAGE("GetField -  iEnd = "<<iEnd<<"; aVtkType = "<<aVtkType);
398     switch(theField.myNbComp) {
399     case 1:
400       for (int i = 0; i < iEnd; i++) 
401         aFloatArray->SetTuple1(k++,anArray[i]);
402       break;
403     case 2:
404       for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*2)
405         aFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+1],0.0);
406       break;
407     case 3:
408       for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*3)
409         aFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+1],anArray[ji+2]);
410       break;
411     default:
412       throw std::runtime_error("GetField - There is algorithm for representation the field !!!");
413     }
414   }
415 }
416
417
418 void VISU_Convertor_impl::FindMesh(const string& theMeshName, VISU::TMesh*& theMesh)
419   throw (std::runtime_error&)
420 {
421   GetMeshMap();
422   if(myMeshMap.find(theMeshName) == myMeshMap.end())
423     throw std::runtime_error("FindMesh >> There is no mesh with the name!!!");
424   theMesh = &myMeshMap[theMeshName];
425 }
426
427
428 void VISU_Convertor_impl::FindMeshOnEntity(const string& theMeshName, VISU::TMesh*& theMesh,
429                                            const VISU::TEntity& theEntity, VISU::TMeshOnEntity*& theMeshOnEntity,
430                                            const string& theFamilyName, VISU::TFamily*& theFamily)
431   throw (std::runtime_error&)
432 {
433   FindMesh(theMeshName,theMesh);
434   VISU::TMeshOnEntityMap& aMeshOnEntityMap = theMesh->myMeshOnEntityMap;
435   if(aMeshOnEntityMap.find(theEntity) == aMeshOnEntityMap.end())
436     throw std::runtime_error("FindMeshOnEntity >> There is no mesh on the entity!!!");
437   theMeshOnEntity = &aMeshOnEntityMap[theEntity];
438   theFamily = VISU::GetFamily(*theMeshOnEntity,theFamilyName);
439 }
440
441
442 vtkIdType VISU_Convertor_impl::GetMeshOnEntitySize(const std::string& theMeshName, 
443                                                    const VISU::TEntity& theEntity,
444                                                    const std::string& theFamilyName)
445   throw (std::runtime_error&)
446 {
447   VISU::TMesh* pMesh = NULL;
448   VISU::TFamily* pFamily = NULL;
449   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
450   FindMeshOnEntity(theMeshName,pMesh,theEntity,pMeshOnEntity,theFamilyName,pFamily);
451   vtkIdType aResult = 3*pMesh->myNbPoints*sizeof(VISU::TMesh::TCoord);
452   vtkIdType aNbCells, aCellsSize;
453   if(!pFamily){
454     aNbCells = pMeshOnEntity->myNbCells;
455     aCellsSize = pMeshOnEntity->myCellsSize;
456   }else{
457     aNbCells = pFamily->myNbCells;
458     aCellsSize = pFamily->myCellsSize;
459   }
460   //that is Connectivity + (Types + Locations + Links)
461   aResult += aCellsSize*sizeof(vtkIdType) + 
462     aNbCells*(sizeof(char) + sizeof(int) + (sizeof(short) + sizeof(vtkIdType)));
463   if(MYDEBUG) cout<<"VISU_Convertor_impl::GetMeshOnEntitySize = "<<aResult<<endl;
464   return aResult;
465 }
466
467
468 void VISU_Convertor_impl::FindMeshOnGroup(const std::string& theMeshName, VISU::TMesh*& theMesh,
469                                           const std::string& theGroupName, VISU::TGroup*& theGroup)
470   throw (std::runtime_error&)
471 {
472   FindMesh(theMeshName,theMesh);
473   VISU::TGroupMap& aGroupMap = theMesh->myGroupMap;
474   VISU::TGroupMap::iterator aGroupMapIter = aGroupMap.find(theGroupName);
475   if(aGroupMapIter == aGroupMap.end())
476     throw std::runtime_error("FindMeshOnGroup >> There is no the group in the mesh!!!");
477   theGroup = &aGroupMapIter->second;
478 }
479
480
481 vtkIdType VISU_Convertor_impl::GetMeshOnGroupSize(const std::string& theMeshName, 
482                                                   const std::string& theGroupName)
483   throw (std::runtime_error&)
484 {
485   VISU::TMesh* pMesh = NULL;
486   VISU::TGroup* pGroup = NULL;
487   FindMeshOnGroup(theMeshName,pMesh,theGroupName,pGroup);
488   vtkIdType aResult = 3*pMesh->myNbPoints*sizeof(VISU::TMesh::TCoord);
489   aResult += pGroup->myCellsSize*sizeof(vtkIdType);
490   if(MYDEBUG) cout<<"VISU_Convertor_impl::GetMeshOnGroupSize = "<<aResult<<endl;
491   return aResult;
492 }
493
494
495 void VISU_Convertor_impl::FindField(const string& theMeshName, VISU::TMesh*& theMesh,
496                                     const VISU::TEntity& theEntity, 
497                                     VISU::TMeshOnEntity*& theMeshOnEntity,
498                                     VISU::TMeshOnEntity*& theVTKMeshOnEntity,
499                                     const string& theFieldName, VISU::TField*& theField)
500   throw (std::runtime_error&)
501 {
502   VISU::TFamily* pFamily = NULL;
503   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
504   FindMeshOnEntity(theMeshName,theMesh,theEntity,pMeshOnEntity,"",pFamily);
505   theMeshOnEntity = pMeshOnEntity;
506   VISU::TMeshOnEntityMap& aMeshOnEntityMap = theMesh->myMeshOnEntityMap;
507   if(theEntity == VISU::NODE_ENTITY){
508     if(aMeshOnEntityMap.find(VISU::CELL_ENTITY) == aMeshOnEntityMap.end())
509       throw std::runtime_error("FindField >> There is no mesh on CELL_ENTITY!!!");
510     pMeshOnEntity = &aMeshOnEntityMap[VISU::CELL_ENTITY];
511   }
512   theVTKMeshOnEntity = pMeshOnEntity;
513   VISU::TFieldMap& aFieldMap = theMeshOnEntity->myFieldMap;
514   if(aFieldMap.find(theFieldName) == aFieldMap.end())
515     throw std::runtime_error("FindField >> There is no field on the mesh!!!");
516   theField = &aFieldMap[theFieldName];
517 }
518
519
520 vtkIdType VISU_Convertor_impl::GetFieldOnMeshSize(const std::string& theMeshName, 
521                                                   const VISU::TEntity& theEntity,
522                                                   const std::string& theFieldName)
523   throw(std::runtime_error&)
524 {
525   VISU::TMesh* pMesh = NULL;
526   VISU::TField* pField = NULL;
527   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
528   FindField(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
529             theFieldName,pField);
530   
531   vtkIdType aResult = GetMeshOnEntitySize(theMeshName,theEntity);
532   aResult += pField->myDataSize*sizeof(float)*pField->myNbValField;
533   if(MYDEBUG) cout<<"VISU_Convertor_impl::GetFieldOnMeshSize = "<<aResult<<endl;
534   return aResult;
535 }
536
537
538 void VISU_Convertor_impl::FindTimeStamp(const std::string& theMeshName, VISU::TMesh*& theMesh,
539                                         const VISU::TEntity& theEntity, 
540                                         VISU::TMeshOnEntity*& theMeshOnEntity,
541                                         VISU::TMeshOnEntity*& theVTKMeshOnEntity,
542                                         const std::string& theFieldName, VISU::TField*& theField,
543                                         int theStampsNum, VISU::TField::TValForTime*& theValForTime)
544   throw (std::runtime_error&)
545 {
546   FindField(theMeshName,theMesh,theEntity,theMeshOnEntity,theVTKMeshOnEntity,theFieldName,theField);
547   VISU::TField::TValField& aValField = theField->myValField;
548   if(aValField.find(theStampsNum) == aValField.end())
549     throw std::runtime_error("FindTimeStamp >> There is no field with the timestamp!!!");
550   theValForTime = &aValField[theStampsNum];
551 }
552
553
554 vtkIdType VISU_Convertor_impl::GetTimeStampSize(const std::string& theMeshName, 
555                                                 const VISU::TEntity& theEntity,
556                                                 const std::string& theFieldName,
557                                                 int theStampsNum)
558   throw (std::runtime_error&)
559 {
560   VISU::TMesh* pMesh = NULL;
561   VISU::TField* pField = NULL;
562   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
563   VISU::TField::TValForTime* pValForTime = NULL;
564   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
565                 theFieldName,pField,theStampsNum,pValForTime);
566
567   vtkIdType aResult = GetMeshOnEntitySize(theMeshName,theEntity);
568   aResult += pField->myDataSize*sizeof(float);
569   if(MYDEBUG) cout<<"VISU_Convertor_impl::GetTimeStampSize = "<<aResult<<endl;
570   return aResult;
571 }
572
573
574 const VISU::TField& 
575 VISU_Convertor_impl::GetField(const string& theMeshName, 
576                               VISU::TEntity theEntity, 
577                               const string& theFieldName) 
578   throw (std::runtime_error&)
579 {
580   VISU::TMesh* pMesh = NULL;
581   VISU::TField* pField = NULL;
582   VISU::TFamily* pFamily = NULL;
583   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
584   FindField(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
585             theFieldName,pField);
586   return *pField;
587 }
588
589
590 const VISU::TField::TValForTime& 
591 VISU_Convertor_impl::GetTimeStamp(const std::string& theMeshName, 
592                                   const VISU::TEntity& theEntity,
593                                   const std::string& theFieldName,
594                                   int theStampsNum)
595   throw (std::runtime_error&)
596 {
597   VISU::TMesh* pMesh = NULL;
598   VISU::TField* pField = NULL;
599   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
600   VISU::TField::TValForTime* pValForTime = NULL;
601   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
602                 theFieldName,pField,theStampsNum,pValForTime);
603   return *pValForTime;
604 }
605
606