CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
KDTreeLinkerTrackEcal Class Reference

#include <KDTreeLinkerTrackEcal.h>

Inheritance diagram for KDTreeLinkerTrackEcal:
KDTreeLinkerBase

Public Member Functions

void buildTree () override
 
void clear () override
 
void insertFieldClusterElt (reco::PFBlockElement *ecalCluster) override
 
void insertTargetElt (reco::PFBlockElement *track) override
 
 KDTreeLinkerTrackEcal ()
 
void searchLinks () override
 
void updatePFBlockEltWithLinks () override
 
 ~KDTreeLinkerTrackEcal () override
 
- Public Member Functions inherited from KDTreeLinkerBase
const reco::PFBlockElement::TypefieldType () const
 
float getCristalPhiEtaMaxSize () const
 
float getCristalXYMaxSize () const
 
float getPhiOffset () const
 
 KDTreeLinkerBase ()
 
virtual void process ()
 
void setCristalPhiEtaMaxSize (float size)
 
void setCristalXYMaxSize (float size)
 
void setDebug (bool isDebug)
 
void setFieldType (const reco::PFBlockElement::Type &fld)
 
void setPhiOffset (double phiOffset)
 
void setTargetType (const reco::PFBlockElement::Type &tgt)
 
const reco::PFBlockElement::TypetargetType () const
 
virtual ~KDTreeLinkerBase ()
 

Private Attributes

BlockEltSet fieldClusterSet_
 
RecHit2BlockEltMap rechit2ClusterLinks_
 
RecHitSet rechitsSet_
 
BlockElt2BlockEltMap target2ClusterLinks_
 
BlockEltSet targetSet_
 
KDTreeLinkerAlgo tree_
 

Additional Inherited Members

- Protected Attributes inherited from KDTreeLinkerBase
reco::PFBlockElement::Type _fieldType
 
reco::PFBlockElement::Type _targetType
 
float cristalPhiEtaMaxSize_
 
float cristalXYMaxSize_
 
bool debug_
 
float phiOffset_
 

Detailed Description

Definition at line 12 of file KDTreeLinkerTrackEcal.h.

Constructor & Destructor Documentation

KDTreeLinkerTrackEcal::KDTreeLinkerTrackEcal ( )

Definition at line 13 of file KDTreeLinkerTrackEcal.cc.

15 {}
KDTreeLinkerTrackEcal::~KDTreeLinkerTrackEcal ( )
override

Definition at line 17 of file KDTreeLinkerTrackEcal.cc.

References clear().

18 {
19  clear();
20 }

Member Function Documentation

void KDTreeLinkerTrackEcal::buildTree ( )
overridevirtual

Implements KDTreeLinkerBase.

Definition at line 63 of file KDTreeLinkerTrackEcal.cc.

References KDTreeLinkerAlgo< DATA, DIM >::build(), RhoEtaPhi::eta(), KDTreeLinkerBase::getPhiOffset(), M_PI, phi, RhoEtaPhi::phi(), phimax, phimin, rechitsSet_, and tree_.

64 {
65  // List of pseudo-rechits that will be used to create the KDTree
66  std::vector<KDTreeNodeInfo> eltList;
67 
68  // Filling of this list
69  for(RecHitSet::const_iterator it = rechitsSet_.begin();
70  it != rechitsSet_.end(); it++) {
71 
72  const reco::PFRecHit::REPPoint &posrep = (*it)->positionREP();
73 
74  KDTreeNodeInfo rh1 (*it, posrep.eta(), posrep.phi());
75  eltList.push_back(rh1);
76 
77  // Here we solve the problem of phi circular set by duplicating some rechits
78  // too close to -Pi (or to Pi) and adding (substracting) to them 2 * Pi.
79  if (rh1.dim2 > (M_PI - getPhiOffset())) {
80  double phi = rh1.dim2 - 2 * M_PI;
81  KDTreeNodeInfo rh2(*it, posrep.eta(), phi);
82  eltList.push_back(rh2);
83  }
84 
85  if (rh1.dim2 < (M_PI * -1.0 + getPhiOffset())) {
86  double phi = rh1.dim2 + 2 * M_PI;
87  KDTreeNodeInfo rh3(*it, posrep.eta(), phi);
88  eltList.push_back(rh3);
89  }
90  }
91 
92  // Here we define the upper/lower bounds of the 2D space (eta/phi).
93  double phimin = -1.0 * M_PI - getPhiOffset();
94  double phimax = M_PI + getPhiOffset();
95 
96  // etamin-etamax, phimin-phimax
97  KDTreeBox region(-3.0, 3.0, phimin, phimax);
98 
99  // We may now build the KDTree
100  tree_.build(eltList, region);
101 }
float phi() const
momentum azimuthal angle
Definition: PtEtaPhiMass.h:53
void build(std::vector< KDTreeNodeInfoT< DATA, DIM > > &eltList, const KDTreeBoxT< DIM > &region)
#define M_PI
float eta() const
momentum pseudorapidity
Definition: PtEtaPhiMass.h:51
float getPhiOffset() const
void KDTreeLinkerTrackEcal::clear ( void  )
overridevirtual

