CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CaloPropagateTrack.h
Go to the documentation of this file.
1 #ifndef CalibrationIsolatedParticlesCaloPropagateTrack_h
2 #define CalibrationIsolatedParticlesCaloPropagateTrack_h
3 
4 #include <cmath>
5 #include <vector>
6 #include <string>
7 
15 
16 //sim track
23 
27 
29 
30 namespace spr{
31 
32  struct propagatedTrack {
33  propagatedTrack() {ok=false;}
34  bool ok;
37  };
38 
40  propagatedTrackID() {ok=false; okECAL=false; okHCAL=false;}
41  bool ok, okECAL, okHCAL;
43  double etaECAL, etaHCAL;
44  double phiECAL, phiHCAL;
45  reco::TrackCollection::const_iterator trkItr;
46  };
47 
49  propagatedTrackDirection() {ok=false; okECAL=false; okHCAL=false;}
50  bool ok, okECAL, okHCAL;
54  reco::TrackCollection::const_iterator trkItr;
55  };
56 
59  ok=okECAL=okHCAL=false;
60  charge=pdgId=0;
61  }
62  bool ok, okECAL, okHCAL;
66  int charge, pdgId;
67  HepMC::GenEvent::particle_const_iterator trkItr;
68  };
69 
72  ok=okECAL=okHCAL=false;
73  charge=pdgId=0;
74  }
75  bool ok, okECAL, okHCAL;
79  int charge, pdgId;
80  reco::GenParticleCollection::const_iterator trkItr;
81  };
82 
83  struct trackAtOrigin {
84  trackAtOrigin() {ok=false;}
85  bool ok;
86  int charge;
89  };
90 
91  // Returns a vector of DetID's of closest cell on the ECAL/HCAL surface of
92  // all the tracks in the collection. Also saves a boolean if extrapolation
93  // is satisfactory
94  std::vector<spr::propagatedTrackID> propagateCosmicCALO(edm::Handle<reco::TrackCollection>& trkCollection, const CaloGeometry* geo, const MagneticField* bField, std::string & theTrackQuality, bool debug=false);
95  std::vector<spr::propagatedTrackID> propagateCALO(edm::Handle<reco::TrackCollection>& trkCollection, const CaloGeometry* geo, const MagneticField* bField, std::string & theTrackQuality, bool debug=false);
96  void propagateCALO(edm::Handle<reco::TrackCollection>& trkCollection, const CaloGeometry* geo, const MagneticField* bField, std::string & theTrackQuality, std::vector<spr::propagatedTrackID>& vdets, bool debug=false);
97  void propagateCALO(edm::Handle<reco::TrackCollection>& trkCollection, const CaloGeometry* geo, const MagneticField* bField, std::string & theTrackQuality, std::vector<spr::propagatedTrackDirection>& trkDir, bool debug=false);
99  std::vector<spr::propagatedGenTrackID> propagateCALO(const HepMC::GenEvent * genEvent, edm::ESHandle<ParticleDataTable>& pdt, const CaloGeometry* geo, const MagneticField* bField, double etaMax=3.0, bool debug=false);
100  std::vector<spr::propagatedGenParticleID> propagateCALO(edm::Handle<reco::GenParticleCollection>& genParticles, edm::ESHandle<ParticleDataTable>& pdt, const CaloGeometry* geo, const MagneticField* bField, double etaMax=3.0, bool debug=false);
102 
103  // Propagate tracks to the ECAL surface and optionally returns the
104  // extrapolated point (and the track direction at point of extrapolation)
107  std::pair<math::XYZPoint,bool> propagateECAL(const reco::Track*, const MagneticField*, bool debug=false);
108  std::pair<math::XYZPoint,bool> propagateECAL(const GlobalPoint& vertex, const GlobalVector& momentum, int charge, const MagneticField*, bool debug=false);
109 
110  // Propagate tracks to the HCAL surface and optionally returns the
111  // extrapolated point (and the track direction at point of extrapolation)
114  std::pair<math::XYZPoint,bool> propagateHCAL(const reco::Track*, const MagneticField*, bool debug=false);
115  std::pair<math::XYZPoint,bool> propagateHCAL(const GlobalPoint& vertex, const GlobalVector& momentum, int charge, const MagneticField*, bool debug=false);
116 
117  // Propagate the track to the end of the tracker and returns the extrapolated
118  // point and optionally the length of the track upto the end
119  std::pair<math::XYZPoint,bool> propagateTracker(const reco::Track*, const MagneticField*, bool debug=false);
120  std::pair<math::XYZPoint,double> propagateTrackerEnd(const reco::Track*, const MagneticField*, bool debug=false);
121 
122  spr::propagatedTrack propagateCalo(const GlobalPoint& vertex, const GlobalVector& momentum, int charge, const MagneticField*, float zdist, float radius, float corner, bool debug=false);
123 
124  // Gives the vertex and momentum of a SimTrack
126 
127 }
128 #endif
std::vector< spr::propagatedTrackID > propagateCALO(edm::Handle< reco::TrackCollection > &trkCollection, const CaloGeometry *geo, const MagneticField *bField, std::string &theTrackQuality, bool debug=false)
reco::GenParticleCollection::const_iterator trkItr
std::pair< math::XYZPoint, bool > propagateHCAL(const reco::Track *, const MagneticField *, bool debug=false)
spr::trackAtOrigin simTrackAtOrigin(unsigned int thisTrk, edm::Handle< edm::SimTrackContainer > &SimTk, edm::Handle< edm::SimVertexContainer > &SimVtx, bool debug=false)
spr::propagatedTrack propagateCalo(const GlobalPoint &vertex, const GlobalVector &momentum, int charge, const MagneticField *, float zdist, float radius, float corner, bool debug=false)
HepMC::GenEvent::particle_const_iterator trkItr
reco::TrackCollection::const_iterator trkItr
tuple genEvent
Definition: MCTruth.py:33
std::pair< math::XYZPoint, bool > propagateTracker(const reco::Track *, const MagneticField *, bool debug=false)
std::pair< math::XYZPoint, bool > propagateECAL(const reco::Track *, const MagneticField *, bool debug=false)
spr::propagatedTrack propagateTrackToECAL(const reco::Track *, const MagneticField *, bool debug=false)
Definition: DetId.h:18
#define debug
Definition: HDRShower.cc:19
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::pair< math::XYZPoint, double > propagateTrackerEnd(const reco::Track *, const MagneticField *, bool debug=false)
std::vector< spr::propagatedTrackID > propagateCosmicCALO(edm::Handle< reco::TrackCollection > &trkCollection, const CaloGeometry *geo, const MagneticField *bField, std::string &theTrackQuality, bool debug=false)
reco::TrackCollection::const_iterator trkItr
spr::propagatedTrack propagateTrackToHCAL(const reco::Track *, const MagneticField *, bool debug=false)