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
 
void process ()
 
void setDebug (bool isDebug)
 
void setFieldType (const reco::PFBlockElement::Type &fld)
 
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< reco::PFRecHit const * > tree_
 

Additional Inherited Members

- Protected Attributes inherited from KDTreeLinkerBase
reco::PFBlockElement::Type _fieldType
 
reco::PFBlockElement::Type _targetType
 
float cristalPhiEtaMaxSize_ = 0.04
 
float cristalXYMaxSize_ = 3.
 
bool debug_ = false
 
float phiOffset_ = 0.25
 

Detailed Description

Definition at line 12 of file KDTreeLinkerTrackEcal.h.

Constructor & Destructor Documentation

KDTreeLinkerTrackEcal::KDTreeLinkerTrackEcal ( )

Definition at line 13 of file KDTreeLinkerTrackEcal.cc.

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 >::build(), RhoEtaPhi::eta(), M_PI, phi, RhoEtaPhi::phi(), phimax, phimin, KDTreeLinkerBase::phiOffset_, rechitsSet_, and tree_.

64 {
65  // List of pseudo-rechits that will be used to create the KDTree
66  std::vector<KDTreeNodeInfo<reco::PFRecHit const*>> 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<reco::PFRecHit const*> 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.dim[1] > (M_PI - phiOffset_)) {
80  double phi = rh1.dim[1] - 2 * M_PI;
81  KDTreeNodeInfo rh2(*it, posrep.eta(), phi);
82  eltList.push_back(rh2);
83  }
84 
85  if (rh1.dim[1] < (M_PI * -1.0 + phiOffset_)) {
86  double phi = rh1.dim[1] + 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 - phiOffset_;
94  double phimax = M_PI + phiOffset_;
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
#define M_PI
void build(std::vector< KDTreeNodeInfo< DATA > > &eltList, const KDTreeBox &region)
KDTreeLinkerAlgo< reco::PFRecHit const * > tree_
float eta() const
momentum pseudorapidity
Definition: PtEtaPhiMass.h:51
void KDTreeLinkerTrackEcal::clear ( void  )
overridevirtual

Implements KDTreeLinkerBase.

Definition at line 237 of file KDTreeLinkerTrackEcal.cc.

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

Referenced by ~KDTreeLinkerTrackEcal().

238 {
239  targetSet_.clear();
240  fieldClusterSet_.clear();
241 
242  rechitsSet_.clear();
243 
244  rechit2ClusterLinks_.clear();
245  target2ClusterLinks_.clear();
246 
247  tree_.clear();
248 }
RecHit2BlockEltMap rechit2ClusterLinks_
BlockElt2BlockEltMap target2ClusterLinks_
KDTreeLinkerAlgo< reco::PFRecHit const * > tree_
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:37
bool isNull() const
Checks for null.
Definition: Ref.h:248
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, KDTreeLinkerBase::cristalPhiEtaMaxSize_, PFLayer::ECAL_BARREL, reco::PFTrajectoryPoint::ECALShowerMax, PVValHelper::eta, reco::PFTrajectoryPoint::isValid(), M_PI, min(), reco::PFTrajectoryPoint::momentum(), phi, reco::PFTrajectoryPoint::position(), reco::PFTrajectoryPoint::positionREP(), rpcPointValidation_cfi::recHit, rechit2ClusterLinks_, KDTreeLinkerAlgo< DATA >::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 = cristalPhiEtaMaxSize_ * (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<reco::PFRecHit const*> 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(auto const& recHit : recHits) {
145 
146  const auto & cornersxyz = recHit->getCornersXYZ();
147  const auto & posxyz = recHit->position();
148  const auto &rhrep = recHit->positionREP();
149  const auto & corners = recHit->getCornersREP();
150 
151  double rhsizeeta = fabs(corners[3].eta() - corners[1].eta());
152  double rhsizephi = fabs(corners[3].phi() - corners[1].phi());
153  if ( rhsizephi > M_PI ) rhsizephi = 2.*M_PI - rhsizephi;
154 
155  double deta = fabs(rhrep.eta() - tracketa);
156  double dphi = fabs(rhrep.phi() - trackphi);
157  if ( dphi > M_PI ) dphi = 2.*M_PI - dphi;
158 
159  // Find all clusters associated to given rechit
160  RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(recHit);
161 
162  for(BlockEltSet::const_iterator clusterIt = ret->second.begin();
163  clusterIt != ret->second.end(); clusterIt++) {
164 
165  reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
166  double clusterz = clusterref->position().z();
167  int fracsNbr = clusterref->recHitFractions().size();
168 
169  if (clusterref->layer() == PFLayer::ECAL_BARREL){ // BARREL
170  // Check if the track is in the barrel
171  if (fabs(trackz) > 300.) continue;
172 
173  double _rhsizeeta = rhsizeeta * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
174  double _rhsizephi = rhsizephi * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
175 
176  // Check if the track and the cluster are linked
177  if(deta < (_rhsizeeta / 2.) && dphi < (_rhsizephi / 2.))
178  target2ClusterLinks_[*it].insert(*clusterIt);
179 
180 
181  } else { // ENDCAP
182 
183  // Check if the track is in the cap
184  if (fabs(trackz) < 300.) continue;
185  if (trackz * clusterz < 0.) continue;
186 
187  double x[5];
188  double y[5];
189  for ( unsigned jc=0; jc<4; ++jc ) {
190  auto cornerposxyz = cornersxyz[jc];
191  x[3-jc] = cornerposxyz.x() + (cornerposxyz.x()-posxyz.x())
192  * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
193  y[3-jc] = cornerposxyz.y() + (cornerposxyz.y()-posxyz.y())
194  * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
195  }
196 
197  x[4] = x[0];
198  y[4] = y[0];
199 
200  bool isinside = TMath::IsInside(trackx,
201  tracky,
202  5,x,y);
203 
204  // Check if the track and the cluster are linked
205  if( isinside )
206  target2ClusterLinks_[*it].insert(*clusterIt);
207  }
208  }
209  }
210  }
211 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
const math::XYZPoint & position() const
cartesian position (x, y, z)
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
T sqrt(T t)
Definition: SSEVec.h:18
void search(const KDTreeBox &searchBox, std::vector< DATA > &resRecHitList)
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 ?
RecHit2BlockEltMap rechit2ClusterLinks_
BlockElt2BlockEltMap target2ClusterLinks_
KDTreeLinkerAlgo< reco::PFRecHit const * > tree_
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 214 of file KDTreeLinkerTrackEcal.cc.

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

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

Definition at line 57 of file KDTreeLinkerTrackEcal.h.

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