CMS 3D CMS Logo

VoronoiAlgorithm.h
Go to the documentation of this file.
1 #ifndef __VoronoiAlgorithm_h__
2 #define __VoronoiAlgorithm_h__
3 
4 #include <cmath>
5 
6 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
7 #include <CGAL/Delaunay_triangulation_2.h>
8 #include <CGAL/Voronoi_diagram_2.h>
9 #include <CGAL/Delaunay_triangulation_adaptation_traits_2.h>
10 #include <CGAL/Delaunay_triangulation_adaptation_policies_2.h>
11 #include <CGAL/Polygon_2.h>
12 
13 #include <boost/multi_array.hpp>
14 
17 
19 private:
20  //
21  typedef CGAL::Delaunay_triangulation_2<
22  CGAL::Exact_predicates_inexact_constructions_kernel>
25  class particle_t {
26  public:
29  double area;
32  std::set<std::vector<particle_t>::iterator> incident;
34  unsigned int i, double a = NAN,
35  double ps = NAN)
36  : momentum(p), reduced_particle_flow_id(i), area(a),
37  momentum_perp_subtracted(ps),
38  momentum_perp_subtracted_unequalized(ps),
39  incident(std::set<std::vector<particle_t>::
40  iterator>())
41  {
42  }
43  inline operator point_2d_t(void) const
44  {
45  return point_2d_t(momentum.Eta(), momentum.Phi());
46  }
47  };
48  typedef std::vector<particle_t> event_t;
49  // Remaining CGAL classes
50  typedef CGAL::Voronoi_diagram_2<
52  CGAL::Delaunay_triangulation_adaptation_traits_2<
53  delaunay_triangulation_t>,
54  CGAL::
55  Delaunay_triangulation_caching_degeneracy_removal_policy_2<
56  delaunay_triangulation_t> > voronoi_diagram_t;
57  typedef CGAL::Polygon_2<
58  CGAL::Exact_predicates_inexact_constructions_kernel>
60 public:
61  static const size_t nreduced_particle_flow_id = 3;
62  static const size_t nfourier = 5;
63 protected:
64  std::vector<double> _edge_pseudorapidity;
65  std::vector<double> _cms_hcal_edge_pseudorapidity;
66  std::vector<double> _cms_ecal_edge_pseudorapidity;
68  std::pair<double, double> _equalization_threshold;
72  event_t _event;
73  boost::multi_array<double, 4> *_perp_fourier;
74  std::vector<double> _feature;
75  std::vector<bool> _active;
76  std::vector<std::pair<size_t, size_t> > _recombine;
77  std::vector<std::vector<size_t> > _recombine_index;
78  std::vector<std::vector<size_t> > _recombine_unsigned;
79  std::vector<double> _recombine_tie;
80  size_t _ncost;
81  std::vector<size_t> _nblock_subtract;
83  void *_lp_problem;
84  // calibrations
85  const UECalibration *ue;
86 private:
87  void initialize_geometry(void);
88  void allocate(void);
89  void deallocate(void);
90  void event_fourier(void);
91  void feature_extract(void);
92  void voronoi_area_incident(void);
93  void subtract_momentum(void);
94  void recombine_link(void);
95  void lp_populate(void *lp_problem);
96  void equalize(void);
97  void remove_nonpositive(void);
98  void subtract_if_necessary(void);
99 public:
101  const UECalibration *ue,
102  const double dr_max,
103  const std::pair<double, double> equalization_threshold =
104  std::pair<double, double>(5.0, 35.0),
105  const bool remove_nonpositive = true);
106  ~VoronoiAlgorithm(void);
116  void push_back_particle(
117  const double perp, const double pseudorapidity,
118  const double azimuth,
119  const unsigned int reduced_particle_flow_id);
123  void clear(void);
129  std::vector<double> subtracted_equalized_perp(void);
130  std::vector<double> subtracted_unequalized_perp(void);
137  std::vector<double> particle_area(void);
146  std::vector<std::set<size_t> > particle_incident(void);
147  std::vector<double> perp_fourier(void);
148  size_t nedge_pseudorapidity(void) const;
149 };
150 
151 #endif
std::vector< bool > _active
std::vector< std::vector< size_t > > _recombine_index
math::PtEtaPhiELorentzVector momentum
std::vector< double > perp_fourier(void)
void subtract_momentum(void)
void event_fourier(void)
CGAL::Voronoi_diagram_2< delaunay_triangulation_t, CGAL::Delaunay_triangulation_adaptation_traits_2< delaunay_triangulation_t >, CGAL::Delaunay_triangulation_caching_degeneracy_removal_policy_2< delaunay_triangulation_t > > voronoi_diagram_t
double _positive_bound_scale
void remove_nonpositive(void)
std::vector< double > _edge_pseudorapidity
void initialize_geometry(void)
std::pair< double, double > Point
Definition: CaloEllipse.h:18
std::vector< double > _recombine_tie
particle_t(math::PtEtaPhiELorentzVector p, unsigned int i, double a=NAN, double ps=NAN)
const UECalibration * ue
delaunay_triangulation_t::Point point_2d_t
std::vector< double > _cms_ecal_edge_pseudorapidity
std::vector< double > particle_area(void)
void lp_populate(void *lp_problem)
void subtract_if_necessary(void)
void push_back_particle(const double perp, const double pseudorapidity, const double azimuth, const unsigned int reduced_particle_flow_id)
std::vector< particle_t > event_t
std::vector< std::pair< size_t, size_t > > _recombine
static const size_t nreduced_particle_flow_id
std::vector< std::set< size_t > > particle_incident(void)
static const size_t nfourier
std::vector< size_t > _nblock_subtract
std::vector< double > _feature
double _radial_distance_square_max
boost::multi_array< double, 4 > * _perp_fourier
PtEtaPhiELorentzVectorD PtEtaPhiELorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:27
std::vector< std::vector< size_t > > _recombine_unsigned
void feature_extract(void)
VoronoiAlgorithm(const UECalibration *ue, const double dr_max, const std::pair< double, double > equalization_threshold=std::pair< double, double >(5.0, 35.0), const bool remove_nonpositive=true)
T perp() const
Magnitude of transverse component.
void voronoi_area_incident(void)
double a
Definition: hdecay.h:121
std::pair< double, double > _equalization_threshold
size_t nedge_pseudorapidity(void) const
std::vector< double > _cms_hcal_edge_pseudorapidity
CGAL::Polygon_2< CGAL::Exact_predicates_inexact_constructions_kernel > polygon_t
std::vector< double > subtracted_unequalized_perp(void)
std::set< std::vector< particle_t >::iterator > incident
std::vector< double > subtracted_equalized_perp(void)
CGAL::Delaunay_triangulation_2< CGAL::Exact_predicates_inexact_constructions_kernel > delaunay_triangulation_t