test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
KDTreeLinkerTrackEcal.cc
Go to the documentation of this file.
2 
4 #include "TMath.h"
5 
6 // the text name is different so that we can easily
7 // construct it when calling the factory
10  "KDTreeTrackAndECALLinker");
11 
12 
15 {}
16 
18 {
19  clear();
20 }
21 
22 void
24 {
25  if( track->trackRefPF()->extrapolatedPoint( reco::PFTrajectoryPoint::ECALShowerMax ).isValid() ) {
26  targetSet_.insert(track);
27  }
28 }
29 
30 
31 void
33 {
34  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 }
61 
62 void
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 }
102 
103 void
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 std::vector< math::XYZPoint >& cornersxyz = rhit->ptr->getCornersXYZ();
148  const math::XYZPoint& posxyz = rhit->ptr->position();
149  const reco::PFRecHit::REPPoint &rhrep = rhit->ptr->positionREP();
150  const std::vector<reco::PFRecHit::REPPoint>& corners = rhit->ptr->getCornersREP();
151  if(corners.size() != 4) continue;
152 
153  double rhsizeEta = fabs(corners[0].Eta() - corners[2].Eta());
154  double rhsizePhi = fabs(corners[0].Phi() - corners[2].Phi());
155  if ( rhsizePhi > M_PI ) rhsizePhi = 2.*M_PI - rhsizePhi;
156 
157  double deta = fabs(rhrep.Eta() - tracketa);
158  double dphi = fabs(rhrep.Phi() - trackphi);
159  if ( dphi > M_PI ) dphi = 2.*M_PI - dphi;
160 
161  // Find all clusters associated to given rechit
162  RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(rhit->ptr);
163 
164  for(BlockEltSet::const_iterator clusterIt = ret->second.begin();
165  clusterIt != ret->second.end(); clusterIt++) {
166 
167  reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
168  double clusterz = clusterref->position().Z();
169  int fracsNbr = clusterref->recHitFractions().size();
170 
171  if (clusterref->layer() == PFLayer::ECAL_BARREL){ // BARREL
172  // Check if the track is in the barrel
173  if (fabs(trackz) > 300.) continue;
174 
175  double _rhsizeEta = rhsizeEta * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
176  double _rhsizePhi = rhsizePhi * (2.00 + 1.0 / (fracsNbr * std::min(1.,trackPt/2.)));
177 
178  // Check if the track and the cluster are linked
179  if(deta < (_rhsizeEta / 2.) && dphi < (_rhsizePhi / 2.))
180  target2ClusterLinks_[*it].insert(*clusterIt);
181 
182 
183  } else { // ENDCAP
184 
185  // Check if the track is in the cap
186  if (fabs(trackz) < 300.) continue;
187  if (trackz * clusterz < 0.) continue;
188 
189  double x[5];
190  double y[5];
191  for ( unsigned jc=0; jc<4; ++jc ) {
192  math::XYZPoint cornerposxyz = cornersxyz[jc];
193  x[jc] = cornerposxyz.X() + (cornerposxyz.X()-posxyz.X())
194  * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
195  y[jc] = cornerposxyz.Y() + (cornerposxyz.Y()-posxyz.Y())
196  * (1.00+0.50/fracsNbr /std::min(1.,trackPt/2.));
197  }
198 
199  x[4] = x[0];
200  y[4] = y[0];
201 
202  bool isinside = TMath::IsInside(trackx,
203  tracky,
204  5,x,y);
205 
206  // Check if the track and the cluster are linked
207  if( isinside )
208  target2ClusterLinks_[*it].insert(*clusterIt);
209  }
210  }
211  }
212  }
213 }
214 
215 void
217 {
218  //TODO YG : Check if cluster positionREP() is valid ?
219 
220  // Here we save in each track the list of phi/eta values of linked clusters.
221  for (BlockElt2BlockEltMap::iterator it = target2ClusterLinks_.begin();
222  it != target2ClusterLinks_.end(); ++it) {
223  reco::PFMultiLinksTC multitracks(true);
224 
225  for (BlockEltSet::iterator jt = it->second.begin();
226  jt != it->second.end(); ++jt) {
227 
228  double clusterPhi = (*jt)->clusterRef()->positionREP().Phi();
229  double clusterEta = (*jt)->clusterRef()->positionREP().Eta();
230 
231  multitracks.linkedClusters.push_back(std::make_pair(clusterPhi, clusterEta));
232  }
233 
234  it->first->setMultilinks(multitracks);
235  }
236 }
237 
238 void
240 {
241  targetSet_.clear();
242  fieldClusterSet_.clear();
243 
244  rechitsSet_.clear();
245 
246  rechit2ClusterLinks_.clear();
247  target2ClusterLinks_.clear();
248 
249  tree_.clear();
250 }
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
Abstract base class for a PFBlock element (track, cluster...)
tuple ret
prodAgent to be discontinued
void build(std::vector< KDTreeNodeInfo > &eltList, const KDTreeBox &region)
tuple trackPt
Definition: listHistos.py:120
virtual const PFClusterRef & clusterRef() const
const math::XYZPoint & position() const
cartesian position (x, y, z)
void search(const KDTreeBox &searchBox, std::vector< KDTreeNodeInfo > &resRecHitList)
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:35
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
T sqrt(T t)
Definition: SSEVec.h:18
PFMultilinksType linkedClusters
void insertFieldClusterElt(reco::PFBlockElement *ecalCluster)
const Fraction< n, m >::type & fract()
Definition: Fraction.h:38
T min(T a, T b)
Definition: MathUtil.h:58
virtual const PFRecTrackRef & trackRefPF() const
bool isNull() const
Checks for null.
Definition: Ref.h:249
#define M_PI
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
Definition: PFRecHit.h:43
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
bool isValid() const
is this point valid ?
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
void insertTargetElt(reco::PFBlockElement *track)
float getCristalPhiEtaMaxSize() const
RecHit2BlockEltMap rechit2ClusterLinks_
BlockElt2BlockEltMap target2ClusterLinks_
#define DEFINE_EDM_PLUGIN(factory, type, name)
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
float getPhiOffset() const