OGS
XmlStnInterface.cpp
Go to the documentation of this file.
1 
15 #include "XmlStnInterface.h"
16 
17 #include <QFile>
18 #include <QtXml/QDomDocument>
19 #include <limits>
20 #include <rapidxml.hpp>
21 
22 #include "BaseLib/DateTools.h"
23 #include "BaseLib/FileTools.h"
24 #include "BaseLib/Logging.h"
25 #include "GeoLib/GEOObjects.h"
26 #include "GeoLib/StationBorehole.h"
27 
28 namespace GeoLib
29 {
30 namespace IO
31 {
33  : XMLQtInterface("OpenGeoSysSTN.xsd"), _geo_objs(geo_objs)
34 {
35 }
36 
37 int XmlStnInterface::readFile(const QString& fileName)
38 {
39  if (XMLQtInterface::readFile(fileName) == 0)
40  {
41  return 0;
42  }
43 
44  QDomDocument doc("OGS-STN-DOM");
45  doc.setContent(getContent());
46  QDomElement docElement =
47  doc.documentElement(); // root element, used for identifying file-type
48  if (docElement.nodeName().compare("OpenGeoSysSTN"))
49  {
50  ERR("XmlStnInterface::readFile(): Unexpected XML root.");
51  return 0;
52  }
53 
54  QDomNodeList lists = docElement.childNodes();
55  for (int i = 0; i < lists.count(); i++)
56  {
57  // read all the station lists
58  QDomNodeList stationList = lists.at(i).childNodes();
59  auto stations = std::make_unique<std::vector<GeoLib::Point*>>();
60  std::string stnName("[NN]");
61 
62  for (int j = 0; j < stationList.count(); j++)
63  {
64  const QDomNode station_node(stationList.at(j));
65  const QString station_type(station_node.nodeName());
66  if (station_type.compare("name") == 0)
67  {
68  stnName = station_node.toElement().text().toStdString();
69  }
70  else if (station_type.compare("stations") == 0)
71  {
72  readStations(station_node, stations.get(),
73  fileName.toStdString());
74  }
75  else if (station_type.compare("boreholes") == 0)
76  {
77  readStations(station_node, stations.get(),
78  fileName.toStdString());
79  }
80  }
81 
82  if (!stations->empty())
83  {
84  _geo_objs.addStationVec(std::move(stations), stnName);
85  }
86  }
87 
88  return 1;
89 }
90 
91 void XmlStnInterface::readStations(const QDomNode& stationsRoot,
92  std::vector<GeoLib::Point*>* stations,
93  const std::string& station_file_name)
94 {
95  QDomElement station = stationsRoot.firstChildElement();
96  while (!station.isNull())
97  {
98  if (station.hasAttribute("id") && station.hasAttribute("x") &&
99  station.hasAttribute("y"))
100  {
101  std::string stationName("[NN]");
102  std::string sensor_data_file_name;
103  std::string boreholeDate("0000-00-00");
104  double boreholeDepth(0.0);
105  double stationValue(0.0);
106 
107  QDomNodeList stationFeatures = station.childNodes();
108  for (int i = 0; i < stationFeatures.count(); i++)
109  {
110  // check for general station features
111  const QDomNode feature_node(stationFeatures.at(i));
112  const QString feature_name(feature_node.nodeName());
113  const QString element_text(feature_node.toElement().text());
114  if (feature_name.compare("name") == 0)
115  {
116  stationName = element_text.toStdString();
117  }
118  if (feature_name.compare("sensordata") == 0)
119  {
120  sensor_data_file_name = element_text.toStdString();
121  /* add other station features here */
122 
123  // check for general borehole features
124  }
125  else if (feature_name.compare("value") == 0)
126  {
127  stationValue = element_text.toDouble();
128  }
129  else if (feature_name.compare("bdepth") == 0)
130  {
131  boreholeDepth = element_text.toDouble();
132  }
133  else if (feature_name.compare("bdate") == 0)
134  {
135  boreholeDate = element_text.toStdString();
136  }
137  /* add other borehole features here */
138  }
139 
140  double zVal = (station.hasAttribute("z"))
141  ? station.attribute("z").toDouble()
142  : 0.0;
143 
144  if (station.nodeName().compare("station") == 0)
145  {
146  GeoLib::Station* s =
147  new GeoLib::Station(station.attribute("x").toDouble(),
148  station.attribute("y").toDouble(),
149  zVal,
150  stationName);
151  s->setStationValue(stationValue);
152  if (!sensor_data_file_name.empty())
153  {
155  sensor_data_file_name, station_file_name));
156  }
157  stations->push_back(s);
158  }
159  else if (station.nodeName().compare("borehole") == 0)
160  {
163  stationName,
164  station.attribute("x").toDouble(),
165  station.attribute("y").toDouble(),
166  zVal,
167  boreholeDepth,
168  boreholeDate);
169  s->setStationValue(stationValue);
170  /* add stratigraphy to the borehole */
171  for (int j = 0; j < stationFeatures.count(); j++)
172  {
173  if (stationFeatures.at(j).nodeName().compare("strat") == 0)
174  {
175  this->readStratigraphy(stationFeatures.at(j), s);
176  }
177  }
178 
179  stations->push_back(s);
180  }
181  }
182  else
183  {
184  WARN(
185  "XmlStnInterface::readStations(): Attribute missing in "
186  "<station> tag.");
187  }
188  station = station.nextSiblingElement();
189  }
190 }
191 
192 void XmlStnInterface::readStratigraphy(const QDomNode& stratRoot,
193  GeoLib::StationBorehole* borehole)
194 {
195  // borehole->addSoilLayer((*borehole)[0], (*borehole)[1], (*borehole)[2],
196  // "");
197  double depth_check((*borehole)[2]);
198  QDomElement horizon = stratRoot.firstChildElement();
199  while (!horizon.isNull())
200  {
201  if (horizon.hasAttribute("id") && horizon.hasAttribute("x") &&
202  horizon.hasAttribute("y") && horizon.hasAttribute("z"))
203  {
204  std::string horizonName("[NN]");
205 
206  QDomNodeList horizonFeatures = horizon.childNodes();
207  for (int i = 0; i < horizonFeatures.count(); i++)
208  {
209  if (horizonFeatures.at(i).nodeName().compare("name") == 0)
210  {
211  horizonName =
212  horizonFeatures.at(i).toElement().text().toStdString();
213  }
214  }
215  /* add other horizon features here */
216 
217  double depth(horizon.attribute("z").toDouble());
218  if (std::abs(depth - depth_check) >
219  std::numeric_limits<double>::
220  epsilon()) // skip soil-layer if its thickness is zero
221  {
222  borehole->addSoilLayer(horizon.attribute("x").toDouble(),
223  horizon.attribute("y").toDouble(),
224  depth,
225  horizonName);
226  depth_check = depth;
227  }
228  else
229  {
230  WARN(
231  "XmlStnInterface::readStratigraphy(): Skipped layer '{:s}' "
232  "in borehole '{:s}' because of thickness 0.0.",
233  horizonName, borehole->getName());
234  }
235  }
236  else
237  {
238  WARN(
239  "XmlStnInterface::readStratigraphy(): Attribute missing in "
240  "<horizon> tag.");
241  }
242  horizon = horizon.nextSiblingElement();
243  }
244 }
245 
247 {
248  if (export_name.empty())
249  {
250  ERR("XmlStnInterface::write(): No station list specified.");
251  return false;
252  }
253 
254  out << "<?xml version=\"1.0\" encoding=\"ISO-8859-1\"?>\n"; // xml
255  // definition
256  out << "<?xml-stylesheet type=\"text/xsl\" "
257  "href=\"OpenGeoSysSTN.xsl\"?>\n\n"; // stylefile definition
258 
259  QDomDocument doc("OGS-STN-DOM");
260  QDomElement root = doc.createElement("OpenGeoSysSTN");
261  root.setAttribute("xmlns:ogs", "http://www.opengeosys.org");
262  root.setAttribute("xmlns:xsi", "http://www.w3.org/2001/XMLSchema-instance");
263 
264  const std::vector<GeoLib::Point*>* stations(
266  bool const is_borehole =
267  dynamic_cast<GeoLib::StationBorehole*>((*stations)[0]);
268 
269  doc.appendChild(root);
270  QDomElement stationListTag = doc.createElement("stationlist");
271  root.appendChild(stationListTag);
272 
273  QDomElement listNameTag = doc.createElement("name");
274  stationListTag.appendChild(listNameTag);
275  QDomText stationListNameText =
276  doc.createTextNode(QString::fromStdString(export_name));
277  listNameTag.appendChild(stationListNameText);
278  QString listType = is_borehole ? "boreholes" : "stations";
279  QDomElement stationsTag = doc.createElement(listType);
280  stationListTag.appendChild(stationsTag);
281 
282  bool useStationValue(false);
283  double sValue =
284  static_cast<GeoLib::Station*>((*stations)[0])->getStationValue();
285  std::size_t nStations(stations->size());
286  for (std::size_t i = 1; i < nStations; i++)
287  {
288  if ((static_cast<GeoLib::Station*>((*stations)[i])->getStationValue() -
289  sValue) < std::numeric_limits<double>::epsilon())
290  {
291  useStationValue = true;
292  break;
293  }
294  }
295 
296  for (std::size_t i = 0; i < nStations; i++)
297  {
298  QString stationType = is_borehole ? "borehole" : "station";
299  QDomElement stationTag = doc.createElement(stationType);
300  stationTag.setAttribute("id", QString::number(i));
301  stationTag.setAttribute(
302  "x", QString::number((*(*stations)[i])[0], 'f',
303  std::numeric_limits<double>::digits10));
304  stationTag.setAttribute(
305  "y", QString::number((*(*stations)[i])[1], 'f',
306  std::numeric_limits<double>::digits10));
307  stationTag.setAttribute(
308  "z", QString::number((*(*stations)[i])[2], 'f',
309  std::numeric_limits<double>::digits10));
310  stationsTag.appendChild(stationTag);
311 
312  QDomElement stationNameTag = doc.createElement("name");
313  stationTag.appendChild(stationNameTag);
314  QDomText stationNameText = doc.createTextNode(QString::fromStdString(
315  static_cast<GeoLib::Station*>((*stations)[i])->getName()));
316  stationNameTag.appendChild(stationNameText);
317 
318  if (useStationValue)
319  {
320  QDomElement stationValueTag = doc.createElement("value");
321  stationTag.appendChild(stationValueTag);
322  QDomText stationValueText = doc.createTextNode(
323  QString::number(static_cast<GeoLib::Station*>((*stations)[i])
324  ->getStationValue()));
325  stationValueTag.appendChild(stationValueText);
326  }
327 
328  if (is_borehole)
329  {
331  doc, stationTag,
332  static_cast<GeoLib::StationBorehole*>((*stations)[i]));
333  }
334  }
335 
336  std::string xml = doc.toString().toStdString();
337  out << xml;
338  return true;
339 }
340 
341 void XmlStnInterface::writeBoreholeData(QDomDocument& doc,
342  QDomElement& boreholeTag,
343  GeoLib::StationBorehole* borehole) const
344 {
345  QDomElement stationDepthTag = doc.createElement("bdepth");
346  boreholeTag.appendChild(stationDepthTag);
347  QDomText stationDepthText =
348  doc.createTextNode(QString::number(borehole->getDepth(), 'f'));
349  stationDepthTag.appendChild(stationDepthText);
350  if (std::abs(borehole->getDate()) > 0)
351  {
352  QDomElement stationDateTag = doc.createElement("bdate");
353  boreholeTag.appendChild(stationDateTag);
354  QDomText stationDateText = doc.createTextNode(
355  QString::fromStdString(BaseLib::date2string(borehole->getDate())));
356  stationDateTag.appendChild(stationDateText);
357  }
358 
359  std::vector<GeoLib::Point*> profile = borehole->getProfile();
360  std::vector<std::string> soilNames = borehole->getSoilNames();
361  std::size_t nHorizons(profile.size());
362 
363  if (nHorizons > 1)
364  {
365  QDomElement stratTag = doc.createElement("strat");
366  boreholeTag.appendChild(stratTag);
367 
368  for (std::size_t j = 1; j < nHorizons;
369  j++)
371  {
372  QDomElement horizonTag = doc.createElement("horizon");
373  horizonTag.setAttribute("id", QString::number(j));
374  horizonTag.setAttribute("x",
375  QString::number((*profile[j])[0], 'f'));
376  horizonTag.setAttribute("y",
377  QString::number((*profile[j])[1], 'f'));
378  horizonTag.setAttribute("z",
379  QString::number((*profile[j])[2], 'f'));
380  stratTag.appendChild(horizonTag);
381  QDomElement horizonNameTag = doc.createElement("name");
382  horizonTag.appendChild(horizonNameTag);
383  QDomText horizonNameText =
384  doc.createTextNode(QString::fromStdString(soilNames[j]));
385  horizonNameTag.appendChild(horizonNameText);
386  }
387  }
388 }
389 } // namespace IO
390 } // namespace GeoLib
Definition of date helper functions.
Filename manipulation routines.
Definition of the GEOObjects class.
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
void WARN(char const *fmt, Args const &... args)
Definition: Logging.h:37
Definition of the StationBorehole class.
std::string getName(std::string const &line)
Returns the name/title from the "Zone"-description.
Definition of the XmlStnInterface class.
std::ostringstream out
The stream to write to.
Definition: Writer.h:46
QByteArray const & getContent() const
Container class for geometric objects.
Definition: GEOObjects.h:61
void addStationVec(std::unique_ptr< std::vector< Point * >> stations, std::string &name)
Adds a vector of stations with the given name and colour to GEOObjects.
Definition: GEOObjects.cpp:122
const std::vector< GeoLib::Point * > * getStationVec(const std::string &name) const
Returns the station vector with the given name.
Definition: GEOObjects.cpp:131
void readStratigraphy(const QDomNode &stratRoot, GeoLib::StationBorehole *borehole)
Reads the stratigraphy of a borehole from an xml-file.
bool write() override
Writes the object to the internal stream. This method must be implemented by a subclass....
XmlStnInterface(GeoLib::GEOObjects &geo_objs)
void readStations(const QDomNode &stationsRoot, std::vector< GeoLib::Point * > *stations, const std::string &station_file_name)
Reads GeoLib::Station- or StationBorehole-objects from an xml-file.
int readFile(const QString &fileName) override
Reads an xml-file containing station object definitions into the GEOObjects used in the constructor (...
void writeBoreholeData(QDomDocument &doc, QDomElement &boreholeTag, GeoLib::StationBorehole *borehole) const
Writes borehole-specific data to a station-xml-file.
GeoLib::GEOObjects & _geo_objs
A borehole as a geometric object.
double getDate() const
Returns the date entry for the borehole.
const std::vector< std::string > & getSoilNames() const
static StationBorehole * createStation(const std::string &name, double x, double y, double z, double depth, const std::string &date="")
Creates a new borehole object based on the given parameters.
void addSoilLayer(double thickness, const std::string &soil_name)
Add a soil layer to the boreholes stratigraphy.
const std::vector< Point * > & getProfile() const
A Station (observation site) is basically a Point with some additional information.
Definition: Station.h:37
void addSensorDataFromCSV(const std::string &file_name)
Allows to add sensor data from a CSV file to the observation site.
Definition: Station.h:77
std::string const & getName() const
Returns the name of the station.
Definition: Station.h:60
void setStationValue(double station_value)
Allows to set a specific value for this station (e.g. for classification)
Definition: Station.h:74
std::string date2string(double ddate)
Definition: DateTools.cpp:65
std::string copyPathToFileName(const std::string &file_name, const std::string &source)
Definition: FileTools.cpp:196
bool readFile(std::string const &file_name, std::vector< std::unique_ptr< MeshLib::Mesh >> &meshes, DataType const export_type)
Reads the specified file and writes data into internal mesh vector.