]> SALOME platform Git repositories - modules/hydro.git/commitdiff
Salome HOME
debug multiple bathymetries on reloaded study
authorPaul RASCLE <paul.rascle@edf.fr>
Thu, 29 Jun 2017 21:25:43 +0000 (23:25 +0200)
committerPaul RASCLE <paul.rascle@edf.fr>
Thu, 29 Jun 2017 21:47:32 +0000 (23:47 +0200)
src/HYDROData/HYDROData_Bathymetry.cxx
src/HYDROData/HYDROData_Bathymetry.h
src/HYDROData/HYDROData_Document.cxx
src/HYDROData/HYDROData_Document.h

index cd6c718f1f1cbffc9e0050b414a136d952b10081..5fa2c82ebbe86e5e19491dfd7f012da035ece308 100644 (file)
@@ -62,12 +62,11 @@ const int BLOCK_SIZE = 1000;
 
 IMPLEMENT_STANDARD_RTTIEXT(HYDROData_Bathymetry, HYDROData_IAltitudeObject)
 
-//HYDROData_QuadtreeNode* HYDROData_Bathymetry::myQuadtree = 0;
-int HYDROData_Bathymetry::myQuadTreeNumber = 0;
+//int HYDROData_Bathymetry::myQuadTreeNumber = 0;
 std::map<int, HYDROData_QuadtreeNode*> HYDROData_Bathymetry::myQuadtrees;
 
 #ifndef LIGHT_MODE
-int HYDROData_Bathymetry::myDelaunayNumber = 0;
+//int HYDROData_Bathymetry::myDelaunayNumber = 0;
 std::map<int, vtkPolyData*> HYDROData_Bathymetry::myDelaunay2D;
 #endif
 
@@ -185,33 +184,19 @@ HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
       if (aLabel.IsNull())
         return 0;
 
-      Handle(TDataStd_RealArray) aCoordsArray;
-      if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
-        return 0;
-
-      myQuadTreeNumber++;
-      TDataStd_Integer::Set( myLab.FindChild( DataTag_Quadtree ), myQuadTreeNumber );
-      DEBTRACE("GetQuadtreeNodes init " << this << " " << myQuadTreeNumber);
-      HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
-      myQuadtrees[myQuadTreeNumber] = aQuadtree;
-
-      Nodes_3D* aListOfNodes = new Nodes_3D();
-
-      int index =0;
-      for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
+      int aQuadTreeNumber = 0;
+      Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
+      if ( ! aDocument.IsNull() )
         {
-          if (i + 3 > n + 1)
-            break;
-
-          double x = aCoordsArray->Value(i++);
-          double y = aCoordsArray->Value(i++);
-          double z = aCoordsArray->Value(i++);
-          gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
-          index++;
-          aListOfNodes->push_back(aPoint);
+          aQuadTreeNumber = aDocument->GetCountQuadtree();
+          DEBTRACE("aQuadTreeNumber " << aQuadTreeNumber);
+          aQuadTreeNumber++;
+          aDocument->SetCountQuadtree(aQuadTreeNumber);
         }
-      DEBTRACE("  GetQuadtreeNodes call setNodesAndCompute");
-      aQuadtree->setNodesAndCompute(aListOfNodes);
+      else
+          DEBTRACE("document.IsNull()");
+      DEBTRACE("compute Quadtree "<< aQuadTreeNumber);
+      HYDROData_QuadtreeNode* aQuadtree = ComputeQuadtreeNodes(aQuadTreeNumber);
       return aQuadtree;
     }
   else
@@ -221,13 +206,53 @@ HYDROData_QuadtreeNode* HYDROData_Bathymetry::GetQuadtreeNodes() const
         {
           if (myQuadtrees.find(aQuadtreeNum->Get()) != myQuadtrees.end())
             return myQuadtrees[aQuadtreeNum->Get()];
-          else DEBTRACE("myQuadtrees["<<aQuadtreeNum->Get()<<"] does not exists!");
+          else
+            {
+              DEBTRACE("recompute Quadtree "<< aQuadtreeNum->Get());
+              HYDROData_QuadtreeNode* aQuadtree = ComputeQuadtreeNodes(aQuadtreeNum->Get());
+              return aQuadtree;
+            }
         }
       else DEBTRACE("no attribute TDataStd_Integer");
     }
   return 0;
 }
 
