CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Classes | Public Member Functions | Static Public Attributes | Protected Attributes | Private Types | Private Member Functions
VoronoiAlgorithm Class Reference

#include <VoronoiAlgorithm.h>

Classes

class  particle_t
 

Public Member Functions

void clear (void)
 
size_t nedge_pseudorapidity (void) const
 
std::vector< double > particle_area (void)
 
std::vector< std::set< size_t > > particle_incident (void)
 
std::vector< double > perp_fourier (void)
 
void push_back_particle (const double perp, const double pseudorapidity, const double azimuth, const unsigned int reduced_particle_flow_id)
 
std::vector< double > subtracted_equalized_perp (void)
 
std::vector< double > subtracted_unequalized_perp (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)
 
 ~VoronoiAlgorithm (void)
 

Static Public Attributes

static const size_t nfourier = 5
 
static const size_t nreduced_particle_flow_id = 3
 

Protected Attributes

std::vector< bool > _active
 
std::vector< double > _cms_ecal_edge_pseudorapidity
 
std::vector< double > _cms_hcal_edge_pseudorapidity
 
std::vector< double > _edge_pseudorapidity
 
std::pair< double, double > _equalization_threshold
 
event_t _event
 
std::vector< double > _feature
 
void * _lp_environment
 
void * _lp_problem
 
std::vector< size_t > _nblock_subtract
 
size_t _ncost
 
boost::multi_array< double, 4 > * _perp_fourier
 
double _positive_bound_scale
 
double _radial_distance_square_max
 
std::vector< std::pair< size_t,
size_t > > 
_recombine
 
std::vector< std::vector
< size_t > > 
_recombine_index
 
std::vector< double > _recombine_tie
 
std::vector< std::vector
< size_t > > 
_recombine_unsigned
 
bool _remove_nonpositive
 
bool _subtracted
 
const UECalibrationue
 

Private Types

typedef
CGAL::Delaunay_triangulation_2
< CGAL::Exact_predicates_inexact_constructions_kernel > 
delaunay_triangulation_t
 
typedef std::vector< particle_tevent_t
 
typedef
delaunay_triangulation_t::Point 
point_2d_t
 
typedef CGAL::Polygon_2
< CGAL::Exact_predicates_inexact_constructions_kernel > 
polygon_t
 
typedef
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
 

Private Member Functions

void allocate (void)
 
void deallocate (void)
 
void equalize (void)
 
void event_fourier (void)
 
void feature_extract (void)
 
void initialize_geometry (void)
 
void lp_populate (void *lp_problem)
 
void recombine_link (void)
 
void remove_nonpositive (void)
 
void subtract_if_necessary (void)
 
void subtract_momentum (void)
 
void voronoi_area_incident (void)
 

Detailed Description

Definition at line 19 of file VoronoiAlgorithm.h.

Member Typedef Documentation

typedef CGAL::Delaunay_triangulation_2< CGAL::Exact_predicates_inexact_constructions_kernel> VoronoiAlgorithm::delaunay_triangulation_t
private

Definition at line 24 of file VoronoiAlgorithm.h.

typedef std::vector<particle_t> VoronoiAlgorithm::event_t
private

Definition at line 49 of file VoronoiAlgorithm.h.

Definition at line 25 of file VoronoiAlgorithm.h.

typedef CGAL::Polygon_2< CGAL::Exact_predicates_inexact_constructions_kernel> VoronoiAlgorithm::polygon_t
private

Definition at line 60 of file VoronoiAlgorithm.h.

typedef 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> > VoronoiAlgorithm::voronoi_diagram_t
private

Definition at line 57 of file VoronoiAlgorithm.h.

Constructor & Destructor Documentation

VoronoiAlgorithm::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 
)

Definition at line 1830 of file VoronoiAlgorithm.cc.

References _edge_pseudorapidity, allocate(), initialize_geometry(), and nedge_pseudorapidity().

1836  _equalization_threshold(equalization_threshold),
1837  _radial_distance_square_max(dr_max * dr_max),
1838  _positive_bound_scale(0.2),
1839  _subtracted(false),
1840  ue(ue_)
1841  {
1843 
1844  static const size_t nedge_pseudorapidity = 15 + 1;
1845  static const double edge_pseudorapidity[nedge_pseudorapidity] = {
1846  -5.191, -2.650, -2.043, -1.740, -1.479, -1.131, -0.783, -0.522, 0.522, 0.783, 1.131, 1.479, 1.740, 2.043, 2.650, 5.191
1847  };
1848 
1849  _edge_pseudorapidity = std::vector<double>(
1850  edge_pseudorapidity,
1851  edge_pseudorapidity + nedge_pseudorapidity);
1852  allocate();
1853  }
double _positive_bound_scale
void remove_nonpositive(void)
std::vector< double > _edge_pseudorapidity
void initialize_geometry(void)
const UECalibration * ue
double _radial_distance_square_max
std::pair< double, double > _equalization_threshold
size_t nedge_pseudorapidity(void) const
VoronoiAlgorithm::~VoronoiAlgorithm ( void  )

Definition at line 1855 of file VoronoiAlgorithm.cc.

References deallocate().

1856  {
1857  deallocate();
1858  }

Member Function Documentation

void VoronoiAlgorithm::allocate ( void  )
private

Definition at line 930 of file VoronoiAlgorithm.cc.

References _edge_pseudorapidity, _perp_fourier, nfourier, and nreduced_particle_flow_id.

Referenced by VoronoiAlgorithm().