Implements KDTreeLinkerBase.

Definition at line 238 of file KDTreeLinkerTrackEcal.cc.

References KDTreeLinkerAlgo< DATA, DIM >::clear(), fieldClusterSet_, rechit2ClusterLinks_, rechitsSet_, target2ClusterLinks_, targetSet_, and tree_.

Referenced by ~KDTreeLinkerTrackEcal().

239 {
240  targetSet_.clear();
241  fieldClusterSet_.clear();
242 
243  rechitsSet_.clear();
244 
245  rechit2ClusterLinks_.clear();
246  target2ClusterLinks_.clear();
247 
248  tree_.clear();
249 }
RecHit2BlockEltMap rechit2ClusterLinks_
BlockElt2BlockEltMap target2ClusterLinks_
void KDTreeLinkerTrackEcal::insertFieldClusterElt ( reco::PFBlockElement ecalCluster)
overridevirtual

Implements KDTreeLinkerBase.

Definition at line 32 of file KDTreeLinkerTrackEcal.cc.

References reco::PFBlockElement::clusterRef(), fieldClusterSet_, funct::fract(), dedxEstimators_cff::fraction, edm::Ref< C, T, F >::isNull(), rechit2ClusterLinks_, and rechitsSet_.

33 {
34  const reco::PFClusterRef& clusterref = ecalCluster->clusterRef();
35 
36  // This test is more or less done in PFBlockAlgo.h. In others cases, it should be switch on.
37  // if (!((clusterref->layer() == PFLayer::ECAL_ENDCAP) ||
38  // (clusterref->layer() == PFLayer::ECAL_BARREL)))
39  // return;
40 
41  const std::vector<reco::PFRecHitFraction> &fraction = clusterref->recHitFractions();
42 
43  // We create a list of ecalCluster
44  fieldClusterSet_.insert(ecalCluster);
45  for(size_t rhit = 0; rhit < fraction.size(); ++rhit) {
46  const reco::PFRecHitRef& rh = fraction[rhit].recHitRef();
47  double fract = fraction[rhit].fraction();
48 
49  if ((rh.isNull()) || (fract < 1E-4))
50  continue;
51 
52  const reco::PFRecHit& rechit = *rh;
53 
54  // We save the links rechit to EcalClusters
55  rechit2ClusterLinks_[&rechit].insert(ecalCluster);
56 
57  // We create a liste of rechits
58  rechitsSet_.insert(&rechit);
59  }
60 }
virtual const PFClusterRef & clusterRef() const
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:32
const Fraction< n, m >::type & fract()
Definition: Fraction.h:38
bool isNull() const
Checks for null.
Definition: Ref.h:250
RecHit2BlockEltMap rechit2ClusterLinks_
void KDTreeLinkerTrackEcal::insertTargetElt ( reco::PFBlockElement track)
overridevirtual

Implements KDTreeLinkerBase.

Definition at line 23 of file KDTreeLinkerTrackEcal.cc.

References reco::PFTrajectoryPoint::ECALShowerMax, targetSet_, and reco::PFBlockElement::trackRefPF().

24 {
25  if( track->trackRefPF()->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() ) {
26  targetSet_.insert(track);
27  }
28 }
virtual const PFRecTrackRef & trackRefPF() const
void KDTreeLinkerTrackEcal::searchLinks ( )
overridevirtual

Implements KDTreeLinkerBase.

Definition at line 104 of file KDTreeLinkerTrackEcal.cc.

