small_gicp
reduction_omp.hpp
Go to the documentation of this file.
1 // SPDX-FileCopyrightText: Copyright 2024 Kenji Koide
2 // SPDX-License-Identifier: MIT
3 #pragma once
4 
5 #include <Eigen/Core>
6 
7 namespace small_gicp {
8 
9 #ifndef _OPENMP
10 #warning "OpenMP is not available. Parallel reduction will be disabled."
11 inline int omp_get_thread_num() {
12  return 0;
13 }
14 #endif
15 
19 
20  template <typename TargetPointCloud, typename SourcePointCloud, typename TargetTree, typename CorrespondenceRejector, typename Factor>
21  std::tuple<Eigen::Matrix<double, 6, 6>, Eigen::Matrix<double, 6, 1>, double> linearize(
22  const TargetPointCloud& target,
23  const SourcePointCloud& source,
24  const TargetTree& target_tree,
25  const CorrespondenceRejector& rejector,
26  const Eigen::Isometry3d& T,
27  std::vector<Factor>& factors) const {
28  std::vector<Eigen::Matrix<double, 6, 6>> Hs(num_threads, Eigen::Matrix<double, 6, 6>::Zero());
29  std::vector<Eigen::Matrix<double, 6, 1>> bs(num_threads, Eigen::Matrix<double, 6, 1>::Zero());
30  std::vector<double> es(num_threads, 0.0);
31 
32 #pragma omp parallel for num_threads(num_threads) schedule(guided, 8)
33  for (std::int64_t i = 0; i < factors.size(); i++) {
34  Eigen::Matrix<double, 6, 6> H;
35  Eigen::Matrix<double, 6, 1> b;
36  double e;
37 
38  if (!factors[i].linearize(target, source, target_tree, T, i, rejector, &H, &b, &e)) {
39  continue;
40  }
41 
42  const int thread_id = omp_get_thread_num();
43  Hs[thread_id] += H;
44  bs[thread_id] += b;
45  es[thread_id] += e;
46  }
47 
48  for (int i = 1; i < num_threads; i++) {
49  Hs[0] += Hs[i];
50  bs[0] += bs[i];
51  es[0] += es[i];
52  }
53 
54  return {Hs[0], bs[0], es[0]};
55  }
56 
57  template <typename TargetPointCloud, typename SourcePointCloud, typename Factor>
58  double error(const TargetPointCloud& target, const SourcePointCloud& source, const Eigen::Isometry3d& T, std::vector<Factor>& factors) const {
59  double sum_e = 0.0;
60 
61 #pragma omp parallel for num_threads(num_threads) schedule(guided, 8) reduction(+ : sum_e)
62  for (std::int64_t i = 0; i < factors.size(); i++) {
63  sum_e += factors[i].error(target, source, T);
64  }
65  return sum_e;
66  }
67 
69 };
70 
71 } // namespace small_gicp
Definition: flat_container.hpp:12
int omp_get_thread_num()
Definition: reduction_omp.hpp:11
Parallel reduction with OpenMP backend.
Definition: reduction_omp.hpp:17
double error(const TargetPointCloud &target, const SourcePointCloud &source, const Eigen::Isometry3d &T, std::vector< Factor > &factors) const
Definition: reduction_omp.hpp:58
int num_threads
Number of threads.
Definition: reduction_omp.hpp:68
ParallelReductionOMP()
Definition: reduction_omp.hpp:18
std::tuple< Eigen::Matrix< double, 6, 6 >, Eigen::Matrix< double, 6, 1 >, double > linearize(const TargetPointCloud &target, const SourcePointCloud &source, const TargetTree &target_tree, const CorrespondenceRejector &rejector, const Eigen::Isometry3d &T, std::vector< Factor > &factors) const
Definition: reduction_omp.hpp:21