CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
TrackDetectorAssociator.h
Go to the documentation of this file.
1 #ifndef TrackAssociator_TrackDetectorAssociator_h
2 #define TrackAssociator_TrackDetectorAssociator_h 1
3 
4 // -*- C++ -*-
5 //
6 // Package: TrackAssociator
7 // Class: TrackDetectorAssociator
8 //
9 /*
10 
11  Description: main class of tools to associate a track to calorimeter and muon detectors
12 
13 */
14 //
15 // Original Author: Dmytro Kovalskyi
16 // Created: Fri Apr 21 10:59:41 PDT 2006
17 //
18 //
19 
26 
32 
36 
41 
44 
46 
48 
49 
51  public:
54 
57 
72  const edm::EventSetup&,
73  const FreeTrajectoryState&,
74  const AssociatorParameters& );
78  const edm::EventSetup& iSetup,
80  const FreeTrajectoryState* innerState,
81  const FreeTrajectoryState* outerState=0);
84  const edm::EventSetup&,
85  const reco::Track&,
86  const AssociatorParameters&,
87  Direction direction = Any );
90  const edm::EventSetup&,
91  const SimTrack&,
92  const SimVertex&,
93  const AssociatorParameters& );
96  const edm::EventSetup&,
97  const GlobalVector&,
98  const GlobalPoint&,
99  const int,
100  const AssociatorParameters& );
103  {return cachedTrajectory_;}
104 
106  void setPropagator( const Propagator* );
107 
109  void useDefaultPropagator();
110 
113  const reco::Track& );
115  const SimTrack&,
116  const SimVertex& );
118  const GlobalVector&,
119  const GlobalPoint&,
120  const int);
121 
122  static bool crossedIP(const reco::Track& track);
123 
124  private:
125  DetIdAssociator::MapRange getMapRange( const std::pair<float,float>& delta,
126  const float dR ) dso_internal;
127 
128  void fillEcal( const edm::Event&,
131 
132  void fillCaloTowers( const edm::Event&,
134  const AssociatorParameters&) dso_internal;
135 
136  void fillHcal( const edm::Event&,
137  TrackDetMatchInfo&,
138  const AssociatorParameters&) dso_internal;
139 
140  void fillHO( const edm::Event&,
141  TrackDetMatchInfo&,
142  const AssociatorParameters&) dso_internal;
143 
144  void fillPreshower( const edm::Event& iEvent,
145  TrackDetMatchInfo& info,
146  const AssociatorParameters&) dso_internal;
147 
148  void fillMuon( const edm::Event&,
149  TrackDetMatchInfo&,
150  const AssociatorParameters&) dso_internal;
151 
152  void fillCaloTruth( const edm::Event&,
153  TrackDetMatchInfo&,
154  const AssociatorParameters&) dso_internal;
155 
157  const RecSegment*,
158  const AssociatorParameters&) dso_internal;
159 
160  void getTAMuonChamberMatches(std::vector<TAMuonChamberMatch>& matches,
161  const AssociatorParameters& parameters) dso_internal;
162 
163  void init( const edm::EventSetup&) dso_internal;
164 
165  math::XYZPoint getPoint( const GlobalPoint& point) dso_internal
166  {
167  return math::XYZPoint(point.x(),point.y(),point.z());
168  }
169 
171  {
172  return math::XYZPoint(point.x(),point.y(),point.z());
173  }
174 
176  {
177  return math::XYZVector(vec.x(),vec.y(),vec.z());
178  }
179 
181  {
182  return math::XYZVector(vec.x(),vec.y(),vec.z());
183  }
184 
189 
196 
199 
201 };
202 #endif
dbl * delta
Definition: mlp_gen.cc:36
dictionary parameters
Definition: Parameters.py:2
static const TGPicture * info(bool iBackgroundIsBlack)
static bool crossedIP(const reco::Track &track)
static FreeTrajectoryState getFreeTrajectoryState(const edm::EventSetup &, const reco::Track &)
get FreeTrajectoryState from different track representations
void useDefaultPropagator()
use the default propagator
edm::ESWatcher< IdealMagneticFieldRecord > theMagneticFeildWatcher_
math::XYZVector getVector(const GlobalVector &vec)
void fillHcal(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
bool addTAMuonSegmentMatch(TAMuonChamberMatch &, const RecSegment *, const AssociatorParameters &)
math::XYZPoint getPoint(const LocalPoint &point)
const CachedTrajectory & getCachedTrajector() const
trajector information
edm::ESHandle< DetIdAssociator > hcalDetIdAssociator_
void fillCaloTowers(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
void fillCaloTruth(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
int iEvent
Definition: GenABIO.cc:230
void setPropagator(const Propagator *)
use a user configured propagator
void getTAMuonChamberMatches(std::vector< TAMuonChamberMatch > &matches, const AssociatorParameters &parameters)
TrackAssociatorParameters AssociatorParameters
edm::ESHandle< CaloGeometry > theCaloGeometry_
edm::ESHandle< DetIdAssociator > ecalDetIdAssociator_
void fillPreshower(const edm::Event &iEvent, TrackDetMatchInfo &info, const AssociatorParameters &)
math::XYZVector getVector(const LocalVector &vec)
math::XYZPoint getPoint(const GlobalPoint &point)
XYZVectorD XYZVector
spatial vector with cartesian internal representation
Definition: Vector3D.h:30
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::ESHandle< DetIdAssociator > hoDetIdAssociator_
void fillMuon(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
string const
Definition: compareJSON.py:14
void init(const edm::EventSetup &)
DetIdAssociator::MapRange getMapRange(const std::pair< float, float > &delta, const float dR)
#define dso_internal
edm::ESHandle< DetIdAssociator > muonDetIdAssociator_
math::XYZVector XYZPoint
edm::ESHandle< DetIdAssociator > preshowerDetIdAssociator_
void fillEcal(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
TrackDetMatchInfo associate(const edm::Event &, const edm::EventSetup &, const FreeTrajectoryState &, const AssociatorParameters &)
edm::ESHandle< DetIdAssociator > caloDetIdAssociator_
void fillHO(const edm::Event &, TrackDetMatchInfo &, const AssociatorParameters &)
edm::ESHandle< GlobalTrackingGeometry > theTrackingGeometry_
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5