References reco::PFTrajectoryPoint::ClosestApproach, PFLayer::ECAL_BARREL, reco::PFTrajectoryPoint::ECALShowerMax, PVValHelper::eta, KDTreeLinkerBase::getCristalPhiEtaMaxSize(), reco::PFTrajectoryPoint::isValid(), M_PI, min(), reco::PFTrajectoryPoint::momentum(), phi, reco::PFTrajectoryPoint::position(), reco::PFTrajectoryPoint::positionREP(), rechit2ClusterLinks_, KDTreeLinkerAlgo< DATA, DIM >::search(), mathSSE::sqrt(), target2ClusterLinks_, targetSet_, listHistos::trackPt, tree_, x, and y.

105 {
106  // Must of the code has been taken from LinkByRecHit.cc
107 
108  // We iterate over the tracks.
109  for(BlockEltSet::iterator it = targetSet_.begin();
110  it != targetSet_.end(); it++) {
111 
112  reco::PFRecTrackRef trackref = (*it)->trackRefPF();
113 
114  // We set the multilinks flag of the track to true. It will allow us to
115  // use in an optimized way our algo results in the recursive linking algo.
116  (*it)->setIsValidMultilinks(true);
117 
118  const reco::PFTrajectoryPoint& atECAL =
119  trackref->extrapolatedPoint(reco::PFTrajectoryPoint::ECALShowerMax);
120 
121  // The track didn't reach ecal
122  if( ! atECAL.isValid() ) continue;
123 
124  const reco::PFTrajectoryPoint& atVertex =
125  trackref->extrapolatedPoint( reco::PFTrajectoryPoint::ClosestApproach );
126 
127  double trackPt = sqrt(atVertex.momentum().Vect().Perp2());
128  double tracketa = atECAL.positionREP().eta();
129  double trackphi = atECAL.positionREP().phi();
130  double trackx = atECAL.position().X();
131  double tracky = atECAL.position().Y();
132  double trackz = atECAL.position().Z();
133 
134  // Estimate the maximal envelope in phi/eta that will be used to find rechit candidates.
135  // Same envelope for cap et barrel rechits.
136  double range = getCristalPhiEtaMaxSize() * (2.0 + 1.0 / std::min(1., trackPt / 2.));
137 
138  // We search for all candidate recHits, ie all recHits contained in the maximal size envelope.
139  std::vector<KDTreeNodeInfo> recHits;
140  KDTreeBox trackBox(tracketa-range, tracketa+range, trackphi-range, trackphi+range);
141  tree_.search(trackBox, recHits);
142 
143  // Here we check all rechit candidates using the non-approximated method.
144  for(std::vector<KDTreeNodeInfo>::const_iterator rhit = recHits.begin();
145  rhit != recHits.end(); ++rhit) {
146 
147  const auto & cornersxyz = rhit->ptr->getCornersXYZ();
148  const auto & posxyz = rhit->ptr->position();
149  const auto &rhrep = rhit->ptr->positionREP();
150  const auto & corners = rhit->ptr->getCornersREP();
151 
152  double rhsizeeta = fabs(corners[3].eta() - corners[1].eta());
153  double rhsizephi = fabs(corners[3].phi() - corners[1].phi());
154  if ( rhsizephi > M_PI ) rhsizephi = 2.*M_PI - rhsizephi;
155 
156  double deta = fabs(rhrep.eta() - tracketa);
157  double dphi = fabs(rhrep.phi() - trackphi);
158  if ( dphi > M_PI ) dphi = 2.*M_PI - dphi;
159 
160  // Find all clusters associated to given rechit
161  RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(rhit->ptr);
162 
163  for(BlockEltSet::const_iterator clusterIt = ret->second.begin();
164  clusterIt != ret->second.end(); clusterIt++) {
165 
166  reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
167  double clusterz = clusterref->position().z();
168  int fracsNbr = clusterref->recHitFractions().size();
169 
170  if (clusterref->layer() == PFLayer::ECAL_BARREL){ // BARREL
171  // Check if the track is in the barrel
172  if (fabs(trackz) > 300.) continue;
173 
174  double _rhsizeeta = rhsizeeta * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
175  double _rhsizephi = rhsizephi * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
176 
177  // Check if the track and the cluster are linked
178  if(deta < (_rhsizeeta / 2.) && dphi < (_rhsizephi / 2.))
179  target2ClusterLinks_[*it].insert(*clusterIt);
180 
181 
182  } else { // ENDCAP
183 
184  // Check if the track is in the cap
185  if (fabs(trackz) < 300.) continue;
186  if (trackz * clusterz < 0.) continue;
187 
188  double x[5];
189  double y[5];
190  for ( unsigned jc=0; jc<4; ++jc ) {
191  auto cornerposxyz = cornersxyz[jc];
192  x[3-jc] = cornerposxyz.x() + (cornerposxyz.x()-posxyz.x())
193  * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
194  y[3-jc] = cornerposxyz.y() + (cornerposxyz.y()-posxyz.y())
195  * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
196  }
197 
198  x[4] = x[0];
199  y[4] = y[0];
200 
201  bool isinside = TMath::IsInside(trackx,
202  tracky,
203  5,x,y);
204 
205  // Check if the track and the cluster are linked
206  if( isinside )
207  target2ClusterLinks_[*it].insert(*clusterIt);
208  }
209  }
210  }
211  }
212 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
const math::XYZPoint & position() const
cartesian position (x, y, z)
void search(const KDTreeBoxT< DIM > &searchBox, std::vector< KDTreeNodeInfoT< DATA, DIM > > &resRecHitList)
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
T sqrt(T t)
Definition: SSEVec.h:18
T min(T a, T b)
Definition: MathUtil.h:58
#define M_PI
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
bool isValid() const
is this point valid ?
float getCristalPhiEtaMaxSize() const
RecHit2BlockEltMap rechit2ClusterLinks_
BlockElt2BlockEltMap target2ClusterLinks_
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
void KDTreeLinkerTrackEcal::updatePFBlockEltWithLinks ( )
overridevirtual

Implements KDTreeLinkerBase.

Definition at line 215 of file KDTreeLinkerTrackEcal.cc.

References reco::PFMultiLinksTC::linkedClusters, and target2ClusterLinks_.

216 {
217  //TODO YG : Check if cluster positionREP() is valid ?
218 
219  // Here we save in each track the list of phi/eta values of linked clusters.
220  for (BlockElt2BlockEltMap::iterator it = target2ClusterLinks_.begin();
221  it != target2ClusterLinks_.end(); ++it) {
222  reco::PFMultiLinksTC multitracks(true);
223 
224  for (BlockEltSet::iterator jt = it->second.begin();
225  jt != it->second.end(); ++jt) {
226 
227  double clusterphi = (*jt)->clusterRef()->positionREP().phi();
228  double clustereta = (*jt)->clusterRef()->positionREP().eta();
229 
230  multitracks.linkedClusters.push_back(std::make_pair(clusterphi, clustereta));
231  }
232 
233  it->first->setMultilinks(multitracks);
234  }
235 }
BlockElt2BlockEltMap target2ClusterLinks_

Member Data Documentation

BlockEltSet KDTreeLinkerTrackEcal::fieldClusterSet_
private

Definition at line 45 of file KDTreeLinkerTrackEcal.h.

Referenced by clear(), and insertFieldClusterElt().

RecHit2BlockEltMap KDTreeLinkerTrackEcal::rechit2ClusterLinks_
private

Definition at line 54 of file KDTreeLinkerTrackEcal.h.

Referenced by clear(), insertFieldClusterElt(), and searchLinks().

RecHitSet KDTreeLinkerTrackEcal::rechitsSet_
private

Definition at line 48 of file KDTreeLinkerTrackEcal.h.

Referenced by buildTree(), clear(), and insertFieldClusterElt().

BlockElt2BlockEltMap KDTreeLinkerTrackEcal::target2ClusterLinks_
private

Definition at line 51 of file KDTreeLinkerTrackEcal.h.

Referenced by clear(), searchLinks(), and updatePFBlockEltWithLinks().

BlockEltSet KDTreeLinkerTrackEcal::targetSet_
private

Definition at line 44 of file KDTreeLinkerTrackEcal.h.

Referenced by clear(), insertTargetElt(), and searchLinks().

KDTreeLinkerAlgo KDTreeLinkerTrackEcal::tree_
private

Definition at line 57 of file KDTreeLinkerTrackEcal.h.

Referenced by buildTree(), clear(), and searchLinks().