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