CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups 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();
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 
101  : AlignmentMonitorBase(cfg, iC, "AlignmentMonitorTracksFromTrajectories"),
102  m_vertexConstraint(cfg.getParameter<bool>("vertexConstraint")),
103  m_beamSpot(cfg.getParameter<edm::InputTag>("beamSpot")) {
104  theMuonServiceProxy = new MuonServiceProxy(cfg.getParameter<edm::ParameterSet>("ServiceParameters"));
105  theMuonUpdatorAtVertex = new MuonUpdatorAtVertex(cfg.getParameter<edm::ParameterSet>("MuonUpdatorAtVertexParameters"),
107 }
108 
109 //
110 // member functions
111 //
112 
114 // book()
116 
118  m_diMuon_Z = book1D("/iterN/", "diMuon_Z", "Di-muon mass (GeV)", 200, 90. - 50., 90. + 50.);
119  m_diMuon_Zforward = book1D("/iterN/", "diMuon_Zforward", "Di-muon mass (GeV) eta > 1.4", 200, 90. - 50., 90. + 50.);
121  book1D("/iterN/", "diMuon_Zbarrel", "Di-muon mass (GeV) -1.4 < eta < 1.4", 200, 90. - 50., 90. + 50.);
123  book1D("/iterN/", "diMuon_Zbackward", "Di-muon mass (GeV) eta < -1.4", 200, 90. - 50., 90. + 50.);
124  m_diMuon_Ups = book1D("/iterN/", "diMuon_Ups", "Di-muon mass (GeV)", 200, 9. - 3., 9. + 3.);
125  m_diMuon_Jpsi = book1D("/iterN/", "diMuon_Jpsi", "Di-muon mass (GeV)", 200, 3. - 1., 3. + 1.);
126  m_diMuon_log = book1D("/iterN/", "diMuon_log", "Di-muon mass (log GeV)", 300, 0, 3);
127  m_chi2_100 = book1D("/iterN/", "m_chi2_100", "Track chi^2", 100, 0, 100);
128  m_chi2_10000 = book1D("/iterN/", "m_chi2_10000", "Track chi^2", 100, 0, 10000);
129  m_chi2_1000000 = book1D("/iterN/", "m_chi2_1000000", "Track chi^2", 100, 1, 1000000);
130  m_chi2_log = book1D("/iterN/", "m_chi2_log", "Log track chi^2", 100, -3, 7);
131  m_chi2DOF_5 = book1D("/iterN/", "m_chi2DOF_5", "Track chi^2/nDOF", 100, 0, 5);
132  m_chi2DOF_100 = book1D("/iterN/", "m_chi2DOF_100", "Track chi^2/nDOF", 100, 0, 100);
133  m_chi2DOF_1000 = book1D("/iterN/", "m_chi2DOF_1000", "Track chi^2/nDOF", 100, 0, 1000);
134  m_chi2DOF_log = book1D("/iterN/", "m_chi2DOF_log", "Log track chi^2/nDOF", 100, -3, 7);
135  m_chi2_improvement = book1D("/iterN/", "m_chi2_improvement", "Track-by-track chi^2/original chi^2", 100, 0., 10.);
137  "/iterN/", "m_chi2DOF_improvement", "Track-by-track (chi^2/DOF)/(original chi^2/original DOF)", 100, 0., 10.);
138  for (int i = 0; i < 36; i++) {
139  char name[100], title[100];
140  snprintf(name, sizeof(name), "m_pt_phi%d", i);
141  snprintf(title, sizeof(title), "Track pt (GeV) in phi bin %d/36", i);
142  m_pt[i] = book1D("/iterN/", name, title, 100, 0, 100);
143  }
144 }
145 
147 // event()
149 
151  const edm::EventSetup &iSetup,
153  theMuonServiceProxy->update(iSetup);
154 
156  iEvent.getByLabel(m_beamSpot, beamSpot);
157 
158  GlobalVector p1, p2;
159  double e1 = 0.;
160  double e2 = 0.;
161 
162  for (ConstTrajTrackPairCollection::const_iterator it = tracks.begin(); it != tracks.end(); ++it) {
163  const Trajectory *traj = it->first;
164  const reco::Track *track = it->second;
165 
166  int nDOF = 0;
167  TrajectoryStateOnSurface closestTSOS;
168  double closest = 10000.;
169 
170  std::vector<TrajectoryMeasurement> measurements = traj->measurements();
171  for (std::vector<TrajectoryMeasurement>::const_iterator im = measurements.begin(); im != measurements.end(); ++im) {
172  const TrajectoryMeasurement meas = *im;
173  const TransientTrackingRecHit *hit = &(*meas.recHit());
174  // const DetId id = hit->geographicalId();
175 
176  nDOF += hit->dimension();
177 
179  GlobalPoint where = TSOS.surface().toGlobal(LocalPoint(0, 0, 0));
180 
181  if (where.mag() < closest) {
182  closest = where.mag();
183  closestTSOS = TSOS;
184  }
185 
186  } // end loop over measurements
187  nDOF -= 5;
188 
189  if (closest != 10000.) {
190  std::pair<bool, FreeTrajectoryState> state;
191 
192  if (m_vertexConstraint) {
193  state = theMuonUpdatorAtVertex->propagateWithUpdate(closestTSOS, *beamSpot);
194  // add in chi^2 contribution from vertex contratint?
195  } else {
196  state = theMuonUpdatorAtVertex->propagate(closestTSOS, *beamSpot);
197  }
198 
199  if (state.first) {
200  double chi2 = traj->chiSquared();
201  double chi2DOF = chi2 / double(nDOF);
202 
203  m_chi2_100->Fill(chi2);
204  m_chi2_10000->Fill(chi2);
205  m_chi2_1000000->Fill(chi2);
206  m_chi2_log->Fill(log10(chi2));
207 
208  m_chi2DOF_5->Fill(chi2DOF);
209  m_chi2DOF_100->Fill(chi2DOF);
210  m_chi2DOF_1000->Fill(chi2DOF);
211  m_chi2DOF_log->Fill(log10(chi2DOF));
212  m_chi2_improvement->Fill(chi2 / track->chi2());
213  m_chi2DOF_improvement->Fill(chi2DOF / track->normalizedChi2());
214 
215  GlobalVector momentum = state.second.momentum();
216  double energy = momentum.mag();
217 
218  if (energy > e1) {
219  e2 = e1;
220  e1 = energy;
221  p2 = p1;
222  p1 = momentum;
223  } else if (energy > e2) {
224  e2 = energy;
225  p2 = momentum;
226  }
227 
228  double pt = momentum.perp();
229  double phi = momentum.phi();
230  while (phi < -M_PI)
231  phi += 2. * M_PI;
232  while (phi > M_PI)
233  phi -= 2. * M_PI;
234 
235  double phibin = (phi + M_PI) / (2. * M_PI) * 36.;
236 
237  for (int i = 0; i < 36; i++) {
238  if (phibin < i + 1) {
239  m_pt[i]->Fill(pt);
240  break;
241  }
242  }
243 
244  } // end if propagate was successful
245  } // end if we have a closest TSOS
246  } // end loop over tracks
247 
248  if (e1 > 0. && e2 > 0.) {
249  double energy_tot = e1 + e2;
250  GlobalVector momentum_tot = p1 + p2;
251  double mass = sqrt(energy_tot * energy_tot - momentum_tot.mag2());
252  double eta = momentum_tot.eta();
253 
254  m_diMuon_Z->Fill(mass);
255  if (eta > 1.4)
256  m_diMuon_Zforward->Fill(mass);
257  else if (eta > -1.4)
258  m_diMuon_Zbarrel->Fill(mass);
259  else
260  m_diMuon_Zbackward->Fill(mass);
261 
262  m_diMuon_Ups->Fill(mass);
263  m_diMuon_Jpsi->Fill(mass);
264  m_diMuon_log->Fill(log10(mass));
265  } // end if we have two tracks
266 }
267 
269 // afterAlignment()
271 
273 
274 //
275 // const member functions
276 //
277 
278 //
279 // static member functions
280 //
281 
282 //
283 // SEAL definitions
284 //
285 
288  "AlignmentMonitorTracksFromTrajectories");
GlobalPoint toGlobal(const Point2DBase< Scalar, LocalTag > lp) const
Definition: Surface.h:79
virtual int dimension() const =0
AlignmentMonitorTracksFromTrajectories(const edm::ParameterSet &cfg, edm::ConsumesCollector iC)
T mag2() const
Definition: PV3DBase.h:63
tuple cfg
Definition: looper.py:296
Point3DBase< Scalar, LocalTag > LocalPoint
Definition: Definitions.h:30
T perp() const
Definition: PV3DBase.h:69
ConstRecHitPointer const & recHit() const
const TString p2
Definition: fwPaths.cc:13
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:593
void book()
Book or retrieve histograms; MUST be reimplemented.
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
std::pair< bool, FreeTrajectoryState > propagate(const TrajectoryStateOnSurface &tsos, const reco::BeamSpot &beamSpot) const
Propagate the state to the 2D-PCA.
auto const & tracks
cannot be loose
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:178
int iEvent
Definition: GenABIO.cc:224
const SurfaceType & surface() const
T mag() const
Definition: PV3DBase.h:64
double chi2() const
chi-squared of the fit
Definition: TrackBase.h:587
T sqrt(T t)
Definition: SSEVec.h:19
const TString p1
Definition: fwPaths.cc:12
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
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:500
#define M_PI
TrajectoryStateOnSurface TSOS
Definition: TestHits.cc:19
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
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:241
T eta() const
Definition: PV3DBase.h:73
void update(const edm::EventSetup &setup, bool duringEvent=true)
update the services each event
#define DEFINE_EDM_PLUGIN(factory, type, name)
TrajectoryStateOnSurface const & backwardPredictedState() const
Access to backward predicted state (from smoother)