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 Member Functions | Private Attributes
ConvBremPFTrackFinder Class Reference

#include <ConvBremPFTrackFinder.h>

Public Member Functions

 ConvBremPFTrackFinder (const TransientTrackBuilder &builder, double mvaBremConvCut, std::string mvaWeightFileConvBrem)
 
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, bool useNuclear, bool useConversions, bool useV0, const reco::PFClusterCollection &theEClus, 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, bool useNuclear, bool useConversions, bool useV0, const reco::PFClusterCollection &theEClus, reco::GsfPFRecTrack gsfpfrectk)
 

Private Attributes

TransientTrackBuilder builder_
 
float detaBremKF
 
float Epout
 
bool found_
 
double mvaBremConvCut_
 
std::string mvaWeightFileConvBrem_
 
float nHITS1
 
std::vector< reco::PFRecTrackRefpfRecTrRef_vec_
 
float ptRatioGsfKF
 
float secPin
 
float secPout
 
float secR
 
float sTIP
 
TMVA::Reader * tmvaReader_
 

Detailed Description

Definition at line 28 of file ConvBremPFTrackFinder.h.

Constructor & Destructor Documentation

ConvBremPFTrackFinder::ConvBremPFTrackFinder ( const TransientTrackBuilder builder,
double  mvaBremConvCut,
std::string  mvaWeightFileConvBrem 
)

Definition at line 22 of file ConvBremPFTrackFinder.cc.

References detaBremKF, Epout, nHITS1, ptRatioGsfKF, secPin, secR, sTIP, and tmvaReader_.

24  :
25  builder_(builder),
26  mvaBremConvCut_(mvaBremConvCut),
27  mvaWeightFileConvBrem_(mvaWeightFileConvBrem)
28 {
29  tmvaReader_ = new TMVA::Reader("!Color:Silent");
30  tmvaReader_->AddVariable("secR",&secR);
31  tmvaReader_->AddVariable("sTIP",&sTIP);
32  tmvaReader_->AddVariable("nHITS1",&nHITS1);
33  tmvaReader_->AddVariable("secPin",&secPin);
34  tmvaReader_->AddVariable("Epout",&Epout);
35  tmvaReader_->AddVariable("detaBremKF",&detaBremKF);
36  tmvaReader_->AddVariable("ptRatioGsfKF",&ptRatioGsfKF);
37  tmvaReader_->BookMVA("BDT",mvaWeightFileConvBrem.c_str());
38 }
TransientTrackBuilder builder_
ConvBremPFTrackFinder::~ConvBremPFTrackFinder ( )

Definition at line 39 of file ConvBremPFTrackFinder.cc.

References tmvaReader_.

39 {delete tmvaReader_;}

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,
bool  useNuclear,
bool  useConversions,
bool  useV0,
const reco::PFClusterCollection theEClus,
reco::GsfPFRecTrack  gsfpfrectk 
)
inline

Definition at line 36 of file ConvBremPFTrackFinder.h.

References found_, and runConvBremFinder().

Referenced by PFElecTkProducer::produce().

46  {
47  found_ = false;
48  runConvBremFinder(thePfRecTrackCol,primaryVertex,
49  pfNuclears,pfConversions,
50  pfV0,useNuclear,
52  theEClus,gsfpfrectk);
53  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, bool useNuclear, bool useConversions, bool useV0, const reco::PFClusterCollection &theEClus, reco::GsfPFRecTrack gsfpfrectk)
const std::vector<reco::PFRecTrackRef>& ConvBremPFTrackFinder::getConvBremPFRecTracks ( )
inline

Definition at line 56 of file ConvBremPFTrackFinder.h.

References pfRecTrRef_vec_.

Referenced by PFElecTkProducer::produce().

56 {return pfRecTrRef_vec_;};
std::vector< reco::PFRecTrackRef > pfRecTrRef_vec_
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,
bool  useNuclear,
bool  useConversions,
bool  useV0,
const reco::PFClusterCollection theEClus,
reco::GsfPFRecTrack  gsfpfrectk 
)
private

Definition at line 42 of file ConvBremPFTrackFinder.cc.

References TrackingRecHit::all, TransientTrackBuilder::build(), builder_, reco::PFCluster::calculatePositionREP(), reco::PFTrack::calculatePositionREP(), gather_cfg::cout, debug, detaBremKF, ExpressReco_HICollisions_FallBack::e, reco::PFTrajectoryPoint::ECALEntrance, PFEnergyCalibration::energyEm(), Epout, eta(), found_, reco::GsfPFRecTrack::gsfTrackRef(), i, edm::Ref< C, T, F >::isNonnull(), reco::PFTrajectoryPoint::isValid(), reco::GsfPFRecTrack::kfPFRecTrackRef(), reco::PFTrajectoryPoint::momentum(), mvaBremConvCut_, nHITS1, notFound, L1TEmulatorMonitor_cff::p, reco::GsfPFRecTrack::PFRecBrem(), pfRecTrRef_vec_, phi, Pi, edm::Handle< T >::product(), ptRatioGsfKF, secPin, secPout, secR, IPTools::signedTransverseImpactParameter(), findQualityFiles::size, mathSSE::sqrt(), sTIP, LinkByRecHit::testTrackAndClusterByRecHit(), tmvaReader_, and TwoPi.

Referenced by foundConvBremPFRecTrack().

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

Member Data Documentation

TransientTrackBuilder ConvBremPFTrackFinder::builder_
private

Definition at line 73 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::detaBremKF
private

Definition at line 78 of file ConvBremPFTrackFinder.h.

Referenced by ConvBremPFTrackFinder(), and runConvBremFinder().

float ConvBremPFTrackFinder::Epout
private

Definition at line 78 of file ConvBremPFTrackFinder.h.

Referenced by ConvBremPFTrackFinder(), and runConvBremFinder().

bool ConvBremPFTrackFinder::found_
private

Definition at line 72 of file ConvBremPFTrackFinder.h.

Referenced by foundConvBremPFRecTrack(), and runConvBremFinder().

double ConvBremPFTrackFinder::mvaBremConvCut_
private

Definition at line 74 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

std::string ConvBremPFTrackFinder::mvaWeightFileConvBrem_
private

Definition at line 75 of file ConvBremPFTrackFinder.h.

float ConvBremPFTrackFinder::nHITS1
private

Definition at line 80 of file ConvBremPFTrackFinder.h.

Referenced by ConvBremPFTrackFinder(), and runConvBremFinder().

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

Definition at line 77 of file ConvBremPFTrackFinder.h.

Referenced by getConvBremPFRecTracks(), and runConvBremFinder().

float ConvBremPFTrackFinder::ptRatioGsfKF
private

Definition at line 78 of file ConvBremPFTrackFinder.h.

Referenced by ConvBremPFTrackFinder(), and runConvBremFinder().

float ConvBremPFTrackFinder::secPin
private

Definition at line 78 of file ConvBremPFTrackFinder.h.

Referenced by ConvBremPFTrackFinder(), and runConvBremFinder().

float ConvBremPFTrackFinder::secPout
private

Definition at line 78 of file ConvBremPFTrackFinder.h.

Referenced by runConvBremFinder().

float ConvBremPFTrackFinder::secR
private

Definition at line 78 of file ConvBremPFTrackFinder.h.

Referenced by ConvBremPFTrackFinder(), and runConvBremFinder().

float ConvBremPFTrackFinder::sTIP
private

Definition at line 78 of file ConvBremPFTrackFinder.h.

Referenced by ConvBremPFTrackFinder(), and runConvBremFinder().

TMVA::Reader* ConvBremPFTrackFinder::tmvaReader_
private