931  {
932  _perp_fourier = new boost::multi_array<double, 4>(
933  boost::extents[_edge_pseudorapidity.size() - 1]
935  }
std::vector< double > _edge_pseudorapidity
static const size_t nreduced_particle_flow_id
static const size_t nfourier
boost::multi_array< double, 4 > * _perp_fourier
void VoronoiAlgorithm::clear ( void  )

Clears the list of unsubtracted particles

Definition at line 1882 of file VoronoiAlgorithm.cc.

References _event, and _subtracted.

Referenced by VoronoiBackgroundProducer::produce().

1883  {
1884  _event.clear();
1885  _subtracted = false;
1886  }
void VoronoiAlgorithm::deallocate ( void  )
private

Definition at line 936 of file VoronoiAlgorithm.cc.

References _perp_fourier.

Referenced by ~VoronoiAlgorithm().

937  {
938  delete _perp_fourier;
939  }
boost::multi_array< double, 4 > * _perp_fourier
void VoronoiAlgorithm::equalize ( void  )
private

Definition at line 1771 of file VoronoiAlgorithm.cc.

References _event, _lp_environment, _nblock_subtract, _ncost, _recombine, plotBeamSpotDB::first, relval_steps::k, lp_populate(), pi, recombine_link(), edm::second(), and x.

Referenced by subtract_if_necessary().

1772  {
1773  bpmpd_problem_t lp_problem = reinterpret_cast<bpmpd_environment_t *>(_lp_environment)->problem();
1774 
1775  recombine_link();
1776  lp_populate(&lp_problem);
1777  lp_problem.optimize();
1778 
1779  int solution_status;
1780  double objective_value;
1781  std::vector<double> x;
1782  std::vector<double> pi;
1783 
1784  lp_problem.solve(solution_status, objective_value,
1785  x, pi);
1786 
1787  for (size_t k = _ncost; k < x.size(); k++) {
1788  if (_event[_recombine[k - _ncost].first].
1789  momentum_perp_subtracted < 0 &&
1791  momentum_perp_subtracted >= 0 && x[k] >= 0) {
1792  _event[_recombine[k - _ncost].first].
1793  momentum_perp_subtracted += x[k];
1794  _event[_recombine[k - _ncost].second].
1795  momentum_perp_subtracted -= x[k];
1796  }
1797  }
1798  for (size_t k = 0; k < _event.size(); k++) {
1799  if (_nblock_subtract[k] != 0 &&
1800  x[_nblock_subtract[k]] >= 0) {
1801  _event[k].momentum_perp_subtracted -=
1802  x[_nblock_subtract[k]];
1803  }
1804  }
1805  }
const Double_t pi
U second(std::pair< T, U > const &p)
void lp_populate(void *lp_problem)
std::vector< std::pair< size_t, size_t > > _recombine
std::vector< size_t > _nblock_subtract
Definition: DDAxes.h:10
void VoronoiAlgorithm::event_fourier ( void  )
private

Definition at line 940 of file VoronoiAlgorithm.cc.

References _edge_pseudorapidity, _event, _perp_fourier, funct::cos(), lumiContext::fill, relval_steps::k, prof2calltree::l, nfourier, and funct::sin().

Referenced by subtract_if_necessary().

941  {
942  std::fill(_perp_fourier->data(),
943  _perp_fourier->data() +
944  _perp_fourier->num_elements(),
945  0);
946 
947  for (std::vector<particle_t>::const_iterator iterator =
948  _event.begin();
949  iterator != _event.end(); iterator++) {
950  const unsigned int reduced_id =
951  iterator->reduced_particle_flow_id;
952 
953  for (size_t k = 1; k < _edge_pseudorapidity.size();
954  k++) {
955  if (iterator->momentum.Eta() >=
956  _edge_pseudorapidity[k - 1] &&
957  iterator->momentum.Eta() <
959  const double azimuth =
960  iterator->momentum.Phi();
961 
962  for (size_t l = 0; l < nfourier; l++) {
963  (*_perp_fourier)[k - 1][reduced_id]
964  [l][0] +=
965  iterator->momentum.Pt() *
966  cos(l * azimuth);
967  (*_perp_fourier)[k - 1][reduced_id]
968  [l][1] +=
969  iterator->momentum.Pt() *
970  sin(l * azimuth);
971  }
972  }
973  }
974  }
975  }
string fill
Definition: lumiContext.py:319
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< double > _edge_pseudorapidity
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static const size_t nfourier
boost::multi_array< double, 4 > * _perp_fourier
void VoronoiAlgorithm::feature_extract ( void  )
private

Definition at line 976 of file VoronoiAlgorithm.cc.

References _edge_pseudorapidity, _feature, j, relval_steps::k, nfourier, nreduced_particle_flow_id, pileupReCalc_HLTpaths::scale, and mathSSE::sqrt().

Referenced by subtract_if_necessary().

977  {
978  const size_t nfeature = 2 * nfourier - 1;
979 
980  _feature.resize(nfeature);
981 
982  // Scale factor to get 95% of the coefficient below 1.0
983  // (where however one order of magnitude tolerance is
984  // acceptable). This is valid for nfourier < 18 (where
985  // interference behavior with the HF geometry starts to
986  // appear)
987 
988  std::vector<double> scale(nfourier, 1.0 / 200.0);
989 
990  if (nfourier >= 1) {
991  scale[0] = 1.0 / 5400.0;
992  }
993  if (nfourier >= 2) {
994  scale[1] = 1.0 / 130.0;
995  }
996  if (nfourier >= 3) {
997  scale[2] = 1.0 / 220.0;
998  }
999 
1000  const size_t index_edge_end =
1001  _edge_pseudorapidity.size() - 2;
1002 
1003  _feature[0] = 0;
1004  for (size_t j = 0; j < nreduced_particle_flow_id; j++) {
1005  _feature[0] += scale[0] *
1006  ((*_perp_fourier)[0 ][j][0][0] +
1007  (*_perp_fourier)[index_edge_end][j][0][0]);
1008  }
1009 
1010  for (size_t k = 1; k < nfourier; k++) {
1011  _feature[2 * k - 1] = 0;
1012  for (size_t j = 0; j < nreduced_particle_flow_id; j++) {
1013  _feature[2 * k - 1] += scale[k] *
1014  ((*_perp_fourier)[0 ][j][k][0] +
1015  (*_perp_fourier)[index_edge_end][j][k][0]);
1016  }
1017  _feature[2 * k] = 0;
1018  for (size_t j = 0; j < nreduced_particle_flow_id; j++) {
1019  _feature[2 * k] += scale[k] *
1020  ((*_perp_fourier)[0 ][j][k][1] +
1021  (*_perp_fourier)[index_edge_end][j][k][1]);
1022  }
1023  }
1024 
1025 
1026 #if 0
1027  const double event_plane = atan2(_feature[4], _feature[3]);
1028  const double v2 =
1029  sqrt(_feature[3] * _feature[3] +
1030  _feature[4] * _feature[4]) / _feature[0];
1031 #endif
1032  }
std::vector< double > _edge_pseudorapidity
T sqrt(T t)
Definition: SSEVec.h:48
static const size_t nreduced_particle_flow_id
static const size_t nfourier
int j
Definition: DBlmapReader.cc:9
std::vector< double > _feature
void VoronoiAlgorithm::initialize_geometry ( void  )
private

Definition at line 896 of file VoronoiAlgorithm.cc.

References _cms_ecal_edge_pseudorapidity, _cms_hcal_edge_pseudorapidity, and i.

Referenced by VoronoiAlgorithm().

897  {
898  static const size_t ncms_hcal_edge_pseudorapidity = 82 + 1;
899  static const double cms_hcal_edge_pseudorapidity[
900  ncms_hcal_edge_pseudorapidity] = {
901  -5.191, -4.889, -4.716, -4.538, -4.363, -4.191, -4.013,
902  -3.839, -3.664, -3.489, -3.314, -3.139, -2.964, -2.853,
903  -2.650, -2.500, -2.322, -2.172, -2.043, -1.930, -1.830,
904  -1.740, -1.653, -1.566, -1.479, -1.392, -1.305, -1.218,
905  -1.131, -1.044, -0.957, -0.879, -0.783, -0.696, -0.609,
906  -0.522, -0.435, -0.348, -0.261, -0.174, -0.087,
907  0.000,
908  0.087, 0.174, 0.261, 0.348, 0.435, 0.522, 0.609,
909  0.696, 0.783, 0.879, 0.957, 1.044, 1.131, 1.218,
910  1.305, 1.392, 1.479, 1.566, 1.653, 1.740, 1.830,
911  1.930, 2.043, 2.172, 2.322, 2.500, 2.650, 2.853,
912  2.964, 3.139, 3.314, 3.489, 3.664, 3.839, 4.013,
913  4.191, 4.363, 4.538, 4.716, 4.889, 5.191
914  };
915 
916  _cms_hcal_edge_pseudorapidity = std::vector<double>(
917  cms_hcal_edge_pseudorapidity,
918  cms_hcal_edge_pseudorapidity +
919  ncms_hcal_edge_pseudorapidity);
920 
921  static const size_t ncms_ecal_edge_pseudorapidity = 344 + 1;
922 
923  for (size_t i = 0; i < ncms_ecal_edge_pseudorapidity; i++) {
925  i * (2 * 2.9928 /
926  (ncms_ecal_edge_pseudorapidity - 1)) -
927  2.9928);
928  }
929  }
int i
Definition: DBlmapReader.cc:9
std::vector< double > _cms_ecal_edge_pseudorapidity
std::vector< double > _cms_hcal_edge_pseudorapidity
void VoronoiAlgorithm::lp_populate ( void *  lp_problem)
private

