]> SALOME platform Git repositories - modules/visu.git/commitdiff
Salome HOME
To introduce ID's mapping
authorapo <apo@opencascade.com>
Mon, 29 Aug 2005 06:15:12 +0000 (06:15 +0000)
committerapo <apo@opencascade.com>
Mon, 29 Aug 2005 06:15:12 +0000 (06:15 +0000)
src/CONVERTOR/VISU_Convertor.hxx
src/CONVERTOR/VISU_Convertor_impl.cxx
src/CONVERTOR/VISU_Convertor_impl.hxx
src/CONVERTOR/VISU_MedConvertor.cxx
src/CONVERTOR/VISU_MedConvertor.hxx

index e7d04964b62a791dc7d0c21a8d3f2dde41934c52..718e18ffb8f37db97f2f6ceaeaa02a93c047fbbe 100644 (file)
@@ -137,11 +137,11 @@ namespace VISU
   {
     virtual 
     vtkIdType 
-    GetNodeObjId(int theVtkI) const = 0;
+    GetNodeObjID(int theVtkI) const = 0;
 
     virtual 
     vtkIdType 
-    GetElemObjId(int theVtkI) const = 0;
+    GetElemObjID(int theVtkI) const = 0;
   };
 
 
@@ -174,7 +174,7 @@ namespace VISU
   {
     virtual 
     TGaussPointID 
-    GetObjId(int theVtkI) const = 0;
+    GetObjID(int theVtkI) const = 0;
   };
 
 
index ffda7dc5707d246f3a04e08349aa6361ed2a4b99..65973f9f43cec3908befed78e451b8f63841ddef 100644 (file)
@@ -153,10 +153,12 @@ namespace VISU
   void
   TNamedPointCoords
   ::Init(vtkIdType theNbPoints,
-        vtkIdType theDim)
+        vtkIdType theDim,
+        const TVectorID& theVectorID)
   {
     TPointCoords::Init(theNbPoints,theDim);
     myPointsDim.resize(theDim);
+    myVectorID = theVectorID;
   }
 
   std::string&
@@ -173,8 +175,20 @@ namespace VISU
     return myPointsDim[theDim];
   }
 
+  vtkIdType
+  TNamedPointCoords
+  ::GetObjID(vtkIdType theID) const
+  {
+    if(myVectorID.empty())
+      return theID;
+    else
+      return myVectorID[theID];
+  }
+
+
   //---------------------------------------------------------------
   TMeshImpl::TMeshImpl():
