CMS 3D CMS Logo

KDTreeLinkerPSEcal.cc
Go to the documentation of this file.
4 
5 #include "TMath.h"
6 
7 // This class is used to find all links between PreShower clusters and ECAL clusters
8 // using a KDTree algorithm.
9 // It is used in PFBlockAlgo.cc in the function links().
11 public:
13  ~KDTreeLinkerPSEcal() override;
14 
15  // With this method, we create the list of psCluster that we want to link.
16  void insertTargetElt(reco::PFBlockElement *psCluster) override;
17 
18  // Here, we create the list of ecalCluster that we want to link. From ecalCluster
19  // and fraction, we will create a second list of rechits that will be used to
20  // build the KDTree.
21  void insertFieldClusterElt(reco::PFBlockElement *ecalCluster) override;
22 
23  // The KDTree building from rechits list.
24  void buildTree() override;
25 
26  // Here we will iterate over all psCluster. For each one, we will search the closest
27  // rechits in the KDTree, from rechits we will find the ecalClusters and after that
28  // we will check the links between the psCluster and all closest ecalClusters.
29  void searchLinks() override;
30 
31  // Here, we will store all PS/ECAL founded links in the PFBlockElement class
32  // of each psCluster in the PFmultilinks field.
33  void updatePFBlockEltWithLinks() override;
34 
35  // Here we free all allocated structures.
36  void clear() override;
37 
38 private:
39  // This method allows us to build the "tree" from the "rechitsSet".
41 
42 private:
43  // Some const values.
44  const double resPSpitch_;
45  const double resPSlength_;
46  const double ps1ToEcal_; // ratio : zEcal / zPS1
47  const double ps2ToEcal_; // ration : zEcal / zPS2
48 
49  // Data used by the KDTree algorithm : sets of PS and ECAL clusters.
52 
53  // Sets of rechits that compose the ECAL clusters. We differenctiate
54  // the rechits by their Z value.
57 
58  // Map of linked PS/ECAL clusters.
60 
61  // Map of the ECAL clusters associated to a rechit.
63 
64  // KD trees
67 };
68 
69 // the text name is different so that we can easily
70 // construct it when calling the factory
71 DEFINE_EDM_PLUGIN(KDTreeLinkerFactory, KDTreeLinkerPSEcal, "KDTreePreshowerAndECALLinker");
72 
74  : KDTreeLinkerBase(conf), resPSpitch_(0.19), resPSlength_(6.1), ps1ToEcal_(1.072), ps2ToEcal_(1.057) {}
75 
77 
79  // This test is more or less done in PFBlockAlgo.h. In others cases, it should be switch on.
80  //if (!((psCluster->clusterRef()->layer() == PFLayer::PS1) || (psCluster->clusterRef()->layer() == PFLayer::PS2)))
81  // return;
82  targetSet_.insert(psCluster);
83 }
84 
86  const reco::PFClusterRef &clusterref = ecalCluster->clusterRef();
87 
88  if (clusterref->layer() != PFLayer::ECAL_ENDCAP)
89  return;
90 
91  const std::vector<reco::PFRecHitFraction> &fraction = clusterref->recHitFractions();
92 
93  // We create a list of cluster
94  fieldClusterSet_.insert(ecalCluster);
95 
96  double clusterz = clusterref->position().Z();
97  RecHitSet &rechitsSet = (clusterz < 0) ? rechitsNegSet_ : rechitsPosSet_;
98 
99  for (size_t rhit = 0; rhit < fraction.size(); ++rhit) {
100  const reco::PFRecHitRef &rh = fraction[rhit].recHitRef();
101  double fract = fraction[rhit].fraction();
102 
103  if ((rh.isNull()) || (fract < cutOffFrac))
104  continue;
105 
106  const reco::PFRecHit &rechit = *rh;
107 
108  // We save the links rechit to Clusters
109  rechit2ClusterLinks_[&rechit].insert(ecalCluster);
110 
111  // We create a liste of rechits
112  rechitsSet.insert(&rechit);
113  }
114 }
115 
119 }
120 
122  // List of pseudo-rechits that will be used to create the KDTree
123  std::vector<KDTreeNodeInfo<reco::PFRecHit const *, 2>> eltList;
124 
125  // Filling of this eltList
126  for (RecHitSet::const_iterator it = rechitsSet.begin(); it != rechitsSet.end(); it++) {
127  const reco::PFRecHit *rh = *it;
128  const auto &posxyz = rh->position();
129 
130  KDTreeNodeInfo<reco::PFRecHit const *, 2> rhinfo{rh, posxyz.x(), posxyz.y()};
131  eltList.push_back(rhinfo);
132  }
133 
134  // xmin-xmax, ymain-ymax
135  KDTreeBox region{-150.f, 150.f, -150.f, 150.f};
136 
137  // We may now build the KDTree
138  tree.build(eltList, region);
139 }
140 
142  // Most of the code has been taken from LinkByRecHit.cc
143 
144  // We iterate over the PS clusters.
145  for (BlockEltSet::iterator it = targetSet_.begin(); it != targetSet_.end(); it++) {
146  reco::PFClusterRef clusterPSRef = (*it)->clusterRef();
147  const reco::PFCluster &clusterPS = *clusterPSRef;
148 
149  // PS cluster position, extrapolated to ECAL
150  double zPS = clusterPS.position().Z();
151  double xPS = clusterPS.position().X();
152  double yPS = clusterPS.position().Y();
153 
154  double etaPS = fabs(clusterPS.positionREP().eta());
155  double deltaX = 0.;
156  double deltaY = 0.;
157  float xPSonEcal = xPS;
158  float yPSonEcal = yPS;
159 
160  if (clusterPS.layer() == PFLayer::PS1) { // PS1
161 
162  // vertical strips, measure x with pitch precision
163  deltaX = resPSpitch_;
164  deltaY = resPSlength_;
165  xPSonEcal *= ps1ToEcal_;
166  yPSonEcal *= ps1ToEcal_;
167 
168  } else { // PS2
169 
170  // horizontal strips, measure y with pitch precision
171  deltaY = resPSpitch_;
172  deltaX = resPSlength_;
173  xPSonEcal *= ps2ToEcal_;
174  yPSonEcal *= ps2ToEcal_;
175  }
176 
177  // Estimate the maximal envelope in phi/eta that will be used to find rechit candidates.
178  // Same envelope for cap et barrel rechits.
179 
180  double maxEcalRadius = cristalXYMaxSize_ / 2.;
181 
182  // The inflation factor includes the approximate projection from Preshower to ECAL
183  double inflation = 2.4 - (etaPS - 1.6);
184  float rangeX = maxEcalRadius * (1 + (0.05 + 1.0 / maxEcalRadius * deltaX / 2.)) * inflation;
185  float rangeY = maxEcalRadius * (1 + (0.05 + 1.0 / maxEcalRadius * deltaY / 2.)) * inflation;
186 
187  // We search for all candidate recHits, ie all recHits contained in the maximal size envelope.
188  std::vector<reco::PFRecHit const *> recHits;
189  KDTreeBox trackBox(xPSonEcal - rangeX, xPSonEcal + rangeX, yPSonEcal - rangeY, yPSonEcal + rangeY);
190 
191  if (zPS < 0)
192  treeNeg_.search(trackBox, recHits);
193  else
194  treePos_.search(trackBox, recHits);
195 
196  for (auto const &recHit : recHits) {
197  const auto &corners = recHit->getCornersXYZ();
198 
199  // Find all clusters associated to given rechit
200  RecHit2BlockEltMap::iterator ret = rechit2ClusterLinks_.find(recHit);
201 
202  for (BlockEltSet::const_iterator clusterIt = ret->second.begin(); clusterIt != ret->second.end(); clusterIt++) {
203  reco::PFClusterRef clusterref = (*clusterIt)->clusterRef();
204  double clusterz = clusterref->position().z();
205 
206  const auto &posxyz = recHit->position() * zPS / clusterz;
207 
208  double x[5];
209  double y[5];
210  for (unsigned jc = 0; jc < 4; ++jc) {
211  auto cornerpos = corners[jc].basicVector() * zPS / clusterz;
212  x[3 - jc] = cornerpos.x() +
213  (cornerpos.x() - posxyz.x()) * (0.05 + 1.0 / fabs((cornerpos.x() - posxyz.x())) * deltaX / 2.);
214  y[3 - jc] = cornerpos.y() +
215  (cornerpos.y() - posxyz.y()) * (0.05 + 1.0 / fabs((cornerpos.y() - posxyz.y())) * deltaY / 2.);
216  }
217 
218  x[4] = x[0];
219  y[4] = y[0];
220 
221  bool isinside = TMath::IsInside(xPS, yPS, 5, x, y);
222 
223  // Check if the track and the cluster are linked
224  if (isinside)
225  target2ClusterLinks_[*it].insert(*clusterIt);
226  }
227  }
228  }
229 }
230 
232  //TODO YG : Check if cluster positionREP() is valid ?
233 
234  // Here we save in each PS the list of phi/eta values of linked ECAL clusters.
235  for (BlockElt2BlockEltMap::iterator it = target2ClusterLinks_.begin(); it != target2ClusterLinks_.end(); ++it) {
236  const auto &psElt = it->first;
237  const auto &ecalEltSet = it->second;
238  reco::PFMultiLinksTC multitracks(true);
239 
240  for (const auto &ecalElt : ecalEltSet) {
241  double clusterphi = ecalElt->clusterRef()->positionREP().phi();
242  double clustereta = ecalElt->clusterRef()->positionREP().eta();
243 
244  multitracks.linkedClusters.push_back(std::make_pair(clusterphi, clustereta));
245 
246  // We set the multilinks flag of the ECAL element (for links to PS) to true. It will allow us to
247  // use it in an optimized way in prefilter
248  ecalElt->setIsValidMultilinks(true, _targetType);
249  }
250 
251  // We set multilinks of the PS element (for links to ECAL)
252  psElt->setMultilinks(multitracks, _fieldType);
253  }
254 }
255 
257  targetSet_.clear();
258  fieldClusterSet_.clear();
259 
260  rechitsNegSet_.clear();
261  rechitsPosSet_.clear();
262 
263  rechit2ClusterLinks_.clear();
264  target2ClusterLinks_.clear();
265 
266  treeNeg_.clear();
267  treePos_.clear();
268 }
KDTreeLinkerAlgo.h
runTheMatrix.ret
ret
prodAgent to be discontinued
Definition: runTheMatrix.py:367
KDTreeLinkerPSEcal::insertFieldClusterElt
void insertFieldClusterElt(reco::PFBlockElement *ecalCluster) override
Definition: KDTreeLinkerPSEcal.cc:85
KDTreeLinkerPSEcal::updatePFBlockEltWithLinks
void updatePFBlockEltWithLinks() override
Definition: KDTreeLinkerPSEcal.cc:231
DDAxes::y
KDTreeLinkerPSEcal::searchLinks
void searchLinks() override
Definition: KDTreeLinkerPSEcal.cc:141
KDTreeLinkerBase::cristalXYMaxSize_
float cristalXYMaxSize_
Definition: KDTreeLinkerBase.h:78
KDTreeLinkerPSEcal::KDTreeLinkerPSEcal
KDTreeLinkerPSEcal(const edm::ParameterSet &conf)
Definition: KDTreeLinkerPSEcal.cc:73
reco::PFCluster::layer
PFLayer::Layer layer() const
cluster layer, see PFLayer.h in this directory
Definition: PFCluster.cc:56
edm::Ref::isNull
bool isNull() const
Checks for null.
Definition: Ref.h:235
tree
Definition: tree.py:1
KDTreeLinkerPSEcal::rechitsPosSet_
RecHitSet rechitsPosSet_
Definition: KDTreeLinkerPSEcal.cc:56
KDTreeLinkerPSEcal::~KDTreeLinkerPSEcal
~KDTreeLinkerPSEcal() override
Definition: KDTreeLinkerPSEcal.cc:76
DDAxes::x
rpcPointValidation_cfi.recHit
recHit
Definition: rpcPointValidation_cfi.py:7
KDTreeLinkerBase::cutOffFrac
const float cutOffFrac
Definition: KDTreeLinkerBase.h:87
edm::Ref< PFClusterCollection >
KDTreeLinkerPSEcal::target2ClusterLinks_
BlockElt2BlockEltMap target2ClusterLinks_
Definition: KDTreeLinkerPSEcal.cc:59
PFLayer::PS1
Definition: PFLayer.h:31
KDTreeLinkerAlgo::search
void search(const KDTreeBox< DIM > &searchBox, std::vector< DATA > &resRecHitList)
Definition: KDTreeLinkerAlgo.h:203
KDTreeLinkerBase::_targetType
reco::PFBlockElement::Type _targetType
Definition: KDTreeLinkerBase.h:75
KDTreeLinkerPSEcal::treePos_
KDTreeLinkerAlgo< reco::PFRecHit const * > treePos_
Definition: KDTreeLinkerPSEcal.cc:66
KDTreeLinkerBase.h
HLT_FULL_cff.fraction
fraction
Definition: HLT_FULL_cff.py:52858
KDTreeLinkerPSEcal::treeNeg_
KDTreeLinkerAlgo< reco::PFRecHit const * > treeNeg_
Definition: KDTreeLinkerPSEcal.cc:65
KDTreeBox
Definition: KDTreeLinkerAlgo.h:14
PFCluster.h
KDTreeLinkerAlgo< reco::PFRecHit const * >
reco::PFMultiLinksTC
Definition: PFMultilinksTC.h:14
reco::PFMultiLinksTC::linkedClusters
PFMultilinksType linkedClusters
Definition: PFMultilinksTC.h:17
reco::PFCluster::positionREP
const REPPoint & positionREP() const
cluster position: rho, eta, phi
Definition: PFCluster.h:92
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
KDTreeLinkerPSEcal::targetSet_
BlockEltSet targetSet_
Definition: KDTreeLinkerPSEcal.cc:50
FastTrackerRecHitMaskProducer_cfi.recHits
recHits
Definition: FastTrackerRecHitMaskProducer_cfi.py:8
KDTreeLinkerPSEcal::fieldClusterSet_
BlockEltSet fieldClusterSet_
Definition: KDTreeLinkerPSEcal.cc:51
KDTreeLinkerPSEcal
Definition: KDTreeLinkerPSEcal.cc:10
edm::ParameterSet
Definition: ParameterSet.h:47
RecHit2BlockEltMap
std::map< const reco::PFRecHit *, BlockEltSet > RecHit2BlockEltMap
Definition: KDTreeLinkerBase.h:16
edmplugin::PluginFactory
Definition: PluginFactory.h:34
KDTreeLinkerPSEcal::rechit2ClusterLinks_
RecHit2BlockEltMap rechit2ClusterLinks_
Definition: KDTreeLinkerPSEcal.cc:62
HLT_FULL_cff.region
region
Definition: HLT_FULL_cff.py:84949
KDTreeLinkerPSEcal::insertTargetElt
void insertTargetElt(reco::PFBlockElement *psCluster) override
Definition: KDTreeLinkerPSEcal.cc:78
funct::fract
const Fraction< n, m >::type & fract()
Definition: Fraction.h:36
KDTreeNodeInfo
Definition: KDTreeLinkerAlgo.h:34
reco::PFBlockElement
Abstract base class for a PFBlock element (track, cluster...)
Definition: PFBlockElement.h:26
reco::CaloCluster::position
const math::XYZPoint & position() const
cluster centroid position
Definition: CaloCluster.h:154
KDTreeLinkerPSEcal::clear
void clear() override
Definition: KDTreeLinkerPSEcal.cc:256
BlockElt2BlockEltMap
std::map< reco::PFBlockElement *, BlockEltSet > BlockElt2BlockEltMap
Definition: KDTreeLinkerBase.h:17
KDTreeLinkerPSEcal::ps1ToEcal_
const double ps1ToEcal_
Definition: KDTreeLinkerPSEcal.cc:46
KDTreeLinkerAlgo::clear
void clear()
Definition: KDTreeLinkerAlgo.h:122
RecHitSet
std::set< const reco::PFRecHit * > RecHitSet
Definition: KDTreeLinkerBase.h:14
reco::PFCluster
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
BlockEltSet
std::set< reco::PFBlockElement * > BlockEltSet
Definition: KDTreeLinkerBase.h:13
KDTreeLinkerBase::_fieldType
reco::PFBlockElement::Type _fieldType
Definition: KDTreeLinkerBase.h:75
reco::PFRecHit::position
PositionType const & position() const
rechit cell centre x, y, z
Definition: PFRecHit.h:117
reco::PFRecHit
Particle flow rechit (rechit + geometry and topology information). See clustering algorithm in PFClus...
Definition: PFRecHit.h:31
KDTreeLinkerBase
Definition: KDTreeLinkerBase.h:19
KDTreeLinkerPSEcal::resPSlength_
const double resPSlength_
Definition: KDTreeLinkerPSEcal.cc:45
KDTreeLinkerPSEcal::ps2ToEcal_
const double ps2ToEcal_
Definition: KDTreeLinkerPSEcal.cc:47
KDTreeLinkerPSEcal::resPSpitch_
const double resPSpitch_
Definition: KDTreeLinkerPSEcal.cc:44
KDTreeLinkerPSEcal::buildTree
void buildTree() override
Definition: KDTreeLinkerPSEcal.cc:116
PFLayer::ECAL_ENDCAP
Definition: PFLayer.h:32
KDTreeLinkerPSEcal::rechitsNegSet_
RecHitSet rechitsNegSet_
Definition: KDTreeLinkerPSEcal.cc:55
reco::PFBlockElement::clusterRef
virtual const PFClusterRef & clusterRef() const
Definition: PFBlockElement.h:90