]> SALOME platform Git repositories - modules/visu.git/blob - src/CONVERTOR/VISU_MeshValue.cxx
Salome HOME
Merge from BR_V5_DEV 16Feb09
[modules/visu.git] / src / CONVERTOR / VISU_MeshValue.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:
24 //  Author:
25 //  Module : VISU
26 //
27 #include "VISU_MeshValue.hxx"
28 #include "VISU_ElnoMeshValue.hxx"
29 #include "VISU_Structures_impl.hxx"
30 #include "VISU_ConvertorUtils.hxx"
31
32 #include "VISU_PointCoords.hxx"
33 #include "VISU_VTKTypeList.hxx"
34
35 #include <vtkUnstructuredGrid.h>
36 #include <vtkPolyData.h>
37
38 #include <vtkPointData.h>
39 #include <vtkCellData.h>
40
41 #include <vtkCharArray.h>
42 #include <vtkUnsignedCharArray.h>
43 #include <vtkShortArray.h>
44 #include <vtkUnsignedShortArray.h>
45 #include <vtkIntArray.h>
46 #include <vtkUnsignedIntArray.h>
47 #include <vtkLongArray.h>
48 #include <vtkUnsignedLongArray.h> 
49 #include <vtkFloatArray.h> 
50 #include <vtkDoubleArray.h> 
51
52 #include <string>
53 #include <algorithm>
54
55 #ifdef _DEBUG_
56 static int MYDEBUG = 0;
57 #else
58 static int MYDEBUG = 0;
59 #endif
60
61 namespace VISU
62 {
63   //---------------------------------------------------------------
64   std::string
65   GenerateFieldName(const PFieldImpl& theField,
66                     const PValForTimeImpl& theValForTime)
67   {
68     const VISU::TTime& aTime = theValForTime->myTime;
69     std::string aFieldName = theField->myMeshName + ", " + theField->myName + ": " + 
70       VISU_Convertor::GenerateName(aTime);
71     return aFieldName;
72   }
73
74
75   //---------------------------------------------------------------
76   void
77   TMeshValueBase
78   ::Init(vtkIdType theNbElem,
79          vtkIdType theNbGauss,
80          vtkIdType theNbComp)
81   {
82     myNbElem = theNbElem;
83     myNbGauss = theNbGauss;
84     myNbComp = theNbComp;
85     myStep = theNbComp*theNbGauss;
86   }
87
88   vtkIdType 
89   TMeshValueBase
90   ::GetNbElem() const
91   {
92     return myNbElem;
93   }
94
95   vtkIdType 
96   TMeshValueBase
97   ::GetNbComp() const
98   {
99     return myNbComp;
100   }
101
102   vtkIdType 
103   TMeshValueBase
104   ::GetNbGauss() const
105   {
106     return myNbGauss;
107   }
108
109   size_t 
110   TMeshValueBase
111   ::size() const
112   {
113     return myNbElem * myStep;
114   }
115
116
117   //----------------------------------------------------------------------------
118   template<int EDataType>
119   void 
120   InitTimeStampOnProfile(const PUnstructuredGrid& theSource,
121                          const PFieldImpl& theField, 
122                          const PValForTimeImpl& theValForTime,
123                          const VISU::TEntity& theEntity);
124
125
126   //----------------------------------------------------------------------------
127   void 
128   GetTimeStampOnProfile(const PUnstructuredGrid& theSource,
129                         const PFieldImpl& theField, 
130                         const PValForTimeImpl& theValForTime,
131                         const VISU::TEntity& theEntity)
132   {
133     vtkIdType aDataType = theField->GetDataType();
134     switch(aDataType){
135     case VTK_DOUBLE:
136       InitTimeStampOnProfile<VTK_DOUBLE>(theSource, theField, theValForTime, theEntity);
137       break;
138     case VTK_FLOAT:
139       InitTimeStampOnProfile<VTK_FLOAT>(theSource, theField, theValForTime, theEntity);
140       break;
141     case VTK_INT:
142       InitTimeStampOnProfile<VTK_INT>(theSource, theField, theValForTime, theEntity);
143       break;
144     case VTK_LONG:
145       InitTimeStampOnProfile<VTK_LONG>(theSource, theField, theValForTime, theEntity);
146       break;
147     default:
148       EXCEPTION(std::runtime_error,
149                 "GetTimeStampOnProfile - handling unsupported data type - "<<aDataType);
150     }
151   }
152
153
154   //----------------------------------------------------------------------------
155   template<int EDataType>
156   struct TDataArrayHolder
157   {
158     typedef typename TL::TEnum2VTKArrayType<EDataType>::TResult TVTKDataArray;
159     typedef typename TL::TEnum2VTKBasicType<EDataType>::TResult TVTKBasicType;
160     TVTKDataArray* myDataArray;
161
162     TDataArrayHolder(TVTKDataArray* theDataArray):
163       myDataArray(theDataArray)
164     {}
165     
166     void
167     WritePointer(TVTKDataArray* theDataArray,
168                  vtkIdType theTupleId,
169                  TVTKBasicType* thePointer)
170     {
171       vtkIdType aNumberOfComponents = theDataArray->GetNumberOfComponents();
172       vtkIdType aPosition = theTupleId * aNumberOfComponents;
173       TVTKBasicType *aPtr = theDataArray->WritePointer(aPosition, aNumberOfComponents);
174       for(vtkIdType anId = 0; anId < aNumberOfComponents; anId++)
175         *aPtr++ = *thePointer++;
176     }
177
178     virtual
179     void
180     SetTuple(vtkIdType theTupleId, 
181              TVTKBasicType* thePointer)
182     {
183       this->WritePointer(myDataArray, theTupleId, thePointer);
184     }
185   };
186
187
188   //----------------------------------------------------------------------------
189   template<int EDataType>
190   struct TDataArrayHolder2: TDataArrayHolder<EDataType>
191   {
192     typedef TDataArrayHolder<EDataType> TSuperClass;
193     typedef typename TSuperClass::TVTKDataArray TVTKDataArray;
194     typedef typename TSuperClass::TVTKBasicType TVTKBasicType;
195     TVTKDataArray* myDataArray2;
196
197     TDataArrayHolder2(TVTKDataArray* theDataArray,
198                       TVTKDataArray* theDataArray2):
199       TSuperClass(theDataArray),
200       myDataArray2(theDataArray2)
201     {}
202     
203     virtual
204     void
205     SetTuple(vtkIdType theTupleId, 
206              TVTKBasicType* thePointer)
207     {
208       this->WritePointer(this->myDataArray, theTupleId, thePointer);
209       this->WritePointer(this->myDataArray2, theTupleId, thePointer);
210     }
211   };
212
213
214   //----------------------------------------------------------------------------
215   template<int EDataType>
216   struct TTimeStampOnProfileInitArray
217   {
218     typedef typename TL::TEnum2VTKArrayType<EDataType>::TResult TVTKDataArray;
219     typedef typename TL::TEnum2VTKBasicType<EDataType>::TResult TVTKBasicType;
220     typedef TTMeshValue<TVTKBasicType> TMeshValue;
221     typedef MED::SharedPtr<TMeshValue> TMeshValuePtr;
222
223     typedef TDataArrayHolder<EDataType> TTDataArrayHolder;
224     typedef MED::SharedPtr<TTDataArrayHolder> PDataArrayHolder;
225     PDataArrayHolder myDataArrayHolder;
226
227     TTimeStampOnProfileInitArray(const PDataArrayHolder& theDataArrayHolder):
228       myDataArrayHolder(theDataArrayHolder)
229     {}
230     
231     void
232     Execute(const PFieldImpl& theField,
233             const PValForTimeImpl& theValForTime)
234     {
235       vtkIdType aNbComp = theField->myNbComp;
236       vtkIdType aSize = std::max(vtkIdType(3), aNbComp);
237       TVector<TVTKBasicType> aDataValues(aSize);
238       
239       const TGeom2MeshValue& aGeom2MeshValue = theValForTime->GetGeom2MeshValue();
240       TGeom2MeshValue::const_iterator anIter = aGeom2MeshValue.begin();
241       for(int aTupleId = 0; anIter != aGeom2MeshValue.end(); anIter++){
242         EGeometry aEGeom = anIter->first;
243         const TMeshValuePtr aMeshValue = anIter->second;
244         
245         vtkIdType aNbElem = aMeshValue->GetNbElem();
246         vtkIdType aNbGauss = aMeshValue->GetNbGauss();
247         
248         INITMSG(MYDEBUG,
249                 "- aEGeom = "<<aEGeom<<
250                 "; aNbElem = "<<aNbElem<<
251                 "; aNbGauss = "<<aNbGauss<<
252                 std::endl);
253         
254         for(vtkIdType iElem = 0; iElem < aNbElem; iElem++, aTupleId++){
255           typename TMeshValue::TCValueSliceArr aValueSliceArr = aMeshValue->GetCompValueSliceArr(iElem);
256           for(vtkIdType iComp = 0; iComp < aNbComp; iComp++){
257             const typename TMeshValue::TCValueSlice& aValueSlice = aValueSliceArr[iComp];
258             aDataValues[iComp] = TVTKBasicType();
259             for(vtkIdType iGauss = 0; iGauss < aNbGauss; iGauss++){
260               aDataValues[iComp] += aValueSlice[iGauss];
261             }
262             aDataValues[iComp] /= aNbGauss;
263           }
264           this->myDataArrayHolder->SetTuple(aTupleId, &aDataValues[0]);
265         }
266       }
267     }
268   };
269
270
271   //----------------------------------------------------------------------------
272   template<int EDataType>
273   void 
274   InitTimeStampOnProfile(const PUnstructuredGrid& theSource,
275                          const PFieldImpl& theField, 
276                          const PValForTimeImpl& theValForTime,
277                          const VISU::TEntity& theEntity)
278   {
279     vtkIdType aNbTuples = theField->myDataSize / theField->myNbComp;
280     std::string aFieldName = VISU::GenerateFieldName(theField, theValForTime);
281     
282     vtkDataSetAttributes* aDataSetAttributes;
283     switch ( theEntity ) {
284     case VISU::NODE_ENTITY : 
285       aDataSetAttributes = theSource->GetPointData();
286       break;
287     default: 
288       aDataSetAttributes = theSource->GetCellData();
289     }
290
291     typedef typename TL::TEnum2VTKArrayType<EDataType>::TResult TVTKDataArray;
292     TVTKDataArray *aSelectedDataArray = TVTKDataArray::New();
293     vtkIdType aNbComp = theField->myNbComp;
294
295     switch ( aNbComp ) {
296     case 1:
297       aSelectedDataArray->SetNumberOfComponents( 1 );
298       aDataSetAttributes->SetScalars( aSelectedDataArray );
299       break;
300     default:
301       aSelectedDataArray->SetNumberOfComponents( 3 );
302       aDataSetAttributes->SetVectors( aSelectedDataArray );
303     }
304     aSelectedDataArray->SetNumberOfTuples( aNbTuples );
305     aSelectedDataArray->SetName( aFieldName.c_str() );
306
307     TVTKDataArray *aFullDataArray = TVTKDataArray::New();
308     aFullDataArray->SetNumberOfComponents( aNbComp );
309     aFullDataArray->SetNumberOfTuples( aNbTuples );
310     aFullDataArray->SetName( "VISU_FIELD" );
311     aDataSetAttributes->AddArray( aFullDataArray );
312
313     INITMSG(MYDEBUG,"InitTimeStampOnProfile "<<
314             "- theEntity = "<<theEntity<<
315             "; aNbTuples = "<<aNbTuples<<
316             "; aNbComp = "<<aNbComp<<
317             std::endl);
318
319     TTimerLog aTimerLog(MYDEBUG,"InitTimeStampOnProfile");
320     
321     const TGeom2MeshValue& aGeom2MeshValue = theValForTime->GetGeom2MeshValue();
322     typedef typename TL::TEnum2VTKBasicType< EDataType >::TResult TVTKBasicType;
323     typedef TTMeshValue< TVTKBasicType > TMeshValue;
324     typedef MED::SharedPtr< TMeshValue > TMeshValuePtr;
325
326     typedef TDataArrayHolder< EDataType > TTDataArrayHolder;
327     typedef MED::SharedPtr< TTDataArrayHolder > PDataArrayHolder;
328
329     TMeshValuePtr aMeshValue = theValForTime->GetFirstMeshValue();
330     if ( aGeom2MeshValue.size() == 1 && aMeshValue->GetNbGauss() == 1 ) {
331       aFullDataArray->SetVoidArray(aMeshValue->GetPointer(),
332                                    aMeshValue->size(),
333                                    true);
334       INITMSG(MYDEBUG,"InitTimeStampOnProfile - aFullDataArray->SetVoidArray()"<<std::endl);
335       if ( aNbComp == 1 ) {
336         aSelectedDataArray->SetVoidArray( aMeshValue->GetPointer(),
337                                           aMeshValue->size(),
338                                           true );
339         INITMSG(MYDEBUG,"InitTimeStampOnProfile - aSelectedDataArray->SetVoidArray()"<<std::endl);
340       }else{
341         PDataArrayHolder aDataArrayHolder(new TTDataArrayHolder(aSelectedDataArray));
342         TTimeStampOnProfileInitArray<EDataType>(aDataArrayHolder).Execute(theField, theValForTime);
343       }
344     }else{
345       typedef TDataArrayHolder2<EDataType> TTDataArrayHolder2;
346       PDataArrayHolder aDataArrayHolder(new TTDataArrayHolder2(aSelectedDataArray, aFullDataArray));
347       TTimeStampOnProfileInitArray<EDataType>(aDataArrayHolder).Execute(theField, theValForTime);
348     }
349
350     aSelectedDataArray->Delete();
351     aFullDataArray->Delete();
352
353     // Process the case for ELNO data
354     //-------------------------------
355     if ( theField->myIsELNO ) {
356       // To calculate effective number of components for the VTK compatibel ELNO data representation
357       vtkIdType aEffectNbTuples = 0;
358       TGeom2MeshValue::const_iterator anIter = aGeom2MeshValue.begin();
359       for ( ; anIter != aGeom2MeshValue.end(); anIter++ ) {
360         const PMeshValue& aMeshValue = anIter->second;
361         aEffectNbTuples += aMeshValue->GetNbElem() * aMeshValue->GetNbGauss();
362       }
363
364       vtkIdType anEffectNbComp = ( aEffectNbTuples * aNbComp ) / aNbTuples + 1;
365     
366       // To create corresponding VTK representation for the ELNO data
367       TSetElnoNodeData< EDataType > aSetElnoNodeData( anEffectNbComp,
368                                                       aNbComp,
369                                                       aNbTuples,
370                                                       "ELNO_FIELD",
371                                                       "ELNO_COMPONENT_MAPPER" );
372
373       std::vector< TVTKBasicType > aDataValues( aNbComp ); // To reserve a temproary value holder
374
375       // To initilize these VTK representation for the ELNO data from the MED
376       anIter = aGeom2MeshValue.begin();
377       for ( ; anIter != aGeom2MeshValue.end(); anIter++ ) {
378         EGeometry aEGeom = anIter->first;
379         const TMeshValuePtr aMeshValue = anIter->second;
380         
381         vtkIdType aNbElem = aMeshValue->GetNbElem();
382         vtkIdType aNbGauss = aMeshValue->GetNbGauss();
383         
384         INITMSG(MYDEBUG,
385                 "- aEGeom = "<<aEGeom<<
386                 "; aNbElem = "<<aNbElem<<
387                 "; aNbGauss = "<<aNbGauss<<
388                 std::endl);
389         std::vector<int> med2visu(aNbGauss);
390         InitMed2VisuArray(med2visu,aEGeom);
391         for ( vtkIdType iElem = 0; iElem < aNbElem; iElem++ ) {
392           const typename TMeshValue::TValueSliceArr& aValueSliceArr = aMeshValue->GetGaussValueSliceArr( iElem );          
393
394           for( vtkIdType iGauss = 0; iGauss < aNbGauss; iGauss++ ) {
395             const typename TMeshValue::TCValueSlice& aValueSlice = aValueSliceArr[ med2visu[iGauss] ];
396
397             for( vtkIdType iComp = 0; iComp < aNbComp; iComp++ ) {
398               aDataValues[ iComp ] = aValueSlice[ iComp ];
399             }
400
401             aSetElnoNodeData.AddNextPointData( &aDataValues[ 0 ] );
402           }
403
404           aSetElnoNodeData.InsertNextCellData();
405         }
406       }
407       
408       // Assign the ELNO data on the corresponding VTK data set attribute 
409       aSetElnoNodeData.AddData( aDataSetAttributes );
410     }
411     //-------------------------------
412   }
413
414
415   //----------------------------------------------------------------------------
416   template<int EDataType>
417   void 
418   InitTimeStampOnGaussMesh(const PPolyData& theSource,
419                            const PFieldImpl& theField, 
420                            const PValForTimeImpl& theValForTime);
421
422   void 
423   GetTimeStampOnGaussMesh(const PPolyData& theSource,
424                           const PFieldImpl& theField, 
425                           const PValForTimeImpl& theValForTime)
426   {
427     vtkIdType aDataType = theField->GetDataType();
428     switch(aDataType){
429     case VTK_DOUBLE:
430       InitTimeStampOnGaussMesh<VTK_DOUBLE>(theSource, theField, theValForTime);
431       break;
432     case VTK_FLOAT:
433       InitTimeStampOnGaussMesh<VTK_FLOAT>(theSource, theField, theValForTime);
434       break;
435     case VTK_INT:
436       InitTimeStampOnGaussMesh<VTK_INT>(theSource, theField, theValForTime);
437       break;
438     case VTK_LONG:
439       InitTimeStampOnGaussMesh<VTK_LONG>(theSource, theField, theValForTime);
440       break;
441     default:
442       EXCEPTION(std::runtime_error,
443                 "GetTimeStampOnGaussMesh - handling unsupported data type - "<<aDataType);
444     }
445   }
446
447   //----------------------------------------------------------------------------
448   template<int EDataType>
449   struct TTimeStampOnGaussMeshInitArray
450   {
451     typedef typename TL::TEnum2VTKArrayType<EDataType>::TResult TVTKDataArray;
452     typedef typename TL::TEnum2VTKBasicType<EDataType>::TResult TVTKBasicType;
453     typedef TTMeshValue<TVTKBasicType> TMeshValue;
454     typedef MED::SharedPtr<TMeshValue> TMeshValuePtr;
455
456     typedef TDataArrayHolder<EDataType> TTDataArrayHolder;
457     typedef MED::SharedPtr<TTDataArrayHolder> PDataArrayHolder;
458     PDataArrayHolder myDataArrayHolder;
459
460     TTimeStampOnGaussMeshInitArray(const PDataArrayHolder& theDataArrayHolder):
461       myDataArrayHolder(theDataArrayHolder)
462     {}
463     
464     void
465     Execute(const PFieldImpl& theField,
466             const PValForTimeImpl& theValForTime)
467     {
468       vtkIdType aNbComp = theField->myNbComp;
469       vtkIdType aSize = std::max(vtkIdType(3), aNbComp);
470       TVector<TVTKBasicType> aDataValues(aSize);
471
472       const TGeom2MeshValue& aGeom2MeshValue = theValForTime->GetGeom2MeshValue();
473
474       PGaussMeshImpl aGaussMesh = theValForTime->myGaussMesh;
475       const TGeom2GaussSubMesh& aGeom2GaussSubMesh = aGaussMesh->myGeom2GaussSubMesh;
476       TGeom2GaussSubMesh::const_iterator anIter = aGeom2GaussSubMesh.begin();
477       for(int aTupleId = 0; anIter != aGeom2GaussSubMesh.end(); anIter++){
478         EGeometry aEGeom = anIter->first;
479         
480         PGaussSubMeshImpl aGaussSubMesh = anIter->second;
481         if(!aGaussSubMesh->myIsDone)
482           continue;
483         
484         TGeom2MeshValue::const_iterator anIter2 = aGeom2MeshValue.find(aEGeom);
485         if(anIter2 == aGeom2MeshValue.end()){
486           EXCEPTION(std::runtime_error,
487                     "TTimeStampOnGaussMeshInitArray >> Can't find values for corresponding Gauss Points SubMesh");
488         }
489         TMeshValuePtr aMeshValue = anIter2->second;
490         vtkIdType aNbGauss = aMeshValue->GetNbGauss();
491         vtkIdType aNbElem = aMeshValue->GetNbElem();
492         
493         if(aNbGauss < 1)
494           continue;
495         
496         const TPointCoords& aCoords = aGaussSubMesh->myPointCoords;
497         
498         INITMSG(MYDEBUG,
499                 "- aEGeom = "<<aEGeom<<
500                 "; aNbElem = "<<aNbElem<<
501                 "; aNbGauss = "<<aNbGauss<<
502                 "; aCoords.GetNbPoints() = "<<aCoords.GetNbPoints()<<
503                 std::endl);
504
505         if(aCoords.GetNbPoints() == aNbElem*aNbGauss){
506           for(int iElem = 0; iElem < aNbElem; iElem++){
507             typename TMeshValue::TCValueSliceArr aValueSliceArr = aMeshValue->GetGaussValueSliceArr(iElem);
508             for(int iGauss = 0; iGauss < aNbGauss; iGauss++, aTupleId++){
509               const typename TMeshValue::TCValueSlice& aValueSlice = aValueSliceArr[iGauss];
510               for(int iComp = 0; iComp < aNbComp; iComp++){
511                 aDataValues[iComp] = aValueSlice[iComp];
512               }
513               this->myDataArrayHolder->SetTuple(aTupleId, &aDataValues[0]);
514             }
515           }
516         }else{
517           for(int iElem = 0; iElem < aNbElem; iElem++, aTupleId++){
518             typename TMeshValue::TCValueSliceArr aValueSliceArr = aMeshValue->GetCompValueSliceArr(iElem);
519             for(int iComp = 0; iComp < aNbComp; iComp++){
520               const typename TMeshValue::TCValueSlice& aValueSlice = aValueSliceArr[iComp];
521               aDataValues[iComp] = TVTKBasicType();
522               for(int iGauss = 0; iGauss < aNbGauss; iGauss++){
523                 aDataValues[iComp] += aValueSlice[iGauss];
524               }
525               aDataValues[iComp] /= aNbGauss;
526             }
527             this->myDataArrayHolder->SetTuple(aTupleId, &aDataValues[0]);
528           }
529         }
530       }
531     }
532   };
533
534
535   template<int EDataType>
536   void 
537   InitTimeStampOnGaussMesh(const PPolyData& theSource,
538                            const PFieldImpl& theField, 
539                            const PValForTimeImpl& theValForTime)
540   {
541     vtkIdType aNbTuples = theSource->GetNumberOfPoints();
542     std::string aFieldName = VISU::GenerateFieldName(theField, theValForTime);
543
544     vtkDataSetAttributes* aDataSetAttributes = theSource->GetPointData();
545
546     typedef typename TL::TEnum2VTKArrayType<EDataType>::TResult TVTKDataArray;
547     TVTKDataArray *aSelectedDataArray = TVTKDataArray::New();
548     vtkIdType aNbComp = theField->myNbComp;
549     switch(aNbComp){
550     case 1:
551       aSelectedDataArray->SetNumberOfComponents(1);
552       aDataSetAttributes->SetScalars(aSelectedDataArray);
553       break;
554     default:
555       aSelectedDataArray->SetNumberOfComponents(3);
556       aDataSetAttributes->SetVectors(aSelectedDataArray);
557     }
558     aSelectedDataArray->SetNumberOfTuples(aNbTuples);
559     aSelectedDataArray->SetName(aFieldName.c_str());
560
561     TVTKDataArray *aFullDataArray = TVTKDataArray::New();
562     aFullDataArray->SetNumberOfComponents(aNbComp);
563     aFullDataArray->SetNumberOfTuples(aNbTuples);
564     aFullDataArray->SetName("VISU_FIELD");
565     aDataSetAttributes->AddArray(aFullDataArray);
566
567     INITMSG(MYDEBUG,"InitTimeStampOnGaussMesh "<<
568             "- aNbTuples = "<<aNbTuples<<
569             "; aNbComp = "<<aNbComp<<
570             std::endl);
571     TTimerLog aTimerLog(MYDEBUG,"InitTimeStampOnGaussMesh");
572     
573     const TGeom2MeshValue& aGeom2MeshValue = theValForTime->GetGeom2MeshValue();
574     typedef typename TL::TEnum2VTKBasicType<EDataType>::TResult TVTKBasicType;
575     typedef TTMeshValue<TVTKBasicType> TMeshValue;
576     typedef MED::SharedPtr<TMeshValue> TMeshValuePtr;
577
578     typedef TDataArrayHolder<EDataType> TTDataArrayHolder;
579     typedef MED::SharedPtr<TTDataArrayHolder> PDataArrayHolder;
580
581     TMeshValuePtr aMeshValue = theValForTime->GetFirstMeshValue();
582     if(aGeom2MeshValue.size() == 1){
583       aFullDataArray->SetVoidArray(aMeshValue->GetPointer(),
584                                    aMeshValue->size(),
585                                    true);
586       INITMSG(MYDEBUG,"InitTimeStampOnGaussMesh - aFullDataArray->SetVoidArray()"<<std::endl);
587       if(aNbComp == 1 || aNbComp == 3){
588         aSelectedDataArray->SetVoidArray(aMeshValue->GetPointer(),
589                                          aMeshValue->size(),
590                                          true);
591         INITMSG(MYDEBUG,"InitTimeStampOnGaussMesh - aSelectedDataArray->SetVoidArray()"<<std::endl);
592       }else{
593         PDataArrayHolder aDataArrayHolder(new TTDataArrayHolder(aSelectedDataArray));
594         TTimeStampOnGaussMeshInitArray<EDataType>(aDataArrayHolder).Execute(theField, theValForTime);
595       }
596     }else{
597       typedef TDataArrayHolder2<EDataType> TTDataArrayHolder2;
598       PDataArrayHolder aDataArrayHolder(new TTDataArrayHolder2(aSelectedDataArray, aFullDataArray));
599       TTimeStampOnGaussMeshInitArray<EDataType>(aDataArrayHolder).Execute(theField, theValForTime);
600     }
601
602     aSelectedDataArray->Delete();
603     aFullDataArray->Delete();
604   }
605
606   
607   //----------------------------------------------------------------------------
608
609   void InitMed2VisuArray(std::vector<int>& anArray, EGeometry aEGeom){
610     switch(aEGeom){
611 #if !(defined(VTK_QUADRATIC_EDGE) && defined(VISU_USE_VTK_QUADRATIC)) && defined(VISU_ENABLE_QUADRATIC)
612     case eSEG3:
613       anArray[0] = 0;
614       anArray[2] = 1;  
615       anArray[1] = 2;
616       break;
617 #endif
618
619 #if !(defined(VTK_QUADRATIC_TRIANGLE) && defined(VISU_USE_VTK_QUADRATIC)) && defined(VISU_ENABLE_QUADRATIC)
620     case eTRIA6:
621       anArray[0] = 0;
622       anArray[2] = 1;  
623       anArray[4] = 2;  
624       
625       anArray[1] = 3;
626       anArray[3] = 4;  
627       anArray[5] = 5;
628       break;
629 #endif
630
631 #if !(defined(VTK_QUADRATIC_QUAD) && defined(VISU_USE_VTK_QUADRATIC)) && defined(VISU_ENABLE_QUADRATIC)
632     case eQUAD8:
633       anArray[0] = 0;
634       anArray[2] = 1;  
635       anArray[4] = 2;  
636       anArray[6] = 3;  
637       
638       anArray[1] = 4;
639       anArray[3] = 5;  
640       anArray[5] = 6;  
641       anArray[7] = 7;
642       break;
643 #endif
644     case eTETRA4:
645       anArray[0] = 0;
646       anArray[1] = 2;
647       anArray[2] = 1;  
648       anArray[3] = 3;
649       break;
650     case ePYRA5:
651       anArray[0] = 0;
652       anArray[1] = 3;  
653       anArray[2] = 2;
654       anArray[3] = 1;  
655       anArray[4] = 4;
656       break;
657 #if (defined(VTK_QUADRATIC_TETRA) && defined(VISU_USE_VTK_QUADRATIC)) && defined(VISU_ENABLE_QUADRATIC)
658     case eTETRA10:
659       anArray[0] = 0;
660       anArray[1] = 2;
661       anArray[2] = 1;  
662       anArray[3] = 3;  
663       
664       anArray[4] = 6;
665       anArray[5] = 5;
666       anArray[6] = 4;  
667       
668       anArray[7] = 7;  
669       anArray[8] = 9;  
670       anArray[9] = 8;
671       break;
672 #endif
673
674 #if (defined(VTK_QUADRATIC_PYRAMID) && defined(VISU_USE_VTK_QUADRATIC)) && defined(VISU_ENABLE_QUADRATIC)
675     case ePYRA13:
676       anArray[0] = 0;
677       anArray[1] = 3;
678       anArray[2] = 2;  
679       anArray[3] = 1;  
680       anArray[4] = 4;
681       
682       anArray[5] = 8;
683       anArray[6] = 7;  
684       anArray[7] = 6;  
685       anArray[8] = 5;  
686       
687       anArray[9] = 9;  
688       anArray[10] = 12;  
689       anArray[11] = 11;  
690       anArray[12] = 10;  
691       break;
692 #endif
693     default:
694       for(int i=0;i<anArray.size();i++){
695         anArray[i] = i;
696       }
697       break;
698     }
699   }
700 }