1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365 366 367 368 369 370 371 372 373 374 375 376 377 378 379 380 381 382 383 384 385 386 387 388 389 390 391 392 393 394 395 396 397 398 399 400 401 402 403 404 405 406 407 408 409 410 411 412 413 414 415 416 417 418 419 420 421 422 423 424 425 426 427 428 429 430 431 432 433 434 435 436 437 438 439 440 441 442 443 444 445 446 447 448 449 450 451 452 453 454 455 456 457 458 459
|
#include "config.h"
#include "Scene_points_with_normal_item.h"
#include <CGAL/Three/CGAL_Lab_plugin_helper.h>
#include <CGAL/Three/CGAL_Lab_plugin_interface.h>
#include <CGAL/pca_estimate_normals.h>
#include <CGAL/jet_estimate_normals.h>
#include <CGAL/vcm_estimate_normals.h>
#include <CGAL/mst_orient_normals.h>
#include <CGAL/scanline_orient_normals.h>
#include <CGAL/Timer.h>
#include <CGAL/Memory_sizer.h>
#include <QObject>
#include <QAction>
#include <QMainWindow>
#include <QApplication>
#include <QtPlugin>
#include <QInputDialog>
#include <QMessageBox>
#include <QMultipleInputDialog.h>
#include "run_with_qprogressdialog.h"
#include "ui_Point_set_normal_estimation_plugin.h"
// Concurrency
typedef CGAL::Parallel_if_available_tag Concurrency_tag;
struct PCA_estimate_normals_functor
: public Functor_with_signal_callback
{
Point_set* points;
int neighborhood_size;
unsigned int sharpness_angle;
PCA_estimate_normals_functor (Point_set* points,
int neighborhood_size)
: points (points), neighborhood_size (neighborhood_size)
{ }
void operator()()
{
CGAL::pca_estimate_normals<Concurrency_tag>(points->all_or_selection_if_not_empty(),
neighborhood_size,
points->parameters().
callback (*(this->callback())));
}
};
struct Jet_estimate_normals_functor
: public Functor_with_signal_callback
{
Point_set* points;
int neighborhood_size;
unsigned int sharpness_angle;
Jet_estimate_normals_functor (Point_set* points,
int neighborhood_size)
: points (points), neighborhood_size (neighborhood_size)
{ }
void operator()()
{
CGAL::jet_estimate_normals<Concurrency_tag>(points->all_or_selection_if_not_empty(),
neighborhood_size,
points->parameters().
callback (*(this->callback())));
}
};
struct Vector_to_pmap
{
typedef boost::readable_property_map_tag category;
typedef Point_set::Index key_type;
typedef bool value_type;
typedef value_type reference;
std::vector<bool>* vec;
Vector_to_pmap (std::vector<bool>* vec = nullptr) : vec (vec) { }
friend inline
value_type get(const Vector_to_pmap& map, key_type i)
{
return (*map.vec)[i];
}
};
using namespace CGAL::Three;
class CGAL_Lab_point_set_normal_estimation_plugin :
public QObject,
public CGAL_Lab_plugin_helper
{
Q_OBJECT
Q_INTERFACES(CGAL::Three::CGAL_Lab_plugin_interface)
Q_PLUGIN_METADATA(IID "com.geometryfactory.CGALLab.PluginInterface/1.0")
QAction* actionNormalEstimation;
QAction* actionNormalOrientation;
QAction* actionNormalInversion;
public:
void init(QMainWindow* mainWindow, CGAL::Three::Scene_interface* scene_interface, Messages_interface*) {
scene = scene_interface;
mw = mainWindow;
actionNormalEstimation = new QAction(tr("Normal Estimation"), mainWindow);
actionNormalEstimation->setObjectName("actionNormalEstimation");
actionNormalEstimation->setProperty("subMenuName","Point Set Processing");
actionNormalOrientation = new QAction(tr("Normal Orientation"), mainWindow);
actionNormalOrientation->setObjectName("actionNormalOrientation");
actionNormalOrientation->setProperty("subMenuName","Point Set Processing");
actionNormalInversion = new QAction(tr("Inverse Normal Orientations"), mainWindow);
actionNormalInversion->setObjectName("actionNormalInversion");
actionNormalInversion->setProperty("subMenuName","Point Set Processing");
autoConnectActions();
}
QList<QAction*> actions() const {
return QList<QAction*>() << actionNormalEstimation << actionNormalOrientation << actionNormalInversion;
}
bool applicable(QAction* action) const {
Scene_points_with_normal_item* item = qobject_cast<Scene_points_with_normal_item*>(scene->item(scene->mainSelectionIndex()));
if (action==actionNormalEstimation)
return item;
else
return item && item->has_normals();
}
public Q_SLOTS:
void on_actionNormalEstimation_triggered();
void on_actionNormalOrientation_triggered();
void on_actionNormalInversion_triggered();
}; // end PS_demo_smoothing_plugin
class Point_set_demo_normal_estimation_dialog : public QDialog, private Ui::NormalEstimationDialog
{
Q_OBJECT
public:
Point_set_demo_normal_estimation_dialog(QWidget* /*parent*/ = nullptr)
{
setupUi(this);
m_offset_radius->setMinimum(0.01);
}
int pca_neighbors() const { return m_pca_neighbors->value(); }
int jet_neighbors() const { return m_jet_neighbors->value(); }
unsigned int convolution_neighbors() const { return m_convolution_neighbors->value(); }
double convolution_radius() const { return m_convolution_radius->value(); }
double offset_radius() const { return m_offset_radius->value(); }
unsigned int method () const
{
if (buttonPCA->isChecked ()) return 0;
if (buttonJet->isChecked ()) return 1;
if (buttonVCM->isChecked ()) return 2;
return -1;
}
bool use_convolution_radius () const
{
return buttonRadius->isChecked();
}
};
void CGAL_Lab_point_set_normal_estimation_plugin::on_actionNormalInversion_triggered()
{
const CGAL::Three::Scene_interface::Item_id index = scene->mainSelectionIndex();
Scene_points_with_normal_item* item =
qobject_cast<Scene_points_with_normal_item*>(scene->item(index));
if(item)
{
// Gets point set
Point_set* points = item->point_set();
if(points == nullptr)
return;
for(Point_set::iterator it = points->begin_or_selection_begin(); it != points->end(); ++it){
points->normal(*it) = -1 * points->normal(*it);
}
item->invalidateOpenGLBuffers();
scene->itemChanged(item);
}
}
void CGAL_Lab_point_set_normal_estimation_plugin::on_actionNormalEstimation_triggered()
{
const CGAL::Three::Scene_interface::Item_id index = scene->mainSelectionIndex();
Scene_points_with_normal_item* item =
qobject_cast<Scene_points_with_normal_item*>(scene->item(index));
if(item)
{
// Gets options
Point_set_demo_normal_estimation_dialog dialog;
if(!dialog.exec())
return;
// Gets point set
Point_set* points = item->point_set();
if(points == nullptr)
return;
if (!(points->has_normal_map()))
points->add_normal_map();
QApplication::setOverrideCursor(Qt::BusyCursor);
QApplication::processEvents();
//***************************************
// normal estimation
//***************************************
if (dialog.method() == 0) // PCA
{
CGAL::Timer task_timer; task_timer.start();
std::cerr << "Estimates normal direction by PCA (k=" << dialog.pca_neighbors() <<")...\n";
// Estimates normals direction.
PCA_estimate_normals_functor functor (points, dialog.pca_neighbors());
run_with_qprogressdialog (functor, "Estimating normals by PCA...", mw);
std::size_t memory = CGAL::Memory_sizer().virtual_size();
std::cerr << "Estimates normal direction: " << task_timer.time() << " seconds, "
<< (memory>>20) << " Mb allocated"
<< std::endl;
}
else if (dialog.method() == 1) // Jet
{
CGAL::Timer task_timer; task_timer.start();
std::cerr << "Estimates normal direction by Jet Fitting (k=" << dialog.jet_neighbors() <<")...\n";
// Estimates normals direction.
Jet_estimate_normals_functor functor (points, dialog.jet_neighbors());
run_with_qprogressdialog (functor, "Estimating normals by jet fitting...", mw);
std::size_t memory = CGAL::Memory_sizer().virtual_size();
std::cerr << "Estimates normal direction: " << task_timer.time() << " seconds, "
<< (memory>>20) << " Mb allocated"
<< std::endl;
}
else if (dialog.method() == 2) // VCM
{
CGAL::Timer task_timer; task_timer.start();
// Estimates normals direction.
if (dialog.use_convolution_radius())
{
std::cerr << "Estimates Normals Direction using VCM (R="
<< dialog.offset_radius() << " and r=" << dialog.convolution_radius() << ")...\n";
CGAL::vcm_estimate_normals
(points->all_or_selection_if_not_empty(),
dialog.offset_radius(),
dialog.convolution_radius(),
points->parameters());
}
else
{
std::cerr << "Estimates Normals Direction using VCM (R="
<< dialog.offset_radius() << " and k=" << dialog.convolution_neighbors() << ")...\n";
CGAL::vcm_estimate_normals
(points->all_or_selection_if_not_empty(),
dialog.offset_radius(),
dialog.convolution_neighbors(),
points->parameters());
}
std::size_t memory = CGAL::Memory_sizer().virtual_size();
std::cerr << "Estimates normal direction: " << task_timer.time() << " seconds, "
<< (memory>>20) << " Mb allocated"
<< std::endl;
}
item->resetMenu();
item->setRenderingMode(ShadedPoints);
// Updates scene
item->invalidateOpenGLBuffers();
scene->itemChanged(index);
item->itemChanged();
QApplication::restoreOverrideCursor();
}
}
void CGAL_Lab_point_set_normal_estimation_plugin::on_actionNormalOrientation_triggered()
{
const CGAL::Three::Scene_interface::Item_id index = scene->mainSelectionIndex();
Scene_points_with_normal_item* item =
qobject_cast<Scene_points_with_normal_item*>(scene->item(index));
if(item)
{
// Gets point set
Point_set* points = item->point_set();
if(points == nullptr)
return;
// Chose method
QMultipleInputDialog method ("Normal Orientation", mw);
QRadioButton* mst = method.add<QRadioButton> ("Orient by Minimum Spanning Tree");
mst->setChecked(true);
QRadioButton* scanline = method.add<QRadioButton> ("Orient 2.5D airborne acquired scanlines");
scanline->setChecked(false);
if (!method.exec())
return;
if (mst->isChecked())
{
// Gets options
QMultipleInputDialog dialog ("Normal Orientation", mw);
QSpinBox* neighborhood = dialog.add<QSpinBox> ("Neighborhood Size: ");
neighborhood->setRange (1, 10000000);
neighborhood->setValue (18);
QRadioButton* use_seed_points = nullptr;
QRadioButton* orient_selection = nullptr;
if (points->nb_selected_points() != 0)
{
use_seed_points = dialog.add<QRadioButton> ("Use selection as seed points and orient the unselected points");
use_seed_points->setChecked(true);
orient_selection = dialog.add<QRadioButton> ("Orient selection");
orient_selection->setChecked(false);
}
if(!dialog.exec())
return;
QApplication::setOverrideCursor(Qt::BusyCursor);
QApplication::processEvents();
// First point to delete
Point_set::iterator first_unoriented_point = points->end();
//***************************************
// normal orientation
//***************************************
CGAL::Timer task_timer; task_timer.start();
std::cerr << "Orient normals with a Minimum Spanning Tree (k=" << neighborhood->value() << ")...\n";
// Tries to orient normals
if (points->nb_selected_points() != 0 && use_seed_points->isChecked())
{
std::vector<bool> constrained_map (points->size(), false);
for (Point_set::iterator it = points->first_selected(); it != points->end(); ++ it)
constrained_map[*it] = true;
first_unoriented_point =
CGAL::mst_orient_normals(*points,
std::size_t(neighborhood->value()),
points->parameters().
point_is_constrained_map(Vector_to_pmap(&constrained_map)));
}
else
first_unoriented_point =
CGAL::mst_orient_normals(points->all_or_selection_if_not_empty(),
std::size_t(neighborhood->value()),
points->parameters());
std::size_t nb_unoriented_normals = std::distance(first_unoriented_point, points->end());
std::size_t memory = CGAL::Memory_sizer().virtual_size();
std::cerr << "Orient normals: " << nb_unoriented_normals << " point(s) with an unoriented normal are selected ("
<< task_timer.time() << " seconds, "
<< (memory>>20) << " Mb allocated)"
<< std::endl;
// Selects points with an unoriented normal
points->set_first_selected (first_unoriented_point);
// Updates scene
item->invalidateOpenGLBuffers();
scene->itemChanged(index);
QApplication::restoreOverrideCursor();
// Warns user
if (nb_unoriented_normals > 0)
{
QMessageBox::information(nullptr,
tr("Points with an unoriented normal"),
tr("%1 point(s) with an unoriented normal are selected.\nPlease orient them or remove them before running Poisson reconstruction.")
.arg(nb_unoriented_normals));
}
}
else // scanline method
{
QApplication::setOverrideCursor(Qt::BusyCursor);
QApplication::processEvents();
//***************************************
// normal orientation
//***************************************
CGAL::Timer task_timer; task_timer.start();
std::cerr << "Orient normals with along 2.5D scanlines..." << std::endl;
std::optional<Point_set::Property_map<float>> scan_angle = points->property_map<float>("scan_angle");
std::optional<Point_set::Property_map<unsigned char>> scan_direction_flag = points->property_map<unsigned char>("scan_direction_flag");
if (!scan_angle.has_value() && !scan_direction_flag.has_value())
{
std::cerr << " using no additional properties" << std::endl;
CGAL::scanline_orient_normals(points->all_or_selection_if_not_empty(),
points->parameters());
}
else if (!scan_angle.has_value() && scan_direction_flag.has_value())
{
std::cerr << " using scan direction flag" << std::endl;
CGAL::scanline_orient_normals(points->all_or_selection_if_not_empty(),
points->parameters().
scanline_id_map (scan_direction_flag.value()));
}
else if (scan_angle.has_value() && !scan_direction_flag.has_value())
{
std::cerr << " using scan angle" << std::endl;
CGAL::scanline_orient_normals(points->all_or_selection_if_not_empty(),
points->parameters().
scan_angle_map (scan_angle.value()));
}
else // if (angle_found && flag_found)
{
std::cerr << " using scan angle and direction flag" << std::endl;
CGAL::scanline_orient_normals(points->all_or_selection_if_not_empty(),
points->parameters().
scan_angle_map (scan_angle.value()).
scanline_id_map (scan_direction_flag.value()));
}
std::size_t memory = CGAL::Memory_sizer().virtual_size();
std::cerr << "Orient normals: "
<< task_timer.time() << " seconds, "
<< (memory>>20) << " Mb allocated)"
<< std::endl;
// Updates scene
item->invalidateOpenGLBuffers();
scene->itemChanged(index);
QApplication::restoreOverrideCursor();
}
}
}
#include "Point_set_normal_estimation_plugin.moc"
|