Definition at line 1402 of file VoronoiAlgorithm.cc.

References _active, _equalization_threshold, _event, _nblock_subtract, _ncost, _positive_bound_scale, _recombine_index, _recombine_tie, _recombine_unsigned, alignCSCRings::e, cond::serialization::equal(), i, infinity, M_PI, bookConverter::max, min(), nedge_pseudorapidity(), AlCaHLTBitMon_ParallelJobs::p, jetcorrextractor::sign(), and histoStyle::weight.

Referenced by equalize().

1403  {
1404  bpmpd_problem_t *p = reinterpret_cast<bpmpd_problem_t *>(lp_problem);
1405 
1406  // The minimax problem is transformed into the LP notation
1407  // using the cost variable trick:
1408  //
1409  // Minimize c
1410  // Subject to:
1411  // c + sum_l t_kl + n_k >= 0 for negative cells n_k
1412  // c - sum_k t_kl + p_l >= 0 for positive cells p_l
1413 
1414  // Common LP mistakes during code development and their
1415  // CPLEX errors when running CPLEX in data checking mode:
1416  //
1417  // Error 1201 (column index ... out of range): Bad column
1418  // indexing, usually index_column out of bound for the
1419  // cost variables.
1420  //
1421  // Error 1222 (duplicate entry): Forgetting to increment
1422  // index_row, or index_column out of bound for the cost
1423  // variables.
1424 
1425  p->set_objective_sense(bpmpd_problem_t::minimize);
1426 
1427  // Rows (RHS of the constraints) of the LP problem
1428 
1429  static const size_t nsector_azimuth = 12;
1430 
1431  // Approximatively 2 pi / nsector_azimuth segmentation of
1432  // the CMS HCAL granularity
1433 
1434  static const size_t ncms_hcal_edge_pseudorapidity = 19 + 1;
1435  static const double cms_hcal_edge_pseudorapidity[
1436  ncms_hcal_edge_pseudorapidity] = {
1437  -5.191, -4.538, -4.013,
1438  -3.489, -2.853, -2.322, -1.830, -1.305, -0.783, -0.261,
1439  0.261, 0.783, 1.305, 1.830, 2.322, 2.853, 3.489,
1440  4.013, 4.538, 5.191
1441  };
1442 
1443  size_t nedge_pseudorapidity;
1444  const double *edge_pseudorapidity;
1445 
1446  nedge_pseudorapidity = ncms_hcal_edge_pseudorapidity;
1447  edge_pseudorapidity = cms_hcal_edge_pseudorapidity;
1448 
1449  const size_t nsuperblock = (nedge_pseudorapidity - 2) *
1450  nsector_azimuth;
1451 
1452  size_t index_row = 0;
1453  for (size_t index_pseudorapidity = 0;
1454  index_pseudorapidity < nedge_pseudorapidity - 2;
1455  index_pseudorapidity++) {
1456  for (size_t index_azimuth = 0;
1457  index_azimuth < nsector_azimuth - 1;
1458  index_azimuth++) {
1459  const size_t index_column =
1460  index_pseudorapidity * nsector_azimuth +
1461  index_azimuth;
1462  p->push_back_row(
1463  bpmpd_problem_t::greater_equal, 0);
1464  p->push_back_coefficient(
1465  index_row, index_column, 1);
1466  p->push_back_coefficient(
1467  index_row, nsuperblock + index_column, -1);
1468  index_row++;
1469  p->push_back_row(
1470  bpmpd_problem_t::greater_equal, 0);
1471  p->push_back_coefficient(
1472  index_row, index_column, 1);
1473  p->push_back_coefficient(
1474  index_row, nsuperblock + index_column + 1, -1);
1475  index_row++;
1476  p->push_back_row(
1477  bpmpd_problem_t::greater_equal, 0);
1478  p->push_back_coefficient(
1479  index_row, index_column, 1);
1480  p->push_back_coefficient(
1481  index_row,
1482  nsuperblock + index_column + nsector_azimuth, -1);
1483  index_row++;
1484  p->push_back_row(
1485  bpmpd_problem_t::greater_equal, 0);
1486  p->push_back_coefficient(
1487  index_row, index_column, 1);
1488  p->push_back_coefficient(
1489  index_row,
1490  nsuperblock + index_column + nsector_azimuth + 1,
1491  -1);
1492  index_row++;
1493  }
1494  const size_t index_column =
1495  index_pseudorapidity * nsector_azimuth +
1496  nsector_azimuth - 1;
1497  p->push_back_row(
1498  bpmpd_problem_t::greater_equal, 0);
1499  p->push_back_coefficient(
1500  index_row, index_column, 1);
1501  p->push_back_coefficient(
1502  index_row, nsuperblock + index_column, -1);
1503  index_row++;
1504  p->push_back_row(
1505  bpmpd_problem_t::greater_equal, 0);
1506  p->push_back_coefficient(
1507  index_row, index_column, 1);
1508  p->push_back_coefficient(
1509  index_row,
1510  nsuperblock + index_column - (nsector_azimuth - 1),
1511  -1);
1512  index_row++;
1513  p->push_back_row(
1514  bpmpd_problem_t::greater_equal, 0);
1515  p->push_back_coefficient(
1516  index_row, index_column, 1);
1517  p->push_back_coefficient(
1518  index_row,
1519  nsuperblock + index_column + nsector_azimuth, -1);
1520  index_row++;
1521  p->push_back_row(
1522  bpmpd_problem_t::greater_equal, 0);
1523  p->push_back_coefficient(
1524  index_row, index_column, 1);
1525  p->push_back_coefficient(
1526  index_row,
1527  nsuperblock + index_column + nsector_azimuth -
1528  (nsector_azimuth - 1),
1529  -1);
1530  index_row++;
1531  }
1532 
1533  const size_t nstaggered_block =
1534  (nedge_pseudorapidity - 1) * nsector_azimuth;
1535  const size_t nblock = nsuperblock + 2 * nstaggered_block;
1536 
1537  _nblock_subtract = std::vector<size_t>(_event.size(), 0);
1538 
1539  std::vector<size_t>
1540  positive_index(_event.size(), _event.size());
1541  size_t positive_count = 0;
1542 
1543  for (std::vector<particle_t>::const_iterator iterator =
1544  _event.begin();
1545  iterator != _event.end(); iterator++) {
1546  if (iterator->momentum_perp_subtracted >= 0) {
1547  positive_index[iterator - _event.begin()] =
1548  positive_count;
1549  positive_count++;
1550  }
1551  }
1552 
1553  _ncost = nblock + positive_count;
1554 
1555  const double sum_unequalized_0 = _equalization_threshold.first;
1556  const double sum_unequalized_1 = (2.0 / 3.0) * _equalization_threshold.first + (1.0 / 3.0) * _equalization_threshold.second;
1557  const double sum_unequalized_2 = (1.0 / 3.0) * _equalization_threshold.first + (2.0 / 3.0) * _equalization_threshold.second;
1558  const double sum_unequalized_3 = _equalization_threshold.second;
1559 
1560  std::vector<particle_t>::const_iterator
1561  iterator_particle = _event.begin();
1562  std::vector<bool>::const_iterator iterator_active =
1563  _active.begin();
1564  std::vector<std::vector<size_t> >::const_iterator
1565  iterator_recombine_index_outer =
1566  _recombine_index.begin();
1567  std::vector<std::vector<size_t> >::const_iterator
1568  iterator_recombine_unsigned_outer =
1569  _recombine_unsigned.begin();
1570  size_t index_column_max = _ncost - 1;
1571  for (; iterator_particle != _event.end();
1572  iterator_particle++, iterator_active++,
1573  iterator_recombine_index_outer++,
1574  iterator_recombine_unsigned_outer++) {
1575  if (*iterator_active) {
1576  int index_pseudorapidity = -1;
1577 
1579  for (size_t i = 1; i < nedge_pseudorapidity; i++) {
1580  if (iterator_particle->momentum.Eta() >= edge_pseudorapidity[i - 1] &&
1581  iterator_particle->momentum.Eta() < edge_pseudorapidity[i]) {
1582  index_pseudorapidity = i - 1;
1583  }
1584  }
1585 
1586  const int index_azimuth = floor(
1587  (iterator_particle->momentum.Phi() + M_PI) *
1588  ((nsector_azimuth >> 1) / M_PI));
1589 
1590  if (index_pseudorapidity != -1) {
1591  // p_i - sum t - u = c_i
1592  // or: c_i + u + sum_t = p_i
1593  // n_i + sum t - u <= 0
1594  // or: u - sum_t >= n_i
1595 
1596  // Inequality RHS
1597  p->push_back_row(
1598  iterator_particle->momentum_perp_subtracted >= 0 ?
1600  bpmpd_problem_t::greater_equal,
1601  iterator_particle->momentum_perp_subtracted);
1602 
1603  // Energy transfer coefficients t_kl
1604  const double sign = iterator_particle->momentum_perp_subtracted >= 0 ? 1 : -1;
1605  const size_t index_column_block_subtract =
1606  nsuperblock +
1607  (nedge_pseudorapidity - 1) * nsector_azimuth +
1608  index_pseudorapidity * nsector_azimuth +
1609  index_azimuth;
1610 
1611  _nblock_subtract[iterator_particle - _event.begin()] =
1612  index_column_block_subtract;
1613 
1614  if (iterator_particle->momentum_perp_subtracted >= 0) {
1615  const size_t index_column_cost =
1616  nblock + positive_index[iterator_particle - _event.begin()];
1617 
1618  p->push_back_coefficient(
1619  index_row, index_column_cost, 1);
1620  index_column_max =
1621  std::max(index_column_max, index_column_cost);
1622  }
1623  p->push_back_coefficient(
1624  index_row, index_column_block_subtract, 1);
1625  index_column_max =
1626  std::max(index_column_max, index_column_block_subtract);
1627 
1628  for (std::vector<size_t>::const_iterator
1629  iterator_recombine_index_inner =
1630  iterator_recombine_index_outer->begin();
1631  iterator_recombine_index_inner !=
1632  iterator_recombine_index_outer->end();
1633  iterator_recombine_index_inner++) {
1634  const size_t index_column =
1635  *iterator_recombine_index_inner +
1636  _ncost;
1637 
1638  p->push_back_coefficient(
1639  index_row, index_column, sign);
1640  index_column_max =
1641  std::max(index_column_max, index_column);
1642  }
1643  index_row++;
1644 
1645  const size_t index_column_block =
1646  nsuperblock +
1647  index_pseudorapidity * nsector_azimuth +
1648  index_azimuth;
1649 
1650  // sum_R c_i - o_i >= -d
1651  // or: d + sum_R c_i >= o_i
1652  // sum_R c_i - o_i <= d
1653  // or: d - sum_R c_i >= -o_i
1654 
1655  double sum_unequalized;
1656 
1657  sum_unequalized = 0;
1658  for (std::vector<size_t>::const_iterator
1659  iterator_recombine_unsigned_inner =
1660  iterator_recombine_unsigned_outer->begin();
1661  iterator_recombine_unsigned_inner !=
1662  iterator_recombine_unsigned_outer->end();
1663  iterator_recombine_unsigned_inner++) {
1664  sum_unequalized +=
1665  _event[*iterator_recombine_unsigned_inner].momentum_perp_subtracted;
1666  }
1667  sum_unequalized = std::max(0.0, sum_unequalized);
1668 
1669  if (sum_unequalized >= sum_unequalized_3 ||
1670  (sum_unequalized >= sum_unequalized_2 &&
1671  (iterator_particle - _event.begin()) % 2 == 0) ||
1672  (sum_unequalized >= sum_unequalized_1 &&
1673  (iterator_particle - _event.begin()) % 4 == 0) ||
1674  (sum_unequalized >= sum_unequalized_0 &&
1675  (iterator_particle - _event.begin()) % 8 == 0)) {
1676 
1677  const double weight = sum_unequalized *
1678  std::min(1.0, std::max(1e-3,
1679  iterator_particle->area));
1680 
1681  if (weight > 0) {
1682  p->push_back_row(
1683  bpmpd_problem_t::greater_equal,
1684  sum_unequalized);
1685 
1686  p->push_back_coefficient(
1687  index_row, index_column_block, 1.0 / weight);
1688 
1689  for (std::vector<size_t>::const_iterator
1690  iterator_recombine_unsigned_inner =
1691  iterator_recombine_unsigned_outer->begin();
1692  iterator_recombine_unsigned_inner !=
1693  iterator_recombine_unsigned_outer->end();
1694  iterator_recombine_unsigned_inner++) {
1695  if (_event[*iterator_recombine_unsigned_inner].momentum_perp_subtracted >= 0) {
1696  const size_t index_column_cost =
1697  nblock +
1698  positive_index[*iterator_recombine_unsigned_inner];
1699 
1700  p->push_back_coefficient(
1701  index_row, index_column_cost, 1);
1702  index_column_max =
1703  std::max(index_column_max, index_column_cost);
1704  }
1705  }
1706  index_row++;
1707 
1708  p->push_back_row(
1709  bpmpd_problem_t::greater_equal,
1710  -sum_unequalized);
1711 
1712  p->push_back_coefficient(
1713  index_row, index_column_block, _positive_bound_scale / weight);
1714 
1715  for (std::vector<size_t>::const_iterator iterator_recombine_unsigned_inner = iterator_recombine_unsigned_outer->begin();
1716  iterator_recombine_unsigned_inner != iterator_recombine_unsigned_outer->end();
1717  iterator_recombine_unsigned_inner++) {
1718  if (_event[*iterator_recombine_unsigned_inner].momentum_perp_subtracted >= 0) {
1719  const size_t index_column_cost =
1720  nblock +
1721  positive_index[*iterator_recombine_unsigned_inner];
1722 
1723  p->push_back_coefficient(
1724  index_row, index_column_cost, -1);
1725  index_column_max =
1726  std::max(index_column_max, index_column_cost);
1727  }
1728  }
1729  index_row++;
1730  }
1731 
1732  }
1733 
1734  }
1735  }
1736  }
1737 
1738  // Epsilon that breaks the degeneracy, in the same units
1739  // as the pT of the event (i.e. GeV)
1740  static const double epsilon_degeneracy = 1e-2;
1741 
1742  // Columns (variables and the objective coefficients) of
1743  // the LP problem
1744  //
1745  // Cost variables (objective coefficient 1)
1746  for (size_t i = 0; i < nsuperblock; i++) {
1747  p->push_back_column(
1749  }
1750  for (size_t i = nsuperblock; i < nsuperblock + nstaggered_block; i++) {
1751  p->push_back_column(
1753  }
1754  for (size_t i = nsuperblock + nstaggered_block; i < nsuperblock + 2 * nstaggered_block; i++) {
1755  p->push_back_column(
1757  }
1758  for (size_t i = nsuperblock + 2 * nstaggered_block; i < _ncost; i++) {
1759  p->push_back_column(
1761  }
1762  //fprintf(stderr, "%s:%d: %lu %lu\n", __FILE__, __LINE__, index_column_max, recombine_tie.size());
1763  // Energy transfer coefficients t_kl (objective
1764  // coefficient 0 + epsilon)
1765  for (size_t i = _ncost; i <= index_column_max; i++) {
1766  p->push_back_column(
1767  epsilon_degeneracy * _recombine_tie[i - _ncost],
1769  }
1770  }
std::vector< bool > _active
std::vector< std::vector< size_t > > _recombine_index
int i
Definition: DBlmapReader.cc:9
double sign(double x)
double _positive_bound_scale
std::vector< double > _recombine_tie
bool equal(const T &first, const T &second)
Definition: Equal.h:34
const double infinity
std::vector< size_t > _nblock_subtract
T min(T a, T b)
Definition: MathUtil.h:58
#define M_PI
std::vector< std::vector< size_t > > _recombine_unsigned
std::pair< double, double > _equalization_threshold
size_t nedge_pseudorapidity(void) const
int weight
Definition: histoStyle.py:50
size_t VoronoiAlgorithm::nedge_pseudorapidity ( void  ) const

Definition at line 1980 of file VoronoiAlgorithm.cc.

References _edge_pseudorapidity.

Referenced by lp_populate(), and VoronoiAlgorithm().

1981  {
1982  return _edge_pseudorapidity.size();
1983  }
std::vector< double > _edge_pseudorapidity
std::vector< double > VoronoiAlgorithm::particle_area ( void  )

Returns the area in the Voronoi diagram diagram occupied by a given particle

Returns
vector of area

Definition at line 1926 of file VoronoiAlgorithm.cc.

References _event, run_regression::ret, and subtract_if_necessary().

Referenced by VoronoiBackgroundProducer::produce().

1927  {
1929 
1930  std::vector<double> ret;
1931 
1932  for (std::vector<particle_t>::const_iterator iterator =
1933  _event.begin();
1934  iterator != _event.end(); iterator++) {
1935  ret.push_back(iterator->area);
1936  }
1937 
1938  return ret;
1939  }
void subtract_if_necessary(void)
std::vector< std::set< size_t > > VoronoiAlgorithm::particle_incident ( void  )

Returns the incident particles in the Delaunay diagram (particles that has a given particle as the nearest neighbor)

Returns
vector of sets of incident particles indices, using the original indexing

Definition at line 1948 of file VoronoiAlgorithm.cc.

References _event, alignCSCRings::e, run_regression::ret, and subtract_if_necessary().

1949  {
1951 
1952  std::vector<std::set<size_t> > ret;
1953 
1954  for (std::vector<particle_t>::const_iterator
1955  iterator_outer = _event.begin();
1956  iterator_outer != _event.end(); iterator_outer++) {
1957  std::set<size_t> e;
1958 
1959  for (std::set<std::vector<particle_t>::iterator>::
1960  const_iterator iterator_inner =
1961  iterator_outer->incident.begin();
1962  iterator_inner != iterator_outer->incident.begin();
1963  iterator_inner++) {
1964  e.insert(*iterator_inner - _event.begin());
1965  }
1966  ret.push_back(e);
1967  }
1968 
1969  return ret;
1970  }
void subtract_if_necessary(void)
std::vector< double > VoronoiAlgorithm::perp_fourier ( void  )

Definition at line 1971 of file VoronoiAlgorithm.cc.

References _perp_fourier, and subtract_if_necessary().

Referenced by VoronoiBackgroundProducer::produce().

1972  {
1974 
1975  return std::vector<double>(
1976  _perp_fourier->data(),
1977  _perp_fourier->data() +
1978  _perp_fourier->num_elements());
1979  }
void subtract_if_necessary(void)
boost::multi_array< double, 4 > * _perp_fourier
void VoronoiAlgorithm::push_back_particle ( const double  perp,
const double  pseudorapidity,
const double  azimuth,
const unsigned int  reduced_particle_flow_id 
)

Add a new unsubtracted particle to the current event

Parameters
[in]perptransverse momentum
[in]pseudorapiditypseudorapidity
[in]azimuthazimuth
[in]reduced_particle_flow_idreduced particle flow ID, between 0 and 2 (inclusive)

Definition at line 1869 of file VoronoiAlgorithm.cc.

References _event, and AlCaHLTBitMon_ParallelJobs::p.

Referenced by VoronoiBackgroundProducer::produce().

1873  {
1874  math::PtEtaPhiELorentzVector p(perp, pseudorapidity, azimuth, NAN);
1875 
1876  p.SetE(p.P());
1877  _event.push_back(particle_t(p, reduced_particle_flow_id));
1878  }
PtEtaPhiELorentzVectorD PtEtaPhiELorentzVector
Lorentz vector with cartesian internal representation.
Definition: LorentzVector.h:27
T perp() const
Magnitude of transverse component.
void VoronoiAlgorithm::recombine_link ( void  )
private

Definition at line 1274 of file VoronoiAlgorithm.cc.

References _active, _event, _radial_distance_square_max, _recombine, _recombine_index, _recombine_tie, _recombine_unsigned, begin, i, j, and python.multivaluedict::sort().

Referenced by equalize().

1275  {
1276  boost::multi_array<double, 2> radial_distance_square(
1277  boost::extents[_event.size()][_event.size()]);
1278 
1279  for (std::vector<particle_t>::const_iterator
1280  iterator_outer = _event.begin();
1281  iterator_outer != _event.end(); iterator_outer++) {
1282  radial_distance_square
1283  [iterator_outer - _event.begin()]
1284  [iterator_outer - _event.begin()] = 0;
1285 
1286  for (std::vector<particle_t>::const_iterator
1287  iterator_inner = _event.begin();
1288  iterator_inner != iterator_outer;
1289  iterator_inner++) {
1290  const double deta = iterator_outer->momentum.Eta() -
1291  iterator_inner->momentum.Eta();
1292  const double dphi = angular_range_reduce(
1293  iterator_outer->momentum.Phi() -
1294  iterator_inner->momentum.Phi());
1295 
1296  radial_distance_square
1297  [iterator_outer - _event.begin()]
1298  [iterator_inner - _event.begin()] =
1299  deta * deta + dphi * dphi;
1300  radial_distance_square
1301  [iterator_inner - _event.begin()]
1302  [iterator_outer - _event.begin()] =
1303  radial_distance_square
1304  [iterator_outer - _event.begin()]
1305  [iterator_inner - _event.begin()];
1306  }
1307  }
1308 
1309  _active.clear();
1310 
1311  for (std::vector<particle_t>::const_iterator
1312  iterator_outer = _event.begin();
1313  iterator_outer != _event.end(); iterator_outer++) {
1314  double incident_area_sum = iterator_outer->area;
1315 
1316  for (std::set<std::vector<particle_t>::iterator>::
1317  const_iterator iterator_inner =
1318  iterator_outer->incident.begin();
1319  iterator_inner !=
1320  iterator_outer->incident.end();
1321  iterator_inner++) {
1322  incident_area_sum += (*iterator_inner)->area;
1323  }
1324  _active.push_back(incident_area_sum < 2.0);
1325  }
1326 
1327  _recombine.clear();
1328  _recombine_index = std::vector<std::vector<size_t> >(
1329  _event.size(), std::vector<size_t>());
1330  _recombine_unsigned = std::vector<std::vector<size_t> >(
1331  _event.size(), std::vector<size_t>());
1332  _recombine_tie.clear();
1333 
1334  // 36 cells corresponds to ~ 3 layers, note that for
1335  // hexagonal tiling, cell in proximity = 3 * layer *
1336  // (layer + 1)
1337  static const size_t npair_max = 36;
1338 
1339  for (size_t i = 0; i < _event.size(); i++) {
1340  for (size_t j = 0; j < _event.size(); j++) {
1341  const bool active_i_j = _active[i] && _active[j];
1342  const size_t incident_count =
1343  _event[i].incident.count(_event.begin() + j) +
1344  _event[j].incident.count(_event.begin() + i);
1345 
1346  if (active_i_j &&
1347  (radial_distance_square[i][j] <
1349  incident_count > 0)) {
1350  _recombine_unsigned[i].push_back(j);
1351  }
1352  }
1353 
1354  if (_event[i].momentum_perp_subtracted < 0) {
1355  std::vector<double> radial_distance_square_list;
1356 
1357  for (std::vector<size_t>::const_iterator iterator =
1359  iterator != _recombine_unsigned[i].end();
1360  iterator++) {
1361  const size_t j = *iterator;
1362 
1363  if (_event[j].momentum_perp_subtracted > 0) {
1364  radial_distance_square_list.push_back(
1365  radial_distance_square[i][j]);
1366  }
1367  }
1368 
1369  double radial_distance_square_max_equalization_cut =
1371 
1372  if (radial_distance_square_list.size() >= npair_max) {
1373  std::sort(radial_distance_square_list.begin(),
1374  radial_distance_square_list.end());
1375  radial_distance_square_max_equalization_cut =
1376  radial_distance_square_list[npair_max - 1];
1377  }
1378 
1379  for (std::vector<size_t>::const_iterator iterator =
1381  iterator != _recombine_unsigned[i].end();
1382  iterator++) {
1383  const size_t j = *iterator;
1384 
1385  if (_event[j].momentum_perp_subtracted > 0 &&
1386  radial_distance_square[i][j] <
1387  radial_distance_square_max_equalization_cut) {
1388  _recombine_index[j].push_back(
1389  _recombine.size());
1390  _recombine_index[i].push_back(
1391  _recombine.size());
1392  _recombine.push_back(
1393  std::pair<size_t, size_t>(i, j));
1394  _recombine_tie.push_back(
1395  radial_distance_square[i][j] /
1397  }
1398  }
1399  }
1400  }
1401  }
std::vector< bool > _active
std::vector< std::vector< size_t > > _recombine_index
int i
Definition: DBlmapReader.cc:9
std::vector< double > _recombine_tie
std::vector< std::pair< size_t, size_t > > _recombine
int j
Definition: DBlmapReader.cc:9
double _radial_distance_square_max
std::vector< std::vector< size_t > > _recombine_unsigned
#define begin
Definition: vmac.h:30
void VoronoiAlgorithm::remove_nonpositive ( void  )
private

Definition at line 1806 of file VoronoiAlgorithm.cc.

References _event, and bookConverter::max.

Referenced by subtract_if_necessary().

1807  {
1808  for (std::vector<particle_t>::iterator iterator =
1809  _event.begin();
1810  iterator != _event.end(); iterator++) {
1811  iterator->momentum_perp_subtracted = std::max(
1812  0.0, iterator->momentum_perp_subtracted);
1813  }
1814  }
void VoronoiAlgorithm::subtract_if_necessary ( void  )
private
void VoronoiAlgorithm::subtract_momentum ( void  )
private

Definition at line 1124 of file VoronoiAlgorithm.cc.

References _cms_ecal_edge_pseudorapidity, _cms_hcal_edge_pseudorapidity, _edge_pseudorapidity, _event, _feature, funct::cos(), create_public_lumi_plots::exp, j, prof2calltree::l, visualization-live-secondInstance_cfg::m, M_PI, gen::n, nfourier, nreduced_particle_flow_id, python.connectstrParser::o, AlCaHLTBitMon_ParallelJobs::p, funct::sin(), ue, UECalibration::ue_interpolation_pf0, UECalibration::ue_interpolation_pf1, UECalibration::ue_interpolation_pf2, and UECalibration::ue_predictor_pf.

Referenced by subtract_if_necessary().

1125  {
1126  for (std::vector<particle_t>::iterator iterator =
1127  _event.begin();
1128  iterator != _event.end(); iterator++) {
1129  int predictor_index = -1;
1130  int interpolation_index = -1;
1131  double density = 0;
1132  double pred_0 = 0;
1133 
1134  for (size_t l = 1; l < _edge_pseudorapidity.size(); l++) {
1135  if (iterator->momentum.Eta() >=
1136  _edge_pseudorapidity[l - 1] &&
1137  iterator->momentum.Eta() <
1139  predictor_index = l - 1;
1140  }
1141  }
1142 
1143  for (size_t j = 0; j < nreduced_particle_flow_id; j++) {
1144  if (j == 2) {
1145  // HCAL
1146  for (size_t l = 1;
1147  l < _cms_hcal_edge_pseudorapidity.size(); l++) {
1148  if (iterator->momentum.Eta() >=
1150  iterator->momentum.Eta() <
1152  interpolation_index = l - 1;
1153  }
1154  }
1155  }
1156  else {
1157  // Tracks or ECAL clusters
1158  for (size_t l = 1;
1159  l < _cms_ecal_edge_pseudorapidity.size(); l++) {
1160  if (iterator->momentum.Eta() >=
1162  iterator->momentum.Eta() <
1164  interpolation_index = l - 1;
1165  }
1166  }
1167  }
1168 
1169  if (predictor_index >= 0 && interpolation_index >= 0) {
1170  // Calculate the aggregated prediction and
1171  // interpolation for the pseudorapidity segment
1172 
1173  const double azimuth = iterator->momentum.Phi();
1174  const float (*p)[2][82] =
1175 #ifdef STANDALONE
1176  ue_predictor_pf[j][predictor_index]
1177 #else // STANDALONE
1178  ue->ue_predictor_pf[j][predictor_index]
1179 #endif // STANDALONE
1180  ;
1181  double pred = 0;
1182 
1183  for (size_t l = 0; l < nfourier; l++) {
1184  const size_t norder = l == 0 ? 9 : 1;
1185 
1186  for (size_t m = 0; m < 2; m++) {
1187  float u = p[l][m][0];
1188 
1189  for (size_t n = 0; n < 2 * nfourier - 1; n++) {
1190  if ((l == 0 && n == 0) || (l == 2 && (n == 3 || n == 4))) {
1191  u += p[l][m][9 * n + 1] * _feature[n];
1192  for (size_t o = 2; o < norder + 1; o++) {
1193  u += p[l][m][9 * n + o] *
1194  hermite_h_normalized(
1195  2 * o - 1, _feature[n]) *
1196  exp(-_feature[n] * _feature[n]);
1197  }
1198  }
1199  }
1200 
1201  pred += u * (l == 0 ? 1.0 : 2.0) *
1202  (m == 0 ? cos(l * azimuth) :
1203  sin(l * azimuth));
1204  if (l == 0 && m == 0) {
1205  pred_0 += u /
1206  (2.0 * M_PI *
1207  (_edge_pseudorapidity[predictor_index + 1] -
1208  _edge_pseudorapidity[predictor_index]));
1209  }
1210  }
1211  }
1212 
1213  double interp = 0;
1214 
1215 #ifdef STANDALONE
1216  if (j == 0) {
1217  interp =
1218  ue_interpolation_pf0[predictor_index][
1219  interpolation_index];
1220  }
1221  else if (j == 1) {
1222  interp =
1223  ue_interpolation_pf1[predictor_index][
1224  interpolation_index];
1225  }
1226  else if (j == 2) {
1227  interp =
1228  ue_interpolation_pf2[predictor_index][
1229  interpolation_index];
1230  }
1231 #else // STANDALONE
1232  if (j == 0) {
1233  interp =
1234  ue->ue_interpolation_pf0[predictor_index][
1235  interpolation_index];
1236  }
1237  else if (j == 1) {
1238  interp =
1239  ue->ue_interpolation_pf1[predictor_index][
1240  interpolation_index];
1241  }
1242  else if (j == 2) {
1243  interp =
1244  ue->ue_interpolation_pf2[predictor_index][
1245  interpolation_index];
1246  }
1247 #endif // STANDALONE
1248  // Interpolate down to the finely binned
1249  // pseudorapidity
1250 
1251  density += pred /
1252  (2.0 * M_PI *
1253  (_edge_pseudorapidity[predictor_index + 1] -
1254  _edge_pseudorapidity[predictor_index])) *
1255  interp;
1256  }
1257  }
1258 
1259  if (std::isfinite(iterator->area) && density >= 0) {
1260  // Subtract the PF candidate by density times
1261  // Voronoi cell area
1262  iterator->momentum_perp_subtracted =
1263  iterator->momentum.Pt() -
1264  density * iterator->area;
1265  }
1266  else {
1267  iterator->momentum_perp_subtracted =
1268  iterator->momentum.Pt();
1269  }
1270  iterator->momentum_perp_subtracted_unequalized =
1271  iterator->momentum_perp_subtracted;
1272  }
1273  }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< double > _edge_pseudorapidity
const UECalibration * ue
std::vector< double > _cms_ecal_edge_pseudorapidity
static const size_t nreduced_particle_flow_id
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static const size_t nfourier
int j
Definition: DBlmapReader.cc:9
std::vector< double > _feature
#define M_PI
float ue_interpolation_pf0[15][344]
float ue_predictor_pf[3][15][5][2][82]
std::vector< double > _cms_hcal_edge_pseudorapidity
float ue_interpolation_pf1[15][344]
float ue_interpolation_pf2[15][82]
std::vector< double > VoronoiAlgorithm::subtracted_equalized_perp ( void  )

Returns the transverse momenta of the subtracted particles

Returns
vector of transverse momenta

Definition at line 1892 of file VoronoiAlgorithm.cc.

References _event, run_regression::ret, and subtract_if_necessary().

Referenced by VoronoiBackgroundProducer::produce().

1893  {
1895 
1896  std::vector<double> ret;
1897 
1898  for (std::vector<particle_t>::const_iterator iterator =
1899  _event.begin();
1900  iterator != _event.end(); iterator++) {
1901  ret.push_back(iterator->momentum_perp_subtracted);
1902  }
1903 
1904  return ret;
1905  }
void subtract_if_necessary(void)
std::vector< double > VoronoiAlgorithm::subtracted_unequalized_perp ( void  )

Definition at line 1906 of file VoronoiAlgorithm.cc.

References _event, run_regression::ret, and subtract_if_necessary().

Referenced by VoronoiBackgroundProducer::produce().

1907  {
1909 
1910  std::vector<double> ret;
1911 
1912  for (std::vector<particle_t>::const_iterator iterator =
1913  _event.begin();
1914  iterator != _event.end(); iterator++) {
1915  ret.push_back(iterator->momentum_perp_subtracted_unequalized);
1916  }
1917 
1918  return ret;
1919  }
void subtract_if_necessary(void)
void VoronoiAlgorithm::voronoi_area_incident ( void  )
private

Definition at line 1033 of file VoronoiAlgorithm.cc.

References _event, patZpeak::handle, edm::eventsetup::heterocontainer::insert(), relval_steps::k, M_PI, NULL, AlCaHLTBitMon_ParallelJobs::p, and query::result.

Referenced by subtract_if_necessary().

1034  {
1035  // Make the Voronoi diagram
1036 
1037  voronoi_diagram_t diagram;
1038 
1039  // Reverse Voronoi face lookup
1040 #ifdef HAVE_SPARSEHASH
1041  // The "empty" or default value of the hash table
1042  const voronoi_diagram_t::Face face_empty;
1043  google::dense_hash_map<voronoi_diagram_t::Face_handle,
1044  size_t, hash<voronoi_diagram_t::Face_handle> >
1045  face_index;
1046 
1047  face_index.set_empty_key(face_empty);
1048 #else // HAVE_SPARSEHASH
1049  std::map<voronoi_diagram_t::Face_handle, size_t>
1050  face_index;
1051 #endif // HAVE_SPARSEHASH
1052 
1053  for (std::vector<particle_t>::const_iterator iterator =
1054  _event.begin();
1055  iterator != _event.end(); iterator++) {
1056  // Make two additional replicas with azimuth +/- 2 pi
1057  // (and use only the middle) to mimick the azimuthal
1058  // cyclicity
1059  for (int k = -1; k <= 1; k++) {
1060  const point_2d_t p(
1061  iterator->momentum.Eta(),
1062  iterator->momentum.Phi() +
1063  k * (2 * M_PI));
1064  const voronoi_diagram_t::Face_handle handle =
1065  diagram.insert(p);
1066 
1067  face_index[handle] = iterator - _event.begin();
1068  }
1069  }
1070 
1071  // Extract the Voronoi cells as polygon and calculate the
1072  // area associated with individual particles
1073 
1074  for (std::vector<particle_t>::iterator iterator =
1075  _event.begin();
1076  iterator != _event.end(); iterator++) {
1077  const voronoi_diagram_t::Locate_result result =
1078  diagram.locate(*iterator);
1079  const voronoi_diagram_t::Face_handle *face =
1080  boost::get<voronoi_diagram_t::Face_handle>(
1081  &result);
1082  double polygon_area;
1083 
1084  if (face != NULL) {
1085  voronoi_diagram_t::Ccb_halfedge_circulator
1086  circulator_start = (*face)->outer_ccb();
1087  bool unbounded = false;
1088  polygon_t polygon;
1089 
1090  voronoi_diagram_t::Ccb_halfedge_circulator
1091  circulator = circulator_start;
1092 
1093  // Circle around the edges and extract the polygon
1094  // vertices
1095  do {
1096  if (circulator->has_target()) {
1097  polygon.push_back(
1098  circulator->target()->point());
1099  _event[face_index[*face]].incident.
1100  insert(
1101  _event.begin() +
1102  face_index[circulator->twin()->
1103  face()]);
1104  }
1105  else {
1106  unbounded = true;
1107  break;
1108  }
1109  }
1110  while (++circulator != circulator_start);
1111  if (unbounded) {
1112  polygon_area = INFINITY;
1113  }
1114  else {
1115  polygon_area = polygon.area();
1116  }
1117  }
1118  else {
1119  polygon_area = NAN;
1120  }
1121  iterator->area = fabs(polygon_area);
1122  }
1123  }
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
#define NULL
Definition: scimark2.h:8
delaunay_triangulation_t::Point point_2d_t
tuple result
Definition: query.py:137
tuple handle
Definition: patZpeak.py:22
bool insert(Storage &iStorage, ItemType *iItem, const IdTag &iIdTag)
Definition: HCMethods.h:49
#define M_PI
CGAL::Polygon_2< CGAL::Exact_predicates_inexact_constructions_kernel > polygon_t

Member Data Documentation

std::vector<bool> VoronoiAlgorithm::_active
protected

Definition at line 76 of file VoronoiAlgorithm.h.

Referenced by lp_populate(), and recombine_link().

std::vector<double> VoronoiAlgorithm::_cms_ecal_edge_pseudorapidity
protected

Definition at line 67 of file VoronoiAlgorithm.h.

Referenced by initialize_geometry(), and subtract_momentum().

std::vector<double> VoronoiAlgorithm::_cms_hcal_edge_pseudorapidity
protected

Definition at line 66 of file VoronoiAlgorithm.h.

Referenced by initialize_geometry(), and subtract_momentum().

std::vector<double> VoronoiAlgorithm::_edge_pseudorapidity
protected
std::pair<double, double> VoronoiAlgorithm::_equalization_threshold
protected

Definition at line 69 of file VoronoiAlgorithm.h.

Referenced by lp_populate().

event_t VoronoiAlgorithm::_event
protected
std::vector<double> VoronoiAlgorithm::_feature
protected

Definition at line 75 of file VoronoiAlgorithm.h.

Referenced by feature_extract(), and subtract_momentum().

void* VoronoiAlgorithm::_lp_environment
protected

Definition at line 83 of file VoronoiAlgorithm.h.

Referenced by equalize().

void* VoronoiAlgorithm::_lp_problem
protected

Definition at line 84 of file VoronoiAlgorithm.h.

std::vector<size_t> VoronoiAlgorithm::_nblock_subtract
protected

Definition at line 82 of file VoronoiAlgorithm.h.

Referenced by equalize(), and lp_populate().

size_t VoronoiAlgorithm::_ncost
protected

Definition at line 81 of file VoronoiAlgorithm.h.

Referenced by equalize(), and lp_populate().

boost::multi_array<double, 4>* VoronoiAlgorithm::_perp_fourier
protected

Definition at line 74 of file VoronoiAlgorithm.h.

Referenced by allocate(), deallocate(), event_fourier(), and perp_fourier().

double VoronoiAlgorithm::_positive_bound_scale
protected

Definition at line 71 of file VoronoiAlgorithm.h.

Referenced by lp_populate().

double VoronoiAlgorithm::_radial_distance_square_max
protected

Definition at line 70 of file VoronoiAlgorithm.h.

Referenced by recombine_link().

std::vector<std::pair<size_t, size_t> > VoronoiAlgorithm::_recombine
protected

Definition at line 77 of file VoronoiAlgorithm.h.

Referenced by equalize(), and recombine_link().

std::vector<std::vector<size_t> > VoronoiAlgorithm::_recombine_index
protected

Definition at line 78 of file VoronoiAlgorithm.h.

Referenced by lp_populate(), and recombine_link().

std::vector<double> VoronoiAlgorithm::_recombine_tie
protected

Definition at line 80 of file VoronoiAlgorithm.h.

Referenced by lp_populate(), and recombine_link().

std::vector<std::vector<size_t> > VoronoiAlgorithm::_recombine_unsigned
protected

Definition at line 79 of file VoronoiAlgorithm.h.

Referenced by lp_populate(), and recombine_link().

bool VoronoiAlgorithm::_remove_nonpositive
protected

Definition at line 68 of file VoronoiAlgorithm.h.

Referenced by subtract_if_necessary().

bool VoronoiAlgorithm::_subtracted
protected

Definition at line 72 of file VoronoiAlgorithm.h.

Referenced by clear(), and subtract_if_necessary().

const size_t VoronoiAlgorithm::nfourier = 5
static

Definition at line 63 of file VoronoiAlgorithm.h.

Referenced by allocate(), event_fourier(), feature_extract(), and subtract_momentum().

const size_t VoronoiAlgorithm::nreduced_particle_flow_id = 3
static

Definition at line 62 of file VoronoiAlgorithm.h.

Referenced by allocate(), feature_extract(), and subtract_momentum().

const UECalibration* VoronoiAlgorithm::ue
protected

Definition at line 86 of file VoronoiAlgorithm.h.

Referenced by subtract_momentum().