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