]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
0020518: [CEA 357] VISU representation of uniform cell FIELD
authordmv <dmv@opencascade.com>
Wed, 7 Oct 2009 10:44:28 +0000 (10:44 +0000)
committerdmv <dmv@opencascade.com>
Wed, 7 Oct 2009 10:44:28 +0000 (10:44 +0000)
src/PIPELINE/VISU_Extractor.cxx
src/PIPELINE/VISU_Extractor.hxx

index 8dce35e755628f384a0997bd6ac25a87aee9e5ff..4cdfb5834cabe328fa2130f483c94d836286389d 100644 (file)
@@ -108,11 +108,11 @@ vtkFloatingPointType CutValue (vtkFloatingPointType theValue, int theDecimals)
     aDegree = (long long)log10(abs((long long)v)) + 1;
   aDegree = theDecimals - aDegree;
   //printf("$$$ 1 v = %.20g , aDegree = %lld \n", v, (long long)aDegree);
-  if (aDegree > 0) {
-    aDegree = pow(10, aDegree);
-    v = ((vtkFloatingPointType)((long long)(v * aDegree))) / aDegree;
-    //printf("$$$ 2 v = %.20g , aDegree = %lld \n", v, (long long)aDegree);
-  }
+
+  aDegree = pow(10, aDegree);
+  v = ((vtkFloatingPointType)((long long)(v * aDegree))) / aDegree;
+  //printf("$$$ 2 v = %.20g , aDegree = %lld \n", v, (long long)aDegree);
+
   return v;
 }
 
@@ -123,11 +123,6 @@ Module2Scalars(vtkDataArray *theInputDataArray,
               TValueType* theOutputPtr,
               vtkIdType theNbOfTuples)
 {
-  // jfa try begin
-  SUIT_ResourceMgr* aResourceMgr = SUIT_Session::session()->resourceMgr();
-  int aDecimals = aResourceMgr->integerValue("VISU", "floating_point_precision", 6);
-  // jfa try end
-
   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
   std::vector<vtkFloatingPointType> anArray(aNbComp < 3? 3: aNbComp);
   for(vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++){
@@ -136,10 +131,7 @@ Module2Scalars(vtkDataArray *theInputDataArray,
     vtkFloatingPointType aScalar = sqrt(aVector[0]*aVector[0] + 
                                        aVector[1]*aVector[1] + 
                                        aVector[2]*aVector[2]);
-    //*theOutputPtr = TValueType(aScalar);
-    // jfa try begin
-    *theOutputPtr = TValueType(CutValue(aScalar, aDecimals));
-    // jfa try end
+    *theOutputPtr = TValueType(aScalar);
     theOutputPtr++;
   }
 }
@@ -151,11 +143,6 @@ Module2ScalarsMOD(vtkDataArray *theInputDataArray,
                   vtkIdType theNbOfTuples,
                   VISU::TGaussMetric theGaussMetric)
 {
-  // jfa try begin
-  SUIT_ResourceMgr* aResourceMgr = SUIT_Session::session()->resourceMgr();
-  int aDecimals = aResourceMgr->integerValue("VISU", "floating_point_precision", 6);
-  // jfa try end
-
   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
   if (aNbComp != 3) // Min, Max, Avg
     return;
@@ -167,9 +154,6 @@ Module2ScalarsMOD(vtkDataArray *theInputDataArray,
       case VISU::MAXIMUM_METRIC: *theOutputPtr = TValueType(anArray[1]); break;
       case VISU::AVERAGE_METRIC: *theOutputPtr = TValueType(anArray[2]); break;
     }
-    // jfa try begin
-    *theOutputPtr = TValueType(CutValue(*theOutputPtr, aDecimals));
-    // jfa try end
     theOutputPtr++;
   }
 }
@@ -184,22 +168,31 @@ Component2Scalars(vtkDataArray *theInputDataArray,
                   vtkIdType theNbOfTuples,
                   vtkIdType theComponentId)
 {
-  // jfa try begin
-  SUIT_ResourceMgr* aResourceMgr = SUIT_Session::session()->resourceMgr();
-  int aDecimals = aResourceMgr->integerValue("VISU", "floating_point_precision", 6);
-  // jfa try end
-
   vtkIdType aNbComp = theInputDataArray->GetNumberOfComponents();
   for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
     *theOutputPtr = *(theInputPtr + theComponentId);
-    // jfa try begin
-    *theOutputPtr = TValueType(CutValue(*theOutputPtr, aDecimals));
-    // jfa try end
     theInputPtr += aNbComp;
     theOutputPtr++;
   }
 }
 
