CMS 3D CMS Logo

AlignmentMonitorMuonSystemMap1D.cc
Go to the documentation of this file.
1 /*
2  * Package: CommonAlignmentProducer
3  * Class : AlignmentMonitorMuonSystemMap1D
4  *
5  * Original Author: Jim Pivarski
6  * Created: Mon Nov 12 13:30:14 CST 2007
7  *
8  * $Id: AlignmentMonitorMuonSystemMap1D.cc,v 1.5 2011/04/15 23:09:37 khotilov Exp $
9  */
10 
14 
23 
24 #include "TH1F.h"
25 #include "TH2F.h"
26 
33 
35 public:
38 
39  void book() override;
40 
41  void event(const edm::Event &iEvent,
42  const edm::EventSetup &iSetup,
43  const ConstTrajTrackPairCollection &iTrajTracks) override;
45 
46  void afterAlignment() override;
47 
48 private:
49  // parameters
51  double m_minTrackPt;
52  double m_maxTrackPt;
53  double m_minTrackP;
54  double m_maxTrackP;
55  double m_maxDxy;
63  bool m_doDT;
64  bool m_doCSC;
67 
68  // counter
79 
80  // histogram helper
82  public:
85  int bins,
86  double low,
87  double high,
88  bool xy,
89  bool add_1d);
90 
91  void fill_x_1d(double residx, double chi2, int dof);
92  void fill_x(char charge, double abscissa, double residx, double chi2, int dof);
93  void fill_y(char charge, double abscissa, double residy, double chi2, int dof);
94  void fill_dxdz(char charge, double abscissa, double slopex, double chi2, int dof);
95  void fill_dydz(char charge, double abscissa, double slopey, double chi2, int dof);
96 
97  private:
99  int m_bins;
100  bool m_xy;
101  bool m_1d;
102  TH1F *m_x_1d;
104  };
105 
106  MuonSystemMapPlot1D *m_DTvsz_station[4][14]; // [station][sector]
107  MuonSystemMapPlot1D *m_CSCvsr_me[2][4][36]; // [endcap][station][chamber]
108  MuonSystemMapPlot1D *m_DTvsphi_station[4][5]; // [station][wheel]
109  MuonSystemMapPlot1D *m_CSCvsphi_me[2][4][3]; // [endcap][station][ring]
110 
111  std::vector<MuonSystemMapPlot1D *> m_plots;
112 
113  std::string num02d(int num);
114 
115  // optional debug ntuple
116  TTree *m_cscnt;
117 
118  struct MyCSCDetId {
119  void init(CSCDetId &id) {
120  e = id.endcap();
121  s = id.station();
122  r = id.ring();
123  c = id.chamber();
124  t = id.iChamberType();
125  }
126  Short_t e, s, r, c;
127  Short_t t; // type 1-10: ME1/a,1/b,1/2,1/3,2/1...4/2
128  };
130 
131  struct MyTrack {
132  Int_t q;
133  Float_t pt, pz;
134  };
136 
137  struct MyResidual {
138  Float_t res, slope, rho, phi, z;
139  };
141 
142  UInt_t m_run;
143 };
144 
146  : AlignmentMonitorBase(cfg, "AlignmentMonitorMuonSystemMap1D"),
147  m_muonCollectionTag(cfg.getParameter<edm::InputTag>("muonCollectionTag")),
148  m_minTrackPt(cfg.getParameter<double>("minTrackPt")),
149  m_maxTrackPt(cfg.getParameter<double>("maxTrackPt")),
150  m_minTrackP(cfg.getParameter<double>("minTrackP")),
151  m_maxTrackP(cfg.getParameter<double>("maxTrackP")),
152  m_maxDxy(cfg.getParameter<double>("maxDxy")),
153  m_minTrackerHits(cfg.getParameter<int>("minTrackerHits")),
154  m_maxTrackerRedChi2(cfg.getParameter<double>("maxTrackerRedChi2")),
155  m_allowTIDTEC(cfg.getParameter<bool>("allowTIDTEC")),
156  m_minNCrossedChambers(cfg.getParameter<int>("minNCrossedChambers")),
157  m_minDT13Hits(cfg.getParameter<int>("minDT13Hits")),
158  m_minDT2Hits(cfg.getParameter<int>("minDT2Hits")),
159  m_minCSCHits(cfg.getParameter<int>("minCSCHits")),
160  m_doDT(cfg.getParameter<bool>("doDT")),
161  m_doCSC(cfg.getParameter<bool>("doCSC")),
162  m_useStubPosition(cfg.getParameter<bool>("useStubPosition")),
163  m_createNtuple(cfg.getParameter<bool>("createNtuple")) {
164  if (m_createNtuple) {
166  m_cscnt = fs->make<TTree>("mualNtuple", "mualNtuple");
167  m_cscnt->Branch("id", &m_id.e, "e/S:s:r:c:t");
168  m_cscnt->Branch("tr", &m_tr.q, "q/I:pt/F:pz");
169  m_cscnt->Branch("re", &m_re.res, "res/F:slope:rho:phi:z");
170  m_cscnt->Branch("run", &m_run, "run/i");
171  }
172 }
173 
175  assert(num >= 0 && num < 100);
176  char tmp[4];
177  sprintf(tmp, "%02d", num);
178  return std::string(tmp);
179 }
180 
182  std::string wheel_label[5] = {"A", "B", "C", "D", "E"};
183 
184  for (unsigned char station = 1; station <= 4; station++) {
185  std::string s_station = std::to_string(station);
186 
187  bool do_y = true;
188  if (station == 4)
189  do_y = false;
190 
191  // *** DT ***
192  if (m_doDT)
193  for (int sector = 1; sector <= 14; sector++) {
194  if ((station < 4 && sector <= 12) || station == 4) {
195  m_DTvsz_station[station - 1][sector - 1] = new MuonSystemMapPlot1D(
196  "DTvsz_st" + s_station + "sec" + num02d(sector), this, 60, -660., 660., do_y, false);
197  m_plots.push_back(m_DTvsz_station[station - 1][sector - 1]);
198  }
199  }
200 
201  if (m_doDT)
202  for (int wheel = -2; wheel <= 2; wheel++) {
204  "DTvsphi_st" + s_station + "wh" + wheel_label[wheel + 2], this, 180, -M_PI, M_PI, do_y, false);
205  m_plots.push_back(m_DTvsphi_station[station - 1][wheel + 2]);
206  }
207 
208  // *** CSC ***
209  if (m_doCSC)
210  for (int endcap = 1; endcap <= 2; endcap++) {
211  std::string s_endcap("m");
212  if (endcap == 1)
213  s_endcap = "p";
214 
215  for (int chamber = 1; chamber <= 36; chamber++) {
217  "CSCvsr_me" + s_endcap + s_station + "ch" + num02d(chamber), this, 60, 100., 700., false, false);
218  m_plots.push_back(m_CSCvsr_me[endcap - 1][station - 1][chamber - 1]);
219  }
220 
221  for (int ring = 1; ring <= 3; ring++) // the ME1/a (ring4) is not independent from ME1/b (ring1)
222  {
223  std::string s_ring = std::to_string(ring);
224  if ((station > 1 && ring <= 2) || station == 1) {
225  m_CSCvsphi_me[endcap - 1][station - 1][ring - 1] =
226  new MuonSystemMapPlot1D("CSCvsphi_me" + s_endcap + s_station + s_ring,
227  this,
228  180,
229  -M_PI / 180. * 5.,
230  M_PI * (2. - 5. / 180.),
231  false,
232  true);
233  m_plots.push_back(m_CSCvsphi_me[endcap - 1][station - 1][ring - 1]);
234  }
235  }
236  } // endcaps
237  } // stations
238 
239  m_counter_event = 0;
240  m_counter_track = 0;
242  m_counter_trackdxy = 0;
244  m_counter_dt = 0;
246  m_counter_2numhits = 0;
247  m_counter_csc = 0;
249 }
250 
252  const edm::EventSetup &iSetup,
253  const ConstTrajTrackPairCollection &trajtracks) {
254  m_counter_event++;
255 
257  iSetup.get<GlobalTrackingGeometryRecord>().get(globalGeometry);
258 
260  iEvent.getByLabel(m_beamSpotTag, beamSpot);
261 
262  edm::ESHandle<DetIdAssociator> muonDetIdAssociator_;
263  iSetup.get<DetIdAssociatorRecord>().get("MuonDetIdAssociator", muonDetIdAssociator_);
264 
266  iSetup.get<TrackingComponentsRecord>().get("SteppingHelixPropagatorAny", prop);
267 
270 
271  if (m_muonCollectionTag.label().empty()) // use trajectories
272  {
273  for (ConstTrajTrackPairCollection::const_iterator trajtrack = trajtracks.begin(); trajtrack != trajtracks.end();
274  ++trajtrack) {
275  m_counter_track++;
276  const Trajectory *traj = (*trajtrack).first;
277  const reco::Track *track = (*trajtrack).second;
278 
279  if (m_minTrackPt < track->pt() && track->pt() < m_maxTrackPt && m_minTrackP < track->p() &&
280  track->p() < m_maxTrackP) {
282  if (fabs(track->dxy(beamSpot->position())) < m_maxDxy) {
284 
285  MuonResidualsFromTrack muonResidualsFromTrack(
286  iSetup, magneticField, globalGeometry, muonDetIdAssociator_, prop, traj, track, pNavigator(), 1000.);
287  processMuonResidualsFromTrack(muonResidualsFromTrack, iEvent);
288  }
289  } // end if track has acceptable momentum
290  } // end loop over tracks
291  } else {
293  iEvent.getByLabel(m_muonCollectionTag, muons);
294 
295  for (reco::MuonCollection::const_iterator muon = muons->begin(); muon != muons->end(); ++muon) {
296  if (!(muon->isTrackerMuon() && muon->innerTrack().isNonnull()))
297  continue;
298 
299  m_counter_track++;
300 
301  if (m_minTrackPt < muon->pt() && muon->pt() < m_maxTrackPt && m_minTrackP < muon->p() &&
302  muon->p() < m_maxTrackP) {
304  if (fabs(muon->innerTrack()->dxy(beamSpot->position())) < m_maxDxy) {
306 
307  MuonResidualsFromTrack muonResidualsFromTrack(globalGeometry, &(*muon), pNavigator(), 100.);
308  processMuonResidualsFromTrack(muonResidualsFromTrack, iEvent);
309  }
310  }
311  }
312  }
313 }
314 
316  const edm::Event &iEvent) {
317  if (mrft.trackerNumHits() < m_minTrackerHits)
318  return;
319  if (!m_allowTIDTEC && mrft.contains_TIDTEC())
320  return;
321  if (mrft.normalizedChi2() > m_maxTrackerRedChi2)
322  return;
323 
324  int nMuChambers = 0;
325  std::vector<DetId> chamberIds = mrft.chamberIds();
326  for (unsigned ch = 0; ch < chamberIds.size(); ch++)
327  if (chamberIds[ch].det() == DetId::Muon)
328  nMuChambers++;
329  if (nMuChambers < m_minNCrossedChambers)
330  return;
331 
332  char charge = (mrft.getTrack()->charge() > 0 ? 1 : -1);
333  // double qoverpt = track->charge() / track->pt();
334  // double qoverpz = track->charge() / track->pz();
335 
337 
338  for (std::vector<DetId>::const_iterator chamberId = chamberIds.begin(); chamberId != chamberIds.end(); ++chamberId) {
339  if (chamberId->det() != DetId::Muon)
340  continue;
341 
342  if (m_doDT && chamberId->subdetId() == MuonSubdetId::DT) {
345  DTChamberId id(chamberId->rawId());
346 
347  m_counter_dt++;
348 
349  if (id.station() < 4 && dt13 != nullptr && dt13->numHits() >= m_minDT13Hits && dt2 != nullptr &&
350  dt2->numHits() >= m_minDT2Hits && (dt2->chi2() / double(dt2->ndof())) < 2.0) {
352 
353  double residual = dt13->global_residual();
354  double resslope = dt13->global_resslope();
355  double chi2 = dt13->chi2();
356  int dof = dt13->ndof();
357 
358  align::GlobalPoint gpos;
359  if (m_useStubPosition)
360  gpos = dt13->global_stubpos();
361  else
362  gpos = dt13->global_trackpos();
363  double phi = atan2(gpos.y(), gpos.x());
364  double z = gpos.z();
365 
366  assert(1 <= id.sector() && id.sector() <= 14);
367 
368  m_DTvsz_station[id.station() - 1][id.sector() - 1]->fill_x(charge, z, residual, chi2, dof);
369  m_DTvsz_station[id.station() - 1][id.sector() - 1]->fill_dxdz(charge, z, resslope, chi2, dof);
370  m_DTvsphi_station[id.station() - 1][id.wheel() + 2]->fill_x(charge, phi, residual, chi2, dof);
371  m_DTvsphi_station[id.station() - 1][id.wheel() + 2]->fill_dxdz(charge, phi, resslope, chi2, dof);
372 
374 
375  residual = dt2->global_residual();
376  resslope = dt2->global_resslope();
377  chi2 = dt2->chi2();
378  dof = dt2->ndof();
379 
380  if (m_useStubPosition)
381  gpos = dt2->global_stubpos();
382  else
383  gpos = dt2->global_trackpos();
384  phi = atan2(gpos.y(), gpos.x());
385  z = gpos.z();
386 
387  assert(1 <= id.sector() && id.sector() <= 14);
388 
389  m_DTvsz_station[id.station() - 1][id.sector() - 1]->fill_y(charge, z, residual, chi2, dof);
390  m_DTvsz_station[id.station() - 1][id.sector() - 1]->fill_dydz(charge, z, resslope, chi2, dof);
391  m_DTvsphi_station[id.station() - 1][id.wheel() + 2]->fill_y(charge, phi, residual, chi2, dof);
392  m_DTvsphi_station[id.station() - 1][id.wheel() + 2]->fill_dydz(charge, phi, resslope, chi2, dof);
393  }
394 
395  if (id.station() == 4 && dt13 != nullptr && dt13->numHits() >= m_minDT13Hits) {
397 
398  double residual = dt13->global_residual();
399  double resslope = dt13->global_resslope();
400  double chi2 = dt13->chi2();
401  int dof = dt13->ndof();
402 
403  align::GlobalPoint gpos;
404  if (m_useStubPosition)
405  gpos = dt13->global_stubpos();
406  else
407  gpos = dt13->global_trackpos();
408  double phi = atan2(gpos.y(), gpos.x());
409  double z = gpos.z();
410 
411  assert(1 <= id.sector() && id.sector() <= 14);
412 
413  m_DTvsz_station[id.station() - 1][id.sector() - 1]->fill_x(charge, z, residual, chi2, dof);
414  m_DTvsz_station[id.station() - 1][id.sector() - 1]->fill_dxdz(charge, z, resslope, chi2, dof);
415  m_DTvsphi_station[id.station() - 1][id.wheel() + 2]->fill_x(charge, phi, residual, chi2, dof);
416  m_DTvsphi_station[id.station() - 1][id.wheel() + 2]->fill_dxdz(charge, phi, resslope, chi2, dof);
417  }
418  }
419 
420  else if (m_doCSC && chamberId->subdetId() == MuonSubdetId::CSC) {
422  CSCDetId id(chamberId->rawId());
423 
424  int ring = id.ring();
425  if (id.ring() == 4)
426  ring = 1; // combine ME1/a + ME1/b
427 
428  m_counter_csc++;
429 
430  if (csc != nullptr && csc->numHits() >= m_minCSCHits) {
432 
433  double residual = csc->global_residual();
434  double resslope = csc->global_resslope();
435  double chi2 = csc->chi2();
436  int dof = csc->ndof();
437 
438  align::GlobalPoint gpos;
439  if (m_useStubPosition)
440  gpos = csc->global_stubpos();
441  else
442  gpos = csc->global_trackpos();
443  double phi = atan2(gpos.y(), gpos.x());
444  // start phi from -5deg
445  if (phi < -M_PI / 180. * 5.)
446  phi += 2. * M_PI;
447  double R = sqrt(pow(gpos.x(), 2) + pow(gpos.y(), 2));
448 
449  int chamber = id.chamber() - 1;
450  if (id.station() > 1 && ring == 1)
451  chamber *= 2;
452 
453  assert(1 <= id.endcap() && id.endcap() <= 2 && 0 <= chamber && chamber <= 35);
454 
455  if (R > 0.)
456  m_CSCvsphi_me[id.endcap() - 1][id.station() - 1][ring - 1]->fill_x_1d(residual / R, chi2, dof);
457 
458  m_CSCvsr_me[id.endcap() - 1][id.station() - 1][chamber]->fill_x(charge, R, residual, chi2, dof);
459  m_CSCvsr_me[id.endcap() - 1][id.station() - 1][chamber]->fill_dxdz(charge, R, resslope, chi2, dof);
460  m_CSCvsphi_me[id.endcap() - 1][id.station() - 1][ring - 1]->fill_x(charge, phi, residual, chi2, dof);
461  m_CSCvsphi_me[id.endcap() - 1][id.station() - 1][ring - 1]->fill_dxdz(charge, phi, resslope, chi2, dof);
462 
463  if (m_createNtuple && chi2 > 0.) // && TMath::Prob(chi2, dof) < 0.95)
464  {
465  m_id.init(id);
466  m_tr.q = charge;
467  m_tr.pt = mrft.getTrack()->pt();
468  m_tr.pz = mrft.getTrack()->pz();
469  m_re.res = residual;
470  m_re.slope = resslope;
471  m_re.rho = R;
472  m_re.phi = phi;
473  m_re.z = gpos.z();
474  m_run = iEvent.id().run();
475  m_cscnt->Fill();
476  }
477  }
478  }
479 
480  //else { assert(false); }
481  } // end loop over chambers
482 }
483 
485  std::cout << "AlignmentMonitorMuonSystemMap1D counters:" << std::endl;
486  std::cout << " monitor m_counter_event = " << m_counter_event << std::endl;
487  std::cout << " monitor m_counter_track = " << m_counter_track << std::endl;
488  std::cout << " monitor m_counter_trackppt = " << m_counter_trackmoment << std::endl;
489  std::cout << " monitor m_counter_trackdxy = " << m_counter_trackdxy << std::endl;
490  std::cout << " monitor m_counter_trackokay = " << m_counter_trackokay << std::endl;
491  std::cout << " monitor m_counter_dt = " << m_counter_dt << std::endl;
492  std::cout << " monitor m_counter_13numhits = " << m_counter_13numhits << std::endl;
493  std::cout << " monitor m_counter_2numhits = " << m_counter_2numhits << std::endl;
494  std::cout << " monitor m_counter_csc = " << m_counter_csc << std::endl;
495  std::cout << " monitor m_counter_cscnumhits = " << m_counter_cscnumhits << std::endl;
496 }
497 
499  std::string name, AlignmentMonitorMuonSystemMap1D *module, int bins, double low, double high, bool xy, bool add_1d)
500  : m_name(name), m_bins(bins), m_xy(xy), m_1d(add_1d) {
501  m_x_2d = m_y_2d = m_dxdz_2d = m_dydz_2d = nullptr;
502  std::stringstream name_x_2d, name_y_2d, name_dxdz_2d, name_dydz_2d;
503  name_x_2d << m_name << "_x_2d";
504  name_y_2d << m_name << "_y_2d";
505  name_dxdz_2d << m_name << "_dxdz_2d";
506  name_dydz_2d << m_name << "_dydz_2d";
507 
508  const int nbins = 200;
509  const double window = 100.;
510 
511  m_x_2d = module->book2D("/iterN/", name_x_2d.str(), "", m_bins, low, high, nbins, -window, window);
512  if (m_xy)
513  m_y_2d = module->book2D("/iterN/", name_y_2d.str(), "", m_bins, low, high, nbins, -window, window);
514  m_dxdz_2d = module->book2D("/iterN/", name_dxdz_2d.str(), "", m_bins, low, high, nbins, -window, window);
515  if (m_xy)
516  m_dydz_2d = module->book2D("/iterN/", name_dydz_2d.str(), "", m_bins, low, high, nbins, -window, window);
517 
518  m_x_1d = nullptr;
519  if (m_1d) {
520  std::stringstream name_x_1d; //, name_y_1d, name_dxdz_1d, name_dydz_1d;
521  name_x_1d << m_name << "_x_1d";
522  m_x_1d = module->book1D("/iterN/", name_x_1d.str(), "", nbins, -window, window);
523  }
524 }
525 
527  if (m_1d && chi2 > 0.) {
528  // assume that residx was in radians
529  double residual = residx * 1000.;
530  m_x_1d->Fill(residual);
531  }
532 }
533 
535  char charge, double abscissa, double residx, double chi2, int dof) {
536  if (chi2 > 0.) {
537  double residual = residx * 10.;
538  //double weight = dof / chi2;
539  m_x_2d->Fill(abscissa, residual);
540  }
541 }
542 
544  char charge, double abscissa, double residy, double chi2, int dof) {
545  if (m_xy && chi2 > 0.) {
546  double residual = residy * 10.;
547  //double weight = dof / chi2;
548  m_y_2d->Fill(abscissa, residual);
549  }
550 }
551 
553  char charge, double abscissa, double slopex, double chi2, int dof) {
554  if (chi2 > 0.) {
555  double residual = slopex * 1000.;
556  //double weight = dof / chi2;
557  m_dxdz_2d->Fill(abscissa, residual);
558  }
559 }
560 
562  char charge, double abscissa, double slopey, double chi2, int dof) {
563  if (m_xy && chi2 > 0.) {
564  double residual = slopey * 1000.;
565  //double weight = dof / chi2;
566  m_dydz_2d->Fill(abscissa, residual);
567  }
568 }
569 
AlignmentMonitorMuonSystemMap1D::m_minCSCHits
int m_minCSCHits
Definition: AlignmentMonitorMuonSystemMap1D.cc:62
Propagator.h
AlignmentMonitorMuonSystemMap1D::m_plots
std::vector< MuonSystemMapPlot1D * > m_plots
Definition: AlignmentMonitorMuonSystemMap1D.cc:111
PDWG_BPHSkim_cff.muons
muons
Definition: PDWG_BPHSkim_cff.py:47
AlignmentMonitorMuonSystemMap1D::m_minNCrossedChambers
int m_minNCrossedChambers
Definition: AlignmentMonitorMuonSystemMap1D.cc:59
MuonSubdetId::CSC
static constexpr int CSC
Definition: MuonSubdetId.h:12
AlignmentMonitorMuonSystemMap1D::m_counter_2numhits
long m_counter_2numhits
Definition: AlignmentMonitorMuonSystemMap1D.cc:76
electrons_cff.bool
bool
Definition: electrons_cff.py:393
AlignmentMonitorMuonSystemMap1D::m_counter_trackdxy
long m_counter_trackdxy
Definition: AlignmentMonitorMuonSystemMap1D.cc:72
DOFs::dof
dof
Definition: AlignPCLThresholdsWriter.cc:37
pwdgSkimBPark_cfi.beamSpot
beamSpot
Definition: pwdgSkimBPark_cfi.py:5
AlignmentMonitorMuonSystemMap1D::m_minDT2Hits
int m_minDT2Hits
Definition: AlignmentMonitorMuonSystemMap1D.cc:61
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11713
MessageLogger.h
AlignmentMonitorMuonSystemMap1D::MyTrack::pt
Float_t pt
Definition: AlignmentMonitorMuonSystemMap1D.cc:133
AlignmentMonitorMuonSystemMap1D::MyResidual
Definition: AlignmentMonitorMuonSystemMap1D.cc:137
muon
Definition: MuonCocktails.h:17
PV3DBase::x
T x() const
Definition: PV3DBase.h:59
AlignmentMonitorMuonSystemMap1D::m_counter_13numhits
long m_counter_13numhits
Definition: AlignmentMonitorMuonSystemMap1D.cc:75
DiDispStaMuonMonitor_cfi.pt
pt
Definition: DiDispStaMuonMonitor_cfi.py:39
AlignmentMonitorMuonSystemMap1D::m_minTrackerHits
int m_minTrackerHits
Definition: AlignmentMonitorMuonSystemMap1D.cc:56
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_dydz_2d
TH2F * m_dydz_2d
Definition: AlignmentMonitorMuonSystemMap1D.cc:103
relativeConstraints.station
station
Definition: relativeConstraints.py:67
edm
HLT enums.
Definition: AlignableModifier.h:19
AlignmentMonitorMuonSystemMap1D::m_minTrackPt
double m_minTrackPt
Definition: AlignmentMonitorMuonSystemMap1D.cc:51
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::MuonSystemMapPlot1D
MuonSystemMapPlot1D(std::string name, AlignmentMonitorMuonSystemMap1D *module, int bins, double low, double high, bool xy, bool add_1d)
Definition: AlignmentMonitorMuonSystemMap1D.cc:498
AlignmentMonitorMuonSystemMap1D::m_counter_cscnumhits
long m_counter_cscnumhits
Definition: AlignmentMonitorMuonSystemMap1D.cc:78
gather_cfg.cout
cout
Definition: gather_cfg.py:144
AlignmentMonitorMuonSystemMap1D::m_id
MyCSCDetId m_id
Definition: AlignmentMonitorMuonSystemMap1D.cc:129
HLT_FULL_cff.InputTag
InputTag
Definition: HLT_FULL_cff.py:89287
MuonResidualsFromTrack::chamberResidual
MuonChamberResidual * chamberResidual(DetId chamberId, int type)
Definition: MuonResidualsFromTrack.cc:735
AlignmentMonitorMuonSystemMap1D::m_maxTrackP
double m_maxTrackP
Definition: AlignmentMonitorMuonSystemMap1D.cc:54
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_xy
bool m_xy
Definition: AlignmentMonitorMuonSystemMap1D.cc:100
AlignmentMonitorMuonSystemMap1D::afterAlignment
void afterAlignment() override
Definition: AlignmentMonitorMuonSystemMap1D.cc:484
MuonResidualsFromTrack::trackerNumHits
int trackerNumHits() const
Definition: MuonResidualsFromTrack.h:78
svgfig.window
def window(xmin, xmax, ymin, ymax, x=0, y=0, width=100, height=100, xlogbase=None, ylogbase=None, minusInfinity=-1000, flipx=False, flipy=True)
Definition: svgfig.py:643
AlignmentMonitorMuonSystemMap1D::processMuonResidualsFromTrack
void processMuonResidualsFromTrack(MuonResidualsFromTrack &mrft, const edm::Event &iEvent)
Definition: AlignmentMonitorMuonSystemMap1D.cc:315
cms::cuda::assert
assert(be >=bs)
AlignmentMonitorMuonSystemMap1D::MyCSCDetId::s
Short_t s
Definition: AlignmentMonitorMuonSystemMap1D.cc:126
GlobalTrackingGeometryRecord
Definition: GlobalTrackingGeometryRecord.h:17
AlignmentMonitorMuonSystemMap1D::MyCSCDetId::t
Short_t t
Definition: AlignmentMonitorMuonSystemMap1D.cc:127
AlignmentMonitorPluginFactory
AlignmentMonitorMuonSystemMap1D::m_counter_event
long m_counter_event
Definition: AlignmentMonitorMuonSystemMap1D.cc:69
hltPixelTracks_cff.chi2
chi2
Definition: hltPixelTracks_cff.py:25
makeMuonMisalignmentScenario.endcap
endcap
Definition: makeMuonMisalignmentScenario.py:320
HLT_FULL_cff.magneticField
magneticField
Definition: HLT_FULL_cff.py:348
AlignmentMonitorMuonSystemMap1D::num02d
std::string num02d(int num)
Definition: AlignmentMonitorMuonSystemMap1D.cc:174
createJobs.tmp
tmp
align.sh
Definition: createJobs.py:716
edm::Handle< reco::BeamSpot >
AlignmentMonitorMuonSystemMap1D::MyResidual::rho
Float_t rho
Definition: AlignmentMonitorMuonSystemMap1D.cc:138
MuonChamberResidual::global_stubpos
align::GlobalPoint global_stubpos()
Definition: MuonChamberResidual.cc:45
AlignmentMonitorBase::book2D
TH2F * book2D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
Definition: AlignmentMonitorBase.cc:129
AlignmentMonitorBase::pNavigator
AlignableNavigator * pNavigator()
Definition: AlignmentMonitorBase.h:114
AlignmentMonitorMuonSystemMap1D::m_CSCvsphi_me
MuonSystemMapPlot1D * m_CSCvsphi_me[2][4][3]
Definition: AlignmentMonitorMuonSystemMap1D.cc:109
AlignmentMonitorMuonSystemMap1D::AlignmentMonitorMuonSystemMap1D
AlignmentMonitorMuonSystemMap1D(const edm::ParameterSet &cfg)
Definition: AlignmentMonitorMuonSystemMap1D.cc:145
PV3DBase::z
T z() const
Definition: PV3DBase.h:61
AlignmentMonitorMuonSystemMap1D::m_maxTrackPt
double m_maxTrackPt
Definition: AlignmentMonitorMuonSystemMap1D.cc:52
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_x_1d
TH1F * m_x_1d
Definition: AlignmentMonitorMuonSystemMap1D.cc:102
MuonChamberResidual::global_residual
double global_residual() const
Definition: MuonChamberResidual.cc:49
reco::TrackBase::pt
double pt() const
track transverse momentum
Definition: TrackBase.h:637
AlignmentMonitorMuonSystemMap1D::m_counter_csc
long m_counter_csc
Definition: AlignmentMonitorMuonSystemMap1D.cc:77
IdealMagneticFieldRecord
Definition: IdealMagneticFieldRecord.h:11
edm::InputTag::label
std::string const & label() const
Definition: InputTag.h:36
csc
Definition: L1Track.h:19
AlignmentMonitorMuonSystemMap1D::MyTrack::pz
Float_t pz
Definition: AlignmentMonitorMuonSystemMap1D.cc:133
AlignmentMonitorMuonSystemMap1D::MyTrack
Definition: AlignmentMonitorMuonSystemMap1D.cc:131
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_name
std::string m_name
Definition: AlignmentMonitorMuonSystemMap1D.cc:98
edm::EventSetup::get
T get() const
Definition: EventSetup.h:80
HLT_FULL_cff.muon
muon
Definition: HLT_FULL_cff.py:11710
Service.h
mathSSE::sqrt
T sqrt(T t)
Definition: SSEVec.h:19
DDAxes::z
reco::Track
Definition: Track.h:27
IdealMagneticFieldRecord.h
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D
Definition: AlignmentMonitorMuonSystemMap1D.cc:81
MuonResidualsFromTrack
Definition: MuonResidualsFromTrack.h:52
edm::ESHandle< GlobalTrackingGeometry >
MuonResidualsFromTrack::getTrack
const reco::Track * getTrack()
Definition: MuonResidualsFromTrack.h:75
AlignmentMonitorMuonSystemMap1D::MyResidual::phi
Float_t phi
Definition: AlignmentMonitorMuonSystemMap1D.cc:138
AlignmentMonitorMuonSystemMap1D::event
void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks) override
Called for each event (by "run()"): may be reimplemented.
Definition: AlignmentMonitorMuonSystemMap1D.cc:251
LaserClient_cfi.nbins
nbins
Definition: LaserClient_cfi.py:51
MuonChamberResidual::chi2
double chi2() const
Definition: MuonChamberResidual.h:59
MuonResidualsFromTrack::chamberIds
const std::vector< DetId > chamberIds() const
Definition: MuonResidualsFromTrack.h:86
AlignmentMonitorMuonSystemMap1D::m_counter_trackmoment
long m_counter_trackmoment
Definition: AlignmentMonitorMuonSystemMap1D.cc:71
reco::TrackBase::charge
int charge() const
track electric charge
Definition: TrackBase.h:596
MuonSubdetId::DT
static constexpr int DT
Definition: MuonSubdetId.h:11
Point3DBase< Scalar, GlobalTag >
MuonResidualsFromTrack::normalizedChi2
double normalizedChi2() const
Definition: MuonResidualsFromTrack.cc:729
DEFINE_EDM_PLUGIN
#define DEFINE_EDM_PLUGIN(factory, type, name)
Definition: PluginFactory.h:124
geometryCSVtoXML.xy
xy
Definition: geometryCSVtoXML.py:19
AlignmentMonitorPluginFactory.h
GlobalTrackingGeometryRecord.h
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
AlignmentMonitorMuonSystemMap1D::m_doCSC
bool m_doCSC
Definition: AlignmentMonitorMuonSystemMap1D.cc:64
TFileService.h
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_1d
bool m_1d
Definition: AlignmentMonitorMuonSystemMap1D.cc:101
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
MuonChamberResidual::kDT2
Definition: MuonChamberResidual.h:28
AlignmentMonitorMuonSystemMap1D::m_minDT13Hits
int m_minDT13Hits
Definition: AlignmentMonitorMuonSystemMap1D.cc:60
edm::ParameterSet
Definition: ParameterSet.h:47
MuonChamberResidual::numHits
int numHits() const
Definition: MuonChamberResidual.h:52
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::fill_x
void fill_x(char charge, double abscissa, double residx, double chi2, int dof)
Definition: AlignmentMonitorMuonSystemMap1D.cc:534
AlignmentMonitorMuonSystemMap1D::m_re
MyResidual m_re
Definition: AlignmentMonitorMuonSystemMap1D.cc:140
AlignmentMonitorMuonSystemMap1D::m_doDT
bool m_doDT
Definition: AlignmentMonitorMuonSystemMap1D.cc:63
AlignmentMonitorMuonSystemMap1D::m_cscnt
TTree * m_cscnt
Definition: AlignmentMonitorMuonSystemMap1D.cc:116
AlignmentMonitorMuonSystemMap1D::m_muonCollectionTag
edm::InputTag m_muonCollectionTag
Definition: AlignmentMonitorMuonSystemMap1D.cc:50
AlignmentMonitorBase.h
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_dxdz_2d
TH2F * m_dxdz_2d
Definition: AlignmentMonitorMuonSystemMap1D.cc:103
AlignmentMonitorMuonSystemMap1D::m_CSCvsr_me
MuonSystemMapPlot1D * m_CSCvsr_me[2][4][36]
Definition: AlignmentMonitorMuonSystemMap1D.cc:107
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::fill_dxdz
void fill_dxdz(char charge, double abscissa, double slopex, double chi2, int dof)
Definition: AlignmentMonitorMuonSystemMap1D.cc:552
CSCDetId
Definition: CSCDetId.h:26
AlignmentMonitorMuonSystemMap1D::m_run
UInt_t m_run
Definition: AlignmentMonitorMuonSystemMap1D.cc:142
AlignmentMonitorMuonSystemMap1D::MyResidual::slope
Float_t slope
Definition: AlignmentMonitorMuonSystemMap1D.cc:138
PV3DBase::y
T y() const
Definition: PV3DBase.h:60
AlignmentMonitorBase
Definition: AlignmentMonitorBase.h:42
makeMuonMisalignmentScenario.wheel
wheel
Definition: makeMuonMisalignmentScenario.py:319
edm::Service< TFileService >
createfilelist.int
int
Definition: createfilelist.py:10
MuonChamberResidual
Definition: MuonChamberResidual.h:26
iEvent
int iEvent
Definition: GenABIO.cc:224
AlignmentMonitorMuonSystemMap1D::MyTrack::q
Int_t q
Definition: AlignmentMonitorMuonSystemMap1D.cc:132
AlignmentMonitorMuonSystemMap1D::MyResidual::z
Float_t z
Definition: AlignmentMonitorMuonSystemMap1D.cc:138
AlignmentMonitorMuonSystemMap1D::MyCSCDetId::e
Short_t e
Definition: AlignmentMonitorMuonSystemMap1D.cc:126
AlignmentMonitorMuonSystemMap1D::m_useStubPosition
bool m_useStubPosition
Definition: AlignmentMonitorMuonSystemMap1D.cc:65
M_PI
#define M_PI
Definition: BXVectorInputProducer.cc:49
MuonResidualsFromTrack::contains_TIDTEC
bool contains_TIDTEC() const
Definition: MuonResidualsFromTrack.h:84
AlignmentMonitorMuonSystemMap1D::MyResidual::res
Float_t res
Definition: AlignmentMonitorMuonSystemMap1D.cc:138
AlignmentMonitorBase::m_beamSpotTag
const edm::InputTag m_beamSpotTag
Definition: AlignmentMonitorBase.h:116
DetIdAssociatorRecord.h
AlignmentMonitorMuonSystemMap1D::m_counter_dt
long m_counter_dt
Definition: AlignmentMonitorMuonSystemMap1D.cc:74
MagneticField.h
edm::EventSetup
Definition: EventSetup.h:57
MuonChamberResidual::kCSC
Definition: MuonChamberResidual.h:28
MuonResidualsFromTrack.h
AlignmentMonitorMuonSystemMap1D::m_allowTIDTEC
bool m_allowTIDTEC
Definition: AlignmentMonitorMuonSystemMap1D.cc:58
get
#define get
DetIdAssociatorRecord
Definition: DetIdAssociatorRecord.h:13
EgammaValidation_cff.num
num
Definition: EgammaValidation_cff.py:34
AlignmentMonitorMuonSystemMap1D::m_maxDxy
double m_maxDxy
Definition: AlignmentMonitorMuonSystemMap1D.cc:55
InputTag.h
MuonChamberResidual::global_trackpos
align::GlobalPoint global_trackpos()
Definition: MuonChamberResidual.cc:41
looper.cfg
cfg
Definition: looper.py:297
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_bins
int m_bins
Definition: AlignmentMonitorMuonSystemMap1D.cc:99
AlignmentMonitorBase::book1D
TH1F * book1D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX)
Definition: AlignmentMonitorBase.cc:107
DDAxes::phi
AlignmentMonitorMuonSystemMap1D::m_createNtuple
bool m_createNtuple
Definition: AlignmentMonitorMuonSystemMap1D.cc:66
AlignmentMonitorMuonSystemMap1D::MyCSCDetId::r
Short_t r
Definition: AlignmentMonitorMuonSystemMap1D.cc:126
LaserClient_cfi.high
high
Definition: LaserClient_cfi.py:50
AlignmentMonitorMuonSystemMap1D::m_tr
MyTrack m_tr
Definition: AlignmentMonitorMuonSystemMap1D.cc:135
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_x_2d
TH2F * m_x_2d
Definition: AlignmentMonitorMuonSystemMap1D.cc:103
Trajectory
Definition: Trajectory.h:38
triggerObjects_cff.id
id
Definition: triggerObjects_cff.py:31
DetIdAssociator.h
TrackingComponentsRecord.h
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::fill_x_1d
void fill_x_1d(double residx, double chi2, int dof)
Definition: AlignmentMonitorMuonSystemMap1D.cc:526
AlignmentMonitorMuonSystemMap1D::m_maxTrackerRedChi2
double m_maxTrackerRedChi2
Definition: AlignmentMonitorMuonSystemMap1D.cc:57
relativeConstraints.ring
ring
Definition: relativeConstraints.py:68
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::fill_y
void fill_y(char charge, double abscissa, double residy, double chi2, int dof)
Definition: AlignmentMonitorMuonSystemMap1D.cc:543
relativeConstraints.chamber
chamber
Definition: relativeConstraints.py:53
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
EventSetup.h
AlignmentMonitorMuonSystemMap1D::m_DTvsphi_station
MuonSystemMapPlot1D * m_DTvsphi_station[4][5]
Definition: AlignmentMonitorMuonSystemMap1D.cc:108
AlignmentMonitorMuonSystemMap1D::~AlignmentMonitorMuonSystemMap1D
~AlignmentMonitorMuonSystemMap1D() override
Definition: AlignmentMonitorMuonSystemMap1D.cc:37
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::m_y_2d
TH2F * m_y_2d
Definition: AlignmentMonitorMuonSystemMap1D.cc:103
AlignmentMonitorMuonSystemMap1D::MuonSystemMapPlot1D::fill_dydz
void fill_dydz(char charge, double abscissa, double slopey, double chi2, int dof)
Definition: AlignmentMonitorMuonSystemMap1D.cc:561
funct::pow
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29
reco::TrackBase::pz
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:646
AlignmentMonitorBase::ConstTrajTrackPairCollection
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
Definition: AlignmentMonitorBase.h:45
AlignmentMonitorMuonSystemMap1D::MyCSCDetId::c
Short_t c
Definition: AlignmentMonitorMuonSystemMap1D.cc:126
DetId::Muon
Definition: DetId.h:26
MuonChamberResidual::ndof
int ndof() const
Definition: MuonChamberResidual.h:60
DTChamberId
Definition: DTChamberId.h:14
trigObjTnPSource_cfi.bins
bins
Definition: trigObjTnPSource_cfi.py:20
AlignmentMonitorMuonSystemMap1D::m_counter_track
long m_counter_track
Definition: AlignmentMonitorMuonSystemMap1D.cc:70
AlignmentMonitorMuonSystemMap1D::m_DTvsz_station
MuonSystemMapPlot1D * m_DTvsz_station[4][14]
Definition: AlignmentMonitorMuonSystemMap1D.cc:106
ParameterSet.h
AlignmentMonitorMuonSystemMap1D::MyCSCDetId
Definition: AlignmentMonitorMuonSystemMap1D.cc:118
GlobalTrackingGeometry.h
MuonChamberResidual::kDT13
Definition: MuonChamberResidual.h:28
edm::Event
Definition: Event.h:73
AlignmentMonitorMuonSystemMap1D::MyCSCDetId::init
void init(CSCDetId &id)
Definition: AlignmentMonitorMuonSystemMap1D.cc:119
MuonChamberResidual::global_resslope
double global_resslope() const
Definition: MuonChamberResidual.cc:51
dttmaxenums::R
Definition: DTTMax.h:29
LaserClient_cfi.low
low
Definition: LaserClient_cfi.py:52
edm::InputTag
Definition: InputTag.h:15
AlignmentMonitorMuonSystemMap1D
Definition: AlignmentMonitorMuonSystemMap1D.cc:34
AlignmentMonitorMuonSystemMap1D::m_minTrackP
double m_minTrackP
Definition: AlignmentMonitorMuonSystemMap1D.cc:53
TFileService::make
T * make(const Args &... args) const
make new ROOT object
Definition: TFileService.h:64
AlignmentMonitorMuonSystemMap1D::m_counter_trackokay
long m_counter_trackokay
Definition: AlignmentMonitorMuonSystemMap1D.cc:73
AlignmentMonitorMuonSystemMap1D::book
void book() override
Book or retrieve histograms; MUST be reimplemented.
Definition: AlignmentMonitorMuonSystemMap1D.cc:181
TrackingComponentsRecord
Definition: TrackingComponentsRecord.h:12