Coverage Control Library
Loading...
Searching...
No Matches
world_idf.cpp
Go to the documentation of this file.
1/*
2 * This file is part of the CoverageControl library
3 *
4 * Author: Saurav Agarwal
5 * Contact: sauravag@seas.upenn.edu, agr.saurav1@gmail.com
6 * Repository: https://github.com/KumarRobotics/CoverageControl
7 *
8 * Copyright (c) 2024, Saurav Agarwal
9 *
10 * The CoverageControl library is free software: you can redistribute it and/or
11 * modify it under the terms of the GNU General Public License as published by
12 * the Free Software Foundation, either version 3 of the License, or (at your
13 * option) any later version.
14 *
15 * The CoverageControl library is distributed in the hope that it will be
16 * useful, but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General
18 * Public License for more details.
19 *
20 * You should have received a copy of the GNU General Public License along with
21 * CoverageControl library. If not, see <https://www.gnu.org/licenses/>.
22 */
23
33#ifdef COVERAGECONTROL_WITH_CUDA
35#endif
36
37namespace CoverageControl {
38
42 /* std::cout << "Generating map using CPU" << std::endl; */
43 std::vector<Polygon_2> cgal_polygons;
44 cgal_polygons.reserve(polygon_features_.size());
45 for (auto const &poly : polygon_features_) {
46 Polygon_2 cgal_poly;
47 for (auto const &pt : poly.poly) {
48 cgal_poly.push_back(CGAL_Point2(pt.x(), pt.y()));
49 }
50 cgal_polygons.push_back(cgal_poly);
51 }
52 float max_importance = 0;
53 float res = static_cast<float>(params_.pResolution);
54 for (int i = 0; i < params_.pWorldMapSize; ++i) { // Row (x index)
55 float x1 = res * i; // Left x-coordinate of pixel
56 float x2 = x1 + res; // Right x-coordinate of pixel
57 for (int j = 0; j < params_.pWorldMapSize; ++j) { // Column (y index)
58 float y1 = res * j; // Lower y-coordinate of pixel
59 float y2 = y1 + res; // Upper y-coordinate of pixel
60 Point2f pt1(x1, y1);
61 Point2f pt2(x2, y2);
62 float importance = ComputeImportanceBND<Point2f, float>(pt1, pt2);
63 /* importance += ComputeImportancePoly(pt1, pt2); */
64 CGAL_Point2 mid((x1 + x2) / 2, (y1 + y2) / 2);
65 float importance_poly = 0;
66 for (size_t k = 0; k < cgal_polygons.size(); ++k) {
67 auto &poly = cgal_polygons[k];
68 if (poly.bounded_side(mid) == CGAL::ON_BOUNDED_SIDE) {
69 importance_poly = std::max(importance_poly, polygon_features_[k].imp);
70 }
71 }
72 importance += importance_poly;
73 /* auto importance = ComputeImportanceBND(pt1, pt2); */
74 /* if (std::abs(importance) < kEps) { */
75 /* importance = 0; */
76 /* } */
77 world_map_(i, j) = importance;
78 if (importance > max_importance) {
79 max_importance = importance;
80 }
81 }
82 }
83
84 if (max_importance < kEps) {
85 normalization_factor_ = static_cast<float>(params_.pNorm);
86 } else {
87 normalization_factor_ = static_cast<float>(params_.pNorm) / max_importance;
88 }
89
90 if (not(normalization_factor_ > 1e-5)) { // Parity with CUDA code
91 return;
92 }
93
94 // Normalize the world map
95#pragma omp parallel for
96 for (int i = 0; i < params_.pWorldMapSize; ++i) {
97 for (int j = 0; j < params_.pWorldMapSize; ++j) {
98 world_map_(i, j) *= normalization_factor_;
99 }
100 }
101}
102
103#ifdef COVERAGECONTROL_WITH_CUDA
104void WorldIDF::GenerateMapCuda(float const resolution, float const truncation,
105 int const map_size) {
106 int num_dists = normal_distributions_.size();
107 /* std::cout << "num_dists: " << num_dists << std::endl; */
108
109 /* BND_Cuda *host_dists = (BND_Cuda*) malloc(num_dists * sizeof(BND_Cuda));
110 */
111 BND_Cuda *host_dists = new BND_Cuda[num_dists];
112
113 for (int i = 0; i < num_dists; ++i) {
114 Point2 mean = normal_distributions_[i].GetMean();
115 host_dists[i].mean_x = static_cast<float>(mean.x());
116 host_dists[i].mean_y = static_cast<float>(mean.y());
117 Point2 sigma = normal_distributions_[i].GetSigma();
118 host_dists[i].sigma_x = static_cast<float>(sigma.x());
119 host_dists[i].sigma_y = static_cast<float>(sigma.y());
120 host_dists[i].rho = static_cast<float>(normal_distributions_[i].GetRho());
121 host_dists[i].scale =
122 static_cast<float>(normal_distributions_[i].GetScale());
123 }
124
125 Polygons_Cuda_Host host_polygons;
126 for (auto const &poly_feature : polygon_features_) {
127 std::vector<PointVector> partition_polys;
128 PolygonYMonotonePartition(poly_feature.poly, partition_polys);
129 for (auto const &poly : partition_polys) {
130 Bounds bounds;
131 for (Point2 const &pt : poly) {
132 float x = static_cast<float>(pt.x());
133 float y = static_cast<float>(pt.y());
134 host_polygons.x.push_back(x);
135 host_polygons.y.push_back(y);
136 if (x < bounds.xmin) {
137 bounds.xmin = x;
138 }
139 if (y < bounds.ymin) {
140 bounds.ymin = y;
141 }
142 if (x > bounds.xmax) {
143 bounds.xmax = x;
144 }
145 if (y > bounds.ymax) {
146 bounds.ymax = y;
147 }
148 }
149 host_polygons.imp.push_back(poly_feature.imp);
150 host_polygons.sz.push_back(poly.size());
151 host_polygons.bounds.push_back(bounds);
152 }
153 }
154
155 host_polygons.num_pts = host_polygons.x.size();
156 host_polygons.num_polygons = host_polygons.imp.size();
157
158 float f_norm = 0;
159 generate_world_map_cuda(host_dists, host_polygons, num_dists, map_size,
160 resolution, truncation, params_.pNorm,
161 world_map_.data(), f_norm);
162 normalization_factor_ = static_cast<double>(f_norm);
163}
164#endif
165
166int WorldIDF::WriteDistributions(std::string const &file_name) const {
167 std::ofstream file(file_name);
168 if (!file.is_open()) {
169 std::cerr << "Could not open file: " << file_name << std::endl;
170 return -1;
171 }
172 file << std::setprecision(kMaxPrecision);
173 for (auto const &dist : normal_distributions_) {
174 Point2 sigma = dist.GetSigma();
175 if (sigma.x() == sigma.y()) {
176 file << "CircularBND" << std::endl;
177 file << dist.GetMean().x() << " " << dist.GetMean().y() << " "
178 << sigma.x() << " " << dist.GetScale() << std::endl;
179 } else {
180 file << "BND" << std::endl;
181 file << dist.GetMean().x() << " " << dist.GetMean().y() << " "
182 << dist.GetSigma().x() << " " << dist.GetSigma().y() << " "
183 << dist.GetRho() << " " << dist.GetScale() << std::endl;
184 }
185 }
186 for (auto const &poly : polygon_features_) {
187 file << "Uniform" << std::endl;
188 file << poly.poly.size() << " " << poly.imp << std::endl;
189 for (auto const &pt : poly.poly) {
190 file << pt.x() << " " << pt.y() << std::endl;
191 }
192 }
193 file.close();
194 return 0;
195}
196} // namespace CoverageControl
int WriteDistributions(std::string const &file_name) const
Contains the configuration for the CGAL library.
CGAL::Polygon_2< K > Polygon_2
Definition config.h:57
K::Point_2 CGAL_Point2
Definition config.h:50
Constants for the CoverageControl library.
Contains structs and functions to interface with the CUDA code to generate the world map.
double const kEps
Definition constants.h:48
constexpr auto kMaxPrecision
Definition constants.h:56
Eigen::Vector2d Point2
Definition typedefs.h:44
Eigen::Vector2f Point2f
Definition typedefs.h:45
Namespace for the CoverageControl library.
void PolygonYMonotonePartition(PointVector const &polygon, std::vector< PointVector > &y_monotone_polygons)
Partition a polygon into y-monotone polygons.
void generate_world_map_cuda(BND_Cuda *, Polygons_Cuda_Host const &, int const, int const, float const, float const, float const, float *importance_vec, float &)
Function to generate the world map on the device.
Provides utilities for polygon manipulation using CGAL.
Structure to store the parameters of the Bivariate Normal Distribution.
Contains the class for Importance Density Function (IDF) for the world.