+//----------------------------------------------------------------------------
+template<typename TValueType>
+void
+CutFieldDataTempl(vtkDataArray *theInputDataArray,
+                  TValueType* theDataPtr,
+                  vtkIdType theNbOfTuples)
+{
+  SUIT_ResourceMgr* aResourceMgr = SUIT_Session::session()->resourceMgr();
+  int aDecimals = aResourceMgr->integerValue("VISU", "floating_point_precision", 6);
+
+  for (vtkIdType aTupleId = 0; aTupleId < theNbOfTuples; aTupleId++) {
+    double aValue;
+    theInputDataArray->GetTuple(aTupleId, &aValue);
+    *theDataPtr = TValueType(CutValue(aValue, aDecimals));
+    theDataPtr++;
+  }
+}
 
 //----------------------------------------------------------------------------
 template<typename TDataSetAttributesType> void
@@ -289,6 +282,7 @@ VISU_Extractor
     int aNbElems = anInput->GetNumberOfPoints();
     if ( anInputPointData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputPointData, anOutputPointData );
+    CutFieldData( aNbElems, anOutputPointData );
   }
 
   vtkCellData *anInputCellData = anInput->GetCellData();
@@ -298,7 +292,41 @@ VISU_Extractor
     int aNbElems = anInput->GetNumberOfCells();
     if ( anInputCellData->GetAttribute( vtkDataSetAttributes::VECTORS ) )
       ExecuteScalars( aNbElems, myScalarMode, myGaussMetric, anInputCellData, anOutputCellData );
+    CutFieldData( aNbElems, anOutputCellData );
   }
 
   return 1;
 }
+
+//---------------------------------------------------------------
+
+void VISU_Extractor::CutFieldData(int theNbElems, vtkDataSetAttributes *theData )
+{
+  if (theNbElems < 1) return;
+
+  vtkDataArray* aFieldArray = NULL;
+  switch (myGaussMetric) {
+  case VISU::AVERAGE_METRIC: aFieldArray = theData->GetArray("VISU_FIELD"); break;
+  case VISU::MINIMUM_METRIC: aFieldArray = theData->GetArray("VISU_FIELD_GAUSS_MIN"); break;
+  case VISU::MAXIMUM_METRIC: aFieldArray = theData->GetArray("VISU_FIELD_GAUSS_MAX"); break;
+  }
+  if( !aFieldArray ) return;
+  
+  vtkIdType anInputDataType = aFieldArray->GetDataType();
+  vtkDataArray *aScalars = vtkDataArray::CreateDataArray(anInputDataType);
+  aScalars->SetNumberOfComponents( 1 );
+  aScalars->SetNumberOfTuples( theNbElems );
+  void *aPtr = aScalars->GetVoidPointer(0);
+
+  switch (anInputDataType) {
+    vtkTemplateMacro3(CutFieldDataTempl,
+                      aFieldArray,
+                      (VTK_TT *)(aPtr),
+                      theNbElems);
+  default:
+    break;
+  }
+
+  theData->SetScalars(aScalars);
+  aScalars->Delete();
+}
index 5d6b8874faa76143e5df641ebd6e0c7bcefc74de..23a284ed8453b8d99d6f6cfbf6fde11b66404cf2 100644 (file)
@@ -31,6 +31,7 @@
 #include "VISU_ConvertorDef.hxx"
 #include <vtkDataSetAlgorithm.h>
 
+class vtkDataSetAttributes;
 
 //----------------------------------------------------------------------------
 class VISU_PIPELINE_EXPORT VISU_Extractor : public vtkDataSetAlgorithm
@@ -68,6 +69,8 @@ protected:
   int
   RequestData(vtkInformation *, vtkInformationVector **, vtkInformationVector *);
 
+  void CutFieldData(int theNbElems, vtkDataSetAttributes *theData);
+
   int myScalarMode;
   VISU::TGaussMetric myGaussMetric;
 };