]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
Fix for Bug IPAL11380
authorapo <apo@opencascade.com>
Thu, 9 Feb 2006 13:12:37 +0000 (13:12 +0000)
committerapo <apo@opencascade.com>
Thu, 9 Feb 2006 13:12:37 +0000 (13:12 +0000)
   Position of the parent mesh element is incorrect:

src/CONVERTOR/VISU_Convertor_impl.cxx
src/CONVERTOR/VISU_Convertor_impl.hxx
src/CONVERTOR/VISU_MedConvertor.cxx
src/CONVERTOR/VISU_MedConvertor.hxx

index ba2a8ad5a02b542fc21907c9bb9f915c677b9111..8d1fe336eba9086ae6cb4f4679757a6beab0a9f3 100644 (file)
@@ -326,11 +326,11 @@ namespace VISU
     if(myIsAll)
       return myMeshOnEntity->GetElemObjID(theID);
 
-    vtkIdType anInputID;
+    vtkIdType anInputID, aStartID, anInputDataSetID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
-    vtkIdType aID = anAppendFilter->GetCellInputID(theID,anInputID);
-    PSubProfileImpl aSubProfileImpl = mySubProfileArr[anInputID];
-    return aSubProfileImpl->GetElemObjID(aID);
+    anAppendFilter->GetCellInputID(theID,anInputID,aStartID,anInputDataSetID);
+    PSubProfileImpl aSubProfileImpl = mySubProfileArr[anInputDataSetID];
+    return aSubProfileImpl->GetElemObjID(anInputID);
   }
 
   vtkIdType
@@ -383,12 +383,12 @@ namespace VISU
     if(myIsAll)
       return myMeshOnEntity->GetElemName(theObjID);
 
-    vtkIdType anInputID;
     vtkIdType aVTKId = GetElemVTKID(theObjID);
+    vtkIdType anInputID, aStartID, anInputDataSetID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
-    vtkIdType aSubID = anAppendFilter->GetCellInputID(aVTKId,anInputID);
-    PSubProfileImpl aSubProfileImpl = mySubProfileArr[anInputID];
-    vtkIdType anEntityObjId = aSubProfileImpl->GetElemObjID(aSubID);
+    anAppendFilter->GetCellInputID(aVTKId,anInputID,aStartID,anInputDataSetID);
+    PSubProfileImpl aSubProfileImpl = mySubProfileArr[anInputDataSetID];
+    vtkIdType anEntityObjId = aSubProfileImpl->GetElemObjID(anInputID);
     return myMeshOnEntity->GetElemName(anEntityObjId);
   }
 
@@ -475,9 +475,10 @@ namespace VISU
   
   TGaussPointID
   TGaussSubMeshImpl
