]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_Convertor_impl.cxx
Salome HOME
Join modifications from BR_Dev_For_4_0 tag V4_1_1.
[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.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 //
23 //  File   : VISU_Convertor_impl.cxx
24 //  Author : Alexey PETROV
25 //  Module : VISU
26
27 #include "VISU_Convertor_impl.hxx"
28 #include "VISU_Structures_impl.hxx"
29 #include "VISU_PointCoords.hxx"
30 #include "VISU_MeshValue.hxx"
31
32 #include "VISU_AppendFilter.hxx"
33 #include "VISU_AppendPolyData.hxx"
34 #include "VTKViewer_CellLocationsArray.h"
35 #include "VISU_CommonCellsFilter.hxx"
36
37 #include "VISU_ConvertorUtils.hxx"
38
39 #include <vtkUnstructuredGrid.h>
40 #include <vtkPolyData.h>
41
42 #include <vtkPoints.h>
43 #include <vtkPointData.h>
44 #include <vtkCellData.h>
45
46 #include <vtkIdList.h>
47 #include <vtkCellType.h>
48 #include <vtkCellArray.h>
49 #include <vtkCellLinks.h>
50 #include <vtkUnsignedCharArray.h>
51
52 #include <qstring.h>
53 #include <qfileinfo.h>
54
55 #include <memory>
56
57 static vtkFloatingPointType ERR_SIZE_CALC = 1.00;
58
59 static int MYVTKDEBUG = 0;
60
61 #ifdef _DEBUG_
62 static int MYDEBUG = 0;
63 static int MYDEBUGWITHFILES = 0;
64 //#define _DEXCEPT_
65 #else
66 static int MYDEBUG = 0;
67 static int MYDEBUGWITHFILES = 0;
68 #endif
69
70
71 namespace
72 {
73   //---------------------------------------------------------------
74   template<class T> 
75   std::string 
76   dtos(const std::string& fmt, T val)
77   {
78     static QString aString;
79     aString.sprintf(fmt.c_str(),val);
80     return aString.latin1();
81   }
82
83   //--------------------------------------------------------------
84
85
86   
87   //---------------------------------------------------------------
88   void
89   PrintCells(int& theStartId,
90              vtkCellArray* theConnectivity, 
91              const VISU::TConnect& theVector)
92   {
93     vtkIdList *anIdList = vtkIdList::New();
94     int kEnd = theVector.size();
95     anIdList->SetNumberOfIds(kEnd);
96     for(int k = 0; k < kEnd; k++)
97       anIdList->SetId(k,theVector[k]);
98     theConnectivity->InsertNextCell(anIdList);
99     anIdList->Delete();
100   }
101
102
103   //---------------------------------------------------------------
104   void
105   GetCellsOnSubMesh(const VISU::PUnstructuredGrid& theSource,
106                     const VISU::PMeshOnEntityImpl& theMeshOnEntity, 
107                     const VISU::PSubMeshImpl& theSubMesh,
108                     const vtkIdType theGeom) 
109   {
110     VISU::TTimerLog aTimerLog(MYDEBUG,"GetCellsOnSubMesh");
111     const VISU::TCell2Connect& anArray = theSubMesh->myCell2Connect;
112     vtkIdType aCellsSize = theSubMesh->myCellsSize;
113     vtkIdType aNbCells = theSubMesh->myNbCells;
114     INITMSG(MYDEBUG,"GetCellsOnSubMesh "<<
115             "- theGeom = "<<theGeom<<
116             "; aNbCells = "<<aNbCells<<
117             endl);
118
119
120     vtkCellArray* aConnectivity = vtkCellArray::New();
121     aConnectivity->Allocate(aCellsSize,0);
122     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
123     aCellTypesArray->SetNumberOfComponents(1);
124     aCellTypesArray->SetNumberOfTuples(aNbCells);
125
126     for(vtkIdType anID = 0; anID < aNbCells; anID++){
127       PrintCells(anID,aConnectivity,anArray[anID]);
128       aCellTypesArray->SetValue(anID,(unsigned char)theGeom);
129     }
130
131     {
132       int aNbTuples = aNbCells;
133       int anEntity = int(theMeshOnEntity->myEntity);
134       vtkIntArray *aDataArray = vtkIntArray::New();
135       aDataArray->SetName("VISU_CELLS_MAPPER");
136       aDataArray->SetNumberOfComponents(2);
137       aDataArray->SetNumberOfTuples(aNbTuples);
138       int *aDataArrayPtr = aDataArray->GetPointer(0);
139       for(int aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
140         int anObjID = theSubMesh->GetElemObjID(aTupleId);
141         *aDataArrayPtr++ = anObjID;
142         *aDataArrayPtr++ = anEntity;
143       }
144       theSource->GetCellData()->AddArray(aDataArray);
145       aDataArray->Delete();
146     }
147
148     vtkIdType *pts = 0, npts = 0;
149     VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
150     aCellLocationsArray->SetNumberOfComponents(1);
151     aCellLocationsArray->SetNumberOfTuples(aNbCells);
152     aConnectivity->InitTraversal();
153     for(int i=0; aConnectivity->GetNextCell(npts,pts); i++)
154       aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
155     theSource->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
156
157     if(MYVTKDEBUG) aConnectivity->DebugOn();
158
159     aCellLocationsArray->Delete();
160     aCellTypesArray->Delete();
161     aConnectivity->Delete();
162   } 
163   
164   
165   //---------------------------------------------------------------
166   void
167   GetCellsOnFamily(const VISU::PUnstructuredGrid& theSource,
168                    const VISU::PMeshOnEntityImpl& theMeshOnEntity, 
169                    const VISU::PFamilyImpl& theFamily) 
170   {
171     INITMSG(MYDEBUG,"GetCellsOnFamily"<<endl);
172
173     vtkIdType aNbCells = theFamily->myNbCells;
174     vtkIdType aCellsSize = theFamily->myCellsSize;
175
176     vtkCellArray* aConnectivity = vtkCellArray::New();
177     aConnectivity->Allocate(aCellsSize,0);
178     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
179     aCellTypesArray->SetNumberOfComponents(1);
180     aCellTypesArray->SetNumberOfTuples(aNbCells);
181
182     VISU::TSubMeshID& aMeshID = theFamily->myMeshID;
183     aMeshID.resize(aNbCells);
184
185     vtkIntArray *aDataArray = vtkIntArray::New();
186     int anEntity = int(theMeshOnEntity->myEntity);
187     aDataArray->SetName("VISU_CELLS_MAPPER");
188     aDataArray->SetNumberOfComponents(2);
189     aDataArray->SetNumberOfTuples(aNbCells);
190     int *aDataArrayPtr = aDataArray->GetPointer(0);
191
192     VISU::TID2ID& anElemObj2VTKID = theFamily->myElemObj2VTKID;
193
194     const VISU::TGeom2SubMesh& aGeom2SubMesh = theMeshOnEntity->myGeom2SubMesh;
195     VISU::TGeom2SubMesh::const_iterator anIter = aGeom2SubMesh.begin();
196     for(vtkIdType aCellId = 0; anIter != aGeom2SubMesh.end(); anIter++){
197       VISU::EGeometry aEGeom = anIter->first;
198       vtkIdType aVGeom = VISUGeom2VTK(aEGeom);
199
200       const VISU::TSubMeshImpl& aSubMesh = anIter->second;
201       const VISU::TCell2Connect& anArray = aSubMesh.myCell2Connect;
202
203       const VISU::TGeom2SubMeshID& aGeom2SubMeshID = theFamily->myGeom2SubMeshID;
204       if(aGeom2SubMeshID.empty()) 
205         EXCEPTION(std::runtime_error,"GetCells >> There is no elements on the family !!!");
206
207       VISU::TGeom2SubMeshID::const_iterator aGeom2SubMeshIDIter = aGeom2SubMeshID.find(aEGeom);
208       if(aGeom2SubMeshIDIter == aGeom2SubMeshID.end()) 
209         continue;
210
211       const VISU::TSubMeshID& aSubMeshID = aGeom2SubMeshIDIter->second;
212
213       INITMSG(MYDEBUG,
214               " - aEGeom = "<<aEGeom<<
215               "; aVGeom = "<<aVGeom<<
216               "; aSubMeshID.size() = "<<aSubMeshID.size()<<
217               endl);
218
219       VISU::TSubMeshID::const_iterator aSubMeshIDIter = aSubMeshID.begin();
220       for(; aSubMeshIDIter != aSubMeshID.end(); aSubMeshIDIter++, aCellId++){
221         vtkIdType anID = *aSubMeshIDIter;
222         PrintCells(aCellId, aConnectivity, anArray[anID]);
223         aCellTypesArray->SetValue(aCellId, (unsigned char)aVGeom);
224         vtkIdType anObjID = aSubMesh.GetElemObjID(anID);
225         anElemObj2VTKID[anObjID] = aCellId;
226         aMeshID[aCellId] = anObjID;
227         *aDataArrayPtr++ = anObjID;
228         *aDataArrayPtr++ = anEntity;
229       }
230     }
231
232     theSource->GetCellData()->AddArray(aDataArray);
233     aDataArray->Delete();
234
235     vtkIdType *pts = 0, npts = 0;
236     VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
237     aCellLocationsArray->SetNumberOfComponents(1);
238     aCellLocationsArray->SetNumberOfTuples(aNbCells);
239     aConnectivity->InitTraversal();
240     for(int i=0; aConnectivity->GetNextCell(npts,pts); i++)
241       aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
242     theSource->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
243
244     if(MYVTKDEBUG) aConnectivity->DebugOn();
245
246     aCellLocationsArray->Delete();
247     aCellTypesArray->Delete();
248     aConnectivity->Delete();
249   }
250   
251   
252   //---------------------------------------------------------------
253   void
254   GetCells(const VISU::PUnstructuredGrid& theSource,
255            const VISU::PSubProfileImpl& theSubProfile,
256            const VISU::PProfileImpl& theProfile,
257            const VISU::PMeshOnEntityImpl& theMeshOnEntity)
258   {
259     vtkIdType aNbCells = theSubProfile->myNbCells;
260     vtkIdType aCellsSize = theSubProfile->myCellsSize;
261     VISU::EGeometry aEGeom = theSubProfile->myGeom;
262     vtkIdType aNbNodes = VISUGeom2NbNodes(aEGeom);
263     vtkIdType aVGeom = VISUGeom2VTK(aEGeom);
264
265     INITMSG(MYDEBUG,"GetCells - aVGeom = "<<aVGeom<<endl);
266
267     const VISU::TSubMeshID& aSubMeshID = theSubProfile->mySubMeshID;
268
269     const VISU::TGeom2SubMesh& aGeom2SubMesh = theMeshOnEntity->myGeom2SubMesh;
270     VISU::TGeom2SubMesh::const_iterator anIter = aGeom2SubMesh.find(aEGeom);
271     if(anIter == aGeom2SubMesh.end())
272       EXCEPTION(std::runtime_error,"GetCells >> There is no elements for the GEOM("<<aEGeom<<")");
273     
274     const VISU::TSubMeshImpl& aSubMesh = anIter->second;
275     const VISU::TCell2Connect& aCell2Connect = aSubMesh.myCell2Connect;
276     
277     vtkCellArray* aConnectivity = vtkCellArray::New();
278     aConnectivity->Allocate(aCellsSize,0);
279     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
280     aCellTypesArray->SetNumberOfComponents(1);
281     aCellTypesArray->SetNumberOfTuples(aNbCells);
282     
283     if(theSubProfile->myStatus == VISU::eAddAll){
284       VISU::TCell2Connect::const_iterator anIter = aCell2Connect.begin();
285       for(vtkIdType anId = 0, aConnId = 0; anIter != aCell2Connect.end(); anIter++){
286         const VISU::TConnect& anArray = aCell2Connect[anId];
287         PrintCells(aConnId,aConnectivity,anArray);
288         aCellTypesArray->SetValue(anId,(unsigned char)aVGeom);
289         aConnId += aNbNodes;
290         anId++;
291       }
292     }else{
293       VISU::TSubMeshID::const_iterator anIter = aSubMeshID.begin();
294       for(vtkIdType anId = 0, aConnId = 0; anIter != aSubMeshID.end(); anIter++){
295         vtkIdType aSubId = *anIter;
296         const VISU::TConnect& anArray = aCell2Connect[aSubId];
297         PrintCells(aConnId,aConnectivity,anArray);
298         aCellTypesArray->SetValue(anId,(unsigned char)aVGeom);
299         aConnId += aNbNodes;
300         anId++;
301       }
302     }
303     
304     vtkIdType *pts = 0, npts = 0;
305     VTKViewer_CellLocationsArray* aCellLocationsArray = VTKViewer_CellLocationsArray::New();
306     
307     aCellLocationsArray->SetNumberOfComponents(1);
308     aCellLocationsArray->SetNumberOfTuples(aNbCells);
309     aConnectivity->InitTraversal();
310     for(int i=0; aConnectivity->GetNextCell(npts,pts); i++)
311       aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
312     theSource->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
313     
314     {
315       int aNbTuples = aNbCells;
316       int anEntity = int(theMeshOnEntity->myEntity);
317       vtkIntArray *aDataArray = vtkIntArray::New();
318       aDataArray->SetName("VISU_CELLS_MAPPER");
319       aDataArray->SetNumberOfComponents(2);
320       aDataArray->SetNumberOfTuples(aNbTuples);
321       int *aDataArrayPtr = aDataArray->GetPointer(0);
322       for(int aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
323         int anObjID = theSubProfile->GetElemObjID(aTupleId);
324         *aDataArrayPtr++ = anObjID;
325         *aDataArrayPtr++ = anEntity;
326       }
327       theSource->GetCellData()->AddArray(aDataArray);
328       aDataArray->Delete();
329     }
330     
331     aCellLocationsArray->Delete();
332     aCellTypesArray->Delete();
333     aConnectivity->Delete();
334   }
335   
336   
337   //---------------------------------------------------------------
338   void
339   GetMeshOnSubProfile(const VISU::PMeshImpl& theMesh,
340                       const VISU::PMeshOnEntityImpl& theMeshOnEntity,
341                       const VISU::PProfileImpl& theProfile,
342                       const VISU::PSubProfileImpl& theSubProfile)
343   {
344     INITMSG(MYDEBUG,"GetMeshOnSubProfile - aEGeom = "<<theSubProfile->myGeom<<endl);
345     
346     const VISU::PUnstructuredGrid& aSource = theSubProfile->GetSource();
347     if(theSubProfile->myIsVTKDone)
348       return;
349     
350     aSource->ShallowCopy(theMesh->GetPointSet());
351     INITMSGA(MYDEBUG,0,"GetNumberOfPoints - "<<aSource->GetNumberOfPoints()<<endl);
352     GetCells(aSource,theSubProfile,theProfile,theMeshOnEntity);
353     BEGMSG(MYDEBUG,"GetNumberOfCells - "<<aSource->GetNumberOfCells()<<endl);
354     
355     theSubProfile->myIsVTKDone = true;
356   }
357   
358   
359   //---------------------------------------------------------------
360   bool
361   GetMeshOnProfile(const VISU::PMeshImpl& theMesh,
362                    const VISU::PMeshOnEntityImpl& theMeshOnEntity,
363                    const VISU::PProfileImpl& theProfile)
364   {
365     if(theProfile->myIsVTKDone)
366       return true;
367    
368     if(theProfile->myMeshOnEntity && theProfile->myMeshOnEntity != theMeshOnEntity.get())
369       return false;
370       
371     VISU::TTimerLog aTimerLog(MYDEBUG,"GetMeshOnProfile");
372     INITMSG(MYDEBUG,"GetMeshOnProfile - anEntity = "<<theMeshOnEntity->myEntity<<std::endl);
373
374     const VISU::PAppendFilter& anAppendFilter = theProfile->GetFilter();
375     anAppendFilter->SetSharedPointSet(theMesh->GetPointSet());
376
377     if(theProfile->myIsAll){
378       vtkUnstructuredGrid* aDataSet = theMeshOnEntity->GetUnstructuredGridOutput();
379       anAppendFilter->AddInput(aDataSet);
380     }else{
381       const VISU::TGeom2SubProfile& aGeom2SubProfile = theProfile->myGeom2SubProfile;
382
383       VISU::TID2ID& anElemObj2VTKID = theProfile->myElemObj2VTKID;
384
385       VISU::TSubProfileArr& aSubProfileArr = theProfile->mySubProfileArr;
386       aSubProfileArr.resize(aGeom2SubProfile.size());
387
388       VISU::TGeom2SubProfile::const_iterator anIter = aGeom2SubProfile.begin();
389       for(vtkIdType anInputID = 0, aCellID = 0; anIter != aGeom2SubProfile.end(); anIter++){
390         VISU::PSubProfileImpl aSubProfile = anIter->second;
391         if(aSubProfile->myStatus == VISU::eRemoveAll)
392           continue;
393         
394         GetMeshOnSubProfile(theMesh,
395                             theMeshOnEntity,
396                             theProfile,
397                             aSubProfile);
398         
399         const VISU::PUnstructuredGrid& aSource = aSubProfile->GetSource();
400         anAppendFilter->AddInput(aSource.GetPointer());
401
402         vtkIdType aNbCells = aSource->GetNumberOfCells();
403         for(vtkIdType aCell = 0; aCell < aNbCells; aCell++, aCellID++){
404           vtkIdType anObjID = aSubProfile->GetElemObjID(aCell);
405           anElemObj2VTKID[anObjID] = aCellID;
406         }
407
408         aSubProfileArr[anInputID++] = aSubProfile;
409       }
410     }
411     anAppendFilter->Update(); // Fix on VTK
412     theProfile->myMeshOnEntity = theMeshOnEntity.get();
413     theProfile->myNamedPointCoords = theMesh->myNamedPointCoords;
414     
415     theProfile->myIsVTKDone = true;
416     return true;
417   }
418   
419   
420   //---------------------------------------------------------------
421   void
422   GetGaussSubMeshSource(const VISU::PPolyData& theSource,
423                         const VISU::PGaussSubMeshImpl& theGaussSubMesh,
424                         const VISU::PMeshOnEntityImpl& theMeshOnEntity)
425   {
426     vtkCellArray* aConnectivity = vtkCellArray::New();
427     vtkIdType aCellsSize = theGaussSubMesh->myCellsSize;
428     aConnectivity->Allocate(aCellsSize, 0);
429     
430     vtkIdList *anIdList = vtkIdList::New();
431     anIdList->SetNumberOfIds(1);
432
433     const VISU::TPointCoords& aCoords = theGaussSubMesh->myPointCoords;
434     vtkIdType aNbPoints = aCoords.GetNbPoints();
435     for(vtkIdType aPointId = 0; aPointId < aNbPoints; aPointId++){
436       anIdList->SetId(0, aPointId);
437       aConnectivity->InsertNextCell(anIdList);
438     }
439     anIdList->Delete();
440     
441     const VISU::PPolyData& aSource = theGaussSubMesh->GetSource();
442     aSource->ShallowCopy(aCoords.GetPointSet());
443     aSource->SetVerts(aConnectivity);
444     
445     aConnectivity->Delete();
446
447     {
448       vtkIdType aNbTuples = aNbPoints;
449       vtkIntArray *aDataArray = vtkIntArray::New();
450       aDataArray->SetName("VISU_POINTS_MAPPER");
451       aDataArray->SetNumberOfComponents(2);
452       aDataArray->SetNumberOfTuples(aNbTuples);
453       int *aDataArrayPtr = aDataArray->GetPointer(0);
454       for(vtkIdType aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
455         vtkIdType aGlobalID = theGaussSubMesh->GetGlobalID(aTupleId);
456         *aDataArrayPtr++ = aGlobalID;
457         *aDataArrayPtr++ = 0;
458       }
459       aSource->GetPointData()->AddArray(aDataArray);
460       aDataArray->Delete();
461     }
462
463     {
464       vtkIdType aNbTuples = aNbPoints;
465       vtkIntArray *aDataArray = vtkIntArray::New();
466       aDataArray->SetName("VISU_CELLS_MAPPER");
467       aDataArray->SetNumberOfComponents(2);
468       aDataArray->SetNumberOfTuples(aNbTuples);
469       int *aDataArrayPtr = aDataArray->GetPointer(0);
470       for(vtkIdType aTupleId = 0; aTupleId < aNbTuples; aTupleId++){
471         VISU::TGaussPointID aGaussPointID = theGaussSubMesh->GetObjID(aTupleId);
472         *aDataArrayPtr++ = aGaussPointID.first;
473         *aDataArrayPtr++ = aGaussPointID.second;
474       }
475       aSource->GetCellData()->AddArray(aDataArray);
476       aDataArray->Delete();
477     }
478   }
479   
480   
481   //---------------------------------------------------------------
482   void
483   GetGaussSubMesh(const VISU::PMeshImpl& theMesh,
484                   const VISU::PMeshOnEntityImpl& theMeshOnEntity,
485                   const VISU::PGaussMeshImpl& theGaussMesh,
486                   const VISU::PGaussSubMeshImpl& theGaussSubMesh)
487   {
488     VISU::PGaussImpl aGauss = theGaussSubMesh->myGauss;
489     
490     if(!theGaussSubMesh->myIsDone)
491       return;
492     
493     if(theGaussSubMesh->myIsVTKDone)
494       return;
495     
496     VISU::TTimerLog aTimerLog(MYDEBUG,"GetGaussSubMesh");
497     INITMSG(MYDEBUG,"GetGaussSubMesh - aVGeom = "<<aGauss->myGeom<<endl);
498
499     const VISU::PPolyData& aSource = theGaussSubMesh->GetSource();
500     GetGaussSubMeshSource(aSource, theGaussSubMesh, theMeshOnEntity);
501
502     INITMSGA(MYDEBUG,0,"GetNumberOfPoints - "<<aSource->GetNumberOfPoints()<<endl);
503     BEGMSG(MYDEBUG,"GetNumberOfCells - "<<aSource->GetNumberOfCells()<<endl);
504     
505     theGaussSubMesh->myIsVTKDone = true;
506   }
507   
508
509   //---------------------------------------------------------------
510   void
511   BuildGaussMesh(const VISU::PMeshImpl& theMesh,
512                  const VISU::PMeshOnEntityImpl& theMeshOnEntity,
513                  const VISU::PGaussMeshImpl& theGaussMesh)
514   {
515     if(theGaussMesh->myIsVTKDone)
516       return;
517
518     VISU::TTimerLog aTimerLog(MYDEBUG,"BuildGaussMesh");
519     const VISU::PAppendPolyData& anAppendFilter = theGaussMesh->GetFilter();
520     const VISU::TGeom2GaussSubMesh& aGeom2GaussSubMesh = theGaussMesh->myGeom2GaussSubMesh;
521     VISU::TGeom2GaussSubMesh::const_iterator anIter = aGeom2GaussSubMesh.begin();
522     for(vtkIdType aStartID = 0; anIter != aGeom2GaussSubMesh.end(); anIter++){
523       VISU::PGaussSubMeshImpl aGaussSubMesh = anIter->second;
524       if(aGaussSubMesh->myStatus == VISU::eRemoveAll)
525         continue;
526
527       aGaussSubMesh->myStartID = aStartID;
528
529       GetGaussSubMesh(theMesh,
530                       theMeshOnEntity,
531                       theGaussMesh,
532                       aGaussSubMesh);
533       
534       const VISU::PPolyData& aSource = aGaussSubMesh->GetSource();
535       aStartID += aSource->GetNumberOfCells();
536
537       anAppendFilter->AddInput(aSource.GetPointer());
538     }
539     anAppendFilter->Update(); // Fix on VTK
540
541     theMeshOnEntity->GetOutput()->Update();
542
543     vtkDataSet* aSource = anAppendFilter->GetOutput();
544     INITMSGA(MYDEBUG,0,"aNbPoints - "<<aSource->GetNumberOfPoints()<<endl);
545     BEGMSG(MYDEBUG,"aNbCells - "<<aSource->GetNumberOfCells()<<endl);
546     
547     theGaussMesh->myIsVTKDone = true;
548   }
549
550
551   //---------------------------------------------------------------
552   void
553   PrintMemorySize(vtkUnstructuredGrid* theDataSet)
554   {
555     theDataSet->Update();
556     BEGMSG(1,"GetPoints() = "<<vtkFloatingPointType(theDataSet->GetPoints()->GetActualMemorySize()*1000)<<endl);
557     BEGMSG(1,"GetCells() = "<<vtkFloatingPointType(theDataSet->GetCells()->GetActualMemorySize()*1000)<<endl);
558     BEGMSG(1,"GetCellTypesArray() = "<<vtkFloatingPointType(theDataSet->GetCellTypesArray()->GetActualMemorySize()*1000)<<endl);
559     BEGMSG(1,"GetCellLocationsArray() = "<<vtkFloatingPointType(theDataSet->GetCellLocationsArray()->GetActualMemorySize()*1000)<<endl);
560     theDataSet->BuildLinks();
561     BEGMSG(1,"GetCellLinks() = "<<vtkFloatingPointType(theDataSet->GetCellLinks()->GetActualMemorySize()*1000)<<endl);
562     BEGMSG(1,"GetPointData() = "<<vtkFloatingPointType(theDataSet->GetPointData()->GetActualMemorySize()*1000)<<endl);
563     BEGMSG(1,"GetCellData() = "<<vtkFloatingPointType(theDataSet->GetCellData()->GetActualMemorySize()*1000)<<endl);
564     BEGMSG(1,"GetActualMemorySize() = "<<vtkFloatingPointType(theDataSet->GetActualMemorySize()*1000)<<endl);
565   }
566 }
567
568
569 //---------------------------------------------------------------
570 VISU_Convertor_impl
571 ::VISU_Convertor_impl()
572 {}
573
574
575 //---------------------------------------------------------------
576 VISU_Convertor_impl
577 ::~VISU_Convertor_impl() 
578 {}
579
580
581 //---------------------------------------------------------------
582 VISU_Convertor* 
583 VISU_Convertor_impl
584 ::Build() 
585
586   if(!myIsDone){ 
587     myIsDone = true;  
588     BuildEntities();
589     BuildFields();
590     BuildMinMax();
591     BuildGroups();
592   }
593   return this;
594 }
595
596 VISU_Convertor* 
597 VISU_Convertor_impl
598 ::BuildEntities() 
599
600   return this;
601 }
602
603 VISU_Convertor* 
604 VISU_Convertor_impl
605 ::BuildFields() 
606
607   BuildEntities();
608   return this;
609 }
610
611 VISU_Convertor* 
612 VISU_Convertor_impl
613 ::BuildMinMax() 
614
615   BuildFields();
616   return this;
617 }
618
619 VISU_Convertor* 
620 VISU_Convertor_impl
621 ::BuildGroups() 
622
623   return this;
624 }
625
626
627 //---------------------------------------------------------------
628 VISU::PNamedIDMapper 
629 VISU_Convertor_impl
630 ::GetMeshOnEntity(const std::string& theMeshName, 
631                   const VISU::TEntity& theEntity)
632 {
633   INITMSG(MYDEBUG,"GetMeshOnEntity"<<
634           "; theMeshName = '"<<theMeshName<<"'"<<
635           "; theEntity = "<<theEntity<<
636           endl);
637
638   //Cheching possibility do the query
639   TFindMeshOnEntity aFindMeshOnEntity = 
640     FindMeshOnEntity(theMeshName,theEntity);
641   
642   VISU::PMeshImpl aMesh = boost::get<0>(aFindMeshOnEntity);;
643   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindMeshOnEntity);
644   
645   //Main part of code
646 #ifndef _DEXCEPT_
647   try{
648 #endif
649     if(!aMeshOnEntity->myIsVTKDone){
650       VISU::TTimerLog aTimerLog(MYDEBUG,"VISU_Convertor_impl::GetMeshOnEntity");
651       const VISU::PAppendFilter& anAppendFilter = aMeshOnEntity->GetFilter();
652       if(MYVTKDEBUG) anAppendFilter->DebugOn();
653
654       LoadMeshOnEntity(aMesh,aMeshOnEntity);
655       anAppendFilter->SetSharedPointSet(aMesh->GetPointSet());
656       
657       const VISU::TGeom2SubMesh& aGeom2SubMesh = aMeshOnEntity->myGeom2SubMesh;
658       VISU::TGeom2SubMesh::const_iterator anIter = aGeom2SubMesh.begin();
659
660       VISU::TID2ID& anElemObj2VTKID = aMeshOnEntity->myElemObj2VTKID;
661       VISU::TSubMeshArr& aSubMeshArr = aMeshOnEntity->mySubMeshArr;
662       aSubMeshArr.resize(aGeom2SubMesh.size());
663
664       for(vtkIdType anID = 0, aCellID = 0; anIter != aGeom2SubMesh.end(); anIter++, anID++){
665         VISU::EGeometry aEGeom = anIter->first;
666         vtkIdType aVGeom = VISUGeom2VTK(aEGeom);
667         VISU::PSubMeshImpl aSubMesh = anIter->second;
668
669         //ENK: 23.11.2006 - PAL13176 - EDF228 VISU : Enhancement of structured datas processing
670         aSubMesh->myIsStructured = aMesh->myIsStructured;
671         aSubMesh->myType = aMesh->myType;
672         aSubMesh->myGrilleStructure = aMesh->myGrilleStructure;
673         aSubMesh->myObjID2StructureID = aMesh->myObjID2StructureID;
674         //ENK: 23.11.2006
675
676         const VISU::PUnstructuredGrid& aSource = aSubMesh->GetSource();
677         aSource->ShallowCopy(aMesh->GetPointSet());
678
679         aSubMesh->myStartID = aCellID;
680         GetCellsOnSubMesh(aSource, aMeshOnEntity, aSubMesh, aVGeom);
681         anAppendFilter->AddInput(aSource.GetPointer());
682         
683         vtkIdType aNbCells = aSource->GetNumberOfCells();
684         for(vtkIdType aCell = 0; aCell < aNbCells; aCell++, aCellID++){
685           vtkIdType anObjID = aSubMesh->GetElemObjID(aCell);
686           anElemObj2VTKID[anObjID] = aCellID;
687         }
688
689         aSubMeshArr[anID] = aSubMesh;
690       }
691       aMeshOnEntity->myNamedPointCoords = aMesh->myNamedPointCoords;
692       
693       //ENK: 23.11.2006 - PAL13176 - EDF228 VISU : Enhancement of structured datas processing
694       aMeshOnEntity->myIsStructured = aMesh->myIsStructured;
695       aMeshOnEntity->myType = aMesh->myType;
696       aMeshOnEntity->myGrilleStructure = aMesh->myGrilleStructure;
697       aMeshOnEntity->myObjID2StructureID = aMesh->myObjID2StructureID;
698       //ENK: 23.11.2006
699
700       aMeshOnEntity->myIsVTKDone = true;
701
702       if(MYDEBUGWITHFILES){
703         std::string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
704         std::string aFileName = std::string(getenv("HOME"))+"/"+getenv("USER")+"-";
705         aFileName += aMeshName + dtos("-%d-",int(theEntity)) + "-MeshOnEntity.vtk";
706         VISU::WriteToFile(anAppendFilter->GetOutput(),aFileName);
707       }
708
709       if(MYVTKDEBUG){
710         GetMeshOnEntitySize(theMeshName,theEntity);
711         PrintMemorySize(anAppendFilter->GetOutput());
712       }
713     }
714
715 #ifndef _DEXCEPT_
716   }catch(...){
717     throw;
718   }
719 #endif
720
721   return aMeshOnEntity;
722 }
723
724
725 //---------------------------------------------------------------
726 VISU::PUnstructuredGridIDMapper 
727 VISU_Convertor_impl
728 ::GetFamilyOnEntity(const std::string& theMeshName, 
729                     const VISU::TEntity& theEntity,
730                     const std::string& theFamilyName)
731 {
732   INITMSG(MYDEBUG,"GetFamilyOnEntity"<<
733           "; theMeshName = '"<<theMeshName<<"'"<<
734           "; theEntity = "<<theEntity<<
735           "; theFamilyName = '"<<theFamilyName<<"'"<<
736           endl);
737
738   //Cheching possibility do the query
739   TFindFamilyOnEntity aFindFamilyOnEntity = 
740     FindFamilyOnEntity(theMeshName,theEntity,theFamilyName);
741
742   VISU::PMeshImpl aMesh = boost::get<0>(aFindFamilyOnEntity);;
743   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindFamilyOnEntity);
744   VISU::PFamilyImpl aFamily = boost::get<2>(aFindFamilyOnEntity);
745
746   //Main part of code
747 #ifndef _DEXCEPT_
748   try{
749 #endif
750     if(!aFamily->myIsVTKDone){
751       const VISU::PUnstructuredGrid& aSource = aFamily->GetSource();
752       if(MYVTKDEBUG) aSource->DebugOn();
753
754       GetMeshOnEntity(theMeshName,theEntity);
755
756       LoadFamilyOnEntity(aMesh,aMeshOnEntity,aFamily);
757       aSource->ShallowCopy(aMesh->GetPointSet());
758       GetCellsOnFamily(aSource,aMeshOnEntity,aFamily);
759
760       aFamily->myNamedPointCoords = aMesh->myNamedPointCoords;
761       //ENK: 23.11.2006 - PAL13176 - EDF228 VISU : Enhancement of structured datas processing
762       aFamily->myIsStructured = aMesh->myIsStructured;
763       aFamily->myType = aMesh->myType;
764       aFamily->myGrilleStructure = aMesh->myGrilleStructure;
765       aFamily->myObjID2StructureID = aMesh->myObjID2StructureID;
766       //ENK: 23.11.2006
767       aFamily->myIsVTKDone = true;
768
769       if(MYDEBUGWITHFILES){
770         std::string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
771         std::string aFamilyName = QString(theFamilyName.c_str()).simplifyWhiteSpace().latin1();
772         std::string aFileName = std::string(getenv("HOME"))+"/"+getenv("USER")+"-";
773         aFileName += aMeshName + dtos("-%d-",int(theEntity)) + aFamilyName + "-FamilyOnEntity.vtk";
774         VISU::WriteToFile(aSource.GetPointer(),aFileName);
775       }
776
777       if(MYVTKDEBUG){
778         GetFamilyOnEntitySize(theMeshName,theEntity,theFamilyName);
779         PrintMemorySize(aSource.GetPointer());
780       }
781     }
782
783 #ifndef _DEXCEPT_
784   }catch(...){
785     throw;
786   }
787 #endif
788
789   return aFamily;
790 }
791
792
793 //---------------------------------------------------------------
794 VISU::PUnstructuredGridIDMapper 
795 VISU_Convertor_impl
796 ::GetMeshOnGroup(const std::string& theMeshName, 
797                  const std::string& theGroupName)
798 {
799   INITMSG(MYDEBUG,"GetMeshOnGroup\n");
800   INITMSGA(MYDEBUG,0,
801            "- theMeshName = '"<<theMeshName<<
802            "'; theGroupName = '"<<theGroupName<<"'"<<
803            endl);
804
805   //Cheching possibility do the query
806   TFindMeshOnGroup aFindMeshOnGroup = FindMeshOnGroup(theMeshName,theGroupName);
807   VISU::PMeshImpl aMesh = boost::get<0>(aFindMeshOnGroup);
808   VISU::PGroupImpl aGroup = boost::get<1>(aFindMeshOnGroup);
809
810   //Main part of code
811 #ifndef _DEXCEPT_
812   try{
813 #endif
814     if(!aGroup->myIsVTKDone){
815       const VISU::PAppendFilter& anAppendFilter = aGroup->GetFilter();
816       const VISU::TFamilySet& aFamilySet = aGroup->myFamilySet;
817
818       LoadMeshOnGroup(aMesh,aFamilySet);
819       anAppendFilter->SetSharedPointSet(aMesh->GetPointSet());
820
821       VISU::TFamilySet::const_iterator anIter = aFamilySet.begin();
822
823       VISU::TID2ID& anElemObj2VTKID = aGroup->myElemObj2VTKID;
824       VISU::TFamilyArr& aFamilyArr = aGroup->myFamilyArr;
825       aFamilyArr.resize(aFamilySet.size());
826
827       for(vtkIdType anID = 0; anIter != aFamilySet.end(); anIter++, anID++){
828         VISU::PFamilyImpl aFamily = *anIter;
829         const std::string& aFamilyName = aFamily->myName;
830         const VISU::TEntity& anEntity = aFamily->myEntity;
831
832         VISU::PIDMapper anIDMapper = GetFamilyOnEntity(theMeshName,anEntity,aFamilyName);
833         vtkDataSet* anOutput = anIDMapper->GetOutput();
834         anAppendFilter->AddInput(anOutput);
835
836         vtkIdType aStartID = anElemObj2VTKID.size();
837         vtkIdType aNbCells = anOutput->GetNumberOfCells();
838         for(vtkIdType aCellID = 0; aCellID < aNbCells; aCellID++){
839           anElemObj2VTKID[aFamily->GetElemObjID(aCellID)] = aStartID + aCellID;
840         }
841         aFamilyArr[anID] = aFamily;
842       }
843       aGroup->myNamedPointCoords = aMesh->myNamedPointCoords;
844       //ENK: 23.11.2006 - PAL13176 - EDF228 VISU : Enhancement of structured datas processing
845       aGroup->myIsStructured = aMesh->myIsStructured;
846       aGroup->myType = aMesh->myType;
847       aGroup->myGrilleStructure = aMesh->myGrilleStructure;
848       aGroup->myObjID2StructureID = aMesh->myObjID2StructureID;
849       //ENK: 23.11.2006
850       aGroup->myIsVTKDone = true;
851
852       if(MYDEBUGWITHFILES){
853         std::string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
854         std::string aGroupName = QString(theGroupName.c_str()).simplifyWhiteSpace().latin1();
855         std::string aFileName = std::string(getenv("HOME"))+"/"+getenv("USER")+"-";
856         aFileName += aMeshName + "-" + aGroupName + "-MeshOnGroup.vtk";
857         VISU::WriteToFile(anAppendFilter->GetOutput(),aFileName);
858       }
859     }
860 #ifndef _DEXCEPT_
861   }catch(...){
862     throw;
863   }
864 #endif
865
866   return aGroup;
867 }
868
869
870 //---------------------------------------------------------------
871 vtkUnstructuredGrid*
872 VISU_Convertor_impl
873 ::GetTimeStampOnProfile(const VISU::PMeshImpl& theMesh,
874                         const VISU::PMeshOnEntityImpl& theMeshOnEntity,
875                         const VISU::PFieldImpl& theField,
876                         const VISU::PValForTimeImpl& theValForTime,
877                         const VISU::PUnstructuredGridIDMapperImpl& theUnstructuredGridIDMapper,
878                         const VISU::PProfileImpl& theProfile,
879                         const VISU::TEntity& theEntity)
880 {
881   vtkUnstructuredGrid* anOutput = NULL;
882   
883   LoadMeshOnEntity(theMesh, theMeshOnEntity);
884   GetMeshOnEntity(theMeshOnEntity->myMeshName, theMeshOnEntity->myEntity);
885   GetMeshOnProfile(theMesh, theMeshOnEntity, theProfile);
886   
887   bool isNeedInCells = false;
888   
889   theUnstructuredGridIDMapper->myIDMapper = theProfile;
890   if(theMeshOnEntity->myEntity == VISU::NODE_ENTITY){
891     // add geometry elements to output,
892     // if timestamp on NODE_ENTITY and
893     // on profiles with status eAddPart
894     VISU::TGeom2SubProfile::const_iterator anIter = theProfile->myGeom2SubProfile.begin();
895     for(; anIter != (theProfile->myGeom2SubProfile).end(); anIter++){
896       const VISU::EGeometry aGeom = anIter->first;
897       const VISU::PSubProfileImpl aSubProfile = anIter->second;
898       if(aSubProfile->myStatus == VISU::eAddPart && aGeom == VISU::ePOINT1){
899         isNeedInCells = true;
900         break;
901       }
902     }
903     if(isNeedInCells){
904       theUnstructuredGridIDMapper->myIsSpecialKey = true;
905       GetTimeStampOnNodalProfile(theMesh,theUnstructuredGridIDMapper,theField,theValForTime,theEntity);
906       anOutput = theUnstructuredGridIDMapper->GetUnstructuredGridOutput();
907     }
908   }
909   if(!isNeedInCells){
910     vtkUnstructuredGrid* anOutput = theUnstructuredGridIDMapper->GetUnstructuredGridOutput();
911     const VISU::PUnstructuredGrid& aSource = theUnstructuredGridIDMapper->mySource.GetSource();
912     VISU::GetTimeStampOnProfile(aSource, theField, theValForTime, theEntity);
913   }
914   return anOutput;
915 }
916
917 void
918 VISU_Convertor_impl
919 ::GetTimeStampOnNodalProfile(const VISU::PMeshImpl& theMesh,
920                              const VISU::PUnstructuredGridIDMapperImpl& theIDMapperFilter,
921                              const VISU::PFieldImpl& theField, 
922                              const VISU::PValForTimeImpl& theValForTime,
923                              const VISU::TEntity& theEntity)
924 {
925   const VISU::TMeshOnEntityMap& aMeshOnEntityMap = theMesh->myMeshOnEntityMap;
926   VISU::TMeshOnEntityMap::const_iterator aIter = aMeshOnEntityMap.begin();
927   for(;aIter!=aMeshOnEntityMap.end();aIter++){
928     VISU::TEntity aEntity = aIter->first;
929     if(aEntity != VISU::NODE_ENTITY){
930       VISU::PNamedIDMapper aMapper = GetMeshOnEntity(theMesh->myName,aEntity);
931       if(aMapper)
932         theIDMapperFilter->myMappers[aEntity] = aMapper;
933     }
934   }
935
936   theIDMapperFilter->GetUnstructuredGridOutput();
937   const VISU::PUnstructuredGrid& aSource = theIDMapperFilter->mySource.GetSource();
938   VISU::GetTimeStampOnProfile(aSource,theField,theValForTime,theEntity);
939 }
940
941 //---------------------------------------------------------------
942 VISU::PUnstructuredGridIDMapper 
943 VISU_Convertor_impl
944 ::GetTimeStampOnMesh(const std::string& theMeshName, 
945                      const VISU::TEntity& theEntity,
946                      const std::string& theFieldName,
947                      int theStampsNum)
948 {
949   INITMSG(MYDEBUG,"GetTimeStampOnMesh"<<
950           "; theMeshName = '"<<theMeshName<<"'"<<
951           "; theEntity = "<<theEntity<<
952           "; theFieldName = '"<<theFieldName<<"'"<<
953           "; theStampsNum = "<<theStampsNum<<
954           endl);
955
956   //Cheching possibility do the query
957   TFindTimeStamp aFindTimeStamp = FindTimeStamp(theMeshName,
958                                                 theEntity,
959                                                 theFieldName,
960                                                 theStampsNum);
961
962   VISU::PMeshImpl aMesh = boost::get<0>(aFindTimeStamp);
963   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindTimeStamp);
964   VISU::PMeshOnEntityImpl aVTKMeshOnEntity = boost::get<2>(aFindTimeStamp);
965   VISU::PValForTimeImpl aValForTime = boost::get<4>(aFindTimeStamp);
966   VISU::PFieldImpl aField = boost::get<3>(aFindTimeStamp);
967
968   //Main part of code
969   VISU::PUnstructuredGridIDMapperImpl anUnstructuredGridIDMapper = aValForTime->myUnstructuredGridIDMapper;
970 #ifndef _DEXCEPT_
971   try{
972 #endif
973     if(!anUnstructuredGridIDMapper->myIsVTKDone){
974       VISU::TTimerLog aTimerLog(MYDEBUG,"VISU_Convertor_impl::GetTimeStampOnMesh");
975       LoadValForTimeOnMesh(aMesh, aMeshOnEntity, aField, aValForTime);
976
977       vtkUnstructuredGrid* anOutput = NULL;
978       try{
979         anOutput = GetTimeStampOnProfile(aMesh,
980                                          aVTKMeshOnEntity,
981                                          aField,
982                                          aValForTime,
983                                          anUnstructuredGridIDMapper,
984                                          aValForTime->myProfile,
985                                          aMeshOnEntity->myEntity);
986       }catch(std::exception& exc){
987         MSG(MYDEBUG,"Follow exception was occured :\n"<<exc.what());
988         anOutput = GetTimeStampOnProfile(aMesh,
989                                          aMeshOnEntity,
990                                          aField,
991                                          aValForTime,
992                                          anUnstructuredGridIDMapper,
993                                          aValForTime->myProfile,
994                                          aVTKMeshOnEntity->myEntity);
995       }
996
997       //ENK: 23.11.2006 - PAL13176 - EDF228 VISU : Enhancement of structured datas processing
998       anUnstructuredGridIDMapper->myIsStructured = aMesh->myIsStructured;
999       anUnstructuredGridIDMapper->myType = aMesh->myType;
1000       anUnstructuredGridIDMapper->myGrilleStructure = aMesh->myGrilleStructure;
1001       anUnstructuredGridIDMapper->myObjID2StructureID = aMesh->myObjID2StructureID;
1002       //ENK: 23.11.2006
1003
1004       anUnstructuredGridIDMapper->myIsVTKDone = true;
1005
1006       if(MYDEBUGWITHFILES){
1007         std::string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
1008         std::string aFieldName = QString(theFieldName.c_str()).simplifyWhiteSpace().latin1();
1009         std::string aPrefix = std::string(getenv("HOME"))+"/"+getenv("USER")+"-";
1010         std::string aFileName = aPrefix + aMeshName + dtos("-%d-",int(theEntity)) + 
1011           aFieldName + dtos("-%d", theStampsNum) + "-TimeStampOnMesh.vtk";
1012         VISU::WriteToFile(anOutput,aFileName);
1013       }
1014       if(MYVTKDEBUG){
1015         GetTimeStampSize(theMeshName, theEntity, theFieldName, theStampsNum);
1016         anOutput->Update();
1017         if(theEntity == VISU::NODE_ENTITY)
1018           BEGMSG(MYVTKDEBUG,"GetPointData() = "<<vtkFloatingPointType(anOutput->GetPointData()->GetActualMemorySize()*1000)<<endl);
1019         else
1020           BEGMSG(MYVTKDEBUG,"GetCellData() = "<<vtkFloatingPointType(anOutput->GetCellData()->GetActualMemorySize()*1000)<<endl);
1021         BEGMSG(MYVTKDEBUG,"GetActualMemorySize() = "<<vtkFloatingPointType(anOutput->GetActualMemorySize()*1000)<<endl);
1022       }
1023       }
1024 #ifndef _DEXCEPT_
1025     }catch(std::exception& exc){
1026     throw;
1027   }catch(...){
1028     throw;
1029   }
1030 #endif
1031
1032   return anUnstructuredGridIDMapper;
1033 }
1034
1035
1036 //---------------------------------------------------------------
1037 VISU::PGaussPtsIDMapper 
1038 VISU_Convertor_impl
1039 ::GetTimeStampOnGaussPts(const std::string& theMeshName, 
1040                          const VISU::TEntity& theEntity,
1041                          const std::string& theFieldName,
1042                          int theStampsNum)
1043 {
1044   INITMSG(MYDEBUG,"GetTimeStampOnGaussPts"<<
1045           "; theMeshName = '"<<theMeshName<<"'"<<
1046           "; theEntity = "<<theEntity<<
1047           "; theFieldName = '"<<theFieldName<<"'"<<
1048           "; theStampsNum = "<<theStampsNum<<
1049           endl);
1050
1051   if(theEntity == VISU::NODE_ENTITY)
1052     EXCEPTION(std::runtime_error, "It is impossible to reate Gauss Points on NODE_ENTITY !!!");
1053
1054   //Cheching possibility do the query
1055   TFindTimeStamp aFindTimeStamp = FindTimeStamp(theMeshName,
1056                                                 theEntity,
1057                                                 theFieldName,
1058                                                 theStampsNum);
1059   
1060   VISU::PMeshImpl aMesh = boost::get<0>(aFindTimeStamp);
1061   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindTimeStamp);
1062   VISU::PMeshOnEntityImpl aVTKMeshOnEntity = aMeshOnEntity;
1063   VISU::PValForTimeImpl aValForTime = boost::get<4>(aFindTimeStamp);
1064   VISU::PFieldImpl aField = boost::get<3>(aFindTimeStamp);
1065
1066   //Main part of code
1067   VISU::PGaussPtsIDFilter aGaussPtsIDFilter = aValForTime->myGaussPtsIDFilter;
1068 #ifndef _DEXCEPT_
1069   try{
1070 #endif
1071     if(!aGaussPtsIDFilter->myIsVTKDone){
1072       VISU::TTimerLog aTimerLog(MYDEBUG,"VISU_Convertor_impl::GetTimeStampOnGaussPts");
1073       LoadValForTimeOnGaussPts(aMesh, aMeshOnEntity, aField, aValForTime);
1074       
1075       GetMeshOnEntity(aVTKMeshOnEntity->myMeshName, aVTKMeshOnEntity->myEntity);
1076       
1077       VISU::PProfileImpl aProfile = aValForTime->myProfile;
1078       GetMeshOnProfile(aMesh, aVTKMeshOnEntity, aProfile);
1079
1080       VISU::PGaussMeshImpl aGaussMesh = aValForTime->myGaussMesh;
1081       if(!aGaussMesh->myIsVTKDone){
1082         BuildGaussMesh(aMesh, aVTKMeshOnEntity, aGaussMesh);
1083         aGaussMesh->myParent = aProfile.get();
1084         aGaussMesh->myIsVTKDone = true;
1085       }
1086
1087       aGaussPtsIDFilter->myIDMapper = aGaussMesh;
1088       aGaussPtsIDFilter->myGaussPtsIDMapper = aGaussMesh;
1089       vtkPolyData* anOutput = aGaussPtsIDFilter->GetPolyDataOutput();
1090       const VISU::PPolyData& aSource = aGaussPtsIDFilter->mySource.GetSource();
1091       VISU::GetTimeStampOnGaussMesh(aSource, aField, aValForTime);
1092
1093       aGaussPtsIDFilter->myIsVTKDone = true;
1094
1095       if(MYDEBUGWITHFILES){
1096         std::string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
1097         std::string aFieldName = QString(theFieldName.c_str()).simplifyWhiteSpace().latin1();
1098         std::string aPrefix = std::string(getenv("HOME"))+"/"+getenv("USER")+"-";
1099         std::string aFileName = aPrefix + aMeshName + dtos("-%d-",int(theEntity)) + 
1100           aFieldName + dtos("-%d",theStampsNum) + "-TimeStampOnGaussPts.vtk";
1101         VISU::WriteToFile(anOutput, aFileName);
1102       }
1103       if(MYVTKDEBUG){
1104         GetTimeStampSize(theMeshName, theEntity, theFieldName, theStampsNum);
1105         anOutput->Update();
1106         if(theEntity == VISU::NODE_ENTITY)
1107           BEGMSG(MYVTKDEBUG,"GetPointData() = "<<vtkFloatingPointType(anOutput->GetPointData()->GetActualMemorySize()*1000)<<endl);
1108         else
1109           BEGMSG(MYVTKDEBUG,"GetCellData() = "<<vtkFloatingPointType(anOutput->GetCellData()->GetActualMemorySize()*1000)<<endl);
1110         BEGMSG(MYVTKDEBUG,"GetActualMemorySize() = "<<vtkFloatingPointType(anOutput->GetActualMemorySize()*1000)<<endl);
1111       }
1112     }
1113 #ifndef _DEXCEPT_
1114   }catch(std::exception& exc){
1115     throw;
1116   }catch(...){
1117     throw;
1118   }
1119 #endif
1120
1121   return aGaussPtsIDFilter;
1122 }
1123
1124 //---------------------------------------------------------------
1125 VISU::PMeshImpl 
1126 VISU_Convertor_impl
1127 ::FindMesh(const std::string& theMeshName)
1128 {
1129   GetMeshMap();
1130   VISU::TMeshMap::iterator aMeshMapIter = myMeshMap.find(theMeshName);
1131   if(aMeshMapIter == myMeshMap.end())
1132     EXCEPTION(std::runtime_error,"FindMesh >> There is no mesh with the name - '"<<theMeshName<<"'!!!");
1133
1134   VISU::PMeshImpl aMesh = aMeshMapIter->second;
1135   return aMesh;
1136 }
1137
1138
1139 //---------------------------------------------------------------
1140 VISU_Convertor_impl::TFindMeshOnEntity
1141 VISU_Convertor_impl
1142 ::FindMeshOnEntity(const std::string& theMeshName,
1143                    const VISU::TEntity& theEntity)
1144 {
1145   VISU::PMeshImpl aMesh = FindMesh(theMeshName);
1146   VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh->myMeshOnEntityMap;
1147   VISU::TMeshOnEntityMap::const_iterator aMeshOnEntityMapIter = aMeshOnEntityMap.find(theEntity);
1148   if(aMeshOnEntityMapIter == aMeshOnEntityMap.end())
1149     EXCEPTION(std::runtime_error,"FindMeshOnEntity >> There is no mesh on the entity - "<<theEntity<<"!!!");
1150
1151   VISU::PMeshOnEntityImpl aMeshOnEntity = aMeshOnEntityMapIter->second;
1152   
1153   return TFindMeshOnEntity(aMesh,
1154                            aMeshOnEntity);
1155 }
1156
1157
1158 //---------------------------------------------------------------
1159 VISU_Convertor_impl::TFindFamilyOnEntity
1160 VISU_Convertor_impl
1161 ::FindFamilyOnEntity(const std::string& theMeshName,
1162                      const VISU::TEntity& theEntity,
1163                      const std::string& theFamilyName)
1164 {
1165   if(theFamilyName != ""){
1166     VISU::PMeshImpl aMesh = FindMesh(theMeshName);
1167     VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh->myMeshOnEntityMap;
1168     VISU::TMeshOnEntityMap::const_iterator aMeshOnEntityMapIter = aMeshOnEntityMap.find(theEntity);
1169     if(aMeshOnEntityMapIter == aMeshOnEntityMap.end())
1170       EXCEPTION(std::runtime_error,"FindFamilyOnEntity >> There is no mesh on the entity - "<<theEntity<<"!!!");
1171
1172     VISU::PMeshOnEntityImpl aMeshOnEntity = aMeshOnEntityMapIter->second;
1173
1174     VISU::TFamilyMap& aFamilyMap = aMeshOnEntity->myFamilyMap;
1175     VISU::TFamilyMap::iterator aFamilyMapIter = aFamilyMap.find(theFamilyName);
1176     if(aFamilyMapIter != aFamilyMap.end()){
1177       const VISU::PFamily& aFamily = aFamilyMapIter->second;
1178       return TFindFamilyOnEntity(aMesh,
1179                                  aMeshOnEntity,
1180                                  aFamily);
1181     }
1182   }
1183   return TFindFamilyOnEntity();
1184 }
1185
1186
1187 //---------------------------------------------------------------
1188 size_t
1189 VISU_Convertor_impl
1190 ::GetSize() 
1191 {
1192   size_t aResult = 0;
1193   const VISU::TMeshMap& aMeshMap = GetMeshMap();
1194   VISU::TMeshMap::const_iterator aMeshMapIter = aMeshMap.begin();
1195   for(; aMeshMapIter != aMeshMap.end(); aMeshMapIter++){
1196     const std::string& aMeshName = aMeshMapIter->first;
1197     const VISU::PMesh aMesh = aMeshMapIter->second;
1198     const VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh->myMeshOnEntityMap;
1199     VISU::TMeshOnEntityMap::const_iterator aMeshOnEntityMapIter;
1200     //Import fields
1201     aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
1202     for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
1203       const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
1204       const VISU::PMeshOnEntity aMeshOnEntity = aMeshOnEntityMapIter->second;
1205       const VISU::TFieldMap& aFieldMap = aMeshOnEntity->myFieldMap;
1206       VISU::TFieldMap::const_iterator aFieldMapIter = aFieldMap.begin();
1207       for(; aFieldMapIter != aFieldMap.end(); aFieldMapIter++){
1208         const std::string& aFieldName = aFieldMapIter->first;
1209         const VISU::PField aField = aFieldMapIter->second;
1210         const VISU::TValField& aValField = aField->myValField;
1211         VISU::TValField::const_iterator aValFieldIter = aValField.begin();
1212         for(; aValFieldIter != aValField.end(); aValFieldIter++){
1213           int aTimeStamp = aValFieldIter->first;
1214           aResult += GetTimeStampSize(aMeshName,anEntity,aFieldName,aTimeStamp);
1215         }
1216       }
1217       //Importing groups
1218       const VISU::TGroupMap& aGroupMap = aMesh->myGroupMap;
1219       VISU::TGroupMap::const_iterator aGroupMapIter = aGroupMap.begin();
1220       for(; aGroupMapIter != aGroupMap.end(); aGroupMapIter++){
1221         const std::string& aGroupName = aGroupMapIter->first;
1222         aResult += GetMeshOnGroupSize(aMeshName,aGroupName);
1223       }
1224       //Import families
1225       const VISU::TFamilyMap& aFamilyMap = aMeshOnEntity->myFamilyMap;
1226       VISU::TFamilyMap::const_iterator aFamilyMapIter = aFamilyMap.begin();
1227       for(; aFamilyMapIter != aFamilyMap.end(); aFamilyMapIter++){
1228         const std::string& aFamilyName = aFamilyMapIter->first;
1229         aResult += GetFamilyOnEntitySize(aMeshName,anEntity,aFamilyName);
1230       }
1231       //Import mesh on entity
1232       aResult += GetMeshOnEntitySize(aMeshName,anEntity);
1233     }
1234   }
1235   MSG(MYDEBUG,"GetSize - aResult = "<<vtkFloatingPointType(aResult));
1236   return aResult;
1237 }
1238
1239
1240 //---------------------------------------------------------------
1241 size_t
1242 VISU_Convertor_impl
1243 ::GetMeshOnEntitySize(const std::string& theMeshName, 
1244                       const VISU::TEntity& theEntity)
1245 {
1246   TFindMeshOnEntity aFindMeshOnEntity = 
1247     FindMeshOnEntity(theMeshName, theEntity);
1248
1249   VISU::PMeshImpl aMesh = boost::get<0>(aFindMeshOnEntity);
1250   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindMeshOnEntity);
1251
1252   size_t aPointsSize = 3*aMesh->GetNbPoints()*sizeof(VISU::TCoord);
1253   size_t aNbCells = aMeshOnEntity->myNbCells;
1254   size_t aCellsSize = aMeshOnEntity->myCellsSize;
1255
1256   size_t aConnectivitySize = aCellsSize*sizeof(vtkIdType);
1257   size_t aTypesSize = aNbCells*sizeof(char);
1258   size_t aLocationsSize = aNbCells*sizeof(int);
1259   vtkFloatingPointType aNbCellsPerPoint = aCellsSize / aNbCells - 1;
1260   size_t aLinksSize = aMesh->GetNbPoints() * 
1261     (vtkIdType(sizeof(vtkIdType)*aNbCellsPerPoint) + sizeof(vtkCellLinks::Link));
1262   aLinksSize = 0;
1263   size_t aResult = aPointsSize + aConnectivitySize + aTypesSize + aLocationsSize + aLinksSize;
1264
1265   MSG(MYDEBUG,"GetMeshOnEntitySize "<<
1266       "- aResult = "<<vtkFloatingPointType(aResult)<<
1267       "; theMeshName = '"<<theMeshName<<"'"<<
1268       "; theEntity = "<<theEntity);
1269   if(MYDEBUG){
1270     INITMSG(MYVTKDEBUG,"- aPointsSize = "<<vtkFloatingPointType(aPointsSize)<<"\n");
1271     BEGMSG(MYVTKDEBUG,"- aConnectivitySize = "<<vtkFloatingPointType(aConnectivitySize)<<"\n");
1272     BEGMSG(MYVTKDEBUG,"- aTypesSize = "<<vtkFloatingPointType(aTypesSize)<<"\n");
1273     BEGMSG(MYVTKDEBUG,"- aLocationsSize = "<<vtkFloatingPointType(aLocationsSize)<<"\n");
1274     BEGMSG(MYVTKDEBUG,"- aLinksSize = "<<vtkFloatingPointType(aLinksSize)<<"\n");
1275   }
1276
1277   aResult = size_t(aResult*ERR_SIZE_CALC);
1278   return aResult;
1279 }
1280
1281
1282 //---------------------------------------------------------------
1283 size_t
1284 VISU_Convertor_impl
1285 ::GetFamilyOnEntitySize(const std::string& theMeshName, 
1286                         const VISU::TEntity& theEntity,
1287                         const std::string& theFamilyName)
1288 {
1289   TFindFamilyOnEntity aFindFamilyOnEntity = 
1290     FindFamilyOnEntity(theMeshName,theEntity,theFamilyName);
1291   VISU::PMeshImpl aMesh = boost::get<0>(aFindFamilyOnEntity);
1292   VISU::PFamilyImpl aFamily = boost::get<2>(aFindFamilyOnEntity);
1293   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindFamilyOnEntity);
1294
1295   size_t aPointsSize = 3*aMesh->GetNbPoints()*sizeof(VISU::TCoord);
1296   size_t aNbCells = aFamily->myNbCells;
1297   size_t aCellsSize = aFamily->myCellsSize;
1298
1299   size_t aConnectivitySize = aCellsSize*sizeof(vtkIdType);
1300   size_t aTypesSize = aNbCells*sizeof(char);
1301   size_t aLocationsSize = aNbCells*sizeof(int);
1302   vtkFloatingPointType aNbCellsPerPoint = aCellsSize / aNbCells - 1;
1303   size_t aLinksSize = aMesh->GetNbPoints() * 
1304     (vtkIdType(sizeof(vtkIdType)*aNbCellsPerPoint) + sizeof(vtkCellLinks::Link));
1305   aLinksSize = 0;
1306   size_t aResult = aPointsSize + aConnectivitySize + aTypesSize + aLocationsSize + aLinksSize;
1307
1308   MSG(MYDEBUG,"GetFamilyOnEntitySize "<<
1309       "- aResult = "<<vtkFloatingPointType(aResult)<<
1310       "; theMeshName = '"<<theMeshName<<"'"<<
1311       "; theEntity = "<<theEntity<<
1312       "; theFamilyName = '"<<theFamilyName<<"'");
1313   if(MYDEBUG){
1314     INITMSG(MYVTKDEBUG,"- aPointsSize = "<<vtkFloatingPointType(aPointsSize)<<"\n");
1315     BEGMSG(MYVTKDEBUG,"- aConnectivitySize = "<<vtkFloatingPointType(aConnectivitySize)<<"\n");
1316     BEGMSG(MYVTKDEBUG,"- aTypesSize = "<<vtkFloatingPointType(aTypesSize)<<"\n");
1317     BEGMSG(MYVTKDEBUG,"- aLocationsSize = "<<vtkFloatingPointType(aLocationsSize)<<"\n");
1318     BEGMSG(MYVTKDEBUG,"- aLinksSize = "<<vtkFloatingPointType(aLinksSize)<<"\n");
1319   }
1320
1321   aResult = size_t(aResult*ERR_SIZE_CALC);
1322   return aResult;
1323 }
1324
1325
1326 //---------------------------------------------------------------
1327 VISU_Convertor_impl::TFindMeshOnGroup
1328 VISU_Convertor_impl
1329 ::FindMeshOnGroup(const std::string& theMeshName, 
1330                   const std::string& theGroupName)
1331 {
1332   VISU::PMeshImpl aMesh = FindMesh(theMeshName);
1333   VISU::TGroupMap& aGroupMap = aMesh->myGroupMap;
1334   VISU::TGroupMap::iterator aGroupMapIter = aGroupMap.find(theGroupName);
1335   if(aGroupMapIter == aGroupMap.end())
1336     EXCEPTION(std::runtime_error,"FindMesh >> There is no the group in the mesh!!! - '"<<theGroupName<<"'");
1337
1338   VISU::PGroupImpl aGroup = aGroupMapIter->second;
1339   return TFindMeshOnGroup(aMesh,aGroup);
1340 }
1341
1342
1343 size_t
1344 VISU_Convertor_impl
1345 ::GetMeshOnGroupSize(const std::string& theMeshName, 
1346                      const std::string& theGroupName)
1347 {
1348   TFindMeshOnGroup aFindMeshOnGroup = FindMeshOnGroup(theMeshName,theGroupName);
1349   VISU::PMeshImpl aMesh = boost::get<0>(aFindMeshOnGroup);
1350   VISU::PGroupImpl aGroup = boost::get<1>(aFindMeshOnGroup);
1351
1352   size_t aPointsSize = 3*aMesh->GetNbPoints()*sizeof(VISU::TCoord);
1353   VISU::TNbASizeCells aNbASizeCells = aGroup->GetNbASizeCells();
1354   size_t aNbCells = aNbASizeCells.first;
1355   size_t aCellsSize = aNbASizeCells.second;
1356   size_t aConnectivityAndTypesSize = aCellsSize*sizeof(vtkIdType);
1357   size_t aLocationsSize = aNbCells*sizeof(int);
1358   vtkFloatingPointType aNbCellsPerPoint = aCellsSize / aNbCells - 1;
1359   size_t aLinksSize = aMesh->GetNbPoints() * 
1360     (vtkIdType(sizeof(vtkIdType)*aNbCellsPerPoint) + sizeof(short));
1361   aLinksSize = 0;
1362   size_t aResult = aPointsSize + aConnectivityAndTypesSize + aLocationsSize + aLinksSize;
1363   if(MYDEBUG){
1364     MSG(MYVTKDEBUG,"GetMeshOnGroupSize - aPointsSize = "<<vtkFloatingPointType(aPointsSize));
1365     MSG(MYVTKDEBUG,"GetMeshOnGroupSize - aConnectivityAndTypesSize = "<<vtkFloatingPointType(aConnectivityAndTypesSize));
1366     MSG(MYVTKDEBUG,"GetMeshOnGroupSize - aLocationsSize = "<<vtkFloatingPointType(aLocationsSize));
1367     MSG(MYVTKDEBUG,"GetMeshOnGroupSize - aLinksSize = "<<vtkFloatingPointType(aLinksSize));
1368   }
1369   MSG(MYDEBUG,"GetMeshOnGroupSize - aResult = "<<vtkFloatingPointType(aResult)<<"; theMeshName = '"
1370       <<theMeshName<<"'; theGroupName = '"<<theGroupName<<"'");
1371
1372   aResult = size_t(aResult*ERR_SIZE_CALC);
1373   return aResult;
1374 }
1375
1376
1377 VISU_Convertor_impl::TFindField
1378 VISU_Convertor_impl
1379 ::FindField(const std::string& theMeshName, 
1380             const VISU::TEntity& theEntity, 
1381             const std::string& theFieldName)
1382 {
1383   TFindMeshOnEntity aFindMeshOnEntity = 
1384     FindMeshOnEntity(theMeshName,theEntity);
1385
1386   VISU::PMeshImpl aMesh = boost::get<0>(aFindMeshOnEntity);;
1387   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindMeshOnEntity);
1388
1389   VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh->myMeshOnEntityMap;
1390   VISU::PMeshOnEntityImpl aVTKMeshOnEntity = aMeshOnEntity;
1391   if(theEntity == VISU::NODE_ENTITY){
1392     if(aMeshOnEntityMap.find(VISU::CELL_ENTITY) != aMeshOnEntityMap.end())
1393       aVTKMeshOnEntity = aMeshOnEntityMap[VISU::CELL_ENTITY];
1394     else if(aMeshOnEntityMap.find(VISU::FACE_ENTITY) != aMeshOnEntityMap.end())
1395       aVTKMeshOnEntity = aMeshOnEntityMap[VISU::FACE_ENTITY];
1396     else if(aMeshOnEntityMap.find(VISU::EDGE_ENTITY) != aMeshOnEntityMap.end())
1397       aVTKMeshOnEntity = aMeshOnEntityMap[VISU::EDGE_ENTITY];
1398     else if(aMeshOnEntityMap.find(VISU::NODE_ENTITY) != aMeshOnEntityMap.end())
1399       aVTKMeshOnEntity = aMeshOnEntityMap[VISU::NODE_ENTITY];
1400   }else
1401     aVTKMeshOnEntity = aMeshOnEntity;
1402   
1403   VISU::TFieldMap& aFieldMap = aMeshOnEntity->myFieldMap;
1404   VISU::TFieldMap::const_iterator aFieldIter= aFieldMap.find(theFieldName);
1405   if(aFieldIter == aFieldMap.end())
1406     EXCEPTION(std::runtime_error,"FindField >> There is no field on the mesh!!!");
1407   
1408   VISU::PFieldImpl aField = aFieldIter->second;
1409
1410   return TFindField(aMesh,
1411                     aMeshOnEntity,
1412                     aVTKMeshOnEntity,
1413                     aField);
1414 }
1415
1416
1417 size_t
1418 VISU_Convertor_impl
1419 ::GetFieldOnMeshSize(const std::string& theMeshName, 
1420                      const VISU::TEntity& theEntity,
1421                      const std::string& theFieldName)
1422 {
1423   TFindField aFindField = FindField(theMeshName,theEntity,theFieldName);
1424   VISU::PMeshOnEntityImpl aVTKMeshOnEntity = boost::get<2>(aFindField);
1425   VISU::PFieldImpl aField = boost::get<3>(aFindField);
1426
1427   size_t aMeshSize = GetMeshOnEntitySize(theMeshName,aVTKMeshOnEntity->myEntity);
1428   size_t aFieldOnMeshSize = size_t(aField->myDataSize*sizeof(vtkFloatingPointType)*aField->myValField.size()*ERR_SIZE_CALC);
1429   size_t aResult = aMeshSize + aFieldOnMeshSize;
1430   if(MYDEBUG)
1431     MSG(MYVTKDEBUG,"GetFieldOnMeshSize - aFieldOnMeshSize = "<<vtkFloatingPointType(aFieldOnMeshSize));
1432   MSG(MYDEBUG,"GetFieldOnMeshSize - aResult = "<<vtkFloatingPointType(aResult)<<"; theMeshName = '"<<theMeshName<<
1433       "'; theEntity = "<<theEntity<<"; theFieldName = '"<<theFieldName<<"'");
1434
1435   return aResult;
1436 }
1437
1438
1439 VISU_Convertor_impl::TFindTimeStamp
1440 VISU_Convertor_impl
1441 ::FindTimeStamp(const std::string& theMeshName, 
1442                 const VISU::TEntity& theEntity, 
1443                 const std::string& theFieldName, 
1444                 int theStampsNum)
1445 {
1446   TFindField aFindField = FindField(theMeshName,theEntity,theFieldName);
1447   VISU::PField aField = boost::get<3>(aFindField);
1448
1449   VISU::TValField& aValField = aField->myValField;
1450   VISU::TValField::const_iterator aValFieldIter= aValField.find(theStampsNum);
1451   if(aValFieldIter == aValField.end())
1452     EXCEPTION(std::runtime_error,"FindTimeStamp >> There is no field with the timestamp!!!");
1453   
1454   VISU::PMeshImpl aMesh = boost::get<0>(aFindField);
1455   VISU::PMeshOnEntityImpl aMeshOnEntity = boost::get<1>(aFindField);
1456   VISU::PMeshOnEntityImpl aVTKMeshOnEntity = boost::get<2>(aFindField);
1457   VISU::PValForTimeImpl aValForTime = aValFieldIter->second;
1458
1459   return TFindTimeStamp(aMesh,
1460                         aMeshOnEntity,
1461                         aVTKMeshOnEntity,
1462                         aField,
1463                         aValForTime);
1464 }
1465
1466
1467 size_t
1468 VISU_Convertor_impl
1469 ::GetTimeStampSize(const std::string& theMeshName, 
1470                    const VISU::TEntity& theEntity,
1471                    const std::string& theFieldName,
1472                    int theStampsNum)
1473 {
1474   TFindTimeStamp aFindTimeStamp = 
1475     FindTimeStamp(theMeshName,theEntity,theFieldName,theStampsNum);
1476   VISU::PMeshOnEntityImpl aVTKMeshOnEntity = boost::get<2>(aFindTimeStamp);
1477   VISU::PFieldImpl aField = boost::get<3>(aFindTimeStamp);
1478   
1479   size_t aMeshSize = GetMeshOnEntitySize(theMeshName, aVTKMeshOnEntity->myEntity);
1480   size_t aTimeStampSize = size_t(aField->myDataSize*sizeof(vtkFloatingPointType) * ERR_SIZE_CALC);
1481   size_t aResult = aMeshSize + aTimeStampSize;
1482
1483   MSG(MYDEBUG && MYVTKDEBUG,"GetTimeStampSize - aTimeStampSize = "<<vtkFloatingPointType(aTimeStampSize));
1484   MSG(MYDEBUG,"GetTimeStampSize - aResult = "<<vtkFloatingPointType(aResult)<<
1485       "; theMeshName = '"<<theMeshName<<"'; theEntity = "<<theEntity<<
1486       "; theFieldName = '"<<theFieldName<<"'; theStampsNum = "<<theStampsNum);
1487
1488   return aResult;
1489 }
1490
1491
1492 size_t
1493 VISU_Convertor_impl
1494 ::GetTimeStampOnMeshSize(const std::string& theMeshName, 
1495                          const VISU::TEntity& theEntity,
1496                          const std::string& theFieldName,
1497                          int theTimeStampNumber,
1498                          bool& theIsEstimated)
1499 {
1500   size_t aSize = 0;
1501
1502   //Cheching possibility do the query
1503   TFindTimeStamp aFindTimeStamp = FindTimeStamp(theMeshName,
1504                                                 theEntity,
1505                                                 theFieldName,
1506                                                 theTimeStampNumber);
1507
1508   VISU::PValForTimeImpl aValForTime = boost::get<4>(aFindTimeStamp);
1509   VISU::PUnstructuredGridIDMapperImpl anUnstructuredGridIDMapper = aValForTime->myUnstructuredGridIDMapper;
1510   if(anUnstructuredGridIDMapper->myIsVTKDone){
1511     VISU::PIDMapper anIDMapper = GetTimeStampOnMesh(theMeshName, 
1512                                                     theEntity, 
1513                                                     theFieldName, 
1514                                                     theTimeStampNumber);
1515     anIDMapper->GetOutput();
1516     aSize += anIDMapper->GetMemorySize();
1517   }else
1518     aSize += GetTimeStampSize(theMeshName, theEntity, theFieldName, theTimeStampNumber);
1519
1520   theIsEstimated = !(anUnstructuredGridIDMapper->myIsVTKDone);
1521
1522   //cout<<"VISU_Convertor_impl::GetTimeStampOnMeshSize - "<<aSize<<"; "<<(anIDMapperFilter->myIsVTKDone)<<endl;
1523   return aSize;
1524 }
1525
1526
1527 size_t
1528 VISU_Convertor_impl
1529 ::GetTimeStampOnGaussPtsSize(const std::string& theMeshName, 
1530                              const VISU::TEntity& theEntity,
1531                              const std::string& theFieldName,
1532                              int theTimeStampNumber,
1533                              bool& theIsEstimated)
1534 {
1535   size_t aSize = 0;
1536
1537   //Cheching possibility do the query
1538   TFindTimeStamp aFindTimeStamp = FindTimeStamp(theMeshName,
1539                                                 theEntity,
1540                                                 theFieldName,
1541                                                 theTimeStampNumber);
1542
1543   VISU::PValForTimeImpl aValForTime = boost::get<4>(aFindTimeStamp);
1544   VISU::PGaussPtsIDFilter aGaussPtsIDFilter = aValForTime->myGaussPtsIDFilter;
1545   if(aGaussPtsIDFilter->myIsVTKDone){
1546     VISU::PGaussPtsIDMapper aGaussPtsIDMapper = GetTimeStampOnGaussPts(theMeshName, 
1547                                                                        theEntity, 
1548                                                                        theFieldName, 
1549                                                                        theTimeStampNumber);
1550     aGaussPtsIDMapper->GetOutput();
1551     aSize += aGaussPtsIDMapper->GetMemorySize();
1552   }else
1553     aSize += GetTimeStampSize(theMeshName, theEntity, theFieldName, theTimeStampNumber);
1554
1555   theIsEstimated = !(aGaussPtsIDFilter->myIsVTKDone);
1556
1557   //cout<<"VISU_Convertor_impl::GetTimeStampOnGaussPtsSize - "<<aSize<<"; "<<(aGaussPtsIDFilter->myIsVTKDone)<<endl;
1558   return aSize;
1559 }
1560
1561
1562 const VISU::PField
1563 VISU_Convertor_impl
1564 ::GetField(const std::string& theMeshName, 
1565            VISU::TEntity theEntity, 
1566            const std::string& theFieldName) 
1567 {
1568   TFindField aFindField = FindField(theMeshName,theEntity,theFieldName);
1569   VISU::PField aField = boost::get<3>(aFindField);
1570   return aField;
1571 }
1572
1573
1574 const VISU::PValForTime 
1575 VISU_Convertor_impl
1576 ::GetTimeStamp(const std::string& theMeshName, 
1577                const VISU::TEntity& theEntity,
1578                const std::string& theFieldName,
1579                int theStampsNum)
1580 {
1581   TFindTimeStamp aFindTimeStamp = 
1582     FindTimeStamp(theMeshName,theEntity,theFieldName,theStampsNum);
1583   VISU::PValForTime aValForTime = boost::get<4>(aFindTimeStamp);
1584   return aValForTime;
1585 }