sgsPy
structurally guided sampling
Loading...
Searching...
No Matches
systematic.h
Go to the documentation of this file.
1/******************************************************************************
2 *
3 * Project: sgs
4 * Purpose: C++ implementation of systematic sampling
5 * Author: Joseph Meyer
6 * Date: July, 2025
7 *
8 ******************************************************************************/
9
14
15#include <iostream>
16#include <random>
17
18#include "utils/access.h"
19#include "utils/existing.h"
20#include "utils/helper.h"
21#include "utils/raster.h"
22#include "utils/vector.h"
23
24namespace sgs {
25namespace systematic {
26
46OGRGeometry *getAccessPolygon(vector::GDALVectorWrapper *p_access, std::string layerName, double buffInner, double buffOuter) {
47 //step 1: create multipolygon buffers
48 OGRMultiPolygon *buffInnerPolygons = new OGRMultiPolygon;
49 OGRMultiPolygon *buffOuterPolygons = new OGRMultiPolygon;
50 OGRGeometry *p_polygonMask;
51
52 //step 2: add geometries to access polygon buffers
53 for (const auto& p_feature : *p_access->getLayer(layerName.c_str())) {
54 OGRGeometry *p_geometry = p_feature->GetGeometryRef();
55 OGRwkbGeometryType type = wkbFlatten(p_geometry->getGeometryType());
56
57 switch (type) {
58 case OGRwkbGeometryType::wkbLineString: {
59 OGRGeometry *p_outer = p_geometry->Buffer(buffOuter);
60 buffOuterPolygons->addGeometry(p_outer);
61 OGRGeometryFactory::destroyGeometry(p_outer);
62
63 if (buffInner != 0) {
64 OGRGeometry *p_inner = p_geometry->Buffer(buffInner);
65 buffInnerPolygons->addGeometry(p_inner);
66 OGRGeometryFactory::destroyGeometry(p_inner);
67 }
68 break;
69 }
70 case OGRwkbGeometryType::wkbMultiLineString: {
71 for (const auto& p_lineString : *p_geometry->toMultiLineString()) {
72 OGRGeometry *p_outer = p_geometry->Buffer(buffOuter);
73 buffOuterPolygons->addGeometry(p_outer);
74 OGRGeometryFactory::destroyGeometry(p_outer);
75 if (buffInner != 0) {
76 OGRGeometry *p_inner = p_geometry->Buffer(buffInner);
77 buffInnerPolygons->addGeometry(p_inner);
78 OGRGeometryFactory::destroyGeometry(p_inner);
79 }
80 }
81 break;
82 }
83 default: {
84 throw std::runtime_error("access polygon geometry type must be LineString or MultiLineString");
85 }
86 }
87 }
88
89 //step 3: generate the polygon mask and free no longer used memory
90 if (buffInner == 0) {
91 p_polygonMask = buffOuterPolygons->UnionCascaded();
92 }
93 else {
94 OGRGeometry *buffOuterUnion = buffOuterPolygons->UnionCascaded();
95 OGRGeometry *buffInnerUnion = buffInnerPolygons->UnionCascaded();
96 p_polygonMask = buffOuterUnion->Difference(buffInnerUnion);
97 OGRGeometryFactory::destroyGeometry(buffOuterUnion);
98 OGRGeometryFactory::destroyGeometry(buffInnerUnion);
99 }
100 delete buffOuterPolygons;
101 delete buffInnerPolygons;
102
103 return p_polygonMask;
104}
105
119inline bool
120checkExtent(double x, double y, double xMin, double xMax, double yMin, double yMax) {
121 return (x >= xMin && x <= xMax && y >= yMin && y <= yMax);
122}
123
133inline bool
134checkAccess(OGRPoint *p_point, OGRGeometry *p_geometry) {
135 return !p_geometry || p_point->Within(p_geometry);
136}
137
148inline bool
150 return !existing.used || !existing.containsCoordinates(x, y);
151}
152
173inline bool
174checkNotNan(raster::GDALRasterWrapper *p_raster, double *IGT, double xCoord, double yCoord, bool force) {
175 if (force) {
176 int x = static_cast<int>(IGT[0] + xCoord * IGT[1] + yCoord * IGT[2]);
177 int y = static_cast<int>(IGT[3] + xCoord * IGT[4] + yCoord * IGT[5]);
178
179 for (int i = 0; i < p_raster->getBandCount(); i++) {
180 GDALRasterBand *p_band = p_raster->getRasterBand(i);
181
182 double val;
183 OGRErr err = p_band->RasterIO(GF_Read, x, y, 1, 1, &val, 1, 1, GDT_Float64, 0, 0);
184 if (err) {
185 throw std::runtime_error("could not read value from raster.");
186 }
187
188 if (val == p_band->GetNoDataValue() || std::isnan(val)) {
189 return false;
190 }
191 }
192 }
193
194 return true;
195}
196
242std::tuple<
243 vector::GDALVectorWrapper *, //GDALDataset containing sample points
244 std::vector<std::vector<double>>, //array of samples to plot
245 std::vector<std::vector<std::vector<double>>> //array of grid to plot
246>
249 double cellSize,
250 std::string shape,
251 std::string location,
252 vector::GDALVectorWrapper *p_existing,
254 std::string layerName,
255 double buffInner,
256 double buffOuter,
257 bool force,
258 bool plot,
259 std::string filename)
260{
261 GDALAllRegister();
262
263 double *GT; //geotransform
264 double IGT[6]; //inverse geotransform
265 GT = p_raster->getGeotransform();
266 GDALInvGeoTransform(GT, IGT);
267
268 //determine raster extent
269 double xMin, xMax, yMin, yMax;
270 xMin = p_raster->getXMin();
271 xMax = p_raster->getXMax();
272 yMin = p_raster->getYMin();
273 yMax = p_raster->getYMax();
274
275 //generate random number generator
276 double xDiff = xMax - xMin;
277 double yDiff = yMax - yMin;
278 double rngMax = yDiff * xDiff;
279 std::mt19937::result_type seed = time(nullptr);
280 auto rng = std::bind(
281 std::uniform_real_distribution<double>(0, rngMax),
282 std::mt19937(seed)
283 );
284
285 //determine random rotation angle within extent polygon
286 double rotation = rng() / (rngMax / 180);
287
288 //determine grid creation function
289 std::string gridFunction;
290 if (shape == "square") {
291 gridFunction = "ST_SquareGrid";
292 }
293 else { //shape == "hexagon"
294 gridFunction = "ST_HexagonalGrid";
295 }
296
297 //create sql query using extent polygon and grid function
298 std::string extentPolygon = "'POLYGON (( "
299 + std::to_string(xMin) + " " + std::to_string(yMin) + ", "
300 + std::to_string(xMin) + " " + std::to_string(yMax) + ", "
301 + std::to_string(xMax) + " " + std::to_string(yMax) + ", "
302 + std::to_string(xMax) + " " + std::to_string(yMin) + ", "
303 + std::to_string(xMin) + " " + std::to_string(yMin) + " ))'";
304
305 std::string queryString = "SELECT RotateCoords("
306 + gridFunction
307 + "(RotateCoords("
308 + "ST_GeomFromText("
309 + extentPolygon
310 + "), " + std::to_string(-rotation) + "), "
311 + std::to_string(cellSize) + "), "
312 + std::to_string(rotation) + ")";
313
314 //query to create grid
315 OGRLayer *p_grid = p_raster->getDataset()->ExecuteSQL(queryString.c_str(), nullptr, "SQLITE");
316
317 //create output dataset before anything that might take a long time (in case creation of dataset fails)
318 GDALDriver *p_driver = GetGDALDriverManager()->GetDriverByName("MEM");
319 if(!p_driver) {
320 throw std::runtime_error("unable to create output sample dataset driver.");
321 }
322 GDALDataset *p_sampleDataset = p_driver->Create("", 0, 0, 0, GDT_Unknown, nullptr);
323 if (!p_sampleDataset) {
324 throw std::runtime_error("unable to create output dataset with driver.");
325 }
326
327 vector::GDALVectorWrapper *p_wrapper = new vector::GDALVectorWrapper(p_sampleDataset, std::string(p_raster->getDataset()->GetProjectionRef()));
328 OGRLayer *p_sampleLayer = p_sampleDataset->CreateLayer("samples", p_wrapper->getSRS(), wkbPoint, nullptr);
329 if (!p_sampleLayer) {
330 throw std::runtime_error("unable to create output dataset layer.");
331 }
332
333 //get access polygon if access is given
334 OGRGeometry *access = nullptr;
335 if (p_access) {
336 //occasionally, samples end up being placed just outside the accessible area,
337 //typically when the buffer sizes are multiples of the pixel size.
338 //
339 //Adjusting the buffers minorly fixes this problem.
340 double pixelSize = std::min(p_raster->getPixelHeight(), p_raster->getPixelWidth());
341 buffOuter = buffOuter - pixelSize / 50;
342 buffInner = buffInner == 0 ? 0 : buffInner + pixelSize / 50;
343
344 access = getAccessPolygon(p_access, layerName, buffInner, buffOuter);
345 }
346
347 //coordinate representations the samples as vectors only if PLOT is true, returned to the (Python) caller
348 std::vector<double> xCoords, yCoords;
349
350 //create existing struct
351 existing::Existing existing(p_existing, p_raster, GT, p_raster->getWidth(), p_sampleLayer, plot, xCoords, yCoords);
352 helper::Field fieldExistingFalse("existing", 0);
353
354 //grid represents the grid used to create the sampels only if PLOT is true, and is returned to the (Python) caller
355 std::vector<std::vector<std::vector<double>>> grid;
356
357 //iterate through the polygons in the grid and populate the ruturn data depending on user inputs
358 for (const auto& p_feature : *p_grid) {
359 OGRGeometry *p_geometry = p_feature->GetGeometryRef();
360 for (const auto& p_polygon : *p_geometry->toMultiPolygon()) {
361 //generate sample depending on 'location' parameter
362 OGRPoint point;
363 OGRPoint secondPoint;
364 if (location == "centers") {
365 p_polygon->Centroid(&point);
366
367 double x = point.getX();
368 double y = point.getY();
369 if (checkExtent(x, y, xMin, xMax, yMin, yMax) &&
370 checkAccess(&point, access) &&
371 checkExisting(x, y, existing) &&
372 checkNotNan(p_raster, IGT, x, y, force))
373 {
374 existing.used ?
375 helper::addPoint(&point, p_sampleLayer, &fieldExistingFalse) :
376 helper::addPoint(&point, p_sampleLayer);
377
378 if (plot) {
379 xCoords.push_back(x);
380 yCoords.push_back(y);
381 }
382 }
383 }
384 else if (location == "corners") {
385 (*p_polygon->begin())->getPoint(0, &point);
386 (*p_polygon->begin())->getPoint(1, &secondPoint);
387
388 double x = point.getX();
389 double y = point.getY();
390 if (checkExtent(x, y, xMin, xMax, yMin, yMax) &&
391 checkAccess(&point, access) &&
392 checkExisting(x, y, existing) &&
393 checkNotNan(p_raster, IGT, x, y, force))
394 {
395 existing.used ?
396 helper::addPoint(&point, p_sampleLayer, &fieldExistingFalse) :
397 helper::addPoint(&point, p_sampleLayer);
398
399 if (plot) {
400 xCoords.push_back(x);
401 yCoords.push_back(y);
402 }
403 }
404
405 x = secondPoint.getX();
406 y = secondPoint.getY();
407 if (checkExtent(x, y, xMin, xMax, yMin, yMax) &&
408 checkAccess(&secondPoint, access) &&
409 checkExisting(x, y, existing) &&
410 checkNotNan(p_raster, IGT, x, y, force))
411 {
412 existing.used ?
413 helper::addPoint(&point, p_sampleLayer, &fieldExistingFalse) :
414 helper::addPoint(&secondPoint, p_sampleLayer);
415
416 if (plot) {
417 xCoords.push_back(x);
418 yCoords.push_back(y);
419 }
420 }
421
422 }
423 else { //location == "random"
424 OGREnvelope envelope;
425 p_polygon->getEnvelope(&envelope);
426 double xMinEnv = envelope.MinX;
427 double xMaxEnv = envelope.MaxX;
428 double xDiffEnv = xMaxEnv - xMinEnv;
429 double yMinEnv = envelope.MinY;
430 double yMaxEnv = envelope.MaxY;
431 double yDiffEnv = yMaxEnv - yMinEnv;
432
433 int tries = 0;
434 bool found = false;
435
436 while (tries < 10 && !found) {
437 point.setX(xMinEnv + rng() / (rngMax / xDiffEnv));
438 point.setY(yMinEnv + rng() / (rngMax / yDiffEnv));
439 while (!p_polygon->Contains(&point)) {
440 point.setX(xMinEnv + rng() / (rngMax / xDiffEnv));
441 point.setY(yMinEnv + rng() / (rngMax / yDiffEnv));
442 }
443
444 double x = point.getX();
445 double y = point.getY();
446 if (checkExtent(x, y, xMin, xMax, yMin, yMax) &&
447 checkAccess(&point, access) &&
448 checkExisting(x, y, existing) &&
449 checkNotNan(p_raster, IGT, x, y, force))
450 {
451 found = true;
452 existing.used ?
453 helper::addPoint(&point, p_sampleLayer, &fieldExistingFalse) :
454 helper::addPoint(&point, p_sampleLayer);
455
456 if (plot) {
457 xCoords.push_back(x);
458 yCoords.push_back(y);
459 }
460 }
461
462 tries++;
463 }
464 }
465
466 //set grid vector to be plot
467 if (plot) {
468 grid.push_back({}); //add new polygon to grid
469 grid.back().push_back({}); //add new x vector to polygon grid
470 grid.back().push_back({}); //add new y vector to polygon grid
471 for (const auto& p_linearRing : *p_polygon) {
472 for (const auto& p_point : *p_linearRing) {
473 grid.back()[0].push_back(p_point.getX());
474 grid.back()[1].push_back(p_point.getY());
475 }
476 }
477 }
478 }
479 }
480
481 if (filename != "") {
482 try {
483 p_wrapper->write(filename);
484 }
485 catch (const std::exception& e) {
486 std::cout << "Exception thrown trying to write file: " << e.what() << std::endl;
487 }
488 }
489
490 //clean up memory
491 p_raster->getDataset()->ReleaseResultSet(p_grid);
492 if (p_access) {
493 OGRGeometryFactory::destroyGeometry(access);
494 }
495
496 return {p_wrapper, {xCoords, yCoords}, grid};
497}
498
499} //namespace systematic
500} //namespace sgs
Definition raster.h:57
double getYMax()
Definition raster.h:525
GDALDataset * getDataset()
Definition raster.h:429
int getWidth()
Definition raster.h:467
double * getGeotransform()
Definition raster.h:590
int getBandCount()
Definition raster.h:485
double getXMax()
Definition raster.h:495
GDALRasterBand * getRasterBand(int band)
Definition raster.h:696
double getXMin()
Definition raster.h:510
double getPixelWidth()
Definition raster.h:555
double getYMin()
Definition raster.h:540
double getPixelHeight()
Definition raster.h:565
Definition vector.h:46
OGRSpatialReference * getSRS(void)
Definition vector.h:457
OGRLayer * getLayer(std::string layerName)
Definition vector.h:246
void write(std::string filename)
Definition vector.h:397
void addPoint(OGRPoint *p_point, OGRLayer *p_layer)
Definition helper.h:700
bool checkExisting(double x, double y, existing::Existing &existing)
Definition systematic.h:149
std::tuple< vector::GDALVectorWrapper *, std::vector< std::vector< double > >, std::vector< std::vector< std::vector< double > > > > systematic(raster::GDALRasterWrapper *p_raster, double cellSize, std::string shape, std::string location, vector::GDALVectorWrapper *p_existing, vector::GDALVectorWrapper *p_access, std::string layerName, double buffInner, double buffOuter, bool force, bool plot, std::string filename)
Definition systematic.h:247
bool checkNotNan(raster::GDALRasterWrapper *p_raster, double *IGT, double xCoord, double yCoord, bool force)
Definition systematic.h:174
bool checkExtent(double x, double y, double xMin, double xMax, double yMin, double yMax)
Definition systematic.h:120
OGRGeometry * getAccessPolygon(vector::GDALVectorWrapper *p_access, std::string layerName, double buffInner, double buffOuter)
Definition systematic.h:46
bool checkAccess(OGRPoint *p_point, OGRGeometry *p_geometry)
Definition systematic.h:134
Definition access.h:23
Definition existing.h:27
Definition pca.h:23
Definition existing.h:38
Definition helper.h:46