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"
27
28namespace GeoLib
29{
30namespace IO
31{
33 : XMLQtInterface("OpenGeoSysSTN.xsd"), _geo_objs(geo_objs)
34{
35}
36
37int 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 std::vector<GeoLib::Point*> stations;
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, fileName.toStdString());
73 }
74 else if (station_type.compare("boreholes") == 0)
75 {
76 readStations(station_node, stations, fileName.toStdString());
77 }
78 }
79
80 if (!stations.empty())
81 {
82 _geo_objs.addStationVec(std::move(stations), stnName);
83 }
84 }
85
86 return 1;
87}
88
89void XmlStnInterface::readStations(const QDomNode& stationsRoot,
90 std::vector<GeoLib::Point*>& stations,
91 const std::string& station_file_name)
92{
93 QDomElement station = stationsRoot.firstChildElement();
94 while (!station.isNull())
95 {
96 if (station.hasAttribute("id") && station.hasAttribute("x") &&
97 station.hasAttribute("y"))
98 {
99 std::string stationName("[NN]");
100 std::string sensor_data_file_name;
101 std::string boreholeDate("0000-00-00");
102 double boreholeDepth(0.0);
103 double stationValue(0.0);
104
105 QDomNodeList stationFeatures = station.childNodes();
106 for (int i = 0; i < stationFeatures.count(); i++)
107 {
108 // check for general station features
109 const QDomNode feature_node(stationFeatures.at(i));
110 const QString feature_name(feature_node.nodeName());
111 const QString element_text(feature_node.toElement().text());
112 if (feature_name.compare("name") == 0)
113 {
114 stationName = element_text.toStdString();
115 }
116 if (feature_name.compare("sensordata") == 0)
117 {
118 sensor_data_file_name = element_text.toStdString();
119 /* add other station features here */
120
121 // check for general borehole features
122 }
123 else if (feature_name.compare("value") == 0)
124 {
125 stationValue = element_text.toDouble();
126 }
127 else if (feature_name.compare("bdepth") == 0)
128 {
129 boreholeDepth = element_text.toDouble();
130 }
131 else if (feature_name.compare("bdate") == 0)
132 {
133 boreholeDate = element_text.toStdString();
134 }
135 /* add other borehole features here */
136 }
137
138 double zVal = (station.hasAttribute("z"))
139 ? station.attribute("z").toDouble()
140 : 0.0;
141
142 if (station.nodeName().compare("station") == 0)
143 {
144 GeoLib::Station* s =
145 new GeoLib::Station(station.attribute("x").toDouble(),
146 station.attribute("y").toDouble(),
147 zVal,
148 stationName);
149 s->setStationValue(stationValue);
150 if (!sensor_data_file_name.empty())
151 {
153 sensor_data_file_name, station_file_name));
154 }
155 stations.push_back(s);
156 }
157 else if (station.nodeName().compare("borehole") == 0)
158 {
161 stationName,
162 station.attribute("x").toDouble(),
163 station.attribute("y").toDouble(),
164 zVal,
165 boreholeDepth,
166 boreholeDate);
167 s->setStationValue(stationValue);
168 /* add stratigraphy to the borehole */
169 for (int j = 0; j < stationFeatures.count(); j++)
170 {
171 if (stationFeatures.at(j).nodeName().compare("strat") == 0)
172 {
173 this->readStratigraphy(stationFeatures.at(j), s);
174 }
175 }
176
177 stations.push_back(s);
178 }
179 }
180 else
181 {
182 WARN(
183 "XmlStnInterface::readStations(): Attribute missing in "
184 "<station> tag.");
185 }
186 station = station.nextSiblingElement();
187 }
188}
189
190void XmlStnInterface::readStratigraphy(const QDomNode& stratRoot,
191 GeoLib::StationBorehole* borehole)
192{
193 // borehole->addSoilLayer((*borehole)[0], (*borehole)[1], (*borehole)[2],
194 // "");
195 double depth_check((*borehole)[2]);
196 QDomElement horizon = stratRoot.firstChildElement();
197 while (!horizon.isNull())
198 {
199 if (horizon.hasAttribute("id") && horizon.hasAttribute("x") &&
200 horizon.hasAttribute("y") && horizon.hasAttribute("z"))
201 {
202 std::string horizonName("[NN]");
203
204 QDomNodeList horizonFeatures = horizon.childNodes();
205 for (int i = 0; i < horizonFeatures.count(); i++)
206 {
207 if (horizonFeatures.at(i).nodeName().compare("name") == 0)
208 {
209 horizonName =
210 horizonFeatures.at(i).toElement().text().toStdString();
211 }
212 }
213 /* add other horizon features here */
214
215 double depth(horizon.attribute("z").toDouble());
216 if (std::abs(depth - depth_check) >
217 std::numeric_limits<double>::
218 epsilon()) // skip soil-layer if its thickness is zero
219 {
220 borehole->addSoilLayer(horizon.attribute("x").toDouble(),
221 horizon.attribute("y").toDouble(),
222 depth,
223 horizonName);
224 depth_check = depth;
225 }
226 else
227 {
228 WARN(
229 "XmlStnInterface::readStratigraphy(): Skipped layer '{:s}' "
230 "in borehole '{:s}' because of thickness 0.0.",
231 horizonName, borehole->getName());
232 }
233 }
234 else
235 {
236 WARN(
237 "XmlStnInterface::readStratigraphy(): Attribute missing in "
238 "<horizon> tag.");
239 }
240 horizon = horizon.nextSiblingElement();
241 }
242}
243
245{
246 if (export_name.empty())
247 {
248 ERR("XmlStnInterface::write(): No station list specified.");
249 return false;
250 }
251
252 out << "<?xml version=\"1.0\" encoding=\"ISO-8859-1\"?>\n"; // xml
253 // definition
254 out << "<?xml-stylesheet type=\"text/xsl\" "
255 "href=\"OpenGeoSysSTN.xsl\"?>\n\n"; // stylefile definition
256
257 QDomDocument doc("OGS-STN-DOM");
258 QDomElement root = doc.createElement("OpenGeoSysSTN");
259 root.setAttribute("xmlns:ogs", "http://www.opengeosys.org");
260 root.setAttribute("xmlns:xsi", "http://www.w3.org/2001/XMLSchema-instance");
261
262 const std::vector<GeoLib::Point*>* stations(
264 bool const is_borehole =
265 dynamic_cast<GeoLib::StationBorehole*>((*stations)[0]);
266
267 doc.appendChild(root);
268 QDomElement stationListTag = doc.createElement("stationlist");
269 root.appendChild(stationListTag);
270
271 QDomElement listNameTag = doc.createElement("name");
272 stationListTag.appendChild(listNameTag);
273 QDomText stationListNameText =
274 doc.createTextNode(QString::fromStdString(export_name));
275 listNameTag.appendChild(stationListNameText);
276 QString listType = is_borehole ? "boreholes" : "stations";
277 QDomElement stationsTag = doc.createElement(listType);
278 stationListTag.appendChild(stationsTag);
279
280 bool useStationValue(false);
281 double sValue =
282 static_cast<GeoLib::Station*>((*stations)[0])->getStationValue();
283 std::size_t nStations(stations->size());
284 for (std::size_t i = 1; i < nStations; i++)
285 {
286 if ((static_cast<GeoLib::Station*>((*stations)[i])->getStationValue() -
287 sValue) < std::numeric_limits<double>::epsilon())
288 {
289 useStationValue = true;
290 break;
291 }
292 }
293
294 for (std::size_t i = 0; i < nStations; i++)
295 {
296 QString stationType = is_borehole ? "borehole" : "station";
297 QDomElement stationTag = doc.createElement(stationType);
298 stationTag.setAttribute("id", QString::number(i));
299 stationTag.setAttribute(
300 "x", QString::number((*(*stations)[i])[0], 'f',
301 std::numeric_limits<double>::digits10));
302 stationTag.setAttribute(
303 "y", QString::number((*(*stations)[i])[1], 'f',
304 std::numeric_limits<double>::digits10));
305 stationTag.setAttribute(
306 "z", QString::number((*(*stations)[i])[2], 'f',
307 std::numeric_limits<double>::digits10));
308 stationsTag.appendChild(stationTag);
309
310 QDomElement stationNameTag = doc.createElement("name");
311 stationTag.appendChild(stationNameTag);
312 QDomText stationNameText = doc.createTextNode(QString::fromStdString(
313 static_cast<GeoLib::Station*>((*stations)[i])->getName()));
314 stationNameTag.appendChild(stationNameText);
315
316 if (useStationValue)
317 {
318 QDomElement stationValueTag = doc.createElement("value");
319 stationTag.appendChild(stationValueTag);
320 QDomText stationValueText = doc.createTextNode(
321 QString::number(static_cast<GeoLib::Station*>((*stations)[i])
322 ->getStationValue()));
323 stationValueTag.appendChild(stationValueText);
324 }
325
326 if (is_borehole)
327 {
329 doc, stationTag,
330 static_cast<GeoLib::StationBorehole*>((*stations)[i]));
331 }
332 }
333
334 std::string xml = doc.toString().toStdString();
335 out << xml;
336 return true;
337}
338
340 QDomElement& boreholeTag,
341 GeoLib::StationBorehole* borehole) const
342{
343 QDomElement stationDepthTag = doc.createElement("bdepth");
344 boreholeTag.appendChild(stationDepthTag);
345 QDomText stationDepthText =
346 doc.createTextNode(QString::number(borehole->getDepth(), 'f'));
347 stationDepthTag.appendChild(stationDepthText);
348 if (std::abs(borehole->getDate()) > 0)
349 {
350 QDomElement stationDateTag = doc.createElement("bdate");
351 boreholeTag.appendChild(stationDateTag);
352 QDomText stationDateText = doc.createTextNode(
353 QString::fromStdString(BaseLib::date2string(borehole->getDate())));
354 stationDateTag.appendChild(stationDateText);
355 }
356
357 std::vector<GeoLib::Point*> profile = borehole->getProfile();
358 std::vector<std::string> soilNames = borehole->getSoilNames();
359 std::size_t nHorizons(profile.size());
360
361 if (nHorizons > 1)
362 {
363 QDomElement stratTag = doc.createElement("strat");
364 boreholeTag.appendChild(stratTag);
365
366 for (std::size_t j = 1; j < nHorizons;
367 j++)
369 {
370 QDomElement horizonTag = doc.createElement("horizon");
371 horizonTag.setAttribute("id", QString::number(j));
372 horizonTag.setAttribute("x",
373 QString::number((*profile[j])[0], 'f'));
374 horizonTag.setAttribute("y",
375 QString::number((*profile[j])[1], 'f'));
376 horizonTag.setAttribute("z",
377 QString::number((*profile[j])[2], 'f'));
378 stratTag.appendChild(horizonTag);
379 QDomElement horizonNameTag = doc.createElement("name");
380 horizonTag.appendChild(horizonNameTag);
381 QDomText horizonNameText =
382 doc.createTextNode(QString::fromStdString(soilNames[j]));
383 horizonNameTag.appendChild(horizonNameText);
384 }
385 }
386}
387} // namespace IO
388} // namespace GeoLib
Definition of date helper functions.
Filename manipulation routines.
Definition of the GEOObjects class.
void ERR(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:45
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
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:47
QByteArray const & getContent() const
Container class for geometric objects.
Definition GEOObjects.h:57
void addStationVec(std::vector< Point * > &&stations, std::string &name)
Adds a vector of stations with the given name and colour to GEOObjects.
const std::vector< GeoLib::Point * > * getStationVec(const std::string &name) const
Returns the station vector with the given name.
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....
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.
XmlStnInterface(GeoLib::GEOObjects &geo_objs)
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.
const std::vector< std::string > & getSoilNames() const
double getDate() const
Returns the date entry for the borehole.
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:83
std::string const & getName() const
Returns the name of the station.
Definition Station.h:60
void setStationValue(double station_value)
Definition Station.h:77
std::string date2string(double ddate)
Definition DateTools.cpp:65
std::string copyPathToFileName(const std::string &file_name, const std::string &source)