OGS
Raster2PointCloud.cpp File Reference

Detailed Description

Definition in file Raster2PointCloud.cpp.

#include <tclap/CmdLine.h>
#include <vtkPolyDataAlgorithm.h>
#include <vtkSmartPointer.h>
#include <vtkXMLPolyDataWriter.h>
#include <algorithm>
#include <memory>
#include <string>
#include "Applications/DataExplorer/VtkVis/VtkGeoImageSource.h"
#include "Applications/DataExplorer/VtkVis/VtkImageDataToPointCloudFilter.h"
#include "BaseLib/FileTools.h"
#include "BaseLib/Logging.h"
#include "BaseLib/MPI.h"
#include "BaseLib/TCLAPArguments.h"
#include "GeoLib/IO/AsciiRasterInterface.h"
#include "InfoLib/GitInfo.h"
Include dependency graph for Raster2PointCloud.cpp:

Go to the source code of this file.

Functions

int main (int argc, char *argv[])
 

Function Documentation

◆ main()

int main ( int argc,
char * argv[] )

Definition at line 28 of file Raster2PointCloud.cpp.

29{
30 TCLAP::CmdLine cmd(
31 "Batch conversion of raster files into VTK point cloud data.\n"
32 "A series of raster files is converted into point clouds by creating a "
33 "number of random points based on the intensity value in each pixel "
34 "(the higher the value, the more points. The [x,y]-extend of these "
35 "random points is limited by the extent of the pixel they represent, "
36 "the elevation is limited by the max-elevation parameter. Likewise, "
37 "the "
38 "maximum number of points per pixel is limited by the max-points "
39 "parameter. The increase of the number of points in relation to the "
40 "pixel value is indicated by gamma, where 1 is a linear increase and "
41 "values smaller or larger than 1 indicate a logarithmic or exponential "
42 "increase, respectively.\n\n"
43 "OpenGeoSys-6 software, version " +
45 ".\n"
46 "Copyright (c) 2012-2025, OpenGeoSys Community "
47 "(http://www.opengeosys.org)",
49 TCLAP::ValueArg<double> gamma_arg(
50 "g", "gamma",
51 "Exponental increase coefficient. A value of '1' indicates a linear "
52 "increase, smaller or larger values indicate a logarithmic or "
53 "exponential increase, respectively (default: 1).",
54 false, 0.0, "exponential increase");
55 cmd.add(gamma_arg);
56 TCLAP::ValueArg<std::size_t> max_points_arg(
57 "p", "max-points",
58 "Maximum number of points per pixel (default: 1000), "
59 "(min = 0).",
60 false, 0, "maximum points per pixel");
61 cmd.add(max_points_arg);
62 TCLAP::ValueArg<double> max_height_arg(
63 "e", "max-elevation",
64 "Maximum elevation of randomly created points "
65 "(default: 5000).",
66 false, 0.0, "maximum point elevation");
67 cmd.add(max_height_arg);
68 TCLAP::ValueArg<std::string> output_arg(
69 "o", "output",
70 "Output (.vtp). Path to output polydata files. "
71 "Specifying 'file.vtp' will "
72 "make the programme write a series of files called 'file0.vtp', "
73 "'file1.vtp', etc.",
74 true, "", "OUTPUT_FILE");
75 cmd.add(output_arg);
76 TCLAP::ValueArg<std::string> input_arg(
77 "i", "input",
78 "Input (.asc). Path to input raster files. "
79 "Specifying 'file.asc' will make "
80 "the programme look for a series of files called 'file0.asc', "
81 "'file1.asc', etc.",
82 true, "", "INPUT_FILE");
83 cmd.add(input_arg);
84 auto log_level_arg = BaseLib::makeLogLevelArg();
85 cmd.add(log_level_arg);
86 cmd.parse(argc, argv);
87
88 BaseLib::MPI::Setup mpi_setup(argc, argv);
89 BaseLib::initOGSLogger(log_level_arg.getValue());
90
91 std::string const input_name = input_arg.getValue().c_str();
92 std::string const output_name = output_arg.getValue().c_str();
93
94 double const max_height =
95 (max_height_arg.isSet()) ? max_height_arg.getValue() : 5000;
96 std::size_t const max_points =
97 (max_points_arg.isSet()) ? max_points_arg.getValue() : 1000;
98 double const gamma = (gamma_arg.isSet()) ? gamma_arg.getValue() : 1;
99
100 std::string const ibase_name = BaseLib::dropFileExtension(input_name);
101 std::string const ifile_ext = BaseLib::getFileExtension(input_name);
102
103 std::string const obase_name = BaseLib::dropFileExtension(output_name);
104 std::string const ofile_ext = BaseLib::getFileExtension(output_name);
105
106 // Rainevents using point cloud filter
107 std::cout << "Starting export... ";
108 double global_min = std::numeric_limits<double>::max();
109 double global_max = 0;
110
111 std::size_t n_rasters = 0;
112 for (std::size_t i = 0;; ++i)
113 {
114 std::string const file_name =
115 ibase_name + std::to_string(i) + ifile_ext;
116 std::unique_ptr<GeoLib::Raster> r(
118 if (r == nullptr)
119 {
120 n_rasters = i;
121 break;
122 }
123
124 for (auto it = r->begin(); it != r->end(); ++it)
125 {
126 if (*it != r->getHeader().no_data)
127 {
128 global_min = std::min(*it, global_min);
129 global_max = std::max(*it, global_max);
130 }
131 }
132 }
133
134 for (std::size_t i = 0; i < n_rasters; ++i)
135 {
136 std::string const file_name =
137 ibase_name + std::to_string(i) + ifile_ext;
138 vtkSmartPointer<VtkGeoImageSource> geo_image =
139 vtkSmartPointer<VtkGeoImageSource>::New();
140 geo_image->readImage(QString::fromStdString(file_name));
141 VtkImageDataToPointCloudFilter* const point_cloud_filter =
143 point_cloud_filter->SetInputConnection(geo_image->GetOutputPort());
144 geo_image->Update();
145 point_cloud_filter->SetMinValueRange(global_min);
146 point_cloud_filter->SetMaxValueRange(global_max);
147 point_cloud_filter->SetMinHeight(0);
148 point_cloud_filter->SetMaxHeight(max_height);
149 point_cloud_filter->SetIsLinear(false);
150 point_cloud_filter->SetMinNumberOfPointsPerCell(0);
151 point_cloud_filter->SetMaxNumberOfPointsPerCell(max_points);
152 point_cloud_filter->SetGamma(gamma);
153 point_cloud_filter->Update();
154 vtkSmartPointer<vtkPolyDataAlgorithm> const pd_alg =
155 dynamic_cast<vtkPolyDataAlgorithm*>(point_cloud_filter);
156 if (pd_alg)
157 {
158 std::string const output =
159 obase_name + std::to_string(i) + ofile_ext;
160 vtkSmartPointer<vtkXMLPolyDataWriter> pd_writer =
161 vtkSmartPointer<vtkXMLPolyDataWriter>::New();
162 pd_writer->SetInputData(pd_alg->GetOutputDataObject(0));
163 pd_writer->SetFileName(output.c_str());
164 pd_writer->Write();
165 }
166 }
167 std::cout << "done." << std::endl;
168 return EXIT_SUCCESS;
169}
static GeoLib::Raster * getRasterFromASCFile(std::string const &fname)
Reads an ArcGis ASC raster file.
static VtkImageDataToPointCloudFilter * New()
Create a new objects (required because of VTKs reference counting)
TCLAP::ValueArg< std::string > makeLogLevelArg()
void initOGSLogger(std::string const &log_level)
Definition Logging.cpp:64
std::string getFileExtension(const std::string &path)
std::string dropFileExtension(std::string const &filename)
GITINFOLIB_EXPORT const std::string ogs_version

References BaseLib::dropFileExtension(), BaseLib::getFileExtension(), FileIO::AsciiRasterInterface::getRasterFromASCFile(), BaseLib::initOGSLogger(), BaseLib::makeLogLevelArg(), VtkImageDataToPointCloudFilter::New(), and GitInfoLib::GitInfo::ogs_version.