OTB  10.0.0
Orfeo Toolbox
otbDEMHandler.h
Go to the documentation of this file.
1 /*
2  * Copyright (C) 2005-2024 Centre National d'Etudes Spatiales (CNES)
3  *
4  * This file is part of Orfeo Toolbox
5  *
6  * https://www.orfeo-toolbox.org/
7  *
8  * Licensed under the Apache License, Version 2.0 (the "License");
9  * you may not use this file except in compliance with the License.
10  * You may obtain a copy of the License at
11  *
12  * http://www.apache.org/licenses/LICENSE-2.0
13  *
14  * Unless required by applicable law or agreed to in writing, software
15  * distributed under the License is distributed on an "AS IS" BASIS,
16  * WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
17  * See the License for the specific language governing permissions and
18  * limitations under the License.
19  */
20 
21 #ifndef otbDEMHandler_h
22 #define otbDEMHandler_h
23 
24 #include "otbImage.h"
25 #include "otbGDALDatasetWrapper.h"
26 #include <string>
27 #include <list>
28 #include <vector>
29 #include <memory>
30 
31 namespace otb
32 {
33 
40 {
41 public:
42  virtual void Update() = 0;
43 protected:
44  DEMObserverInterface() = default;
45  ~DEMObserverInterface() = default;
48 };
49 
56 {
57 public:
58  virtual void AttachObserver(DEMObserverInterface *observer) = 0;
59  virtual void DetachObserver(DEMObserverInterface *observer) = 0;
60  virtual void Notify() const = 0;
61 protected:
62  DEMSubjectInterface() = default;
63  ~DEMSubjectInterface() = default;
66 };
67 
68 // Forward declaration of the class at this point. No need to expose the
69 // full class definition.
70 class DEMHandlerTLS;
71 
114 {
115 public:
116  using Self = DEMHandler;
117  using PointType = itk::Point<double, 2>;
118 
121 
125  void OpenDEMDirectory(std::string DEMDirectory);
126 
130  void OpenDEMFile(std::string path);
131 
135  bool IsValidDEMDirectory(const std::string& DEMDirectory) const;
136 
140  bool OpenGeoidFile(std::string geoidFile);
141 
151  double GetHeightAboveEllipsoid(double lon, double lat) const;
152 
153  double GetHeightAboveEllipsoid(const PointType& geoPoint) const;
154 
164  double GetHeightAboveMSL(double lon, double lat) const;
165  double GetHeightAboveMSL(const PointType& geoPoint) const;
167 
173  double GetGeoidHeight(double lon, double lat) const;
174  double GetGeoidHeight(const PointType& geoPoint) const;
176 
178  std::size_t GetDEMCount() const noexcept;
179 
180  double GetDefaultHeightAboveEllipsoid() const noexcept;
181 
182  void SetDefaultHeightAboveEllipsoid(double height);
183 
189  std::string const& GetDEMDirectory(std::size_t idx = 0) const;
190 
191  std::size_t GetNumberOfDEMDirectories() const noexcept
192  { return m_DEMDirectories.size(); }
193 
195  std::string const& GetGeoidFile() const noexcept;
196 
201 
206  void AttachObserver(DEMObserverInterface *observer) override {m_ObserverList.push_back(observer);};
207 
209  void DetachObserver(DEMObserverInterface *observer) override {m_ObserverList.remove(observer);};
210 
212  void Notify() const override;
213 
215  static constexpr char const* DEM_DATASET_PATH = "/vsimem/otb_dem_dataset.vrt";
216  static constexpr char const* DEM_WARPED_DATASET_PATH = "/vsimem/otb_dem_warped_dataset.vrt";
217  static constexpr char const* DEM_SHIFTED_DATASET_PATH = "/vsimem/otb_dem_shifted_dataset.vrt";
218 
238  DEMHandlerTLS const& GetHandlerForCurrentThread() const;
239 
240 protected:
246 
253 
254 private:
255 
259  std::shared_ptr<DEMHandlerTLS> DoFetchOrCreateHandler() const;
260 
270  void RegisterConfigurationInHandler(DEMHandlerTLS & handler) const;
271 
272  DEMHandler(const Self&) = delete;
273  void operator=(const Self&) = delete;
274 
277 
279  std::vector<std::string> m_DEMDirectories;
280 
282  std::vector<otb::GDALDatasetWrapper::Pointer> m_DatasetList;
283 
285  std::string m_GeoidFilename;
286 
288  std::list<DEMObserverInterface *> m_ObserverList;
289 
291  mutable std::vector<std::shared_ptr<DEMHandlerTLS>> m_tlses;
292 };
293 
294 
308 double GetHeightAboveEllipsoid(DEMHandlerTLS const&, double lon, double lat);
309 double GetHeightAboveMSL (DEMHandlerTLS const&, double lon, double lat);
310 double GetGeoidHeight (DEMHandlerTLS const&, double lon, double lat);
311 double GetHeightAboveEllipsoid(DEMHandlerTLS const&, itk::Point<double, 2> geoPoint);
312 double GetHeightAboveMSL (DEMHandlerTLS const&, itk::Point<double, 2> geoPoint);
313 double GetGeoidHeight (DEMHandlerTLS const&, itk::Point<double, 2> geoPoint);
315 
316 
317 }
318 #endif
Single access point for DEM data retrieval.
bool OpenGeoidFile(std::string geoidFile)
std::vector< std::string > m_DEMDirectories
double GetHeightAboveMSL(const PointType &geoPoint) const
std::vector< otb::GDALDatasetWrapper::Pointer > m_DatasetList
DEMHandlerTLS const & GetHandlerForCurrentThread() const
std::vcl_size_t GetDEMCount() const noexcept
double GetGeoidHeight(const PointType &geoPoint) const
static DEMHandler & GetInstance()
static constexpr char const * DEM_WARPED_DATASET_PATH
double GetHeightAboveEllipsoid(const PointType &geoPoint) const
void Notify() const override
std::vector< std::shared_ptr< DEMHandlerTLS > > m_tlses
itk::Point< double, 2 > PointType
std::shared_ptr< DEMHandlerTLS > DoFetchOrCreateHandler() const
void DetachObserver(DEMObserverInterface *observer) override
std::string const & GetDEMDirectory(std::vcl_size_t idx=0) const
bool IsValidDEMDirectory(const std::string &DEMDirectory) const
double GetDefaultHeightAboveEllipsoid() const noexcept
void RegisterConfigurationInHandler(DEMHandlerTLS &handler) const
std::vcl_size_t GetNumberOfDEMDirectories() const noexcept
double m_DefaultHeightAboveEllipsoid
double GetHeightAboveMSL(double lon, double lat) const
static constexpr char const * DEM_DATASET_PATH
double GetGeoidHeight(double lon, double lat) const
void ClearElevationParameters()
double GetHeightAboveEllipsoid(double lon, double lat) const
void AttachObserver(DEMObserverInterface *observer) override
static constexpr char const * DEM_SHIFTED_DATASET_PATH
void OpenDEMDirectory(std::string DEMDirectory)
void OpenDEMFile(std::string path)
std::string m_GeoidFilename
std::list< DEMObserverInterface * > m_ObserverList
void operator=(const Self &)=delete
void SetDefaultHeightAboveEllipsoid(double height)
std::string const & GetGeoidFile() const noexcept
DEMHandler(const Self &)=delete
Observer design pattern to keep track of DEM configuration changes.
Definition: otbDEMHandler.h:40
DEMObserverInterface(DEMObserverInterface const &)=delete
DEMObserverInterface & operator=(DEMObserverInterface const &)=delete
virtual void Update()=0
Observer design pattern to keep track of DEM configuration changes.
Definition: otbDEMHandler.h:56
virtual void DetachObserver(DEMObserverInterface *observer)=0
virtual void Notify() const =0
DEMSubjectInterface(DEMSubjectInterface const &)=delete
virtual void AttachObserver(DEMObserverInterface *observer)=0
DEMSubjectInterface & operator=(DEMSubjectInterface const &)=delete
The "otb" namespace contains all Orfeo Toolbox (OTB) classes.
double GetHeightAboveMSL(DEMHandlerTLS const &, double lon, double lat)
double GetHeightAboveEllipsoid(DEMHandlerTLS const &, double lon, double lat)
double GetGeoidHeight(DEMHandlerTLS const &, double lon, double lat)