CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 ()
 
void clear ()
 
void insertFieldClusterElt (reco::PFBlockElement *ecalCluster)
 
void insertTargetElt (reco::PFBlockElement *track)
 
 KDTreeLinkerTrackEcal ()
 
void searchLinks ()
 
void updatePFBlockEltWithLinks ()
 
 ~KDTreeLinkerTrackEcal ()
 
- Public Member Functions inherited from KDTreeLinkerBase
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 setPhiOffset (double phiOffset)
 
virtual ~KDTreeLinkerBase ()
 

Private Attributes

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

Additional Inherited Members

- Protected Attributes inherited from KDTreeLinkerBase
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 6 of file KDTreeLinkerTrackEcal.cc.

KDTreeLinkerTrackEcal::~KDTreeLinkerTrackEcal ( )

Definition at line 10 of file KDTreeLinkerTrackEcal.cc.

References clear().

11 {
12  clear();
13 }

Member Function Documentation

void KDTreeLinkerTrackEcal::buildTree ( )
virtual

Implements KDTreeLinkerBase.

Definition at line 54 of file KDTreeLinkerTrackEcal.cc.

References KDTreeLinkerAlgo::build(), KDTreeLinkerBase::getPhiOffset(), M_PI, phi, phimax, phimin, rechitsSet_, and tree_.

55 {
56  // List of pseudo-rechits that will be used to create the KDTree
57  std::vector<KDTreeNodeInfo> eltList;
58 
59  // Filling of this list
60  for(RecHitSet::const_iterator it = rechitsSet_.begin();
61  it != rechitsSet_.end(); it++) {
62 
63  const reco::PFRecHit::REPPoint &posrep = (*it)->positionREP();
64 
65  KDTreeNodeInfo rh1 (*it, posrep.Eta(), posrep.Phi());
66  eltList.push_back(rh1);
67 
68  // Here we solve the problem of phi circular set by duplicating some rechits
69  // too close to -Pi (or to Pi) and adding (substracting) to them 2 * Pi.
70  if (rh1.dim2 > (M_PI - getPhiOffset())) {
71  double phi = rh1.dim2 - 2 * M_PI;
72  KDTreeNodeInfo rh2(*it, posrep.Eta(), phi);
73  eltList.push_back(rh2);
74  }
75 
76  if (rh1.dim2 < (M_PI * -1.0 + getPhiOffset())) {
77  double phi = rh1.dim2 + 2 * M_PI;
78  KDTreeNodeInfo rh3(*it, posrep.Eta(), phi);
79  eltList.push_back(rh3);
80  }
81  }
82 
83  // Here we define the upper/lower bounds of the 2D space (eta/phi).
84  double phimin = -1.0 * M_PI - getPhiOffset();
85  double phimax = M_PI + getPhiOffset();
86 
87  // etamin-etamax, phimin-phimax
88  KDTreeBox region(-3.0, 3.0, phimin, phimax);
89 
90  // We may now build the KDTree
91  tree_.build(eltList, region);
92 }
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< Double32_t > > REPPoint
Definition: PFRecHit.h:39
void build(std::vector< KDTreeNodeInfo > &eltList, const KDTreeBox &region)
#define M_PI
Definition: BFit3D.cc:3
float getPhiOffset() const
Definition: DDAxes.h:10
void KDTreeLinkerTrackEcal::clear ( void  )
virtual

Implements KDTreeLinkerBase.

Definition at line 230 of file KDTreeLinkerTrackEcal.cc.

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

Referenced by ~KDTreeLinkerTrackEcal().

231 {
232  targetSet_.clear();
233  fieldClusterSet_.clear();
234 
235  rechitsSet_.clear();
236 
237  rechit2ClusterLinks_.clear();
238  target2ClusterLinks_.clear();
239 
240  tree_.clear();
241 }
RecHit2BlockEltMap rechit2ClusterLinks_
BlockElt2BlockEltMap target2ClusterLinks_
void KDTreeLinkerTrackEcal::insertFieldClusterElt ( reco::PFBlockElement ecalCluster)
virtual

Implements KDTreeLinkerBase.

Definition at line 23 of file KDTreeLinkerTrackEcal.cc.

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

Referenced by PFBlockAlgo::setInput().

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

Implements KDTreeLinkerBase.

Definition at line 16 of file KDTreeLinkerTrackEcal.cc.

References targetSet_.

Referenced by PFBlockAlgo::setInput().

17 {
18  targetSet_.insert(track);
19 }
void KDTreeLinkerTrackEcal::searchLinks ( )
virtual

Implements KDTreeLinkerBase.

Definition at line 95 of file KDTreeLinkerTrackEcal.cc.

References reco::PFTrajectoryPoint::ClosestApproach, PFLayer::ECAL_BARREL, reco::PFTrajectoryPoint::ECALShowerMax, reco::tau::disc::Eta(), KDTreeLinkerBase::getCristalPhiEtaMaxSize(), reco::PFTrajectoryPoint::isValid(), M_PI, min, reco::PFTrajectoryPoint::momentum(), colinearityKinematic::Phi, reco::PFTrajectoryPoint::position(), reco::PFTrajectoryPoint::positionREP(), rechit2ClusterLinks_, run_regression::ret, KDTreeLinkerAlgo::search(), mathSSE::sqrt(), target2ClusterLinks_, targetSet_, tree_, vdt::x, and detailsBasic3DVector::y.

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

Implements KDTreeLinkerBase.

Definition at line 207 of file KDTreeLinkerTrackEcal.cc.

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

208 {
209  //TODO YG : Check if cluster positionREP() is valid ?
210 
211  // Here we save in each track the list of phi/eta values of linked clusters.
212  for (BlockElt2BlockEltMap::iterator it = target2ClusterLinks_.begin();
213  it != target2ClusterLinks_.end(); ++it) {
214  reco::PFMultiLinksTC multitracks(true);
215 
216  for (BlockEltSet::iterator jt = it->second.begin();
217  jt != it->second.end(); ++jt) {
218 
219  double clusterPhi = (*jt)->clusterRef()->positionREP().Phi();
220  double clusterEta = (*jt)->clusterRef()->positionREP().Eta();
221 
222  multitracks.linkedClusters.push_back(std::make_pair(clusterPhi, clusterEta));
223  }
224 
225  it->first->setMultilinks(multitracks);
226  }
227 }
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().