CMS 3D CMS Logo

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

#include <ConvBremPFTrackFinder.h>

Public Member Functions

 ConvBremPFTrackFinder (const TransientTrackBuilder &builder, double mvaBremConvCutBarrelLowPt, double mvaBremConvCutBarrelHighPt, double mvaBremConvCutEndcapsLowPt, double mvaBremConvCutEndcapsHighPt)
 
bool foundConvBremPFRecTrack (const edm::Handle< reco::PFRecTrackCollection > &thePfRecTrackCol, const edm::Handle< reco::VertexCollection > &primaryVertex, const edm::Handle< reco::PFDisplacedTrackerVertexCollection > &pfNuclears, const edm::Handle< reco::PFConversionCollection > &pfConversions, const edm::Handle< reco::PFV0Collection > &pfV0, const convbremhelpers::HeavyObjectCache *cache, bool useNuclear, bool useConversions, bool useV0, const reco::PFClusterCollection &theEClus, const reco::GsfPFRecTrack &gsfpfrectk)
 
const std::vector< reco::PFRecTrackRef > & getConvBremPFRecTracks ()
 
 ~ConvBremPFTrackFinder ()
 

Private Member Functions

void runConvBremFinder (const edm::Handle< reco::PFRecTrackCollection > &thePfRecTrackCol, const edm::Handle< reco::VertexCollection > &primaryVertex, const edm::Handle< reco::PFDisplacedTrackerVertexCollection > &pfNuclears, const edm::Handle< reco::PFConversionCollection > &pfConversions, const edm::Handle< reco::PFV0Collection > &pfV0, const convbremhelpers::HeavyObjectCache *cache, bool useNuclear, bool useConversions, bool useV0, const reco::PFClusterCollection &theEClus, const reco::GsfPFRecTrack &gsfpfrectk)
 

Private Attributes

TransientTrackBuilder builder_
 
float detaBremKF
 
float Epout
 
bool found_
 
double mvaBremConvCutBarrelHighPt_
 
double mvaBremConvCutBarrelLowPt_
 
double mvaBremConvCutEndcapsHighPt_
 
double mvaBremConvCutEndcapsLowPt_
 
float nHITS1
 
std::vector< reco::PFRecTrackRefpfRecTrRef_vec_
 
float ptRatioGsfKF
 
float secPin
 
float secPout
 
float secR
 
float sTIP
 

Detailed Description

Definition at line 34 of file ConvBremPFTrackFinder.h.

Constructor & Destructor Documentation

◆ ConvBremPFTrackFinder()

ConvBremPFTrackFinder::ConvBremPFTrackFinder ( const TransientTrackBuilder builder,
double  mvaBremConvCutBarrelLowPt,
double  mvaBremConvCutBarrelHighPt,
double  mvaBremConvCutEndcapsLowPt,
double  mvaBremConvCutEndcapsHighPt 
)

Definition at line 21 of file ConvBremPFTrackFinder.cc.

26  : builder_(builder),
27  mvaBremConvCutBarrelLowPt_(mvaBremConvCutBarrelLowPt),
28  mvaBremConvCutBarrelHighPt_(mvaBremConvCutBarrelHighPt),
29  mvaBremConvCutEndcapsLowPt_(mvaBremConvCutEndcapsLowPt),
30  mvaBremConvCutEndcapsHighPt_(mvaBremConvCutEndcapsHighPt) {}
TransientTrackBuilder builder_

◆ ~ConvBremPFTrackFinder()

ConvBremPFTrackFinder::~ConvBremPFTrackFinder ( )

Definition at line 31 of file ConvBremPFTrackFinder.cc.

31 {}

Member Function Documentation

◆ foundConvBremPFRecTrack()

bool ConvBremPFTrackFinder::foundConvBremPFRecTrack ( const edm::Handle< reco::PFRecTrackCollection > &  thePfRecTrackCol,
const edm::Handle< reco::VertexCollection > &  primaryVertex,
const edm::Handle< reco::PFDisplacedTrackerVertexCollection > &  pfNuclears,
const edm::Handle< reco::PFConversionCollection > &  pfConversions,
const edm::Handle< reco::PFV0Collection > &  pfV0,
const convbremhelpers::HeavyObjectCache cache,
bool  useNuclear,
bool  useConversions,
bool  useV0,
const reco::PFClusterCollection theEClus,
const reco::GsfPFRecTrack gsfpfrectk 
)
inline