+HYDROData_QuadtreeNode* HYDROData_Bathymetry::ComputeQuadtreeNodes( int key) const
+{
+  TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
+  if (aLabel.IsNull())
+    return 0;
+
+  Handle(TDataStd_RealArray) aCoordsArray;
+  if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
+    return 0;
+
+  TDataStd_Integer::Set( myLab.FindChild( DataTag_Quadtree ), key );
+  DEBTRACE("GetQuadtreeNodes init " << this << " " << key);
+  HYDROData_QuadtreeNode* aQuadtree = new HYDROData_QuadtreeNode(0, 30, 5, 0.);
+  myQuadtrees[key] = aQuadtree;
+
+  Nodes_3D* aListOfNodes = new Nodes_3D();
+
+  int index =0;
+  for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
+    {
+      if (i + 3 > n + 1)
+        break;
+
+      double x = aCoordsArray->Value(i++);
+      double y = aCoordsArray->Value(i++);
+      double z = aCoordsArray->Value(i++);
+      gpi_XYZ* aPoint = new gpi_XYZ(x, y, z, index);
+      index++;
+      aListOfNodes->push_back(aPoint);
+    }
+  DEBTRACE("  GetQuadtreeNodes call setNodesAndCompute");
+  aQuadtree->setNodesAndCompute(aListOfNodes);
+  return aQuadtree;
+}
+
 #ifndef LIGHT_MODE
 vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
 {
@@ -238,35 +263,19 @@ vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
       if (aLabel.IsNull())
         return 0;
 
-      Handle(TDataStd_RealArray) aCoordsArray;
-      if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
-        return 0;
-
-      myDelaunayNumber++;
-      TDataStd_Integer::Set( myLab.FindChild( DataTag_Delaunay ), myDelaunayNumber );
-      DEBTRACE("GetVtkDelaunay2D init " << this << " " << myDelaunayNumber);
-      vtkPoints *points = vtkPoints::New();
-      points->Allocate(aCoordsArray->Upper() +1);
-      for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
+      int aDelaunayNumber = 0;
+      Handle(HYDROData_Document) aDocument = HYDROData_Document::Document( myLab );
+      if ( ! aDocument.IsNull() )
         {
-          if (i + 3 > n + 1)
-            break;
-          double x = aCoordsArray->Value(i++);
-          double y = aCoordsArray->Value(i++);
-          double z = aCoordsArray->Value(i++);
-          vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
-          //DEBTRACE("  " << index);
+          aDelaunayNumber = aDocument->GetCountDelaunay();
+          DEBTRACE("aDelaunayNumber " << aDelaunayNumber);
+          aDelaunayNumber++;
+          aDocument->SetCountDelaunay(aDelaunayNumber);
         }
-      vtkPolyData* profile = vtkPolyData::New();
-      profile->SetPoints(points);
-      DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
-
-      vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
-      delaunay2D->SetInputData(profile);
-      delaunay2D->Update();
-      vtkPolyData* data = delaunay2D->GetOutput();
-      data->BuildLinks();
-      myDelaunay2D[myDelaunayNumber] = data;
+      else
+          DEBTRACE("document.IsNull()");
+      DEBTRACE("compute Delaunay "<< aDelaunayNumber);
+      vtkPolyData* data = ComputeVtkDelaunay2D(aDelaunayNumber);
       return data;
     }
   else
