CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentMonitorTracksFromTrajectories.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CommonAlignmentProducer
4 // Class : AlignmentMonitorTracksFromTrajectories
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: Jim Pivarski
10 // Created: Thu Jun 28 01:38:33 CDT 2007
11 //
12 
13 // system include files
19 #include "TH1.h"
20 #include "TObject.h"
22 
26 
27 #include <fstream>
28 
29 // user include files
30 
31 //
32 // class definition
33 //
34 
36  public:
39 
40  void book();
41  void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& iTrajTracks);
42  void afterAlignment(const edm::EventSetup &iSetup);
43 
44  private:
49 
50  TH1F *m_diMuon_Z;
54  TH1F *m_diMuon_Ups;
56  TH1F *m_diMuon_log;
57  TH1F *m_chi2_100;
58  TH1F *m_chi2_10000;
60  TH1F *m_chi2_log;
61  TH1F *m_chi2DOF_5;
67  TH1F *m_pt[36];
68 };
69 
70 //
71 // constants, enums and typedefs
72 //
73 
74 //
75 // static data member definitions
76 //
77 
78 //
79 // constructors and destructor
80 //
81 
82 // AlignmentMonitorTracksFromTrajectories::AlignmentMonitorTracksFromTrajectories(const AlignmentMonitorTracksFromTrajectories& rhs)
83 // {
84 // // do actual copying here;
85 // }
86 
87 //
88 // assignment operators
89 //
90 // const AlignmentMonitorTracksFromTrajectories& AlignmentMonitorTracksFromTrajectories::operator=(const AlignmentMonitorTracksFromTrajectories& rhs)
91 // {
92 // //An exception safe implementation is
93 // AlignmentMonitorTracksFromTrajectories temp(rhs);
94 // swap(rhs);
95 //
96 // return *this;
97 // }
98 
100  : AlignmentMonitorBase(cfg, "AlignmentMonitorTracksFromTrajectories")
101  , m_vertexConstraint(cfg.getParameter<bool>("vertexConstraint"))
102  , m_beamSpot(cfg.getParameter<edm::InputTag>("beamSpot"))
103 {
104  theMuonServiceProxy = new MuonServiceProxy(cfg.getParameter<edm::ParameterSet>("ServiceParameters"));
106 }
107 
108 //
109 // member functions
110 //
111 
113 // book()
115 
117  m_diMuon_Z = book1D("/iterN/", "diMuon_Z", "Di-muon mass (GeV)", 200, 90. - 50., 90. + 50.);
118  m_diMuon_Zforward = book1D("/iterN/", "diMuon_Zforward", "Di-muon mass (GeV) eta > 1.4", 200, 90. - 50., 90. + 50.);
119  m_diMuon_Zbarrel = book1D("/iterN/", "diMuon_Zbarrel", "Di-muon mass (GeV) -1.4 < eta < 1.4", 200, 90. - 50., 90. + 50.);
120  m_diMuon_Zbackward = book1D("/iterN/", "diMuon_Zbackward", "Di-muon mass (GeV) eta < -1.4", 200, 90. - 50., 90. + 50.);
121  m_diMuon_Ups = book1D("/iterN/", "diMuon_Ups", "Di-muon mass (GeV)", 200, 9. - 3., 9. + 3.);
122  m_diMuon_Jpsi = book1D("/iterN/", "diMuon_Jpsi", "Di-muon mass (GeV)", 200, 3. - 1., 3. + 1.);
123  m_diMuon_log = book1D("/iterN/", "diMuon_log", "Di-muon mass (log GeV)", 300, 0, 3);
124  m_chi2_100 = book1D("/iterN/", "m_chi2_100", "Track chi^2", 100, 0, 100);
125  m_chi2_10000 = book1D("/iterN/", "m_chi2_10000", "Track chi^2", 100, 0, 10000);
126  m_chi2_1000000 = book1D("/iterN/", "m_chi2_1000000", "Track chi^2", 100, 1, 1000000);
127  m_chi2_log = book1D("/iterN/", "m_chi2_log", "Log track chi^2", 100, -3, 7);
128  m_chi2DOF_5 = book1D("/iterN/", "m_chi2DOF_5", "Track chi^2/nDOF", 100, 0, 5);
129  m_chi2DOF_100 = book1D("/iterN/", "m_chi2DOF_100", "Track chi^2/nDOF", 100, 0, 100);
130  m_chi2DOF_1000 = book1D("/iterN/", "m_chi2DOF_1000", "Track chi^2/nDOF", 100, 0, 1000);
131  m_chi2DOF_log = book1D("/iterN/", "m_chi2DOF_log", "Log track chi^2/nDOF", 100, -3, 7);
132  m_chi2_improvement = book1D("/iterN/", "m_chi2_improvement", "Track-by-track chi^2/original chi^2", 100, 0., 10.);
133  m_chi2DOF_improvement = book1D("/iterN/", "m_chi2DOF_improvement", "Track-by-track (chi^2/DOF)/(original chi^2/original DOF)", 100, 0., 10.);
134  for (int i = 0; i < 36; i++) {
135  char name[100], title[100];
136  snprintf(name, sizeof(name), "m_pt_phi%d", i);
137  snprintf(title, sizeof(title), "Track pt (GeV) in phi bin %d/36", i);
138  m_pt[i] = book1D("/iterN/", name, title, 100, 0, 100);
139  }
140 }
141 
143 // event()
145 
147  theMuonServiceProxy->update(iSetup);
148 
150  iEvent.getByLabel(m_beamSpot, beamSpot);
151 
152  GlobalVector p1, p2;
153  double e1 = 0.;
154  double e2 = 0.;
155 
156  for (ConstTrajTrackPairCollection::const_iterator it = tracks.begin(); it != tracks.end(); ++it) {
157  const Trajectory *traj = it->first;
158  const reco::Track *track = it->second;
159 
160  int nDOF = 0;
161  TrajectoryStateOnSurface closestTSOS;
162  double closest = 10000.;
163 
164  std::vector<TrajectoryMeasurement> measurements = traj->measurements();
165  for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end(); ++im) {
166  const TrajectoryMeasurement meas = *im;
167  const TransientTrackingRecHit* hit = &(*meas.recHit());
168 // const DetId id = hit->geographicalId();
169 
170  nDOF += hit->dimension();
171 
173  GlobalPoint where = TSOS.surface().toGlobal(LocalPoint(0,0,0));
174 
175  if (where.mag() < closest) {
176  closest = where.mag();
177  closestTSOS = TSOS;
178  }
179 
180  } // end loop over measurements
181  nDOF -= 5;
182 
183  if (closest != 10000.) {
184  std::pair<bool, FreeTrajectoryState> state;
185 
186  if (m_vertexConstraint) {
187  state = theMuonUpdatorAtVertex->propagateWithUpdate(closestTSOS, *beamSpot);
188  // add in chi^2 contribution from vertex contratint?
189  }
190  else {
191  state = theMuonUpdatorAtVertex->propagate(closestTSOS, *beamSpot);
192  }
193 
194  if (state.first) {
195  double chi2 = traj->chiSquared();
196  double chi2DOF = chi2 / double(nDOF);
197 
198  m_chi2_100->Fill(chi2);
199  m_chi2_10000->Fill(chi2);
200  m_chi2_1000000->Fill(chi2);
201  m_chi2_log->Fill(log10(chi2));
202 
203  m_chi2DOF_5->Fill(chi2DOF);
204  m_chi2DOF_100->Fill(chi2DOF);
205  m_chi2DOF_1000->Fill(chi2DOF);
206  m_chi2DOF_log->Fill(log10(chi2DOF));
207  m_chi2_improvement->Fill(chi2 / track->chi2());
208  m_chi2DOF_improvement->Fill(chi2DOF / track->normalizedChi2());
209 
210  GlobalVector momentum = state.second.momentum();
211  double energy = momentum.mag();
212 
213  if (energy > e1) {
214  e2 = e1;
215  e1 = energy;
216  p2 = p1;
217  p1 = momentum;
218  }
219  else if (energy > e2) {
220  e2 = energy;
221  p2 = momentum;
222  }
223 
224  double pt = momentum.perp();
225  double phi = momentum.phi();
226  while (phi < -M_PI) phi += 2.* M_PI;
227  while (phi > M_PI) phi -= 2.* M_PI;
228 
229  double phibin = (phi + M_PI) / (2. * M_PI) * 36.;
230 
231  for (int i = 0; i < 36; i++) {
232  if (phibin < i+1) {
233  m_pt[i]->Fill(pt);
234  break;
235  }
236  }
237 
238  } // end if propagate was successful
239  } // end if we have a closest TSOS
240  } // end loop over tracks
241 
242  if (e1 > 0. && e2 > 0.) {
243  double energy_tot = e1 + e2;
244  GlobalVector momentum_tot = p1 + p2;
245  double mass = sqrt(energy_tot*energy_tot - momentum_tot.mag2());
246  double eta = momentum_tot.eta();
247 
248  m_diMuon_Z->Fill(mass);
249  if (eta > 1.4) m_diMuon_Zforward->Fill(mass);
250  else if (eta > -1.4) m_diMuon_Zbarrel->Fill(mass);
251  else m_diMuon_Zbackward->Fill(mass);
252 
253  m_diMuon_Ups->Fill(mass);
254  m_diMuon_Jpsi->Fill(mass);
255  m_diMuon_log->Fill(log10(mass));
256  } // end if we have two tracks
257 }
258 
260 // afterAlignment()
262 
264 }
265 
266 //
267 // const member functions
268 //
269 
270 //
271 // static member functions
272 //
273 
274 //
275 // SEAL definitions
276 //
277 
void update(const edm::EventSetup &setup)
update the services each event
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:114
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual int dimension() const =0
T mag2() const
Definition: PV3DBase.h:66
T perp() const
Definition: PV3DBase.h:72
ConstRecHitPointer const & recHit() const
double normalizedChi2() const
chi-squared divided by n.d.o.f. (or chi-squared * 1e6 if n.d.o.f. is zero)
Definition: TrackBase.h:609
void book()
Book or retrieve histograms; MUST be reimplemented.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
std::pair< bool, FreeTrajectoryState > propagate(const TrajectoryStateOnSurface &tsos, const reco::BeamSpot &beamSpot) const
Propagate the state to the 2D-PCA.
T eta() const
void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks)
Called for each event (by &quot;run()&quot;): may be reimplemented.
DataContainer const & measurements() const
Definition: Trajectory.h:203
int iEvent
Definition: GenABIO.cc:230
const SurfaceType & surface() const
T mag() const
Definition: PV3DBase.h:67
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:597
T sqrt(T t)
Definition: SSEVec.h:48
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
double p2[4]
Definition: TauolaWrapper.h:90
TH1F * book1D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:402
#define M_PI
TrajectoryStateOnSurface TSOS
Definition: TestHits.cc:19
tuple tracks
Definition: testEve_cfg.py:39
std::pair< bool, FreeTrajectoryState > propagateWithUpdate(const TrajectoryStateOnSurface &tsos, const reco::BeamSpot &beamSpot) const
Propagate to the 2D-PCA and apply the vertex constraint.
float chiSquared() const
Definition: Trajectory.h:252
T eta() const
Definition: PV3DBase.h:76
double p1[4]
Definition: TauolaWrapper.h:89
Local3DPoint LocalPoint
Definition: LocalPoint.h:11
#define DEFINE_EDM_PLUGIN(factory, type, name)
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)
Definition: DDAxes.h:10