Definition at line 43 of file ConvBremPFTrackFinder.h.

References utilities::cache(), found_, pfConversions_cfi::pfConversions, pfV0_cfi::pfV0, BeamMonitor_cff::primaryVertex, runConvBremFinder(), RecoParticleFlow_conversion_cff::useConversions, RecoParticleFlow_cff::useNuclear, and pfTrackElec_cfi::useV0.

53  {
54  found_ = false;
55  runConvBremFinder(thePfRecTrackCol,
57  pfNuclears,
59  pfV0,
60  cache,
61  useNuclear,
63  useV0,
64  theEClus,
65  gsfpfrectk);
66  return found_;
67  };
void runConvBremFinder(const edm::Handle< reco::PFRecTrackCollection > &thePfRecTrackCol, const edm::Handle< reco::VertexCollection > &primaryVertex, const edm::Handle< reco::PFDisplacedTrackerVertexCollection > &pfNuclears, const edm::Handle< reco::PFConversionCollection > &pfConversions, const edm::Handle< reco::PFV0Collection > &pfV0, const convbremhelpers::HeavyObjectCache *cache, bool useNuclear, bool useConversions, bool useV0, const reco::PFClusterCollection &theEClus, const reco::GsfPFRecTrack &gsfpfrectk)
def cache(function)
Definition: utilities.py:3
primaryVertex
hltOfflineBeamSpot for HLTMON

◆ getConvBremPFRecTracks()

const std::vector<reco::PFRecTrackRef>& ConvBremPFTrackFinder::getConvBremPFRecTracks ( )
inline

Definition at line 69 of file ConvBremPFTrackFinder.h.

References pfRecTrRef_vec_.

69 { return pfRecTrRef_vec_; };
std::vector< reco::PFRecTrackRef > pfRecTrRef_vec_

◆ runConvBremFinder()

void ConvBremPFTrackFinder::runConvBremFinder ( const edm::Handle< reco::PFRecTrackCollection > &  thePfRecTrackCol,
const edm::Handle< reco::VertexCollection > &  primaryVertex,
const edm::Handle< reco::PFDisplacedTrackerVertexCollection > &  pfNuclears,
const edm::Handle< reco::PFConversionCollection > &  pfConversions,
const edm::Handle< reco::PFV0Collection > &  pfV0,
const convbremhelpers::HeavyObjectCache cache,
bool  useNuclear,
bool  useConversions,
bool  useV0,
const reco::PFClusterCollection theEClus,
const reco::GsfPFRecTrack gsfpfrectk 
)
private

Definition at line 33 of file ConvBremPFTrackFinder.cc.

References TrackingRecHit::all, TransientTrackBuilder::build(), builder_, utilities::cache(), reco::PFCluster::calculatePositionREP(), gather_cfg::cout, debug, detaBremKF, MillePedeFileConverter_cfg::e, reco::PFTrajectoryPoint::ECALEntrance, reco::PFTrajectoryPoint::ECALShowerMax, Epout, PVValHelper::eta, reco::PFTrack::extrapolatedPoint(), found_, reco::GsfPFRecTrack::gsfTrackRef(), mps_fire::i, edm::Ref< C, T, F >::isNonnull(), reco::PFTrajectoryPoint::isValid(), reco::GsfPFRecTrack::kfPFRecTrackRef(), reco::PFTrajectoryPoint::momentum(), mvaBremConvCutBarrelHighPt_, mvaBremConvCutBarrelLowPt_, mvaBremConvCutEndcapsHighPt_, mvaBremConvCutEndcapsLowPt_, nHITS1, notFound, AlCaHLTBitMon_ParallelJobs::p, pfConversions_cfi::pfConversions, reco::GsfPFRecTrack::PFRecBrem(), pfRecTrRef_vec_, pfV0_cfi::pfV0, phi, Pi, BeamMonitor_cff::primaryVertex, edm::Handle< T >::product(), ptRatioGsfKF, MetAnalyzer::pv(), secPin, secPout, secR, IPTools::signedTransverseImpactParameter(), findQualityFiles::size, mathSSE::sqrt(), sTIP, LinkByRecHit::testTrackAndClusterByRecHit(), TwoPi, RecoParticleFlow_conversion_cff::useConversions, RecoParticleFlow_cff::useNuclear, pfTrackElec_cfi::useV0, and HltBtagValidation_cff::Vertex.

Referenced by foundConvBremPFRecTrack().

43  {
44  found_ = false;
45  bool debug = false;
46  bool debugRef = false;
47 
48  if (debug)
49  cout << "runConvBremFinder:: Entering " << endl;
50 
51  const reco::GsfTrackRef& refGsf = gsfpfrectk.gsfTrackRef();
52  const reco::PFRecTrackRef& pfTrackRef = gsfpfrectk.kfPFRecTrackRef();
53  vector<PFBrem> primPFBrem = gsfpfrectk.PFRecBrem();
54 
55  const PFRecTrackCollection& PfRTkColl = *(thePfRecTrackCol.product());
56  reco::PFRecTrackCollection::const_iterator pft = PfRTkColl.begin();
57  reco::PFRecTrackCollection::const_iterator pftend = PfRTkColl.end();
58 
59  vector<PFRecTrackRef> AllPFRecTracks;
60  AllPFRecTracks.clear();
61  unsigned int ipft = 0;
62 
63  for (; pft != pftend; ++pft, ipft++) {
64  // do not consider the kf track already associated to the seed
65  if (pfTrackRef.isNonnull())
66  if (pfTrackRef->trackRef() == pft->trackRef())
67  continue;
68 
69  PFRecTrackRef pfRecTrRef(thePfRecTrackCol, ipft);
70  TrackRef trackRef = pfRecTrRef->trackRef();
71  reco::TrackBaseRef selTrackBaseRef(trackRef);
72 
73  if (debug)
74  cout << "runConvBremFinder:: pushing_back High Purity " << pft->trackRef()->pt() << " eta,phi "
75  << pft->trackRef()->eta() << ", " << pft->trackRef()->phi() << " Memory Address Ref " << &*trackRef
76  << " Memory Address BaseRef " << &*selTrackBaseRef << endl;
77  AllPFRecTracks.push_back(pfRecTrRef);
78  }
79 
80  if (useConversions) {
81  const PFConversionCollection& PfConvColl = *(pfConversions.product());
82  for (unsigned i = 0; i < PfConvColl.size(); i++) {
84 
85  unsigned int trackSize = (convRef->pfTracks()).size();
86  if (convRef->pfTracks().size() < 2)
87  continue;
88  for (unsigned iTk = 0; iTk < trackSize; iTk++) {
89  PFRecTrackRef compPFTkRef = convRef->pfTracks()[iTk];
90  reco::TrackBaseRef newTrackBaseRef(compPFTkRef->trackRef());
91  // do not consider the kf track already associated to the seed
92  if (pfTrackRef.isNonnull()) {
93  reco::TrackBaseRef primaryTrackBaseRef(pfTrackRef->trackRef());
94  if (primaryTrackBaseRef == newTrackBaseRef)
95  continue;
96  }
97  bool notFound = true;
98  for (unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
99  reco::TrackBaseRef selTrackBaseRef(AllPFRecTracks[iPF]->trackRef());
100 
101  if (debugRef)
102  cout << "## Track 1 HP pt " << AllPFRecTracks[iPF]->trackRef()->pt() << " eta, phi "
103  << AllPFRecTracks[iPF]->trackRef()->eta() << ", " << AllPFRecTracks[iPF]->trackRef()->phi()
104  << " Memory Address Ref " << &(*AllPFRecTracks[iPF]->trackRef()) << " Memory Address BaseRef "
105  << &*selTrackBaseRef << endl;
106  if (debugRef)
107  cout << "** Track 2 CONV pt " << compPFTkRef->trackRef()->pt() << " eta, phi "
108  << compPFTkRef->trackRef()->eta() << ", " << compPFTkRef->trackRef()->phi() << " Memory Address Ref "
109  << &*compPFTkRef->trackRef() << " Memory Address BaseRef " << &*newTrackBaseRef << endl;
110  //if(selTrackBaseRef == newTrackBaseRef || AllPFRecTracks[iPF]->trackRef()== compPFTkRef->trackRef()) {
111  if (AllPFRecTracks[iPF]->trackRef() == compPFTkRef->trackRef()) {
112  if (debugRef)
113  cout << " SAME BREM REF " << endl;
114  notFound = false;
115  }
116  }
117  if (notFound) {
118  if (debug)
119  cout << "runConvBremFinder:: pushing_back Conversions " << compPFTkRef->trackRef()->pt() << " eta,phi "
120  << compPFTkRef->trackRef()->eta() << " phi " << compPFTkRef->trackRef()->phi() << endl;
121  AllPFRecTracks.push_back(compPFTkRef);
122  }
123  }
124  }
125  }
126 
127  if (useNuclear) {
128  const PFDisplacedTrackerVertexCollection& PfNuclColl = *(pfNuclears.product());
129  for (unsigned i = 0; i < PfNuclColl.size(); i++) {
130  const reco::PFDisplacedTrackerVertexRef dispacedVertexRef(pfNuclears, i);
131  unsigned int trackSize = dispacedVertexRef->pfRecTracks().size();
132  for (unsigned iTk = 0; iTk < trackSize; iTk++) {
133  reco::PFRecTrackRef newPFRecTrackRef = dispacedVertexRef->pfRecTracks()[iTk];
134  reco::TrackBaseRef newTrackBaseRef(newPFRecTrackRef->trackRef());
135  // do not consider the kf track already associated to the seed
136  if (pfTrackRef.isNonnull()) {
137  reco::TrackBaseRef primaryTrackBaseRef(pfTrackRef->trackRef());
138  if (primaryTrackBaseRef == newTrackBaseRef)
139  continue;
140  }
141  bool notFound = true;
142  for (unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
143  reco::TrackBaseRef selTrackBaseRef(AllPFRecTracks[iPF]->trackRef());
144  if (selTrackBaseRef == newTrackBaseRef)
145  notFound = false;
146  }
147  if (notFound) {
148  if (debug)
149  cout << "runConvBremFinder:: pushing_back displaced Vertex pt " << newPFRecTrackRef->trackRef()->pt()
150  << " eta,phi " << newPFRecTrackRef->trackRef()->eta() << ", " << newPFRecTrackRef->trackRef()->phi()
151  << endl;
152  AllPFRecTracks.push_back(newPFRecTrackRef);
153  }
154  }
155  }
156  }
157 
158  if (useV0) {
159  const PFV0Collection& PfV0Coll = *(pfV0.product());
160  for (unsigned i = 0; i < PfV0Coll.size(); i++) {
161  reco::PFV0Ref v0Ref(pfV0, i);
162  unsigned int trackSize = (v0Ref->pfTracks()).size();
163  for (unsigned iTk = 0; iTk < trackSize; iTk++) {
164  reco::PFRecTrackRef newPFRecTrackRef = (v0Ref->pfTracks())[iTk];
165  reco::TrackBaseRef newTrackBaseRef(newPFRecTrackRef->trackRef());
166  // do not consider the kf track already associated to the seed
167  if (pfTrackRef.isNonnull()) {
168  reco::TrackBaseRef primaryTrackBaseRef(pfTrackRef->trackRef());
169  if (primaryTrackBaseRef == newTrackBaseRef)
170  continue;
171  }
172  bool notFound = true;
173  for (unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
174  reco::TrackBaseRef selTrackBaseRef(AllPFRecTracks[iPF]->trackRef());
175  if (selTrackBaseRef == newTrackBaseRef)
176  notFound = false;
177  }
178  if (notFound) {
179  if (debug)
180  cout << "runConvBremFinder:: pushing_back V0 " << newPFRecTrackRef->trackRef()->pt() << " eta,phi "
181  << newPFRecTrackRef->trackRef()->eta() << ", " << newPFRecTrackRef->trackRef()->phi() << endl;
182  AllPFRecTracks.push_back(newPFRecTrackRef);
183  }
184  }
185  }
186  }
187 
188  pfRecTrRef_vec_.clear();
189 
190  for (unsigned iPF = 0; iPF < AllPFRecTracks.size(); iPF++) {
191  double dphi = fabs(AllPFRecTracks[iPF]->trackRef()->phi() - refGsf->phi());
192  if (dphi > TMath::Pi())
193  dphi -= TMath::TwoPi();
194  double deta = fabs(AllPFRecTracks[iPF]->trackRef()->eta() - refGsf->eta());
195 
196  // limiting the phase space (just for saving cpu-time)
197  if (fabs(dphi) > 1.0 || fabs(deta) > 0.4)
198  continue;
199 
200  double minDEtaBremKF = 1000.;
201  double minDPhiBremKF = 1000.;
202  double minDRBremKF = 1000.;
203  double minDEtaBremKFPos = 1000.;
204  double minDPhiBremKFPos = 1000.;
205  double minDRBremKFPos = 1000.;
206  reco::TrackRef trkRef = AllPFRecTracks[iPF]->trackRef();
207 
208  double secEta = trkRef->innerMomentum().eta();
209  double secPhi = trkRef->innerMomentum().phi();
210 
211  for (unsigned ipbrem = 0; ipbrem < primPFBrem.size(); ipbrem++) {
212  if (primPFBrem[ipbrem].indTrajPoint() == 99)
213  continue;
214  const reco::PFTrajectoryPoint& atPrimECAL =
215  primPFBrem[ipbrem].extrapolatedPoint(reco::PFTrajectoryPoint::ECALEntrance);
216  if (!atPrimECAL.isValid())
217  continue;
218  double bremEta = atPrimECAL.momentum().Eta();
219  double bremPhi = atPrimECAL.momentum().Phi();
220 
221  double deta = fabs(bremEta - secEta);
222  double dphi = fabs(bremPhi - secPhi);
223  if (dphi > TMath::Pi())
224  dphi -= TMath::TwoPi();
225  double DR = sqrt(deta * deta + dphi * dphi);
226 
227  double detaPos = fabs(bremEta - trkRef->innerPosition().eta());
228  double dphiPos = fabs(bremPhi - trkRef->innerPosition().phi());
229  if (dphiPos > TMath::Pi())
230  dphiPos -= TMath::TwoPi();
231  double DRPos = sqrt(detaPos * detaPos + dphiPos * dphiPos);
232 
233  // find the closest track tangent
234  if (DR < minDRBremKF) {
235  minDRBremKF = DR;
236  minDEtaBremKF = deta;
237  minDPhiBremKF = fabs(dphi);
238  }
239 
240  if (DRPos < minDRBremKFPos) {
241  minDRBremKFPos = DR;
242  minDEtaBremKFPos = detaPos;
243  minDPhiBremKFPos = fabs(dphiPos);
244  }
245  }
246 
247  //gsfR
248  float gsfR = sqrt(refGsf->innerPosition().x() * refGsf->innerPosition().x() +
249  refGsf->innerPosition().y() * refGsf->innerPosition().y());
250 
251  // secR
252  secR = sqrt(trkRef->innerPosition().x() * trkRef->innerPosition().x() +
253  trkRef->innerPosition().y() * trkRef->innerPosition().y());
254 
255  // apply loose selection (to be parallel) between the secondary track and brem-tangents.
256  // Moreover if the secR is internal with respect to the GSF track by two pixel layers discard it.
257  if ((minDPhiBremKF < 0.1 || minDPhiBremKFPos < 0.1) && (minDEtaBremKF < 0.02 || minDEtaBremKFPos < 0.02) &&
258  secR > (gsfR - 8)) {
259  if (debug)
260  cout << "runConvBremFinder:: OK Find track and BREM close "
261  << " MinDphi " << minDPhiBremKF << " MinDeta " << minDEtaBremKF << endl;
262 
263  float MinDist = 100000.;
264  float EE_calib = 0.;
265  PFRecTrack pfrectrack = *AllPFRecTracks[iPF];
266  // Find and ECAL associated cluster
267  for (PFClusterCollection::const_iterator clus = theEClus.begin(); clus != theEClus.end(); clus++) {
268  // Removed unusd variable, left this in case it has side effects
269  clus->position();
270  double dist = -1.;
271  PFCluster clust = *clus;
272  clust.calculatePositionREP();
274  ? LinkByRecHit::testTrackAndClusterByRecHit(pfrectrack, clust)
275  : -1.;
276 
277  if (dist > 0. && dist < MinDist) {
278  MinDist = dist;
279  EE_calib = cache->pfcalib_->energyEm(*clus, 0.0, 0.0, false);
280  }
281  }
282  if (MinDist > 0. && MinDist < 100000.) {
283  // compute all the input variables for conv brem selection
284 
285  secPout = sqrt(trkRef->outerMomentum().x() * trkRef->outerMomentum().x() +
286  trkRef->outerMomentum().y() * trkRef->outerMomentum().y() +
287  trkRef->outerMomentum().z() * trkRef->outerMomentum().z());
288 
289  secPin = sqrt(trkRef->innerMomentum().x() * trkRef->innerMomentum().x() +
290  trkRef->innerMomentum().y() * trkRef->innerMomentum().y() +
291  trkRef->innerMomentum().z() * trkRef->innerMomentum().z());
292 
293  // maybe put innter momentum pt?
294  ptRatioGsfKF = trkRef->pt() / (refGsf->ptMode());
295 
296  Vertex dummy;
297  const Vertex* pv = &dummy;
299  if (!primaryVertex->empty()) {
300  pv = &*primaryVertex->begin();
301  // we always use the first vertex (at the moment)
303  } else { // create a dummy PV
305  e(0, 0) = 0.0015 * 0.0015;
306  e(1, 1) = 0.0015 * 0.0015;
307  e(2, 2) = 15. * 15.;
308  Vertex::Point p(0, 0, 0);
309  dummy = Vertex(p, e, 0, 0, 0);
310  }
311 
312  // direction of the Gsf track
313  GlobalVector direction(refGsf->innerMomentum().x(), refGsf->innerMomentum().y(), refGsf->innerMomentum().z());
314 
315  TransientTrack transientTrack = builder_.build(*trkRef);
316  sTIP = IPTools::signedTransverseImpactParameter(transientTrack, direction, *pv).second.significance();
317 
318  Epout = EE_calib / secPout;
319 
320  // eta distance brem-secondary kf track
321  detaBremKF = minDEtaBremKF;
322 
323  // Number of commont hits
324  unsigned int tmpSh = 0;
325  //uint ish=0;
326  int kfhitcounter = 0;
327  for (auto const& nhit : refGsf->recHits())
328  if (nhit->isValid()) {
329  kfhitcounter = 0;
330  for (auto const& ihit : trkRef->recHits())
331  if (ihit->isValid()) {
332  // method 1
333  if (nhit->sharesInput(ihit, TrackingRecHit::all))
334  tmpSh++;
335  kfhitcounter++;
336  // method 2 to switch in case of problem with rechit collections
337  // if((ihit->geographicalId()==nhit->geographicalId())&&
338  // ((nhit->localPosition()-ihit->localPosition()).mag()<0.01)) ish++;
339  }
340  }
341 
342  nHITS1 = tmpSh;
343 
344  TString weightfilepath = "";
345  double mvaValue = -10;
346  double cutvalue = -10;
347 
348  float vars[6] = {secR, sTIP, nHITS1, Epout, detaBremKF, ptRatioGsfKF};
349  if (refGsf->pt() < 20 && fabs(refGsf->eta()) < 1.479) {
350  mvaValue = cache->gbrBarrelLowPt_->GetClassifier(vars);
351  cutvalue = mvaBremConvCutBarrelLowPt_;
352  }
353  if (refGsf->pt() > 20 && fabs(refGsf->eta()) < 1.479) {
354  mvaValue = cache->gbrBarrelHighPt_->GetClassifier(vars);
355  cutvalue = mvaBremConvCutBarrelHighPt_;
356  }
357  if (refGsf->pt() < 20 && fabs(refGsf->eta()) > 1.479) {
358  mvaValue = cache->gbrEndcapsLowPt_->GetClassifier(vars);
359  cutvalue = mvaBremConvCutEndcapsLowPt_;
360  }
361  if (refGsf->pt() > 20 && fabs(refGsf->eta()) > 1.479) {
362  mvaValue = cache->gbrEndcapsHighPt_->GetClassifier(vars);
363  cutvalue = mvaBremConvCutEndcapsHighPt_;
364  }
365 
366  if (debug)
367  cout << "Gsf track Pt, Eta " << refGsf->pt() << " " << refGsf->eta() << endl;
368  if (debug)
369  cout << "Cutvalue " << cutvalue << endl;
370 
371  if ((kfhitcounter - nHITS1) <= 3 && nHITS1 > 3)
372  mvaValue = 2; // this means duplicates tracks, very likely not physical
373 
374  if (debug)
375  cout << " The input variables for conv brem tracks identification " << endl
376  << " secR " << secR << " gsfR " << gsfR << endl
377  << " N shared hits " << nHITS1 << endl
378  << " sTIP " << sTIP << endl
379  << " detaBremKF " << detaBremKF << endl
380  << " E/pout " << Epout << endl
381  << " ptRatioKFGsf " << ptRatioGsfKF << endl
382  << " ***** MVA ***** " << mvaValue << endl;
383 
384  if (mvaValue > cutvalue) {
385  found_ = true;
386  pfRecTrRef_vec_.push_back(AllPFRecTracks[iPF]);
387  }
388  } // end MinDIST
389  } // end selection kf - brem tangents
390  } // loop on the kf tracks
391 }
size
Write out results.
const double TwoPi
const double Pi
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:20
Particle flow cluster, see clustering algorithm in PFClusterAlgo.
Definition: PFCluster.h:42
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
std::pair< bool, Measurement1D > signedTransverseImpactParameter(const reco::TransientTrack &track, const GlobalVector &direction, const reco::Vertex &vertex)
Definition: IPTools.cc:57
T const * product() const
Definition: Handle.h:70
math::Error< dimension >::type Error
covariance error matrix (3x3)
Definition: Vertex.h:44
bool isNonnull() const
Checks for non-null.
Definition: Ref.h:238
std::vector< PFConversion > PFConversionCollection
collection of PFConversion objects
reco::TransientTrack build(const reco::Track *p) const
bool isValid() const
is this point valid ?
void calculatePositionREP()
computes posrep_ once and for all
Definition: PFCluster.h:95
const edm::Ref< std::vector< PFRecTrack > > & kfPFRecTrackRef() const
Definition: GsfPFRecTrack.h:37
T sqrt(T t)
Definition: SSEVec.h:19
TransientTrackBuilder builder_
const std::vector< reco::PFBrem > & PFRecBrem() const
Definition: GsfPFRecTrack.h:45
def pv(vc)
Definition: MetAnalyzer.py:7
std::vector< reco::PFRecTrackRef > pfRecTrRef_vec_
math::XYZPoint Point
point in the space
Definition: Vertex.h:40
std::vector< PFV0 > PFV0Collection
collection of PFV0 objects
Definition: PFV0Fwd.h:9
const reco::PFTrajectoryPoint & extrapolatedPoint(unsigned layerid) const
Definition: PFTrack.cc:46
#define debug
Definition: HDRShower.cc:19
static const GlobalPoint notFound(0, 0, 0)
const reco::GsfTrackRef & gsfTrackRef() const
Definition: GsfPFRecTrack.h:34
def cache(function)
Definition: utilities.py:3
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
vars
Definition: DeepTauId.cc:30
primaryVertex
hltOfflineBeamSpot for HLTMON
static double testTrackAndClusterByRecHit(const reco::PFRecTrack &track, const reco::PFCluster &cluster, bool isBrem=false, bool debug=false)
Definition: LinkByRecHit.cc:18
std::vector< PFRecTrack > PFRecTrackCollection
collection of PFRecTrack objects
Definition: PFRecTrackFwd.h:9
std::vector< PFDisplacedTrackerVertex > PFDisplacedTrackerVertexCollection
collection of DisplacedTrackerVertexs

Member Data Documentation

◆ builder_

TransientTrackBuilder ConvBremPFTrackFinder::builder_
private

Definition at line 85 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ detaBremKF

float ConvBremPFTrackFinder::detaBremKF
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ Epout

float ConvBremPFTrackFinder::Epout
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ found_

bool ConvBremPFTrackFinder::found_
private

Definition at line 84 of file ConvBremPFTrackFinder.h.

Referenced by foundConvBremPFRecTrack(), and runConvBremFinder().

◆ mvaBremConvCutBarrelHighPt_

double ConvBremPFTrackFinder::mvaBremConvCutBarrelHighPt_
private

Definition at line 86 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ mvaBremConvCutBarrelLowPt_

double ConvBremPFTrackFinder::mvaBremConvCutBarrelLowPt_
private

Definition at line 86 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ mvaBremConvCutEndcapsHighPt_

double ConvBremPFTrackFinder::mvaBremConvCutEndcapsHighPt_
private

Definition at line 86 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ mvaBremConvCutEndcapsLowPt_

double ConvBremPFTrackFinder::mvaBremConvCutEndcapsLowPt_
private

Definition at line 86 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ nHITS1

float ConvBremPFTrackFinder::nHITS1
private

Definition at line 92 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ pfRecTrRef_vec_

std::vector<reco::PFRecTrackRef> ConvBremPFTrackFinder::pfRecTrRef_vec_
private

Definition at line 89 of file ConvBremPFTrackFinder.h.

Referenced by getConvBremPFRecTracks(), and runConvBremFinder().

◆ ptRatioGsfKF

float ConvBremPFTrackFinder::ptRatioGsfKF
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ secPin

float ConvBremPFTrackFinder::secPin
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ secPout

float ConvBremPFTrackFinder::secPout
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ secR

float ConvBremPFTrackFinder::secR
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

◆ sTIP

float ConvBremPFTrackFinder::sTIP
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().