@@ -276,13 +285,58 @@ vtkPolyData* HYDROData_Bathymetry::GetVtkDelaunay2D() const
         {
           if (myDelaunay2D.find(aDelaunayNum->Get()) != myDelaunay2D.end())
             return myDelaunay2D[aDelaunayNum->Get()];
-          else DEBTRACE("myDelaunay2D["<<aDelaunayNum->Get()<<"] does not exists!");
+          else
+            {
+              DEBTRACE("recompute Delaunay "<< aDelaunayNum->Get());
+              vtkPolyData* data = ComputeVtkDelaunay2D(aDelaunayNum->Get());
+              return data;
+            }
         }
       else DEBTRACE("no attribute TDataStd_Integer");
     }
   return 0;
 }
+
+vtkPolyData* HYDROData_Bathymetry::ComputeVtkDelaunay2D(int key) const
+{
+  TDF_Label aLabel = myLab.FindChild(DataTag_AltitudePoints, false);
+  if (aLabel.IsNull())
+    return 0;
+
+  Handle(TDataStd_RealArray) aCoordsArray;
+  if (!aLabel.FindAttribute(TDataStd_RealArray::GetID(), aCoordsArray))
+    return 0;
+
+  TDataStd_Integer::Set( myLab.FindChild( DataTag_Delaunay ), key );
+  DEBTRACE("GetVtkDelaunay2D init " << this << " " << key);
+  vtkPoints *points = vtkPoints::New();
+  points->Allocate(aCoordsArray->Upper() +1);
+  for (int i = aCoordsArray->Lower(), n = aCoordsArray->Upper(); i <= n;)
+    {
+      if (i + 3 > n + 1)
+        break;
+      double x = aCoordsArray->Value(i++);
+      double y = aCoordsArray->Value(i++);
+      double z = aCoordsArray->Value(i++);
+      vtkIdType index = points->InsertNextPoint(x, y, z); // same index than in GetQuadtreeNodes
+      //DEBTRACE("  " << index);
+    }
+  vtkPolyData* profile = vtkPolyData::New();
+  profile->SetPoints(points);
+  DEBTRACE("Number of Points: "<< points->GetNumberOfPoints());
+
+  vtkDelaunay2D* delaunay2D = vtkDelaunay2D::New();
+  delaunay2D->SetInputData(profile);
+  delaunay2D->Update();
+  vtkPolyData* data = delaunay2D->GetOutput();
+  data->BuildLinks();
+  myDelaunay2D[key] = data;
+  return data;
+}
+
 #endif
+
+
 void HYDROData_Bathymetry::RemoveAltitudePoints()
 {
   TDF_Label aLabel = myLab.FindChild( DataTag_AltitudePoints, false );
index 5a7620e871a276472699da54815cb0d596100281..67a5868c816d078bd57fb859e85108e8fd1f9b8b 100644 (file)
@@ -96,10 +96,11 @@ public:
    */
   HYDRODATA_EXPORT virtual AltitudePoints   GetAltitudePoints(bool IsConvertToGlobal = false) const;
   HYDRODATA_EXPORT virtual HYDROData_QuadtreeNode* GetQuadtreeNodes() const;
+  HYDRODATA_EXPORT virtual HYDROData_QuadtreeNode* ComputeQuadtreeNodes( int key) const;
 
 #ifndef LIGHT_MODE
   HYDRODATA_EXPORT virtual vtkPolyData* GetVtkDelaunay2D() const;
-  //HYDRODATA_EXPORT bool interpolZtriangle(const gp_XY& point, vtkIdList* triangle, double* z);
+  HYDRODATA_EXPORT virtual vtkPolyData* ComputeVtkDelaunay2D(int key) const;
 #endif
 
   /**
@@ -167,10 +168,10 @@ private:
    */
   bool                                      importFromXYZFile( QFile&          theFile,
                                                                AltitudePoints& thePoints ) const;
-  static int myQuadTreeNumber;
+  //static int myQuadTreeNumber;
   static std::map<int, HYDROData_QuadtreeNode*> myQuadtrees;
 #ifndef LIGHT_MODE
-  static int myDelaunayNumber;
+  //static int myDelaunayNumber;
   static std::map<int, vtkPolyData*> myDelaunay2D;
 #endif
   bool                                      importFromASCFile( QFile&          theFile,
index 2fee31c361d416f89f499c056b6fbf66805f815f..b59211f6e3cf6e0703b456cdca13852c54645f9b 100644 (file)
@@ -49,6 +49,8 @@ static const int TAG_OBJECTS = 2; // tag of the objects sub-tree
 static const int TAG_HISTORY = 3; // tag of the history sub-tree (Root for History)
 static const int TAG_LOCAL_CS = 4; // tag of local coordinate system information
 static const int TAG_DEF_STRICKLER_COEFF = 5; // tag of default strickler coefficient
+static const int TAG_COUNT_QUADTREE = 6; // tag of number of quadtrees created so far
+static const int TAG_COUNT_DELAUNAY = 7; // tag of number of Delaunay triangulations created so far
 static const gp_Pnt2d DEFAULT_LOCAL_CS( 0, 0 );
 
 using namespace std;
@@ -231,6 +233,56 @@ void HYDROData_Document::SetDefaultStricklerCoefficient( double theCoeff ) const
     }
 }
 