+    myNamedPointCoords(new TNamedPointCoords()),
     myPoints(vtkPoints::New()),
     myNbPoints(0) 
   {
@@ -192,14 +206,14 @@ namespace VISU
 
   vtkIdType
   TSubProfileImpl
-  ::GetNodeObjId(vtkIdType theID) const
+  ::GetNodeObjID(vtkIdType theID) const
   {
     return theID;
   }
 
   vtkIdType
   TSubProfileImpl
-  ::GetElemObjId(vtkIdType theID) const
+  ::GetElemObjID(vtkIdType theID) const
   {
     return theID;
   }
@@ -225,24 +239,24 @@ namespace VISU
 
   vtkIdType
   TProfileImpl
-  ::GetNodeObjId(vtkIdType theID) const
+  ::GetNodeObjID(vtkIdType theID) const
   {
     vtkIdType anInputID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
     vtkIdType aID = anAppendFilter->GetNodeObjId(theID,anInputID);
     const TSubProfileImpl& aSubProfileImpl = mySubProfileArr[anInputID];
-    return aSubProfileImpl.GetNodeObjId(aID);
+    return aSubProfileImpl.GetNodeObjID(aID);
   }
 
   vtkIdType
   TProfileImpl
-  ::GetElemObjId(vtkIdType theID) const
+  ::GetElemObjID(vtkIdType theID) const
   {
     vtkIdType anInputID;
     const TVTKAppendFilter& anAppendFilter = GetFilter();
     vtkIdType aID = anAppendFilter->GetElemObjId(theID,anInputID);
     const TSubProfileImpl& aSubProfileImpl = mySubProfileArr[anInputID];
-    return aSubProfileImpl.GetElemObjId(aID);
+    return aSubProfileImpl.GetElemObjID(aID);
   }
 
 
@@ -287,7 +301,7 @@ namespace VISU
   //---------------------------------------------------------------
   TGaussPointID 
   TGaussMeshImpl
-  ::GetObjId(int theVtkI) const
+  ::GetObjID(int theVtkI) const
   {
     TGaussPointID aRetVID;
     int aID, anIndexDS;
index 0d95088c819fe65d052613dd19376c49efea90f2..64b33f0b115de3695a5aaa5d55d3f359e38f6c51 100644 (file)
@@ -113,32 +113,41 @@ namespace VISU
     const TVTKPoints&
     GetPoints() const { return myPoints;}
   };
+  typedef SharedPtr<TPointCoords> PPointCoords;
 
 
   //---------------------------------------------------------------
+  typedef TVector<vtkIdType> TVectorID;
   class TNamedPointCoords: public virtual TPointCoords
   {
     typedef TVector<std::string> TPointsDim;
     TPointsDim myPointsDim;
+    TVectorID myVectorID;
 
   public:
 
     void
     Init(vtkIdType theNbPoints,
-        vtkIdType theDim);
+        vtkIdType theDim,
+        const TVectorID& theVectorID = TVectorID());
     
     std::string&
     GetName(vtkIdType theDim);
     
     const std::string&
     GetName(vtkIdType theDim) const;
+
+    virtual
+    vtkIdType
+    GetObjID(vtkIdType theID) const;
   };
