Salome HOME
DCQ : Merge with Ecole Ete a6.
[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 namespace{
75   void GetPoints(VISU::TVTKSource& theStorage, VISU::TMesh& theMesh) 
76   {
77     vtkPoints* aPoints = theMesh.myPoints.GetPointer();
78     if(!aPoints){
79       aPoints = vtkPoints::New();
80       if(MYVTKDEBUG) aPoints->DebugOn();
81       const VISU::TMesh::TPointsCoord& anArray = theMesh.myPointsCoord;
82       vtkIdType iEnd = theMesh.myPointsCoord.size();
83       vtkIdType aNbPoints = iEnd / theMesh.myDim;
84       aPoints->SetNumberOfPoints(aNbPoints);
85       if(MYDEBUG) 
86         MESSAGE("GetPoints - aNbPoints = "<<aNbPoints<<"; myDim = "<<theMesh.myDim);
87       switch(theMesh.myDim) {
88       case 1:
89         for (vtkIdType i = 0, j = 0; i < iEnd; i += theMesh.myDim, j++) 
90           aPoints->SetPoint(j,anArray[i],0.0,0.0);
91         break;
92       case 2:
93         for (vtkIdType i = 0, j = 0; i < iEnd; i += theMesh.myDim, j++) 
94           aPoints->SetPoint(j,anArray[i],anArray[i+1],0.0);
95         break;
96       case 3: 
97         for (vtkIdType i = 0, j = 0; i < iEnd; i += theMesh.myDim, j++) 
98           aPoints->SetPoint(j,anArray[i],anArray[i+1],anArray[i+2]);
99         break;
100       }
101       theMesh.myPoints = aPoints;
102     }
103     theStorage->SetPoints(aPoints);
104   }
105   
106   
107   inline void PrintCells(int& theStartId,
108                          vtkCellArray* theConnectivity, 
109                          const VISU::TMeshOnEntity::TConnect& theVector)
110   {
111     vtkIdList *anIdList = vtkIdList::New();
112     int kEnd = theVector.size();
113     anIdList->SetNumberOfIds(kEnd);
114     for(int k = 0; k < kEnd; k++)
115       anIdList->SetId(k,theVector[k]);
116     theConnectivity->InsertNextCell(anIdList);
117     anIdList->Delete();
118   }
119
120   void GetCellsOnEntity(VISU::TVTKSource& theStorage,
121                         const VISU::TMeshOnEntity& theMeshOnEntity, 
122                         const string& theFamilyName) 
123   {
124     //Check on existing family
125     const VISU::TFamily* pFamily = VISU::GetFamily(theMeshOnEntity,theFamilyName);
126     bool isFamilyPresent = (pFamily != NULL);
127     const VISU::TFamily& aFamily = *pFamily;
128     //Main part of code
129     pair<int,int> aCellsDim = theMeshOnEntity.GetCellsDims(theFamilyName);
130     int aNbCells = aCellsDim.first, aCellsSize = aCellsDim.second;
131     vtkCellArray* aConnectivity = vtkCellArray::New();
132     aConnectivity->Allocate(aCellsSize,0);
133     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
134     aCellTypesArray->SetNumberOfComponents(1);
135     aCellTypesArray->SetNumberOfTuples(aNbCells);
136     if(MYDEBUG) MESSAGE("GetCellsOnEntity - isFamilyPresent = "<<isFamilyPresent);
137     const VISU::TMeshOnEntity::TCellsConn &aCellsConn = theMeshOnEntity.myCellsConn;
138     VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
139     for(int i = 0, j = 0; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
140       const VISU::TMeshOnEntity::TConnForCellType& anArray = aCellsConnIter->second;
141       int aVtkType = aCellsConnIter->first;
142       if(MYDEBUG) 
143         MESSAGE("GetCellsOnEntity - aVtkType = "<<aVtkType<<"; anArray.size() = "<<anArray.size());
144       if(!isFamilyPresent)
145         for(int k = 0, kEnd = anArray.size(); k < kEnd; k++, i++){
146           PrintCells(i,aConnectivity,anArray[k]);
147           aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
148         }
149       else{
150         const VISU::TFamily::TSubMesh& aSubMesh = aFamily.mySubMesh;
151         if(aSubMesh.empty()) 
152           throw std::runtime_error(EXCEPTION("GetCells >> There is no elements on the family !!!"));
153         VISU::TFamily::TSubMesh::const_iterator aSubMeshIter = aSubMesh.find(aVtkType);
154         if(aSubMeshIter == aSubMesh.end()) continue;
155         const VISU::TFamily::TSubMeshOnCellType& aSubMeshOnCellType = aSubMeshIter->second;
156         if(MYDEBUG) 
157           MESSAGE("GetCellsOnEntity - aSubMeshOnCellType.size() = "<<aSubMeshOnCellType.size());
158         VISU::TFamily::TSubMeshOnCellType::const_iterator aSubMeshOnCellTypeIter = aSubMeshOnCellType.begin();
159         for(; aSubMeshOnCellTypeIter != aSubMeshOnCellType.end(); aSubMeshOnCellTypeIter++, i++){
160           PrintCells(i,aConnectivity,anArray[*aSubMeshOnCellTypeIter]);
161           aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
162         }
163       }
164     }
165     vtkIdType *pts = 0, npts = 0;
166     vtkIntArray* aCellLocationsArray = vtkIntArray::New();
167     aCellLocationsArray->SetNumberOfComponents(1);
168     aCellLocationsArray->SetNumberOfTuples(aNbCells);
169     aConnectivity->InitTraversal();
170     for(int i=0; aConnectivity->GetNextCell(npts,pts); i++)
171       aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
172     theStorage->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
173     if(MYVTKDEBUG) aConnectivity->DebugOn();
174     aCellLocationsArray->Delete();
175     aCellTypesArray->Delete();
176     aConnectivity->Delete();
177   } 
178   
179   
180   void GetCellsOnGroup(VISU::TVTKSource& theStorage,
181                        const VISU::TMesh& theMesh,
182                        const VISU::TFamilyAndEntitySet& theFamilyAndEntitySet) 
183   {
184     //Calculate dimentions of the group
185     int aNbCells = 0, aCellsSize = 0;
186     VISU::TFamilyAndEntitySet::const_iterator aFamilyAndEntitySetIter = theFamilyAndEntitySet.begin();
187     for(; aFamilyAndEntitySetIter != theFamilyAndEntitySet.end(); aFamilyAndEntitySetIter++){
188       const VISU::TFamilyAndEntity& aFamilyAndEntity = *aFamilyAndEntitySetIter;
189       const string& aFamilyName = aFamilyAndEntity.first;
190       const VISU::TEntity& anEntity = aFamilyAndEntity.second;
191       const VISU::TMeshOnEntity& aMeshOnEntity = theMesh.myMeshOnEntityMap.find(anEntity)->second;
192       pair<int,int> aCellsDim = aMeshOnEntity.GetCellsDims(aFamilyName);
193       aNbCells += aCellsDim.first;
194       aCellsSize += aCellsDim.second;
195     }
196     vtkCellArray* aConnectivity = vtkCellArray::New();
197     aConnectivity->Allocate(aCellsSize,0);
198     vtkUnsignedCharArray* aCellTypesArray = vtkUnsignedCharArray::New();
199     aCellTypesArray->SetNumberOfComponents(1);
200     aCellTypesArray->SetNumberOfTuples(aNbCells);
201     aFamilyAndEntitySetIter = theFamilyAndEntitySet.begin();
202     for(int i = 0, j = 0; aFamilyAndEntitySetIter != theFamilyAndEntitySet.end(); aFamilyAndEntitySetIter++){
203       const VISU::TFamilyAndEntity& aFamilyAndEntity = *aFamilyAndEntitySetIter;
204       const string& aFamilyName = aFamilyAndEntity.first;
205       const VISU::TEntity& anEntity = aFamilyAndEntity.second;
206       const VISU::TMeshOnEntity& aMeshOnEntity = theMesh.myMeshOnEntityMap.find(anEntity)->second;
207       const VISU::TFamily& aFamily = *(VISU::GetFamily(aMeshOnEntity,aFamilyName));
208       const VISU::TMeshOnEntity::TCellsConn &aCellsConn = aMeshOnEntity.myCellsConn;
209       VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
210       for(; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
211         const VISU::TMeshOnEntity::TConnForCellType& anArray = aCellsConnIter->second;
212         int aVtkType = aCellsConnIter->first;
213         if(MYDEBUG) 
214           MESSAGE("GetCellsOnGroup - aVtkType = "<<aVtkType<<"; anArray.size() = "<<anArray.size());
215         const VISU::TFamily::TSubMesh& aSubMesh = aFamily.mySubMesh;
216         if(aSubMesh.empty()) 
217           throw std::runtime_error(EXCEPTION("GetCells >> There is no elements on the family !!!"));
218         VISU::TFamily::TSubMesh::const_iterator aSubMeshIter = aSubMesh.find(aVtkType);
219         if(aSubMeshIter == aSubMesh.end()) continue;
220         const VISU::TFamily::TSubMeshOnCellType& aSubMeshOnCellType = aSubMeshIter->second;
221         if(MYDEBUG) 
222           MESSAGE("GetCellsOnGroup - aSubMeshOnCellType.size() = "<<aSubMeshOnCellType.size());
223         VISU::TFamily::TSubMeshOnCellType::const_iterator aSubMeshOnCellTypeIter = aSubMeshOnCellType.begin();
224         for(; aSubMeshOnCellTypeIter != aSubMeshOnCellType.end(); aSubMeshOnCellTypeIter++, i++){
225           PrintCells(i,aConnectivity,anArray[*aSubMeshOnCellTypeIter]);
226           aCellTypesArray->SetValue(j++,(unsigned char)aVtkType);
227         }
228       }
229     }
230     vtkIdType *pts = 0, npts = 0;
231     vtkIntArray* aCellLocationsArray = vtkIntArray::New();
232     aCellLocationsArray->SetNumberOfComponents(1);
233     aCellLocationsArray->SetNumberOfTuples(aNbCells);
234     aConnectivity->InitTraversal();
235     for(int i = 0; aConnectivity->GetNextCell(npts,pts); i++)
236       aCellLocationsArray->SetValue(i,aConnectivity->GetTraversalLocation(npts));
237     theStorage->SetCells(aCellTypesArray,aCellLocationsArray,aConnectivity);
238     aCellLocationsArray->Delete();
239     aCellTypesArray->Delete();
240     aConnectivity->Delete();
241   } 
242   
243   
244   void InitProfile(VISU::TVTKExtractFilter& theFilter,
245                    const VISU::TMeshOnEntity& theMeshOnEntity, 
246                    const VISU::TField::TValForTime& theValForTime)
247   {
248     const VISU::TField::TValForCells& aValForCells = theValForTime.myValForCells;
249     const VISU::TMeshOnEntity::TCellsConn &aCellsConn = theMeshOnEntity.myCellsConn;
250     VISU::TMeshOnEntity::TCellsConn::const_iterator aCellsConnIter = aCellsConn.begin();
251     for(; aCellsConnIter != aCellsConn.end(); aCellsConnIter++){
252       const vtkIdType& aCellType = aCellsConnIter->first;
253       if(aValForCells.find(aCellType) == aValForCells.end())
254         theFilter->RemoveCellsWithType(aCellType);
255     }
256   }
257
258
259   void GetValsOnTimeStamp(vtkFloatArray *theFloatArray, 
260                           const vtkIdType& theNumberOfTuples,
261                           const std::string& theFieldName,
262                           const VISU::TField& theField,
263                           const VISU::TField::TValForTime& theValForTime)
264   {
265     //theFloatArray->DebugOn();
266     theFloatArray->SetNumberOfTuples(theNumberOfTuples);
267     theFloatArray->SetName(theFieldName.c_str());
268     if(MYDEBUG) MESSAGE("GetValsOnTimeStamp - theNumberOfTuples = "<<theNumberOfTuples);
269     const VISU::TField::TValForCells& aValForCells = theValForTime.myValForCells;
270     VISU::TField::TValForCells::const_iterator aValForCellsIter = aValForCells.begin();
271     for(int k = 0; aValForCellsIter != aValForCells.end(); aValForCellsIter++) {
272       const VISU::TField::TValForCellsWithType& anArray = aValForCellsIter->second;
273       int iEnd = anArray.size()/theField.myNbComp;
274       int aVtkType = aValForCellsIter->first;
275       if(MYDEBUG) MESSAGE("GetValsOnTimeStamp -  iEnd = "<<iEnd<<"; aVtkType = "<<aVtkType);
276       switch(theField.myNbComp) {
277       case 1:
278         for (int i = 0; i < iEnd; i++) 
279           theFloatArray->SetTuple1(k++,anArray[i]);
280         break;
281       case 2:
282         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*2)
283           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+1],0.0);
284         break;
285       case 3:
286         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*3)
287           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+1],anArray[ji+2]);
288         break;
289       case 4:
290         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*4)
291           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+2],0.0);
292         break;
293       case 6:
294         for (int i = 0, ji = 0; i < iEnd; ++i, ji = i*4)
295           theFloatArray->SetTuple3(k++,anArray[ji],anArray[ji+2],anArray[ji+5]);
296         break;
297       default:
298         throw std::runtime_error(EXCEPTION("GetValsOnTimeStamp - There is no an algorithm for representation of the field !!!"));
299       }
300     }
301   }
302
303   string GenerateFieldName(const VISU::TField& theField,
304                            const VISU::TField::TValForTime& theValForTime)
305   {
306     const VISU::TField::TTime& aTime = theValForTime.myTime;
307     string aFieldName = theField.myMeshName + ", " + theField.myName + ": " + 
308       VISU_Convertor::GenerateName(aTime);
309     return aFieldName;
310   }
311
312   void GetTimeStamp(VISU::TVTKSource& theStorage,
313                     const VISU::TMesh& theMesh,
314                     const VISU::TMeshOnEntity& theMeshOnEntity, 
315                     const VISU::TField& theField, 
316                     const VISU::TField::TValForTime& theValForTime)
317   {
318     int aNumberOfTuples = theField.myDataSize/theField.myNbComp;
319     string aFieldName = GenerateFieldName(theField,theValForTime);
320     if(MYDEBUG) 
321       MESSAGE("GetTimeStamp(TVTKSource) - aFieldName = "<<aFieldName<<
322               "; aNumberOfTuples = "<<aNumberOfTuples);
323
324     vtkDataSetAttributes* aDataSetAttributes;
325     switch(theField.myEntity){
326     case VISU::NODE_ENTITY : 
327       aDataSetAttributes = theStorage->GetPointData();
328       break;
329     default: 
330       aDataSetAttributes = theStorage->GetCellData();
331     }
332
333     vtkFloatArray *aFloatArray = vtkFloatArray::New();
334     switch(theField.myNbComp) {
335     case 1:
336       aFloatArray->SetNumberOfComponents(1);
337       aDataSetAttributes->SetScalars(aFloatArray);
338       break;
339     default:
340       aFloatArray->SetNumberOfComponents(3);
341       aDataSetAttributes->SetVectors(aFloatArray);
342     }
343
344     GetValsOnTimeStamp(aFloatArray,aNumberOfTuples,aFieldName,theField,theValForTime);
345   }
346
347   void GetTimeStamp(VISU::TVTKAttribyteFilter& theAttribyteFilter,
348                     VISU::TVTKMergetFilter& theMergeFilter,
349                     VISU::TVTKExtractFilter& theExtractFilter,
350                     const VISU::TMesh& theMesh,
351                     const VISU::TMeshOnEntity& theMeshOnEntity, 
352                     const VISU::TField& theField, 
353                     const VISU::TField::TValForTime& theValForTime)
354   {
355     int aNumberOfTuples = theField.myDataSize/theField.myNbComp;
356     string aFieldName = GenerateFieldName(theField,theValForTime);
357     if(MYDEBUG) 
358       MESSAGE("GetTimeStamp(TVTKAttribyteFilter) - aFieldName = "<<aFieldName<<
359               "; aNumberOfTuples = "<<aNumberOfTuples);
360
361     vtkDataObject* aDataObject = vtkDataObject::New();
362     theMergeFilter->SetDataObject(aDataObject);
363     aDataObject->Delete();
364
365     theMergeFilter->SetInput(theExtractFilter->GetOutput());
366     theAttribyteFilter->SetInput(theMergeFilter->GetOutput());
367
368     switch(theField.myEntity){
369     case VISU::NODE_ENTITY : 
370       theMergeFilter->SetOutputFieldToPointDataField();
371       theAttribyteFilter->SetInputFieldToPointDataField();
372       theAttribyteFilter->SetOutputAttributeDataToPointData();
373       break;
374     default: 
375       theMergeFilter->SetOutputFieldToCellDataField();
376       theAttribyteFilter->SetInputFieldToCellDataField();
377       theAttribyteFilter->SetOutputAttributeDataToCellData();
378     }
379
380     vtkFloatArray *aFloatArray = vtkFloatArray::New();
381     switch(theField.myNbComp) {
382     case 1:
383       aFloatArray->SetNumberOfComponents(1);
384       theAttribyteFilter->SetScalarComponent(0,aFieldName.c_str(),0);
385       break;
386     default:
387       aFloatArray->SetNumberOfComponents(3);
388       theAttribyteFilter->SetVectorComponent(0,aFieldName.c_str(),0);
389       theAttribyteFilter->SetVectorComponent(1,aFieldName.c_str(),1);
390       theAttribyteFilter->SetVectorComponent(2,aFieldName.c_str(),2);
391     }
392
393     vtkFieldData* aFieldData = aDataObject->GetFieldData();
394     aFieldData->AddArray(aFloatArray);
395     aFloatArray->Delete();
396
397     GetValsOnTimeStamp(aFloatArray,aNumberOfTuples,aFieldName,theField,theValForTime);
398   }
399 }
400
401 VISU_Convertor_impl::VISU_Convertor_impl() {
402   myIsDone = false;
403 }
404
405 VISU_Convertor_impl::~VISU_Convertor_impl() {}
406
407 VISU_Convertor::TOutput* 
408 VISU_Convertor_impl::GetMeshOnEntity(const string& theMeshName, 
409                                      const VISU::TEntity& theEntity,
410                                      const string& theFamilyName)
411 {
412   if(MYDEBUG) 
413     MESSAGE("GetMeshOnEntity - theMeshName = '"<<theMeshName<<
414             "'; theEntity = "<<theEntity<<"; theFamilyName = '"<<theFamilyName<<"'");
415   //Cheching possibility do the query
416   VISU::TMesh* pMesh = NULL;
417   VISU::TFamily* pFamily = NULL;
418   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
419   FindMeshOnEntity(theMeshName,pMesh,theEntity,pMeshOnEntity,theFamilyName,pFamily);
420   VISU::TMesh& aMesh = *pMesh;
421   VISU::TMeshOnEntity& aMeshOnEntity = *pMeshOnEntity;
422   VISU::TVTKSource* pSource;
423   if(pFamily != NULL)
424     pSource = &(pFamily->myStorage);
425   else
426     pSource = &(aMeshOnEntity.myStorage);
427   VISU::TVTKSource& aSource = *pSource;
428   //Main part of code
429   try{
430     if(aSource.GetPointer() == NULL){
431       aSource = TOutput::New();
432       aSource->Delete();
433       if(MYVTKDEBUG) aSource->DebugOn();
434       LoadMeshOnEntity(aMeshOnEntity,theFamilyName);
435       GetPoints(aSource,aMesh);
436       GetCellsOnEntity(aSource,aMeshOnEntity,theFamilyName);
437       if(MYDEBUGWITHFILES){
438         string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
439         string aFamilyName = QString(theFamilyName.c_str()).simplifyWhiteSpace().latin1();
440         string aFileName = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
441         aFileName += aMeshName + dtos("-%d-",int(theEntity)) + aFamilyName + "-Conv.vtk";
442         VISU::WriteToFile(aSource.GetPointer(),aFileName);
443       }
444     }
445     if(MYVTKDEBUG){
446       GetMeshOnEntitySize(theMeshName,theEntity,theFamilyName);
447       vtkUnstructuredGrid* aDataSet = aSource.GetPointer();
448       aDataSet->Update();
449       MESSAGE("GetMeshOnEntity - GetPoints() = "<<float(aDataSet->GetPoints()->GetActualMemorySize()*1000));
450       MESSAGE("GetMeshOnEntity - GetCells() = "<<float(aDataSet->GetCells()->GetActualMemorySize()*1000));
451       MESSAGE("GetMeshOnEntity - GetCellTypesArray() = "<<float(aDataSet->GetCellTypesArray()->GetActualMemorySize()*1000));
452       MESSAGE("GetMeshOnEntity - GetCellLocationsArray() = "<<float(aDataSet->GetCellLocationsArray()->GetActualMemorySize()*1000));
453       aDataSet->BuildLinks();
454       MESSAGE("GetMeshOnEntity - GetCellLinks() = "<<float(aDataSet->GetCellLinks()->GetActualMemorySize()*1000));
455       MESSAGE("GetMeshOnEntity - GetActualMemorySize() = "<<float(aDataSet->GetActualMemorySize()*1000));
456     }
457   }catch(...){
458     aSource = vtkSmartPointerBase();
459     throw;
460   }
461   return aSource.GetPointer();
462 }
463
464 VISU_Convertor::TOutput* 
465 VISU_Convertor_impl::GetMeshOnGroup(const string& theMeshName, 
466                                     const string& theGroupName)
467 {
468   if(MYDEBUG) MESSAGE("GetMeshOnGroup - theMeshName = '"<<theMeshName<<
469                       "'; theGroupName = '"<<theGroupName<<"'");
470   //Cheching possibility do the query
471   VISU::TMesh* pMesh = NULL;
472   VISU::TGroup* pGroup = NULL;
473   FindMeshOnGroup(theMeshName,pMesh,theGroupName,pGroup);
474   VISU::TMesh& aMesh = *pMesh;
475   VISU::TGroup& aGroup = *pGroup;
476   const VISU::TFamilyAndEntitySet& aFamilyAndEntitySet = aGroup.myFamilyAndEntitySet;
477   VISU::TVTKSource& aSource = aGroup.myStorage;
478   //Main part of code
479   try{
480     if(aSource.GetPointer() == NULL){
481       aSource = TOutput::New();
482       aSource->Delete();
483       LoadMeshOnGroup(aMesh,aFamilyAndEntitySet);
484       GetPoints(aSource,aMesh);
485       GetCellsOnGroup(aSource,aMesh,aFamilyAndEntitySet);
486       if(MYDEBUGWITHFILES){
487         string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
488         string aGroupName = QString(theGroupName.c_str()).simplifyWhiteSpace().latin1();
489         string aFileName = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
490         aFileName += aMeshName + "-" + aGroupName + "-Conv.vtk";
491         VISU::WriteToFile(aSource.GetPointer(),aFileName);
492       }
493     }
494   }catch(...){
495     aSource = vtkSmartPointerBase();
496     throw;
497   }
498   return aSource.GetPointer();
499 }
500
501 VISU_Convertor::TOutput* 
502 VISU_Convertor_impl::GetTimeStampOnMesh(const string& theMeshName, 
503                                         const VISU::TEntity& theEntity,
504                                         const string& theFieldName,
505                                         int theStampsNum)
506 {
507   if(MYDEBUG){
508     MESSAGE("GetTimeStampOnMesh - theMeshName = '"<<theMeshName<<"; theEntity = "<<theEntity);
509     MESSAGE("GetTimeStampOnMesh - theFieldName = '"<<theFieldName<<"'; theStampsNum = "<<theStampsNum);
510   }
511   //Cheching possibility do the query
512   VISU::TMesh* pMesh = NULL;
513   VISU::TField* pField = NULL;
514   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
515   VISU::TField::TValForTime* pValForTime = NULL;
516   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
517                 theFieldName,pField,theStampsNum,pValForTime);
518   VISU::TMesh& aMesh = *pMesh;
519   VISU::TMeshOnEntity& aMeshOnEntity = *pMeshOnEntity;
520   VISU::TField& aField = *pField;
521   VISU::TField::TValForTime& aValForTime = *pValForTime;
522   VISU::TVTKAttribyteFilter& anAttribyteFilter = aValForTime.myAttribyteFilter;
523   VISU::TVTKSource& aSource = aValForTime.myStorage;
524   TOutput* anOutput = NULL;
525   //Main part of code
526   try{
527     if(aSource.GetPointer())
528       return aSource.GetPointer();
529     else if(anAttribyteFilter.GetPointer())
530       return anAttribyteFilter->GetUnstructuredGridOutput();
531     else{
532       LoadFieldOnMesh(aMesh,aMeshOnEntity,aField,aValForTime);
533
534       VISU::TVTKExtractFilter& anExtractFilter = aField.myExtractFilter;
535       if(anExtractFilter.GetPointer() == NULL){
536         anExtractFilter = VISU_ExtractUnstructuredGrid::New();
537         anExtractFilter->Delete();
538         //anExtractFilter->DebugOn();
539         try{
540           LoadMeshOnEntity(*pVTKMeshOnEntity);
541         }catch(std::runtime_error& exc){
542           pVTKMeshOnEntity = pMeshOnEntity;
543           MESSAGE("Follow exception was accured :\n"<<exc.what());
544         }catch(...){
545           pVTKMeshOnEntity = pMeshOnEntity;
546           MESSAGE("Unknown exception was accured!");
547         }
548         GetMeshOnEntity(pVTKMeshOnEntity->myMeshName,pVTKMeshOnEntity->myEntity);
549         
550         anExtractFilter->SetInput(pVTKMeshOnEntity->myStorage.GetPointer());
551         ::InitProfile(anExtractFilter,*pMeshOnEntity,aValForTime);
552       }      
553       if(!anExtractFilter->IsRemoving()){
554         aSource = TOutput::New();
555         aSource->Delete();
556         aSource->ShallowCopy(pVTKMeshOnEntity->myStorage.GetPointer());
557         ::GetTimeStamp(aSource,aMesh,*pVTKMeshOnEntity,aField,aValForTime);
558         anOutput = aSource.GetPointer();
559       }else{
560         anAttribyteFilter = vtkFieldDataToAttributeDataFilter::New();
561         anAttribyteFilter->Delete();
562         //anAttribyteFilter->DebugOn();
563         
564         VISU::TVTKMergetFilter& aMergeFilter = aValForTime.myMergeFilter;
565         aMergeFilter = vtkMergeDataObjectFilter::New();
566         aMergeFilter->Delete();
567         //aMergeFilter->DebugOn();
568
569         ::GetTimeStamp(anAttribyteFilter,aMergeFilter,anExtractFilter,
570                        aMesh,*pVTKMeshOnEntity,aField,aValForTime);
571         anOutput = anAttribyteFilter->GetUnstructuredGridOutput();
572       }
573       if(MYDEBUGWITHFILES){
574         string aMeshName = QString(theMeshName.c_str()).simplifyWhiteSpace().latin1();
575         string aFieldName = QString(theFieldName.c_str()).simplifyWhiteSpace().latin1();
576         string aPrefix = string("/users/")+getenv("USER")+"/"+getenv("USER")+"-";
577         string aFileName = aPrefix + aMeshName + dtos("-%d-",int(theEntity)) + 
578           aFieldName + dtos("-%d",theStampsNum) + "-Conv.vtk";
579         VISU::WriteToFile(anOutput,aFileName);
580       }
581       if(MYVTKDEBUG){
582         GetTimeStampSize(theMeshName,theEntity,theFieldName,theStampsNum);
583         vtkUnstructuredGrid *aDataSet = anAttribyteFilter->GetUnstructuredGridOutput();
584         aDataSet->Update();
585         if(theEntity == VISU::NODE_ENTITY)
586           MESSAGE("GetTimeStampOnMesh - GetData() = "<<float(aDataSet->GetPointData()->GetActualMemorySize()*1000));
587         else
588           MESSAGE("GetMeshOnEntity - GetData() = "<<float(aDataSet->GetCellData()->GetActualMemorySize()*1000));
589         MESSAGE("GetTimeStampOnMesh - GetActualMemorySize() = "<<float(aDataSet->GetActualMemorySize()*1000));
590       }
591     }
592   }catch(...){
593     aSource = vtkSmartPointerBase();
594     anAttribyteFilter = vtkSmartPointerBase();
595     throw;
596   }
597   return anOutput;
598 }
599
600 void VISU_Convertor_impl::FindMesh(const string& theMeshName, VISU::TMesh*& theMesh)
601 {
602   GetMeshMap();
603   if(myMeshMap.find(theMeshName) == myMeshMap.end())
604     throw std::runtime_error(EXCEPT("FindMesh >> There is no mesh with the name - '%1'!!!").
605                              arg(theMeshName.c_str()).latin1());
606   theMesh = &myMeshMap[theMeshName];
607 }
608
609
610 void VISU_Convertor_impl::FindMeshOnEntity(const string& theMeshName, 
611                                            VISU::TMesh*& theMesh,
612                                            const VISU::TEntity& theEntity, 
613                                            VISU::TMeshOnEntity*& theMeshOnEntity,
614                                            const string& theFamilyName, 
615                                            VISU::TFamily*& theFamily)
616 {
617   FindMesh(theMeshName,theMesh);
618   VISU::TMeshOnEntityMap& aMeshOnEntityMap = theMesh->myMeshOnEntityMap;
619   if(aMeshOnEntityMap.find(theEntity) == aMeshOnEntityMap.end())
620     throw std::runtime_error(EXCEPT("FindMeshOnEntity >> There is no mesh on the entity - %1!!!").
621                              arg(int(theEntity)).latin1());
622   theMeshOnEntity = &aMeshOnEntityMap[theEntity];
623   theFamily = VISU::GetFamily(*theMeshOnEntity,theFamilyName);
624 }
625
626
627 float VISU_Convertor_impl::GetSize() {
628   float aResult = 0.0;
629   const VISU::TMeshMap& aMeshMap = GetMeshMap();
630   VISU::TMeshMap::const_iterator aMeshMapIter = aMeshMap.begin();
631   for(; aMeshMapIter != aMeshMap.end(); aMeshMapIter++){
632     const string& aMeshName = aMeshMapIter->first;
633     const VISU::TMesh& aMesh = aMeshMapIter->second;
634     const VISU::TMeshOnEntityMap& aMeshOnEntityMap = aMesh.myMeshOnEntityMap;
635     VISU::TMeshOnEntityMap::const_iterator aMeshOnEntityMapIter;
636     //Import fields
637     aMeshOnEntityMapIter = aMeshOnEntityMap.begin();
638     for(; aMeshOnEntityMapIter != aMeshOnEntityMap.end(); aMeshOnEntityMapIter++){
639       const VISU::TEntity& anEntity = aMeshOnEntityMapIter->first;
640       const VISU::TMeshOnEntity& aMeshOnEntity = aMeshOnEntityMapIter->second;
641       const VISU::TFieldMap& aFieldMap = aMeshOnEntity.myFieldMap;
642       VISU::TFieldMap::const_iterator aFieldMapIter = aFieldMap.begin();
643       for(; aFieldMapIter != aFieldMap.end(); aFieldMapIter++){
644         const string& aFieldName = aFieldMapIter->first;
645         const VISU::TField& aField = aFieldMapIter->second;
646         const VISU::TField::TValField& aValField = aField.myValField;
647         VISU::TField::TValField::const_iterator aValFieldIter = aValField.begin();
648         for(; aValFieldIter != aValField.end(); aValFieldIter++){
649           int aTimeStamp = aValFieldIter->first;
650           aResult += GetTimeStampSize(aMeshName,anEntity,aFieldName,aTimeStamp);
651         }
652       }
653       //Importing groups
654       const VISU::TGroupMap& aGroupMap = aMesh.myGroupMap;
655       VISU::TGroupMap::const_iterator aGroupMapIter = aGroupMap.begin();
656       for(; aGroupMapIter != aGroupMap.end(); aGroupMapIter++){
657         const string& aGroupName = aGroupMapIter->first;
658         aResult += GetMeshOnGroupSize(aMeshName,aGroupName);
659       }
660       //Import families
661       const VISU::TFamilyMap& aFamilyMap = aMeshOnEntity.myFamilyMap;
662       VISU::TFamilyMap::const_iterator aFamilyMapIter = aFamilyMap.begin();
663       for(; aFamilyMapIter != aFamilyMap.end(); aFamilyMapIter++){
664         const string& aFamilyName = aFamilyMapIter->first;
665         aResult += GetMeshOnEntitySize(aMeshName,anEntity,aFamilyName);
666       }
667       //Import mesh on entity
668       aResult += GetMeshOnEntitySize(aMeshName,anEntity);
669     }
670   }
671   if(MYDEBUG)
672     MESSAGE("GetSize - aResult = "<<float(aResult));
673   return aResult;
674 }
675
676
677 float VISU_Convertor_impl::GetMeshOnEntitySize(const std::string& theMeshName, 
678                                                const VISU::TEntity& theEntity,
679                                                const std::string& theFamilyName)
680 {
681   VISU::TMesh* pMesh = NULL;
682   VISU::TFamily* pFamily = NULL;
683   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
684   FindMeshOnEntity(theMeshName,pMesh,theEntity,pMeshOnEntity,theFamilyName,pFamily);
685   vtkIdType aPointsSize = 3*pMesh->myNbPoints*sizeof(VISU::TMesh::TCoord);
686   vtkIdType aNbCells, aCellsSize;
687   if(!pFamily){
688     aNbCells = pMeshOnEntity->myNbCells;
689     aCellsSize = pMeshOnEntity->myCellsSize;
690   }else{
691     aNbCells = pFamily->myNbCells;
692     aCellsSize = pFamily->myCellsSize;
693   }
694   vtkIdType aConnectivitySize = aCellsSize*sizeof(vtkIdType);
695   vtkIdType aTypesSize = aNbCells*sizeof(char);
696   vtkIdType aLocationsSize = aNbCells*sizeof(int);
697   float aNbCellsPerPoint = aCellsSize / aNbCells - 1;
698   vtkIdType aLinksSize = pMesh->myNbPoints * 
699     (vtkIdType(sizeof(vtkIdType)*aNbCellsPerPoint) + sizeof(vtkCellLinks::Link));
700   aLinksSize = 0;
701   vtkIdType aResult = aPointsSize + aConnectivitySize + aTypesSize + aLocationsSize + aLinksSize;
702   if(MYDEBUG){
703     if(MYVTKDEBUG){
704       MESSAGE("GetMeshOnEntitySize - aPointsSize = "<<float(aPointsSize));
705       MESSAGE("GetMeshOnEntitySize - aConnectivitySize = "<<float(aConnectivitySize));
706       MESSAGE("GetMeshOnEntitySize - aTypesSize = "<<float(aTypesSize));
707       MESSAGE("GetMeshOnEntitySize - aLocationsSize = "<<float(aLocationsSize));
708       MESSAGE("GetMeshOnEntitySize - aLinksSize = "<<float(aLinksSize));
709     }
710     MESSAGE("GetMeshOnEntitySize - aResult = "<<float(aResult)<<"; theMeshName = '"<<theMeshName<<
711             "'; theEntity = "<<theEntity<<"; theFamilyName = '"<<theFamilyName<<"'");
712   }
713   aResult = vtkIdType(aResult*ERR_SIZE_CALC);
714   return aResult;
715 }
716
717
718 void VISU_Convertor_impl::FindMeshOnGroup(const std::string& theMeshName, VISU::TMesh*& theMesh,
719                                           const std::string& theGroupName, VISU::TGroup*& theGroup)
720 {
721   FindMesh(theMeshName,theMesh);
722   VISU::TGroupMap& aGroupMap = theMesh->myGroupMap;
723   VISU::TGroupMap::iterator aGroupMapIter = aGroupMap.find(theGroupName);
724   if(aGroupMapIter == aGroupMap.end())
725     throw std::runtime_error(EXCEPT("FindMesh >> There is no the group in the mesh!!! - '%1'").arg(theGroupName.c_str()).latin1());
726   theGroup = &aGroupMapIter->second;
727 }
728
729
730 float VISU_Convertor_impl::GetMeshOnGroupSize(const std::string& theMeshName, 
731                                               const std::string& theGroupName)
732 {
733   VISU::TMesh* pMesh = NULL;
734   VISU::TGroup* pGroup = NULL;
735   FindMeshOnGroup(theMeshName,pMesh,theGroupName,pGroup);
736   vtkIdType aPointsSize = 3*pMesh->myNbPoints*sizeof(VISU::TMesh::TCoord);
737   vtkIdType aNbCells = pGroup->myNbCells, aCellsSize = pGroup->myCellsSize;
738   vtkIdType aConnectivityAndTypesSize = aCellsSize*sizeof(vtkIdType);
739   vtkIdType aLocationsSize = aNbCells*sizeof(int);
740   float aNbCellsPerPoint = aCellsSize / aNbCells - 1;
741   vtkIdType aLinksSize = pMesh->myNbPoints * 
742     (vtkIdType(sizeof(vtkIdType)*aNbCellsPerPoint) + sizeof(short));
743   aLinksSize = 0;
744   vtkIdType aResult = aPointsSize + aConnectivityAndTypesSize + aLocationsSize + aLinksSize;
745   if(MYDEBUG){
746     if(MYVTKDEBUG){
747       MESSAGE("GetMeshOnGroupSize - aPointsSize = "<<float(aPointsSize));
748       MESSAGE("GetMeshOnGroupSize - aConnectivityAndTypesSize = "<<float(aConnectivityAndTypesSize));
749       MESSAGE("GetMeshOnGroupSize - aLocationsSize = "<<float(aLocationsSize));
750       MESSAGE("GetMeshOnGroupSize - aLinksSize = "<<float(aLinksSize));
751     }
752     MESSAGE("GetMeshOnGroupSize - aResult = "<<float(aResult)<<"; theMeshName = '"
753             <<theMeshName<<"'; theGroupName = '"<<theGroupName<<"'");
754   }
755   aResult = vtkIdType(aResult*ERR_SIZE_CALC);
756   return aResult;
757 }
758
759
760 void VISU_Convertor_impl::FindField(const string& theMeshName, VISU::TMesh*& theMesh,
761                                     const VISU::TEntity& theEntity, 
762                                     VISU::TMeshOnEntity*& theMeshOnEntity,
763                                     VISU::TMeshOnEntity*& theVTKMeshOnEntity,
764                                     const string& theFieldName, VISU::TField*& theField)
765 {
766   VISU::TFamily* pFamily = NULL;
767   VISU::TMeshOnEntity* pMeshOnEntity = NULL;
768   FindMeshOnEntity(theMeshName,theMesh,theEntity,pMeshOnEntity,"",pFamily);
769   theMeshOnEntity = pMeshOnEntity;
770   VISU::TMeshOnEntityMap& aMeshOnEntityMap = theMesh->myMeshOnEntityMap;
771   if(theEntity == VISU::NODE_ENTITY){
772     if(aMeshOnEntityMap.find(VISU::CELL_ENTITY) != aMeshOnEntityMap.end())
773       pMeshOnEntity = &aMeshOnEntityMap[VISU::CELL_ENTITY];
774     else if(aMeshOnEntityMap.find(VISU::FACE_ENTITY) != aMeshOnEntityMap.end())
775       pMeshOnEntity = &aMeshOnEntityMap[VISU::FACE_ENTITY];
776     else if(aMeshOnEntityMap.find(VISU::NODE_ENTITY) != aMeshOnEntityMap.end())
777       pMeshOnEntity = &aMeshOnEntityMap[VISU::EDGE_ENTITY];
778   }
779   theVTKMeshOnEntity = pMeshOnEntity;
780   VISU::TFieldMap& aFieldMap = theMeshOnEntity->myFieldMap;
781   if(aFieldMap.find(theFieldName) == aFieldMap.end())
782     throw std::runtime_error(EXCEPTION("FindField >> There is no field on the mesh!!!"));
783   theField = &aFieldMap[theFieldName];
784 }
785
786
787 float VISU_Convertor_impl::GetFieldOnMeshSize(const std::string& theMeshName, 
788                                               const VISU::TEntity& theEntity,
789                                               const std::string& theFieldName)
790 {
791   VISU::TMesh* pMesh = NULL;
792   VISU::TField* pField = NULL;
793   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
794   FindField(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
795             theFieldName,pField);
796   float aMeshSize = GetMeshOnEntitySize(theMeshName,pVTKMeshOnEntity->myEntity);
797   float aFieldOnMeshSize = float(pField->myDataSize*sizeof(float)*pField->myNbValField * ERR_SIZE_CALC);
798   float aResult = aMeshSize + aFieldOnMeshSize;
799   if(MYDEBUG){
800     if(MYVTKDEBUG)
801       MESSAGE("GetFieldOnMeshSize - aFieldOnMeshSize = "<<float(aFieldOnMeshSize));
802     MESSAGE("GetFieldOnMeshSize - aResult = "<<float(aResult)<<"; theMeshName = '"<<theMeshName<<
803             "'; theEntity = "<<theEntity<<"; theFieldName = '"<<theFieldName<<"'");
804   }
805   return aResult;
806 }
807
808
809 void VISU_Convertor_impl::FindTimeStamp(const std::string& theMeshName, VISU::TMesh*& theMesh,
810                                         const VISU::TEntity& theEntity, 
811                                         VISU::TMeshOnEntity*& theMeshOnEntity,
812                                         VISU::TMeshOnEntity*& theVTKMeshOnEntity,
813                                         const std::string& theFieldName, VISU::TField*& theField,
814                                         int theStampsNum, VISU::TField::TValForTime*& theValForTime)
815 {
816   FindField(theMeshName,theMesh,theEntity,theMeshOnEntity,theVTKMeshOnEntity,theFieldName,theField);
817   VISU::TField::TValField& aValField = theField->myValField;
818   if(aValField.find(theStampsNum) == aValField.end())
819     throw std::runtime_error(EXCEPTION("FindTimeStamp >> There is no field with the timestamp!!!"));
820   theValForTime = &aValField[theStampsNum];
821 }
822
823
824 float VISU_Convertor_impl::GetTimeStampSize(const std::string& theMeshName, 
825                                             const VISU::TEntity& theEntity,
826                                             const std::string& theFieldName,
827                                             int theStampsNum)
828 {
829   VISU::TMesh* pMesh = NULL;
830   VISU::TField* pField = NULL;
831   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
832   VISU::TField::TValForTime* pValForTime = NULL;
833   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
834                 theFieldName,pField,theStampsNum,pValForTime);
835   float aMeshSize = GetMeshOnEntitySize(theMeshName,pVTKMeshOnEntity->myEntity);
836   float aTimeStampSize = float(pField->myDataSize*sizeof(float) * ERR_SIZE_CALC);
837   float aResult = aMeshSize + aTimeStampSize;
838   if(MYDEBUG){
839     if(MYVTKDEBUG)
840       MESSAGE("GetTimeStampSize - aTimeStampSize = "<<float(aTimeStampSize));
841     MESSAGE("GetTimeStampSize - aResult = "<<float(aResult)<<"; theMeshName = '"<<theMeshName<<"'; theEntity = "<<theEntity<<
842             "; theFieldName = '"<<theFieldName<<"'; theStampsNum = "<<theStampsNum);
843   }
844   return aResult;
845 }
846
847
848 const VISU::TField& 
849 VISU_Convertor_impl::GetField(const string& theMeshName, 
850                               VISU::TEntity theEntity, 
851                               const string& theFieldName) 
852 {
853   VISU::TMesh* pMesh = NULL;
854   VISU::TField* pField = NULL;
855   VISU::TFamily* pFamily = NULL;
856   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
857   FindField(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
858             theFieldName,pField);
859   return *pField;
860 }
861
862
863 const VISU::TField::TValForTime& 
864 VISU_Convertor_impl::GetTimeStamp(const std::string& theMeshName, 
865                                   const VISU::TEntity& theEntity,
866                                   const std::string& theFieldName,
867                                   int theStampsNum)
868 {
869   VISU::TMesh* pMesh = NULL;
870   VISU::TField* pField = NULL;
871   VISU::TMeshOnEntity *pMeshOnEntity = NULL, *pVTKMeshOnEntity = NULL;
872   VISU::TField::TValForTime* pValForTime = NULL;
873   FindTimeStamp(theMeshName,pMesh,theEntity,pMeshOnEntity,pVTKMeshOnEntity,
874                 theFieldName,pField,theStampsNum,pValForTime);
875   return *pValForTime;
876 }
877
878