-  ::GetObjID(vtkIdType theID) const
+  ::GetObjID(vtkIdType theID,
+            vtkIdType theStartID) const
   {
-    TCellID aCellID = theID / myGauss->myNbPoints;
+    TCellID aCellID = theStartID + theID / myGauss->myNbPoints;
     TLocalPntID aLocalPntID = theID % myGauss->myNbPoints;
     
     return TGaussPointID(aCellID,aLocalPntID);
@@ -514,12 +515,12 @@ namespace VISU
   TGaussMeshImpl
   ::GetObjID(vtkIdType theID) const
   {
+    vtkIdType anInputID, aStartId, anInputDataSetID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
-    vtkIdType anInputDataSetID;
-    vtkIdType anInputID = anAppendFilter->GetCellInputID(theID,anInputDataSetID);
+    anAppendFilter->GetCellInputID(theID,anInputID,aStartId,anInputDataSetID);
     const TGaussSubMeshImpl& aSubMeshImpl = myGaussSubMeshArr[anInputDataSetID];
 
-    return aSubMeshImpl.GetObjID(anInputID);
+    return aSubMeshImpl.GetObjID(anInputID,aStartId);
   }
   
   TVTKOutput* 
@@ -601,11 +602,11 @@ namespace VISU
   TMeshOnEntityImpl
   ::GetElemObjID(vtkIdType theID) const
   {
-    vtkIdType anInputID;
+    vtkIdType anInputID, aStartId, anInputDataSetID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
-    vtkIdType aID = anAppendFilter->GetCellInputID(theID,anInputID);
-    const PSubMeshImpl& aSubMesh = mySubMeshArr[anInputID];
-    return aSubMesh->GetElemObjID(aID);
+    anAppendFilter->GetCellInputID(theID,anInputID,aStartId,anInputDataSetID);
+    const PSubMeshImpl& aSubMesh = mySubMeshArr[anInputDataSetID];
+    return aSubMesh->GetElemObjID(anInputID);
   }
 
   std::string 
@@ -620,11 +621,11 @@ namespace VISU
   ::GetElemName(vtkIdType theObjID) const
   {
     vtkIdType aVTKId = GetElemVTKID(theObjID);
-    vtkIdType anInputID;
+    vtkIdType anInputID, aStartId, anInputDataSetID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
-    vtkIdType aSubID = anAppendFilter->GetCellInputID(aVTKId,anInputID);
-    const PSubMeshImpl& aSubMesh = mySubMeshArr[anInputID];
-    return aSubMesh->GetElemName(aSubID);
+    anAppendFilter->GetCellInputID(aVTKId,anInputID,aStartId,anInputDataSetID);
+    const PSubMeshImpl& aSubMesh = mySubMeshArr[anInputDataSetID];
+    return aSubMesh->GetElemName(anInputID);
   }
 
   //---------------------------------------------------------------
@@ -704,11 +705,11 @@ namespace VISU
   TGroupImpl
   ::GetElemObjID(vtkIdType theID) const
   {
-    vtkIdType anInputID;
+    vtkIdType anInputID, aStartId, anInputDataSetID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
-    vtkIdType anID = anAppendFilter->GetCellInputID(theID,anInputID);
-    const PFamilyImpl& aFamily = myFamilyArr[anInputID];
-    return aFamily->GetElemObjID(anID);
+    anAppendFilter->GetCellInputID(theID,anInputID,aStartId,anInputDataSetID);
+    const PFamilyImpl& aFamily = myFamilyArr[anInputDataSetID];
+    return aFamily->GetElemObjID(anInputID);
   }
 
   vtkIdType 
index c4efc7d25c80c2471f5f83b6b4651e4869cdcbad..20889db78fd36c8067144dc002bdd53b636bda29 100644 (file)
@@ -456,7 +456,8 @@ namespace VISU
     //! To implement the TGaussPtsIDMapper::GetObjID
     virtual
     TGaussPointID
-    GetObjID(vtkIdType theID) const;
+    GetObjID(vtkIdType theID,
+            vtkIdType theStartID) const;
     
     PGaussImpl myGauss; //<! Keep reference to corresponding TGauss structure
 
index 2eaddcfe6d7449aa0b81033270e69e4bc8fdec3b..c5e948621e77a2b588710d97d5a373bda8a587bf 100644 (file)
@@ -894,14 +894,16 @@ namespace
   //---------------------------------------------------------------
   TGaussPointID
   TMEDGaussSubMesh
-  ::GetObjID(vtkIdType theID) const
+  ::GetObjID(vtkIdType theID,
+            vtkIdType theStartID) const
   {
-    TInt aNbPoints = myGauss->myNbPoints;
-    TCellID aCellID = theID / aNbPoints;
+    TCellID aCellID = theID / myGauss->myNbPoints;
+    TLocalPntID aLocalPntID = theID % myGauss->myNbPoints;
+    
     if(myIsElemNum)
       aCellID = myElemNum[aCellID];
-
-    TLocalPntID aLocalPntID = theID % aNbPoints;
+    else
+      aCellID += theStartID;
 
     return TGaussPointID(aCellID,aLocalPntID);
   }
@@ -913,7 +915,10 @@ namespace
   ::Init(const MED::PElemInfo& theElemInfo)
   {
     myIsElemNum = theElemInfo->IsElemNum();
-    myElemNum = theElemInfo->myElemNum;
+
+    if(myIsElemNum)
+      myElemNum = theElemInfo->myElemNum;
+
     if(theElemInfo->IsElemNames())
       myElemInfo = theElemInfo;
   }
index abcb6ad70549e9b5ea9ef8515fa5eb325d978c60..0ce3377aa8bda5b49b3ef9a46e3ea43c6ef4dba2 100644 (file)
@@ -99,7 +99,8 @@ namespace VISU
 
     virtual
     TGaussPointID
-    GetObjID(vtkIdType theID) const;
+    GetObjID(vtkIdType theID,
+            vtkIdType theStartID) const;
   };
   typedef SharedPtr<TMEDGaussSubMesh> PMEDGaussSubMesh;