+int HYDROData_Document::GetCountQuadtree() const
+{
+  int nbQuad = 0;
+  TDF_Label aLabel = myDoc->Main().FindChild(TAG_COUNT_QUADTREE, Standard_False);
+  if ( !aLabel.IsNull() )
+  {
+      Handle(TDataStd_Integer) anAttr;
+      if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anAttr ) )
+          nbQuad = anAttr->Get();
+  }
+  return nbQuad;
+}
+
+void HYDROData_Document::SetCountQuadtree( int nbQuad) const
+{
+  TDF_Label aLabel = myDoc->Main().FindChild(TAG_COUNT_QUADTREE);
+  if ( !aLabel.IsNull() )
+  {
+      Handle(TDataStd_Integer) anAttr;
+      if ( !aLabel.FindAttribute( TDataStd_Integer::GetID(), anAttr ) )
+          aLabel.AddAttribute( anAttr = new TDataStd_Integer() );
+      anAttr->Set( nbQuad );
+  }
+}
+
+int HYDROData_Document::GetCountDelaunay() const
+{
+  int nbDelaunay = 0;
+  TDF_Label aLabel = myDoc->Main().FindChild(TAG_COUNT_DELAUNAY, Standard_False);
+  if ( !aLabel.IsNull() )
+  {
+      Handle(TDataStd_Integer) anAttr;
+      if ( aLabel.FindAttribute( TDataStd_Integer::GetID(), anAttr ) )
+        nbDelaunay = anAttr->Get();
+  }
+  return nbDelaunay;
+}
+
+void HYDROData_Document::SetCountDelaunay( int nbDelaunay) const
+{
+  TDF_Label aLabel = myDoc->Main().FindChild(TAG_COUNT_DELAUNAY);
+  if ( !aLabel.IsNull() )
+  {
+      Handle(TDataStd_Integer) anAttr;
+      if ( !aLabel.FindAttribute( TDataStd_Integer::GetID(), anAttr ) )
+          aLabel.AddAttribute( anAttr = new TDataStd_Integer() );
+      anAttr->Set( nbDelaunay );
+  }
+}
+
 bool HYDROData_Document::DumpToPython( const QString& thePyScriptPath,
                                        const bool     theIsMultiFile ) const
 {
index c342ec3ba9413b5ab0596b992ea793ff12a2b77d..4c0368f81d141b2ff4c450c15d747e02d658c30a 100644 (file)
@@ -228,6 +228,11 @@ public:
   //! Sets default strickler coefficient
   HYDRODATA_EXPORT void SetDefaultStricklerCoefficient( double ) const;
 
+  HYDRODATA_EXPORT int GetCountQuadtree() const;
+  HYDRODATA_EXPORT void SetCountQuadtree( int ) const;
+  HYDRODATA_EXPORT int GetCountDelaunay() const;
+  HYDRODATA_EXPORT void SetCountDelaunay( int ) const;
+
   HYDRODATA_EXPORT QColor GetAssociatedColor( const QString& theStricklerType, const Handle(HYDROData_StricklerTable)& theTable ) const;
 
 protected: