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 35 of file ConvBremPFTrackFinder.h.

Constructor & Destructor Documentation

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

Definition at line 22 of file ConvBremPFTrackFinder.cc.

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

Definition at line 33 of file ConvBremPFTrackFinder.cc.

33 { }

Member Function Documentation

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 45 of file ConvBremPFTrackFinder.h.

References found_, and runConvBremFinder().

56  {
57  found_ = false;
58  runConvBremFinder(thePfRecTrackCol,primaryVertex,
59  pfNuclears,pfConversions,
60  pfV0,cache,useNuclear,
62  theEClus,gsfpfrectk);
63  return found_;};
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)
const std::vector<reco::PFRecTrackRef>& ConvBremPFTrackFinder::getConvBremPFRecTracks ( )
inline
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 36 of file ConvBremPFTrackFinder.cc.

References TrackingRecHit::all, particleFlowSuperClusterECAL_cfi::applyCrackCorrections, TransientTrackBuilder::build(), builder_, reco::PFTrack::calculatePositionREP(), reco::PFCluster::calculatePositionREP(), gather_cfg::cout, debug, detaBremKF, MillePedeFileConverter_cfg::e, reco::PFTrajectoryPoint::ECALEntrance, reco::PFTrajectoryPoint::ECALShowerMax, Epout, PVValHelper::eta, reco::PFTrack::extrapolatedPoint(), found_, convbremhelpers::HeavyObjectCache::gbrBarrelHighPt_, convbremhelpers::HeavyObjectCache::gbrBarrelLowPt_, convbremhelpers::HeavyObjectCache::gbrEndcapsHighPt_, convbremhelpers::HeavyObjectCache::gbrEndcapsLowPt_, 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, convbremhelpers::HeavyObjectCache::pfcalib_, reco::GsfPFRecTrack::PFRecBrem(), pfRecTrRef_vec_, phi, Pi, impactParameterTagInfos_cfi::primaryVertex, edm::Handle< T >::product(), ptRatioGsfKF, MetAnalyzer::pv(), secPin, secPout, secR, IPTools::signedTransverseImpactParameter(), findQualityFiles::size, mathSSE::sqrt(), sTIP, LinkByRecHit::testTrackAndClusterByRecHit(), and TwoPi.

Referenced by foundConvBremPFRecTrack(), and getConvBremPFRecTracks().

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

Member Data Documentation

TransientTrackBuilder ConvBremPFTrackFinder::builder_
private

Definition at line 84 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::detaBremKF
private

Definition at line 88 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::Epout
private

Definition at line 88 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

bool ConvBremPFTrackFinder::found_
private

Definition at line 83 of file ConvBremPFTrackFinder.h.

Referenced by foundConvBremPFRecTrack(), and runConvBremFinder().

double ConvBremPFTrackFinder::mvaBremConvCutBarrelHighPt_
private

Definition at line 85 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

double ConvBremPFTrackFinder::mvaBremConvCutBarrelLowPt_
private

Definition at line 85 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

double ConvBremPFTrackFinder::mvaBremConvCutEndcapsHighPt_
private

Definition at line 85 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

double ConvBremPFTrackFinder::mvaBremConvCutEndcapsLowPt_
private

Definition at line 85 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::nHITS1
private

Definition at line 90 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

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

Definition at line 87 of file ConvBremPFTrackFinder.h.

Referenced by getConvBremPFRecTracks(), and runConvBremFinder().

float ConvBremPFTrackFinder::ptRatioGsfKF
private

Definition at line 88 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::secPin
private

Definition at line 88 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::secPout
private

Definition at line 88 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::secR
private

Definition at line 88 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::sTIP
private

Definition at line 88 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().