]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_Convertor_impl.cxx
Salome HOME
1e3aaa6119ae8ba1dbd3760001c9b6d5e2351b99
[modules/visu.git] / src / CONVERTOR / VISU_Convertor_impl.cxx
1 //  VISU OBJECT : interactive object for VISU entities implementation
2 //
3 //  Copyright (C) 2003  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 //  CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS 
5 // 
6 //  This library is free software; you can redistribute it and/or 
7 //  modify it under the terms of the GNU Lesser General Public 
8 //  License as published by the Free Software Foundation; either 
9 //  version 2.1 of the License. 
10 // 
11 //  This library is distributed in the hope that it will be useful, 
12 //  but WITHOUT ANY WARRANTY; without even the implied warranty of 
13 //  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU 
14 //  Lesser General Public License for more details. 
15 // 
16 //  You should have received a copy of the GNU Lesser General Public 
17 //  License along with this library; if not, write to the Free Software 
18 //  Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA 
19 // 
20 //  See http://www.opencascade.org/SALOME/ or email : webmaster.salome@opencascade.org 
21 //
22 //
23 //  File   : VISU_Convertor_impl.cxx
24 //  Author : Alexey PETROV
25 //  Module : VISU
26
27 #include "VISU_Convertor_impl.hxx"
28 #include "VISU_ConvertorUtils.hxx"
29
30 #include <vtkIdList.h>
31 #include <vtkCellType.h>
32 #include <vtkIntArray.h>
33 #include <vtkCellArray.h>
34 #include <vtkFloatArray.h>
35 #include <vtkUnsignedCharArray.h>
36 #include <vtkPointData.h>
37 #include <vtkCellData.h>
38 #include <vtkCellLinks.h>
39
40 #include <vtkMergeDataObjectFilter.h>
41
42 #include <vtkUnstructuredGridWriter.h>
43
44 #include <qstring.h>
45 #include <qfileinfo.h>
46
47 #include <valarray>     
48 #include <memory>
49
50 using namespace std;
51
52 static float ERR_SIZE_CALC = 1.00;
53
54 static int MYVTKDEBUG = 0;
55
56 #ifdef _DEBUG_
57 static int MYDEBUG = 0;
58 static int MYDEBUGWITHFILES = 0;
59 #else
60 static int MYDEBUG = 0;
61 static int MYDEBUGWITHFILES = 0;
62 #endif
63
64 void VISU::WriteToFile(vtkUnstructuredGrid* theDataSet, const string& theFileName){
65   //vtkUnstructuredGridWriter* aWriter = vtkUnstructuredGridWriter::New();
66   //if(MYVTKDEBUG) aWriter->DebugOn();
67   ////aWriter->SetFileType(VTK_BINARY);
68   //aWriter->SetFileName(theFileName.c_str());
69   //aWriter->SetInput(theDataSet);
70   //aWriter->Write();
71   //aWriter->Delete();
72 }
73
74 enum ECoordName{eX, eY, eZ, eNone};
75 typedef VISU::TMesh::TCoord (*TGetCoord)(const VISU::TMesh::TPointsCoord&, int);
76
77 template<ECoordName TheCoordId>
78 VISU::TMesh::TCoord GetCoord(const VISU::TMesh::TPointsCoord& thePointsCoord, 
79                              int theStartPos)
80 {
81   return thePointsCoord[theStartPos+TheCoordId];
82 }
83
84 template<>
85 VISU::TMesh::TCoord GetCoord<eNone>(const VISU::TMesh::TPointsCoord& thePointsCoord, 
86                                     int theStartPos)
87 {
88   return 0.0;
89 }
90
91
92 static TGetCoord aXYZGetCoord[3] = {
93   &GetCoord<eX>, 
94   &GetCoord<eY>, 
95   &GetCoord<eZ>
96 };
97
98
99 static TGetCoord aXYGetCoord[3] = {
100   &GetCoord<eX>, 
101   &GetCoord<eY>, 
102   &GetCoord<eNone>
103 };
104
105 static TGetCoord aYZGetCoord[3] = {
106   &GetCoord<eNone>,
107   &GetCoord<eX>, 
108   &GetCoord<eY>
109 };
110
111 static TGetCoord aXZGetCoord[3] = {
112   &GetCoord<eX>, 
113   &GetCoord<eNone>,
114   &GetCoord<eY>
115 };
116
117
118 static TGetCoord aXGetCoord[3] = {
119   &GetCoord<eX>, 
120   &GetCoord<eNone>,
121   &GetCoord<eNone>
122 };
123
124 static TGetCoord aYGetCoord[3] = {
125   &GetCoord<eNone>,
126   &GetCoord<eX>, 
127   &GetCoord<eNone>
128 };
129
130 static TGetCoord aZGetCoord[3] = {
131   &GetCoord<eNone>,
132   &GetCoord<eNone>,
133   &GetCoord<eX>
134 };
135
136
137 class TCoordHelper{
138   const VISU::TMesh::TPointsCoord& myPointsCoord;
139   TGetCoord* myGetCoord;
140 public:
141   TCoordHelper(const VISU::TMesh::TPointsCoord& thePointsCoord,
142                TGetCoord* theGetCoord):
143     myPointsCoord(thePointsCoord),
144     myGetCoord(theGetCoord)
145   {}
146   virtual ~TCoordHelper(){}
147   VISU::TMesh::TCoord GetCoord(int theStartPos, int theCoodId){
148     return (*myGetCoord[theCoodId])(myPointsCoord,theStartPos);
149   }
150 };
151 typedef std::auto_ptr<TCoordHelper> TCoordHelperPtr;
152
153
154 namespace{
155   void GetPoints(VISU::TVTKSource& theStorage, VISU::TMesh& theMesh) 
156   {
157     vtkPoints* aPoints = theMesh.myPoints.GetPointer();
158     if(!aPoints){
159       aPoints = vtkPoints::New();
160       TCoordHelperPtr aCoordHelperPtr;
161       const VISU::TMesh::TPointsCoord& anArray = theMesh.myPointsCoord;
162       {
163         int aMeshDimension = theMesh.myDim;
164         bool anIsDimPresent[3] = {false, false, false};
165         for(int iDim = 0; iDim < aMeshDimension; iDim++){
166           string aDimName = theMesh.myPointsDim[iDim];
167           if(aDimName == "x" || aDimName == "X")
168             anIsDimPresent[eX] = true;
169           else if(aDimName == "y" || aDimName == "Y")
170             anIsDimPresent[eY] = true;
171           else if(aDimName == "z" || aDimName == "Z")
172             anIsDimPresent[eZ] = true;
173         }
174
175         switch(aMeshDimension){
176         case 3:
177           aCoordHelperPtr.reset(new TCoordHelper(anArray,aXYZGetCoord));
178           break;
179         case 2:
180           if(anIsDimPresent[eY] && anIsDimPresent[eZ])
181             aCoordHelperPtr.reset(new TCoordHelper(anArray,aYZGetCoord));
182           else if(anIsDimPresent[eX] && anIsDimPresent[eZ])
183             aCoordHelperPtr.reset(new TCoordHelper(anArray,aXZGetCoord));
184           else
185             aCoordHelperPtr.reset(new TCoordHelper(anArray,aXYGetCoord));
186           break;
187         case 1:
188           if(anIsDimPresent[eY])
189             aCoordHelperPtr.reset(new TCoordHelper(anArray,aYGetCoord));
190           else if(anIsDimPresent[eZ])
191             aCoordHelperPtr.reset(new TCoordHelper(anArray,aZGetCoord));
192           else
193             aCoordHelperPtr.reset(new TCoordHelper(anArray,aXGetCoord));
194           break;
195         }
196       }
197
198       if(MYVTKDEBUG) aPoints->DebugOn();
199       vtkIdType iEnd = theMesh.myPointsCoord.size();
200       vtkIdType aNbPoints = iEnd / theMesh.myDim;
201       aPoints->SetNumberOfPoints(aNbPoints);
202       if(MYDEBUG) 
203         MESSAGE("GetPoints - aNbPoints = "<<aNbPoints<<"; myDim = "<<theMesh.myDim);
204       for (vtkIdType i = 0, j = 0; i < iEnd; i += theMesh.myDim, j++) 
205         aPoints->SetPoint(j,
206                           aCoordHelperPtr->GetCoord(i,eX),
207                           aCoordHelperPtr->GetCoord(i,eY),
208                           aCoordHelperPtr->GetCoord(i,eZ));
209       theMesh.myPoints = aPoints;
210     }
211     theStorage->SetPoints(aPoints);
212   }
213   
214   
215   inline void PrintCells(int& theStartId,
216                          vtkCellArray* theConnectivity, 
217                          const VISU::TMeshOnEntity::TConnect& theVector)
218   {
219     vtkIdList *anIdList = vtkIdList::New();
220     int kEnd = theVector.size();
221     anIdList->SetNumberOfIds(kEnd);
222     for(int k = 0; k < kEnd; k++)
223       anIdList->SetId(k,theVector[k]);
224     theConnectivity->InsertNextCell(anIdList);
225     anIdList->Delete();
226   }
227
228   void GetCellsOnEntity(VISU::TVTKSource& theStorage,
229                         const VISU::TMeshOnEntity& theMeshOnEntity, 
230                         const string& theFamilyName) 
231   {
232     //Check on existing family
233     const VISU::TFamily* pFamily = VISU::GetFamily(theMeshOnEntity,theFamilyName);
234     bool isFamilyPresent = (pFamily != NULL);
235     const VISU::TFamily& aFamily = *pFamily;
236     //Main part of code
237     pair<int,int> aCellsDim = theMeshOnEntity.GetCellsDims(theFamilyName);
238     int aNbCells = aCellsDim.first, aCellsSize = aCellsDim.second;
239     vtkCellArray* aConnectivity = vtkCellArray::New();
240     aConnectivity->Allocate(aCellsSize,0);
241     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
242     aCellTypesArray->SetNumberOfComponents(1);
243     aCellTypesArray->SetNumberOfTuples(aNbCells);
244     if(MYDEBUG) MESSAGE("GetCellsOnEntity - isFamilyPresent = "<<isFamilyPresent);
245     const VISU::TMeshOnEntity::TCellsConn &aCellsConn = theMeshOnEntity.myCellsConn;
246     VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
247     for(int i = 0, j = 0; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
248       const VISU::TMeshOnEntity::TConnForCellType& anArray = aCellsConnIter->second;
249       int aVtkType = aCellsConnIter->first;
250       if(MYDEBUG) 
251         MESSAGE("GetCellsOnEntity - aVtkType = "<<aVtkType<<"; anArray.size() = "<<anArray.size());
252       if(!isFamilyPresent)
253         for(int k = 0, kEnd = anArray.size(); k < kEnd; k++, i++){
254           PrintCells(i,aConnectivity,anArray[k]);
255           aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
256         }
257       else{
258         const VISU::TFamily::TSubMesh& aSubMesh = aFamily.mySubMesh;
259         if(aSubMesh.empty()) 
260           throw std::runtime_error(EXCEPTION("GetCells >> There is no elements on the family !!!"));
261         VISU::TFamily::TSubMesh::const_iterator aSubMeshIter = aSubMesh.find(aVtkType);
262         if(aSubMeshIter == aSubMesh.end()) continue;
263         const VISU::TFamily::TSubMeshOnCellType& aSubMeshOnCellType = aSubMeshIter->second;
264         if(MYDEBUG) 
265           MESSAGE("GetCellsOnEntity - aSubMeshOnCellType.size() = "<<aSubMeshOnCellType.size());
266         VISU::TFamily::TSubMeshOnCellType::const_iterator aSubMeshOnCellTypeIter = aSubMeshOnCellType.begin();
267         for(; aSubMeshOnCellTypeIter != aSubMeshOnCellType.end(); aSubMeshOnCellTypeIter++, i++){
268           PrintCells(i,aConnectivity,anArray[*aSubMeshOnCellTypeIter]);
269           aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
270         }
271       }
272     }
273     vtkIdType *pts = 0, npts = 0;
274     vtkIntArray* aCellLocationsArray = vtkIntArray::New();
275     aCellLocationsArray->SetNumberOfComponents(1);
276     aCellLocationsArray->SetNumberOfTuples(aNbCells);
277     aConnectivity->InitTraversal();
278     for(int i=0; aConnectivity->GetNextCell(npts,pts); i++)
279       aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
280     theStorage->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
281     if(MYVTKDEBUG) aConnectivity->DebugOn();
282     aCellLocationsArray->Delete();
283     aCellTypesArray->Delete();
284     aConnectivity->Delete();
285   } 
286   
287   
288   void GetCellsOnGroup(VISU::TVTKSource& theStorage,
289                        const VISU::TMesh& theMesh,
290                        const VISU::TFamilyAndEntitySet& theFamilyAndEntitySet) 
291   {
292     //Calculate dimentions of the group
293     int aNbCells = 0, aCellsSize = 0;
294     VISU::TFamilyAndEntitySet::const_iterator aFamilyAndEntitySetIter = theFamilyAndEntitySet.begin();
295     for(; aFamilyAndEntitySetIter != theFamilyAndEntitySet.end(); aFamilyAndEntitySetIter++){
296       const VISU::TFamilyAndEntity& aFamilyAndEntity = *aFamilyAndEntitySetIter;
297       const string& aFamilyName = aFamilyAndEntity.first;
298       const VISU::TEntity& anEntity = aFamilyAndEntity.second;
299       const VISU::TMeshOnEntity& aMeshOnEntity = theMesh.myMeshOnEntityMap.find(anEntity)->second;
300       pair<int,int> aCellsDim = aMeshOnEntity.GetCellsDims(aFamilyName);
301       aNbCells += aCellsDim.first;
302       aCellsSize += aCellsDim.second;
303     }
304     vtkCellArray* aConnectivity = vtkCellArray::New();
305     aConnectivity->Allocate(aCellsSize,0);
306     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
307     aCellTypesArray->SetNumberOfComponents(1);
308     aCellTypesArray->SetNumberOfTuples(aNbCells);
309     aFamilyAndEntitySetIter = theFamilyAndEntitySet.begin();
310     for(int i = 0, j = 0; aFamilyAndEntitySetIter != theFamilyAndEntitySet.end(); aFamilyAndEntitySetIter++){
311       const VISU::TFamilyAndEntity& aFamilyAndEntity = *aFamilyAndEntitySetIter;
312       const string& aFamilyName = aFamilyAndEntity.first;
313       const VISU::TEntity& anEntity = aFamilyAndEntity.second;
314       const VISU::TMeshOnEntity& aMeshOnEntity = theMesh.myMeshOnEntityMap.find(anEntity)->second;
315       const VISU::TFamily& aFamily = *(VISU::GetFamily(aMeshOnEntity,aFamilyName));
316       const VISU::TMeshOnEntity::TCellsConn &aCellsConn = aMeshOnEntity.myCellsConn;
317       VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
318       for(; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
319         const VISU::TMeshOnEntity::TConnForCellType& anArray = aCellsConnIter->second;
320         int aVtkType = aCellsConnIter->first;
321         if(MYDEBUG) 
322           MESSAGE("GetCellsOnGroup - aVtkType = "<<aVtkType<<"; anArray.size() = "<<anArray.size());
323         const VISU::TFamily::TSubMesh& aSubMesh = aFamily.mySubMesh;
324         if(aSubMesh.empty()) 
325           throw std::runtime_error(EXCEPTION("GetCells >> There is no elements on the family !!!"));
326         VISU::TFamily::TSubMesh::const_iterator aSubMeshIter = aSubMesh.find(aVtkType);
327         if(aSubMeshIter == aSubMesh.end()) continue;
328         const VISU::TFamily::TSubMeshOnCellType& aSubMeshOnCellType = aSubMeshIter->second;
329         if(MYDEBUG) 
330           MESSAGE("GetCellsOnGroup - aSubMeshOnCellType.size() = "<<aSubMeshOnCellType.size());
331         VISU::TFamily::TSubMeshOnCellType::const_iterator aSubMeshOnCellTypeIter = aSubMeshOnCellType.begin();
332         for(; aSubMeshOnCellTypeIter != aSubMeshOnCellType.end(); aSubMeshOnCellTypeIter++, i++){
333           PrintCells(i,aConnectivity,anArray[*aSubMeshOnCellTypeIter]);
334           aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
335         }
336       }
337     }
338     vtkIdType *pts = 0, npts = 0;
339     vtkIntArray* aCellLocationsArray = vtkIntArray::New();
340     aCellLocationsArray->SetNumberOfComponents(1);
341     aCellLocationsArray->SetNumberOfTuples(aNbCells);
342     aConnectivity->InitTraversal();
343     for(int i = 0; aConnectivity->GetNextCell(npts,pts); i++)
344       aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
345     theStorage->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
346     aCellLocationsArray->Delete();
347     aCellTypesArray->Delete();
348     aConnectivity->Delete();
349   } 
350   
351   
352   void InitProfile(VISU::TVTKExtractFilter& theFilter,
353                    const VISU::TMeshOnEntity& theMeshOnEntity, 
354                    const VISU::TField::TValForTime& theValForTime)
355   {
356     const VISU::TField::TValForCells& aValForCells = theValForTime.myValForCells;
357     const VISU::TMeshOnEntity::TCellsConn &aCellsConn = theMeshOnEntity.myCellsConn;
358     VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
359     for(; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
360       const vtkIdType& aCellType = aCellsConnIter->first;
361       if(aValForCells.find(aCellType) == aValForCells.end())
362         theFilter->RemoveCellsWithType(aCellType);
363     }
364   }
365
366
367   void GetValsOnTimeStamp(vtkFloatArray *theFloatArray, 
368                           const vtkIdType& theNumberOfTuples,
369                           const std::string& theFieldName,
370                           const VISU::TField& theField,
371                           const VISU::TField::TValForTime& theValForTime)
372   {
373     //theFloatArray->DebugOn();
374     theFloatArray->SetNumberOfTuples(theNumberOfTuples);
375     theFloatArray->SetName(theFieldName.c_str());
376     if(MYDEBUG) MESSAGE("GetValsOnTimeStamp - theNumberOfTuples = "<<theNumberOfTuples);
377     const VISU::TField::TValForCells& aValForCells = theValForTime.myValForCells;
378     VISU::TField::TValForCells::const_iterator aValForCellsIter = aValForCells.begin();
379     for(int k = 0; aValForCellsIter != aValForCells.end(); aValForCellsIter++) {
380       const VISU::TField::TValForCellsWithType& anArray = aValForCellsIter->second;
381       int iEnd = anArray.size()/theField.myNbComp;
382       int aVtkType = aValForCellsIter->first;
383       if(MYDEBUG) MESSAGE("GetValsOnTimeStamp -  iEnd = "<<iEnd<<"; aVtkType = "<<aVtkType);
384       switch(theField.myNbComp) {
385       case 1:
386         for (int i = 0; i < iEnd; i++) 
387           theFloatArray->SetTuple1(k++,anArray[i]);
388         break;
389       case 2:
390         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*2)
391           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+1],0.0);
392         break;
393       case 3:
394         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*3)
395           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+1],anArray[ji+2]);
396         break;
397       case 4:
398         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*4)
399           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+2],0.0);
400         break;
401       case 6:
402         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*4)
403           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+2],anArray[ji+5]);
404         break;
405       default:
406         throw std::runtime_error(EXCEPTION("GetValsOnTimeStamp - There is no an algorithm for representation of the field !!!"));
407       }
408     }
409   }
410
411   string GenerateFieldName(const VISU::TField& theField,
412                            const VISU::TField::TValForTime& theValForTime)
413   {
414     const VISU::TField::TTime& aTime = theValForTime.myTime;
415     string aFieldName = theField.myMeshName + ", " + theField.myName + ": " + 
416       VISU_Convertor::GenerateName(aTime);
417     return aFieldName;
418   }
419
420   void GetTimeStamp(VISU::TVTKSource& theStorage,
421                     const VISU::TMesh& theMesh,
422                     const VISU::TMeshOnEntity& theMeshOnEntity, 
423                     const VISU::TField& theField, 
424                     const VISU::TField::TValForTime& theValForTime)
425   {
426     int aNumberOfTuples = theField.myDataSize/theField.myNbComp;
427     string aFieldName = GenerateFieldName(theField,theValForTime);
428     if(MYDEBUG) 
429       MESSAGE("GetTimeStamp(TVTKSource) - aFieldName = "<<aFieldName<<
430               "; aNumberOfTuples = "<<aNumberOfTuples);
431
432     vtkDataSetAttributes* aDataSetAttributes;
433     switch(theField.myEntity){
434     case VISU::NODE_ENTITY : 
435       aDataSetAttributes = theStorage->GetPointData();
436       break;
437     default: 
438       aDataSetAttributes = theStorage->GetCellData();
439     }
440
441     vtkFloatArray *aFloatArray = vtkFloatArray::New();
442     switch(theField.myNbComp) {
443     case 1:
444       aFloatArray->SetNumberOfComponents(1);
445       aDataSetAttributes->SetScalars(aFloatArray);
446       break;
447     default:
448       aFloatArray->SetNumberOfComponents(3);
449       aDataSetAttributes->SetVectors(aFloatArray);
450     }
451
452     GetValsOnTimeStamp(aFloatArray,aNumberOfTuples,aFieldName,theField,theValForTime);
453   }
454
455   void GetTimeStamp(VISU::TVTKAttribyteFilter& theAttribyteFilter,
456                     VISU::TVTKMergetFilter& theMergeFilter,
457                     VISU::TVTKExtractFilter& theExtractFilter,
458                     const VISU::TMesh& theMesh,
459                     const VISU::TMeshOnEntity& theMeshOnEntity, 
460                     const VISU::TField& theField, 
461                     const VISU::TField::TValForTime& theValForTime)
462   {
463     int aNumberOfTuples = theField.myDataSize/theField.myNbComp;
464     string aFieldName = GenerateFieldName(theField,theValForTime);
465     if(MYDEBUG) 
466       MESSAGE("GetTimeStamp(TVTKAttribyteFilter) - aFieldName = "<<aFieldName<<
467               "; aNumberOfTuples = "<<aNumberOfTuples);
468
469     vtkDataObject* aDataObject = vtkDataObject::New();
470     theMergeFilter->SetDataObject(aDataObject);
471     aDataObject->Delete();
472
473     theMergeFilter->SetInput(theExtractFilter->GetOutput());
474     theAttribyteFilter->SetInput(theMergeFilter->GetOutput());
475
476     switch(theField.myEntity){
477     case VISU::NODE_ENTITY : 
478       theMergeFilter->SetOutputFieldToPointDataField();
479       theAttribyteFilter->SetInputFieldToPointDataField();
480       theAttribyteFilter->SetOutputAttributeDataToPointData();
481       break;
482     default: 
483       theMergeFilter->SetOutputFieldToCellDataField();
484       theAttribyteFilter->SetInputFieldToCellDataField();
485       theAttribyteFilter->SetOutputAttributeDataToCellData();
486     }
487
488     vtkFloatArray *aFloatArray = vtkFloatArray::New();
489     switch(theField.myNbComp) {
490     case 1:
491       aFloatArray->SetNumberOfComponents(1);
492       theAttribyteFilter->SetScalarComponent(0,aFieldName.c_str(),0);
493       break;
494     default:
495       aFloatArray->SetNumberOfComponents(3);
496       theAttribyteFilter->SetVectorComponent(0,aFieldName.c_str(),0);
497       theAttribyteFilter->SetVectorComponent(1,aFieldName.c_str(),1);
498       theAttribyteFilter->SetVectorComponent(2,aFieldName.c_str(),2);
499     }
500
501     vtkFieldData* aFieldData = aDataObject->GetFieldData();
502     aFieldData->AddArray(aFloatArray);
503     aFloatArray->Delete();
504
505     GetValsOnTimeStamp(aFloatArray,aNumberOfTuples,aFieldName,theField,theValForTime);
506   }
507 }
508
509 VISU_Convertor_impl::VISU_Convertor_impl() {
510   myIsDone = false;
511 }
512
513 VISU_Convertor_impl::~VISU_Convertor_impl() {}
514
515 VISU_Convertor::TOutput* 
516 VISU_Convertor_impl::GetMeshOnEntity(const string& theMeshName, 
517                                      const VISU::TEntity& theEntity,
518                                      const string& theFamilyName)
519 {
520   if(MYDEBUG) 
521     MESSAGE("GetMeshOnEntity - theMeshName = '"<<theMeshName<<
522             "'; theEntity = "<<theEntity<<"; theFamilyName = '"<<theFamilyName<<"'");
523   //Cheching possibility do the query
524   VISU::TMesh* pMesh = NULL;
525   VISU::TFamily* pFamily = NULL;
526   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
527   FindMeshOnEntity(theMeshName,pMesh,theEntity,pMeshOnEntity,theFamilyName,pFamily);
528   VISU::TMesh& aMesh = *pMesh;
529   VISU::TMeshOnEntity& aMeshOnEntity = *pMeshOnEntity;
530   VISU::TVTKSource* pSource;
531   if(pFamily != NULL)
532     pSource = &(pFamily->myStorage);
533   else
534     pSource = &(aMeshOnEntity.myStorage);
535   VISU::TVTKSource& aSource = *pSource;
536   //Main part of code
537   try{
538     if(aSource.GetPointer() == NULL){
539       aSource = TOutput::New();
540       aSource->Delete();
541       if(MYVTKDEBUG) aSource->DebugOn();
542       LoadMeshOnEntity(aMeshOnEntity,theFamilyName);
543       GetPoints(aSource,aMesh);
544       GetCellsOnEntity(aSource,aMeshOnEntity,theFamilyName);
545       if(MYDEBUGWITHFILES){
546         string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
547         string aFamilyName = QString(theFamilyName.c_str()).simplifyWhiteSpace().latin1();
548         string aFileName = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
549         aFileName += aMeshName + dtos("-%d-",int(theEntity)) + aFamilyName + "-Conv.vtk";
550         VISU::WriteToFile(aSource.GetPointer(),aFileName);
551       }
552     }
553     if(MYVTKDEBUG){
554       GetMeshOnEntitySize(theMeshName,theEntity,theFamilyName);
555       vtkUnstructuredGrid* aDataSet = aSource.GetPointer();
556       aDataSet->Update();
557       MESSAGE("GetMeshOnEntity - GetPoints() = "<<float(aDataSet->GetPoints()->GetActualMemorySize()*1000));
558       MESSAGE("GetMeshOnEntity - GetCells() = "<<float(aDataSet->GetCells()->GetActualMemorySize()*1000));
559       MESSAGE("GetMeshOnEntity - GetCellTypesArray() = "<<float(aDataSet->GetCellTypesArray()->GetActualMemorySize()*1000));
560       MESSAGE("GetMeshOnEntity - GetCellLocationsArray() = "<<float(aDataSet->GetCellLocationsArray()->GetActualMemorySize()*1000));
561       aDataSet->BuildLinks();
562       MESSAGE("GetMeshOnEntity - GetCellLinks() = "<<float(aDataSet->GetCellLinks()->GetActualMemorySize()*1000));
563       MESSAGE("GetMeshOnEntity - GetActualMemorySize() = "<<float(aDataSet->GetActualMemorySize()*1000));
564     }
565   }catch(...){
566     aSource = vtkSmartPointerBase();
567     throw;
568   }
569   return aSource.GetPointer();
570 }
571
572 VISU_Convertor::TOutput* 
573 VISU_Convertor_impl::GetMeshOnGroup(const string& theMeshName, 
574                                     const string& theGroupName)
575 {
576   if(MYDEBUG) MESSAGE("GetMeshOnGroup - theMeshName = '"<<theMeshName<<
577                       "'; theGroupName = '"<<theGroupName<<"'");
578   //Cheching possibility do the query
579   VISU::TMesh* pMesh = NULL;
580   VISU::TGroup* pGroup = NULL;
581   FindMeshOnGroup(theMeshName,pMesh,theGroupName,pGroup);
582   VISU::TMesh& aMesh = *pMesh;
583   VISU::TGroup& aGroup = *pGroup;
584   const VISU::TFamilyAndEntitySet& aFamilyAndEntitySet = aGroup.myFamilyAndEntitySet;
585   VISU::TVTKSource& aSource = aGroup.myStorage;
586   //Main part of code
587   try{
588     if(aSource.GetPointer() == NULL){
589       aSource = TOutput::New();
590       aSource->Delete();
591       LoadMeshOnGroup(aMesh,aFamilyAndEntitySet);
592       GetPoints(aSource,aMesh);
593       GetCellsOnGroup(aSource,aMesh,aFamilyAndEntitySet);
594       if(MYDEBUGWITHFILES){
595         string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
596         string aGroupName = QString(theGroupName.c_str()).simplifyWhiteSpace().latin1();
597         string aFileName = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
598         aFileName += aMeshName + "-" + aGroupName + "-Conv.vtk";
599         VISU::WriteToFile(aSource.GetPointer(),aFileName);
600       }
601     }
602   }catch(...){
603     aSource = vtkSmartPointerBase();
604     throw;
605   }
606   return aSource.GetPointer();
607 }
608
609 VISU_Convertor::TOutput* 
610 VISU_Convertor_impl::GetTimeStampOnMesh(const string& theMeshName, 
611                                         const VISU::TEntity& theEntity,
612                                         const string& theFieldName,
613                                         int theStampsNum)
614 {
615   if(MYDEBUG){
616     MESSAGE("GetTimeStampOnMesh - theMeshName = '"<<theMeshName<<"; theEntity = "<<theEntity);
617     MESSAGE("GetTimeStampOnMesh - theFieldName = '"<<theFieldName<<"'; theStampsNum = "<<theStampsNum);
618   }
619   //Cheching possibility do the query
620   VISU::TMesh* pMesh = NULL;
621   VISU::TField* pField = NULL;
622   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
623   VISU::TField::TValForTime* pValForTime = NULL;
624   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
625                 theFieldName,pField,theStampsNum,pValForTime);
626   VISU::TMesh& aMesh = *pMesh;
627   VISU::TMeshOnEntity& aMeshOnEntity = *pMeshOnEntity;
628   VISU::TField& aField = *pField;
629   VISU::TField::TValForTime& aValForTime = *pValForTime;
630   VISU::TVTKAttribyteFilter& anAttribyteFilter = aValForTime.myAttribyteFilter;
631   VISU::TVTKSource& aSource = aValForTime.myStorage;
632   TOutput* anOutput = NULL;
633   //Main part of code
634   try{
635     if(aSource.GetPointer())
636       return aSource.GetPointer();
637     else if(anAttribyteFilter.GetPointer())
638       return anAttribyteFilter->GetUnstructuredGridOutput();
639     else{
640       LoadFieldOnMesh(aMesh,aMeshOnEntity,aField,aValForTime);
641
642       VISU::TVTKExtractFilter& anExtractFilter = aField.myExtractFilter;
643       if(anExtractFilter.GetPointer() == NULL){
644         anExtractFilter = VISU_ExtractUnstructuredGrid::New();
645         anExtractFilter->Delete();
646         //anExtractFilter->DebugOn();
647         try{
648           LoadMeshOnEntity(*pVTKMeshOnEntity);
649         }catch(std::exception& exc){
650           pVTKMeshOnEntity = pMeshOnEntity;
651           MESSAGE("Follow exception was occured :\n"<<exc.what());
652         }catch(...){
653           pVTKMeshOnEntity = pMeshOnEntity;
654           MESSAGE("Unknown exception was occured!");
655         }
656         GetMeshOnEntity(pVTKMeshOnEntity->myMeshName,pVTKMeshOnEntity->myEntity);
657         
658         anExtractFilter->SetInput(pVTKMeshOnEntity->myStorage.GetPointer());
659         ::InitProfile(anExtractFilter,*pMeshOnEntity,aValForTime);
660       }      
661       if(!anExtractFilter->IsRemoving()){
662         aSource = TOutput::New();
663         aSource->Delete();
664         aSource->ShallowCopy(pVTKMeshOnEntity->myStorage.GetPointer());
665         ::GetTimeStamp(aSource,aMesh,*pVTKMeshOnEntity,aField,aValForTime);
666         anOutput = aSource.GetPointer();
667       }else{
668         anAttribyteFilter = vtkFieldDataToAttributeDataFilter::New();
669         anAttribyteFilter->Delete();
670         //anAttribyteFilter->DebugOn();
671         
672         VISU::TVTKMergetFilter& aMergeFilter = aValForTime.myMergeFilter;
673         aMergeFilter = vtkMergeDataObjectFilter::New();
674         aMergeFilter->Delete();
675         //aMergeFilter->DebugOn();
676
677         ::GetTimeStamp(anAttribyteFilter,aMergeFilter,anExtractFilter,
678                        aMesh,*pVTKMeshOnEntity,aField,aValForTime);
679         anOutput = anAttribyteFilter->GetUnstructuredGridOutput();
680       }
681       if(MYDEBUGWITHFILES){
682         string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
683         string aFieldName = QString(theFieldName.c_str()).simplifyWhiteSpace().latin1();
684         string aPrefix = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
685         string aFileName = aPrefix + aMeshName + dtos("-%d-",int(theEntity)) + 
686           aFieldName + dtos("-%d",theStampsNum) + "-Conv.vtk";
687         VISU::WriteToFile(anOutput,aFileName);
688       }
689       if(MYVTKDEBUG){
690         GetTimeStampSize(theMeshName,theEntity,theFieldName,theStampsNum);
691         vtkUnstructuredGrid *aDataSet = anAttribyteFilter->GetUnstructuredGridOutput();
692         aDataSet->Update();
693         if(theEntity == VISU::NODE_ENTITY)
694           MESSAGE("GetTimeStampOnMesh - GetData() = "<<float(aDataSet->GetPointData()->GetActualMemorySize()*1000));
695         else
696           MESSAGE("GetMeshOnEntity - GetData() = "<<float(aDataSet->GetCellData()->GetActualMemorySize()*1000));
697         MESSAGE("GetTimeStampOnMesh - GetActualMemorySize() = "<<float(aDataSet->GetActualMemorySize()*1000));
698       }
699     }
700   }catch(...){
701     aSource = vtkSmartPointerBase();
702     anAttribyteFilter = vtkSmartPointerBase();
703     throw;
704   }
705   return anOutput;
706 }
707
708 void VISU_Convertor_impl::FindMesh(const string& theMeshName, VISU::TMesh*& theMesh)
709 {
710   GetMeshMap();
711   if(myMeshMap.find(theMeshName) == myMeshMap.end())
712     throw std::runtime_error(EXCEPT("FindMesh >> There is no mesh with the name - '%1'!!!").
713                              arg(theMeshName.c_str()).latin1());
714   theMesh = &myMeshMap[theMeshName];
715 }
716
717
718 void VISU_Convertor_impl::FindMeshOnEntity(const string& theMeshName, 
719                                            VISU::TMesh*& theMesh,
720                                            const VISU::TEntity& theEntity, 
721                                            VISU::TMeshOnEntity*& theMeshOnEntity,
722                                            const string& theFamilyName, 
723                                            VISU::TFamily*& theFamily)
724 {
725   FindMesh(theMeshName,theMesh);
726   VISU::TMeshOnEntityMap& aMeshOnEntityMap = theMesh->myMeshOnEntityMap;
727   if(aMeshOnEntityMap.find(theEntity) == aMeshOnEntityMap.end())
728     throw std::runtime_error(EXCEPT("FindMeshOnEntity >> There is no mesh on the entity - %1!!!").
729                              arg(int(theEntity)).latin1());
730   theMeshOnEntity = &aMeshOnEntityMap[theEntity];
731   theFamily = VISU::GetFamily(*theMeshOnEntity,theFamilyName);
732 }
733
734
735 float VISU_Convertor_impl::GetSize() {
736   float aResult = 0.0;
737   const VISU::TMeshMap& aMeshMap = GetMeshMap();
738   VISU::TMeshMap::const_iterator aMeshMapIter = aMeshMap.begin();
739   for(; aMeshMapIter != aMeshMap.end(); aMeshMapIter++){
740     const string& aMeshName = aMeshMapIter->first;
741     const VISU::TMesh& aMesh = aMeshMapIter->second;
742     const VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh.myMeshOnEntityMap;
743     VISU::TMeshOnEntityMap::const_iterator aMeshOnEntityMapIter;
744     //Import fields
745     aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
746     for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
747       const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
748       const VISU::TMeshOnEntity& aMeshOnEntity = aMeshOnEntityMapIter->second;
749       const VISU::TFieldMap& aFieldMap = aMeshOnEntity.myFieldMap;
750       VISU::TFieldMap::const_iterator aFieldMapIter = aFieldMap.begin();
751       for(; aFieldMapIter != aFieldMap.end(); aFieldMapIter++){
752         const string& aFieldName = aFieldMapIter->first;
753         const VISU::TField& aField = aFieldMapIter->second;
754         const VISU::TField::TValField& aValField = aField.myValField;
755         VISU::TField::TValField::const_iterator aValFieldIter = aValField.begin();
756         for(; aValFieldIter != aValField.end(); aValFieldIter++){
757           int aTimeStamp = aValFieldIter->first;
758           aResult += GetTimeStampSize(aMeshName,anEntity,aFieldName,aTimeStamp);
759         }
760       }
761       //Importing groups
762       const VISU::TGroupMap& aGroupMap = aMesh.myGroupMap;
763       VISU::TGroupMap::const_iterator aGroupMapIter = aGroupMap.begin();
764       for(; aGroupMapIter != aGroupMap.end(); aGroupMapIter++){
765         const string& aGroupName = aGroupMapIter->first;
766         aResult += GetMeshOnGroupSize(aMeshName,aGroupName);
767       }
768       //Import families
769       const VISU::TFamilyMap& aFamilyMap = aMeshOnEntity.myFamilyMap;
770       VISU::TFamilyMap::const_iterator aFamilyMapIter = aFamilyMap.begin();
771       for(; aFamilyMapIter != aFamilyMap.end(); aFamilyMapIter++){
772         const string& aFamilyName = aFamilyMapIter->first;
773         aResult += GetMeshOnEntitySize(aMeshName,anEntity,aFamilyName);
774       }
775       //Import mesh on entity
776       aResult += GetMeshOnEntitySize(aMeshName,anEntity);
777     }
778   }
779   if(MYDEBUG)
780     MESSAGE("GetSize - aResult = "<<float(aResult));
781   return aResult;
782 }
783
784
785 float VISU_Convertor_impl::GetMeshOnEntitySize(const std::string& theMeshName, 
786                                                const VISU::TEntity& theEntity,
787                                                const std::string& theFamilyName)
788 {
789   VISU::TMesh* pMesh = NULL;
790   VISU::TFamily* pFamily = NULL;
791   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
792   FindMeshOnEntity(theMeshName,pMesh,theEntity,pMeshOnEntity,theFamilyName,pFamily);
793   vtkIdType aPointsSize = 3*pMesh->myNbPoints*sizeof(VISU::TMesh::TCoord);
794   vtkIdType aNbCells, aCellsSize;
795   if(!pFamily){
796     aNbCells = pMeshOnEntity->myNbCells;
797     aCellsSize = pMeshOnEntity->myCellsSize;
798   }else{
799     aNbCells = pFamily->myNbCells;
800     aCellsSize = pFamily->myCellsSize;
801   }
802   vtkIdType aConnectivitySize = aCellsSize*sizeof(vtkIdType);
803   vtkIdType aTypesSize = aNbCells*sizeof(char);
804   vtkIdType aLocationsSize = aNbCells*sizeof(int);
805   float aNbCellsPerPoint = aCellsSize / aNbCells - 1;
806   vtkIdType aLinksSize = pMesh->myNbPoints * 
807     (vtkIdType(sizeof(vtkIdType)*aNbCellsPerPoint) + sizeof(vtkCellLinks::Link));
808   aLinksSize = 0;
809   vtkIdType aResult = aPointsSize + aConnectivitySize + aTypesSize + aLocationsSize + aLinksSize;
810   if(MYDEBUG){
811     if(MYVTKDEBUG){
812       MESSAGE("GetMeshOnEntitySize - aPointsSize = "<<float(aPointsSize));
813       MESSAGE("GetMeshOnEntitySize - aConnectivitySize = "<<float(aConnectivitySize));
814       MESSAGE("GetMeshOnEntitySize - aTypesSize = "<<float(aTypesSize));
815       MESSAGE("GetMeshOnEntitySize - aLocationsSize = "<<float(aLocationsSize));
816       MESSAGE("GetMeshOnEntitySize - aLinksSize = "<<float(aLinksSize));
817     }
818     MESSAGE("GetMeshOnEntitySize - aResult = "<<float(aResult)<<"; theMeshName = '"<<theMeshName<<
819             "'; theEntity = "<<theEntity<<"; theFamilyName = '"<<theFamilyName<<"'");
820   }
821   aResult = vtkIdType(aResult*ERR_SIZE_CALC);
822   return aResult;
823 }
824
825
826 void VISU_Convertor_impl::FindMeshOnGroup(const std::string& theMeshName, VISU::TMesh*& theMesh,
827                                           const std::string& theGroupName, VISU::TGroup*& theGroup)
828 {
829   FindMesh(theMeshName,theMesh);
830   VISU::TGroupMap& aGroupMap = theMesh->myGroupMap;
831   VISU::TGroupMap::iterator aGroupMapIter = aGroupMap.find(theGroupName);
832   if(aGroupMapIter == aGroupMap.end())
833     throw std::runtime_error(EXCEPT("FindMesh >> There is no the group in the mesh!!! - '%1'").arg(theGroupName.c_str()).latin1());
834   theGroup = &aGroupMapIter->second;
835 }
836
837
838 float VISU_Convertor_impl::GetMeshOnGroupSize(const std::string& theMeshName, 
839                                               const std::string& theGroupName)
840 {
841   VISU::TMesh* pMesh = NULL;
842   VISU::TGroup* pGroup = NULL;
843   FindMeshOnGroup(theMeshName,pMesh,theGroupName,pGroup);
844   vtkIdType aPointsSize = 3*pMesh->myNbPoints*sizeof(VISU::TMesh::TCoord);
845   vtkIdType aNbCells = pGroup->myNbCells, aCellsSize = pGroup->myCellsSize;
846   vtkIdType aConnectivityAndTypesSize = aCellsSize*sizeof(vtkIdType);
847   vtkIdType aLocationsSize = aNbCells*sizeof(int);
848   float aNbCellsPerPoint = aCellsSize / aNbCells - 1;
849   vtkIdType aLinksSize = pMesh->myNbPoints * 
850     (vtkIdType(sizeof(vtkIdType)*aNbCellsPerPoint) + sizeof(short));
851   aLinksSize = 0;
852   vtkIdType aResult = aPointsSize + aConnectivityAndTypesSize + aLocationsSize + aLinksSize;
853   if(MYDEBUG){
854     if(MYVTKDEBUG){
855       MESSAGE("GetMeshOnGroupSize - aPointsSize = "<<float(aPointsSize));
856       MESSAGE("GetMeshOnGroupSize - aConnectivityAndTypesSize = "<<float(aConnectivityAndTypesSize));
857       MESSAGE("GetMeshOnGroupSize - aLocationsSize = "<<float(aLocationsSize));
858       MESSAGE("GetMeshOnGroupSize - aLinksSize = "<<float(aLinksSize));
859     }
860     MESSAGE("GetMeshOnGroupSize - aResult = "<<float(aResult)<<"; theMeshName = '"
861             <<theMeshName<<"'; theGroupName = '"<<theGroupName<<"'");
862   }
863   aResult = vtkIdType(aResult*ERR_SIZE_CALC);
864   return aResult;
865 }
866
867
868 void VISU_Convertor_impl::FindField(const string& theMeshName, VISU::TMesh*& theMesh,
869                                     const VISU::TEntity& theEntity, 
870                                     VISU::TMeshOnEntity*& theMeshOnEntity,
871                                     VISU::TMeshOnEntity*& theVTKMeshOnEntity,
872                                     const string& theFieldName, VISU::TField*& theField)
873 {
874   VISU::TFamily* pFamily = NULL;
875   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
876   FindMeshOnEntity(theMeshName,theMesh,theEntity,pMeshOnEntity,"",pFamily);
877   theMeshOnEntity = pMeshOnEntity;
878   VISU::TMeshOnEntityMap& aMeshOnEntityMap = theMesh->myMeshOnEntityMap;
879   if(theEntity == VISU::NODE_ENTITY){
880     if(aMeshOnEntityMap.find(VISU::CELL_ENTITY) != aMeshOnEntityMap.end())
881       pMeshOnEntity = &aMeshOnEntityMap[VISU::CELL_ENTITY];
882     else if(aMeshOnEntityMap.find(VISU::FACE_ENTITY) != aMeshOnEntityMap.end())
883       pMeshOnEntity = &aMeshOnEntityMap[VISU::FACE_ENTITY];
884     else if(aMeshOnEntityMap.find(VISU::NODE_ENTITY) != aMeshOnEntityMap.end())
885       pMeshOnEntity = &aMeshOnEntityMap[VISU::EDGE_ENTITY];
886   }
887   theVTKMeshOnEntity = pMeshOnEntity;
888   VISU::TFieldMap& aFieldMap = theMeshOnEntity->myFieldMap;
889   if(aFieldMap.find(theFieldName) == aFieldMap.end())
890     throw std::runtime_error(EXCEPTION("FindField >> There is no field on the mesh!!!"));
891   theField = &aFieldMap[theFieldName];
892 }
893
894
895 float VISU_Convertor_impl::GetFieldOnMeshSize(const std::string& theMeshName, 
896                                               const VISU::TEntity& theEntity,
897                                               const std::string& theFieldName)
898 {
899   VISU::TMesh* pMesh = NULL;
900   VISU::TField* pField = NULL;
901   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
902   FindField(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
903             theFieldName,pField);
904   float aMeshSize = GetMeshOnEntitySize(theMeshName,pVTKMeshOnEntity->myEntity);
905   float aFieldOnMeshSize = float(pField->myDataSize*sizeof(float)*pField->myNbValField * ERR_SIZE_CALC);
906   float aResult = aMeshSize + aFieldOnMeshSize;
907   if(MYDEBUG){
908     if(MYVTKDEBUG)
909       MESSAGE("GetFieldOnMeshSize - aFieldOnMeshSize = "<<float(aFieldOnMeshSize));
910     MESSAGE("GetFieldOnMeshSize - aResult = "<<float(aResult)<<"; theMeshName = '"<<theMeshName<<
911             "'; theEntity = "<<theEntity<<"; theFieldName = '"<<theFieldName<<"'");
912   }
913   return aResult;
914 }
915
916
917 void VISU_Convertor_impl::FindTimeStamp(const std::string& theMeshName, VISU::TMesh*& theMesh,
918                                         const VISU::TEntity& theEntity, 
919                                         VISU::TMeshOnEntity*& theMeshOnEntity,
920                                         VISU::TMeshOnEntity*& theVTKMeshOnEntity,
921                                         const std::string& theFieldName, VISU::TField*& theField,
922                                         int theStampsNum, VISU::TField::TValForTime*& theValForTime)
923 {
924   FindField(theMeshName,theMesh,theEntity,theMeshOnEntity,theVTKMeshOnEntity,theFieldName,theField);
925   VISU::TField::TValField& aValField = theField->myValField;
926   if(aValField.find(theStampsNum) == aValField.end())
927     throw std::runtime_error(EXCEPTION("FindTimeStamp >> There is no field with the timestamp!!!"));
928   theValForTime = &aValField[theStampsNum];
929 }
930
931
932 float VISU_Convertor_impl::GetTimeStampSize(const std::string& theMeshName, 
933                                             const VISU::TEntity& theEntity,
934                                             const std::string& theFieldName,
935                                             int theStampsNum)
936 {
937   VISU::TMesh* pMesh = NULL;
938   VISU::TField* pField = NULL;
939   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
940   VISU::TField::TValForTime* pValForTime = NULL;
941   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
942                 theFieldName,pField,theStampsNum,pValForTime);
943   float aMeshSize = GetMeshOnEntitySize(theMeshName,pVTKMeshOnEntity->myEntity);
944   float aTimeStampSize = float(pField->myDataSize*sizeof(float) * ERR_SIZE_CALC);
945   float aResult = aMeshSize + aTimeStampSize;
946   if(MYDEBUG){
947     if(MYVTKDEBUG)
948       MESSAGE("GetTimeStampSize - aTimeStampSize = "<<float(aTimeStampSize));
949     MESSAGE("GetTimeStampSize - aResult = "<<float(aResult)<<"; theMeshName = '"<<theMeshName<<"'; theEntity = "<<theEntity<<
950             "; theFieldName = '"<<theFieldName<<"'; theStampsNum = "<<theStampsNum);
951   }
952   return aResult;
953 }
954
955
956 const VISU::TField& 
957 VISU_Convertor_impl::GetField(const string& theMeshName, 
958                               VISU::TEntity theEntity, 
959                               const string& theFieldName) 
960 {
961   VISU::TMesh* pMesh = NULL;
962   VISU::TField* pField = NULL;
963   VISU::TFamily* pFamily = NULL;
964   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
965   FindField(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
966             theFieldName,pField);
967   return *pField;
968 }
969
970
971 const VISU::TField::TValForTime& 
972 VISU_Convertor_impl::GetTimeStamp(const std::string& theMeshName, 
973                                   const VISU::TEntity& theEntity,
974                                   const std::string& theFieldName,
975                                   int theStampsNum)
976 {
977   VISU::TMesh* pMesh = NULL;
978   VISU::TField* pField = NULL;
979   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
980   VISU::TField::TValForTime* pValForTime = NULL;
981   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
982                 theFieldName,pField,theStampsNum,pValForTime);
983   return *pValForTime;
984 }
985
986