+  typedef SharedPtr<TNamedPointCoords> PNamedPointCoords;
 
 
   //---------------------------------------------------------------
   struct TMeshImpl: virtual TMesh, virtual TIsVTKDone
   {
-    TNamedPointCoords myNamedPointCoords;
+    PNamedPointCoords myNamedPointCoords;
 
     TVTKPoints myPoints;
     vtkIdType myNbPoints;
@@ -162,11 +171,11 @@ namespace VISU
 
     virtual 
     vtkIdType 
-    GetNodeObjId(int theVtkI) const;
+    GetNodeObjID(int theVtkI) const;
 
     virtual 
     vtkIdType 
-    GetElemObjId(int theVtkI) const;
+    GetElemObjID(int theVtkI) const;
 
     ESubMeshStatus myStatus;
     TSubMeshID mySubMeshID;
@@ -185,11 +194,13 @@ namespace VISU
 
     virtual 
     vtkIdType 
-    GetNodeObjId(int theVtkI) const;
+    GetNodeObjID(int theVtkI) const;
 
     virtual 
     vtkIdType 
-    GetElemObjId(int theVtkI) const;
+    GetElemObjID(int theVtkI) const;
+
+    PNamedPointCoords myNamedPointCoords;
 
     TSubProfileArr mySubProfileArr;
     TGeom2SubProfile myGeom2SubProfile;
@@ -233,7 +244,7 @@ namespace VISU
 
     virtual
     TGaussPointID
-    GetObjId(int theVtkI) const;
+    GetObjID(int theVtkI) const;
 
     TGaussSubMeshArr myGaussSubMeshArr;
     TGeom2GaussSubMesh myGeom2GaussSubMesh;
index 0023e3b8e2ce0ddbef70a2ff992686b3eee14b25..330547231a8df48a87a1caeff78935e9177905d6 100644 (file)
@@ -738,16 +738,19 @@ namespace
   //---------------------------------------------------------------
   vtkIdType
   TMEDSubProfile
-  ::GetNodeObjId(vtkIdType theID) const
+  ::GetNodeObjID(vtkIdType theID) const
   {
-    return theID;
+    return myNamedPointCoords->GetObjID(theID);
   }
   
   vtkIdType
   TMEDSubProfile
-  ::GetElemObjId(vtkIdType theID) const
+  ::GetElemObjID(vtkIdType theID) const
   {
-    return theID;
+    if(myIsElemNum)
+      return myElemNum[theID];
+    else
+      return theID;
   }
   
 
@@ -985,7 +988,10 @@ VISU_MedConvertor
       TInt aDim = theMesh->myDim;
 
       TNamedPointCoords& aCoords = theMesh->myNamedPointCoords;
-      aCoords.Init(aNbElem,aDim);
+      if(EBooleen anIsNodeNum = aNodeInfo->IsElemNum())
+       aCoords.Init(aNbElem,aDim,aNodeInfo->myElemNum);
+      else
+       aCoords.Init(aNbElem,aDim);
 
       for(int iDim = 0; iDim < aDim; iDim++)
        aCoords.GetName(iDim) = aNodeInfo->GetCoordName(iDim);
@@ -1231,7 +1237,9 @@ VISU_MedConvertor
 
 //---------------------------------------------------------------
 void
-LoadProfile(MED::TTimeStampVal& theTimeStampVal,
+LoadProfile(const MED::PWrapper& theMed,
+           VISU::PMEDMesh theMesh,
+           MED::TTimeStampVal& theTimeStampVal,
            VISU::TMEDValForTime& theValForTime,
            VISU::TMEDMeshOnEntity& theMeshOnEntity)
 {
@@ -1268,7 +1276,29 @@ LoadProfile(MED::TTimeStampVal& theTimeStampVal,
              endl);
     }
   }
+  {
+    const MED::PMeshInfo& aMeshInfo = theMesh->myMeshInfo;
+    MED::PNodeInfo aNodeInfo = theMed->GetPNodeInfo(aMeshInfo);
+    
+    TEntity aVEntity = theMeshOnEntity.myEntity;
+    MED::EEntiteMaillage aMEntity = VTKEntityToMED(aVEntity);
+    
+    const TGeom2SubProfile& aGeom2SubProfile = aProfile->myGeom2SubProfile;
+    TGeom2SubProfile::const_iterator anIter = aGeom2SubProfile.begin();
+    for(; anIter != aGeom2SubProfile.end(); anIter++){
+      const PMEDSubProfile& aSubProfile = anIter->second;
+      MED::EGeometrieElement aMGeom = aSubProfile->myMGeom;
+      MED::PCellInfo aCellInfo = theMed->GetPCellInfo(aMeshInfo,
+                                                     aMEntity, 
+                                                     aMGeom);
+      
+      aSubProfile->myIsElemNum = aCellInfo->IsElemNum();
+      if(aSubProfile->myIsElemNum)
+       aSubProfile->myElemNum = aCellInfo->myElemNum;
 
+      aSubProfile->myNamedPointCoords = theMesh->myNamedPointCoords;
+    }
+  }
   aProfile->myIsDone = true;
 }
 
@@ -1438,7 +1468,9 @@ LoadValForTime(const MED::PWrapper& theMed,
                             aMKey2Profile,
                             aKey2Gauss);
 
-  LoadProfile(aTimeStampVal,
+  LoadProfile(theMed,
+             theMesh,
+             aTimeStampVal,
              theValForTime,
              theMeshOnEntity);
   
index 11421033aea467faa9eb0699020b465151e0acaa..fd51f28fd5cba609caae5467010ed3dd41b5f455 100644 (file)
@@ -31,16 +31,18 @@ namespace VISU
   struct TMEDSubProfile: virtual TSubProfileImpl
   {
     MED::EGeometrieElement myMGeom;
+
     MED::EBooleen myIsElemNum;
     MED::TElemNum myElemNum;
+    PNamedPointCoords myNamedPointCoords;
 
     virtual 
     vtkIdType 
-    GetNodeObjId(vtkIdType theID) const;
+    GetNodeObjID(vtkIdType theID) const;
 
     virtual 
     vtkIdType 
-    GetElemObjId(vtkIdType theID) const;
+    GetElemObjID(vtkIdType theID) const;
 
   };
   typedef SharedPtr<TMEDSubProfile> PMEDSubProfile;