CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
AlignmentMonitorMuonSystemMap1D.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CommonAlignmentProducer
4 // Class : AlignmentMonitorMuonSystemMap1D
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: Jim Pivarski
10 // Created: Mon Nov 12 13:30:14 CST 2007
11 //
12 
14 
15 //
16 // constants, enums and typedefs
17 //
18 
19 //
20 // static data member definitions
21 //
22 
23 //
24 // member functions
25 //
26 
28  : AlignmentMonitorBase(cfg, "AlignmentMonitorMuonSystemMap1D")
29  , m_minTrackPt(cfg.getParameter<double>("minTrackPt"))
30  , m_maxTrackPt(cfg.getParameter<double>("maxTrackPt"))
31  , m_minTrackerHits(cfg.getParameter<int>("minTrackerHits"))
32  , m_maxTrackerRedChi2(cfg.getParameter<double>("maxTrackerRedChi2"))
33  , m_allowTIDTEC(cfg.getParameter<bool>("allowTIDTEC"))
34  , m_minDT13Hits(cfg.getParameter<int>("minDT13Hits"))
35  , m_minDT2Hits(cfg.getParameter<int>("minDT2Hits"))
36  , m_minCSCHits(cfg.getParameter<int>("minCSCHits"))
37 {}
38 
40  int tens = num / 10;
41  int ones = num % 10;
42 
43  std::string s_tens, s_ones;
44  if (tens == 0) s_tens = std::string("0");
45  if (tens == 1) s_tens = std::string("1");
46  if (tens == 2) s_tens = std::string("2");
47  if (tens == 3) s_tens = std::string("3");
48  if (tens == 4) s_tens = std::string("4");
49  if (tens == 5) s_tens = std::string("5");
50  if (tens == 6) s_tens = std::string("6");
51  if (tens == 7) s_tens = std::string("7");
52  if (tens == 8) s_tens = std::string("8");
53  if (tens == 9) s_tens = std::string("9");
54 
55  if (ones == 0) s_ones = std::string("0");
56  if (ones == 1) s_ones = std::string("1");
57  if (ones == 2) s_ones = std::string("2");
58  if (ones == 3) s_ones = std::string("3");
59  if (ones == 4) s_ones = std::string("4");
60  if (ones == 5) s_ones = std::string("5");
61  if (ones == 6) s_ones = std::string("6");
62  if (ones == 7) s_ones = std::string("7");
63  if (ones == 8) s_ones = std::string("8");
64  if (ones == 9) s_ones = std::string("9");
65 
66  return s_tens + s_ones;
67 }
68 
70  for (int sector = 1; sector <= 14; sector++) {
71  if (sector <= 12) {
72  m_DTvsz_station1[sector-1] = new MuonSystemMapPlot1D(std::string("DTvsz_st1sec") + num02d(sector), this, 60, -660., 660., true); m_plots.push_back(m_DTvsz_station1[sector-1]);
73  m_DTvsz_station2[sector-1] = new MuonSystemMapPlot1D(std::string("DTvsz_st2sec") + num02d(sector), this, 60, -660., 660., true); m_plots.push_back(m_DTvsz_station2[sector-1]);
74  m_DTvsz_station3[sector-1] = new MuonSystemMapPlot1D(std::string("DTvsz_st3sec") + num02d(sector), this, 60, -660., 660., true); m_plots.push_back(m_DTvsz_station3[sector-1]);
75  }
76  m_DTvsz_station4[sector-1] = new MuonSystemMapPlot1D(std::string("DTvsz_st4sec") + num02d(sector), this, 60, -660., 660., false); m_plots.push_back(m_DTvsz_station4[sector-1]);
77  }
78 
79  for (int endcap = 1; endcap <= 2; endcap++) {
80  for (int chamber = 1; chamber <= 36; chamber++) {
81  m_CSCvsr_me1[endcap-1][chamber-1] = new MuonSystemMapPlot1D(std::string("CSCvsr_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("1ch") + num02d(chamber), this, 60, 100., 700., false); m_plots.push_back(m_CSCvsr_me1[endcap-1][chamber-1]);
82  m_CSCvsr_me2[endcap-1][chamber-1] = new MuonSystemMapPlot1D(std::string("CSCvsr_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("2ch") + num02d(chamber), this, 60, 100., 700., false); m_plots.push_back(m_CSCvsr_me2[endcap-1][chamber-1]);
83  m_CSCvsr_me3[endcap-1][chamber-1] = new MuonSystemMapPlot1D(std::string("CSCvsr_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("3ch") + num02d(chamber), this, 60, 100., 700., false); m_plots.push_back(m_CSCvsr_me3[endcap-1][chamber-1]);
84  m_CSCvsr_me4[endcap-1][chamber-1] = new MuonSystemMapPlot1D(std::string("CSCvsr_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("4ch") + num02d(chamber), this, 60, 100., 700., false); m_plots.push_back(m_CSCvsr_me4[endcap-1][chamber-1]);
85  }
86  }
87 
88  for (int wheel = -2; wheel <= 2; wheel++) {
89  std::string s_wheel;
90  if (wheel == -2) s_wheel = std::string("A");
91  else if (wheel == -1) s_wheel = std::string("B");
92  else if (wheel == 0) s_wheel = std::string("C");
93  else if (wheel == +1) s_wheel = std::string("D");
94  else if (wheel == +2) s_wheel = std::string("E");
95 
96  m_DTvsphi_station1[wheel+2] = new MuonSystemMapPlot1D(std::string("DTvsphi_st1wh") + s_wheel, this, 180, -M_PI, M_PI, true); m_plots.push_back(m_DTvsphi_station1[wheel+2]);
97  m_DTvsphi_station2[wheel+2] = new MuonSystemMapPlot1D(std::string("DTvsphi_st2wh") + s_wheel, this, 180, -M_PI, M_PI, true); m_plots.push_back(m_DTvsphi_station2[wheel+2]);
98  m_DTvsphi_station3[wheel+2] = new MuonSystemMapPlot1D(std::string("DTvsphi_st3wh") + s_wheel, this, 180, -M_PI, M_PI, true); m_plots.push_back(m_DTvsphi_station3[wheel+2]);
99  m_DTvsphi_station4[wheel+2] = new MuonSystemMapPlot1D(std::string("DTvsphi_st4wh") + s_wheel, this, 180, -M_PI, M_PI, false); m_plots.push_back(m_DTvsphi_station4[wheel+2]);
100  }
101 
102  for (int endcap = 1; endcap <= 2; endcap++) {
103  m_CSCvsphi_me11[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("11"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me11[endcap-1]);
104  m_CSCvsphi_me12[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("12"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me12[endcap-1]);
105  m_CSCvsphi_me13[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("13"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me13[endcap-1]);
106  m_CSCvsphi_me14[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("14"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me14[endcap-1]);
107  m_CSCvsphi_me21[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("21"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me21[endcap-1]);
108  m_CSCvsphi_me22[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("22"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me22[endcap-1]);
109  m_CSCvsphi_me31[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("31"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me31[endcap-1]);
110  m_CSCvsphi_me32[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("32"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me32[endcap-1]);
111  m_CSCvsphi_me41[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("41"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me41[endcap-1]);
112  m_CSCvsphi_me42[endcap-1] = new MuonSystemMapPlot1D(std::string("CSCvsphi_me") + (endcap == 1 ? std::string("p") : std::string("m")) + std::string("42"), this, 180, -M_PI, M_PI, false); m_plots.push_back(m_CSCvsphi_me42[endcap-1]);
113  }
114 
115  m_counter_event = 0;
116  m_counter_track = 0;
117  m_counter_trackpt = 0;
119  m_counter_dt = 0;
121  m_counter_2numhits = 0;
122  m_counter_csc = 0;
124 }
125 
127  m_counter_event++;
128 
130  iSetup.get<GlobalTrackingGeometryRecord>().get(globalGeometry);
131 
132  for (ConstTrajTrackPairCollection::const_iterator trajtrack = trajtracks.begin(); trajtrack != trajtracks.end(); ++trajtrack) {
133  const Trajectory* traj = (*trajtrack).first;
134  const reco::Track* track = (*trajtrack).second;
135 
136  m_counter_track++;
137 
138  if (m_minTrackPt < track->pt() && track->pt() < m_maxTrackPt) {
139  char charge = (track->charge() > 0 ? 1 : -1);
140  // double qoverpt = track->charge() / track->pt();
141  // double qoverpz = track->charge() / track->pz();
142  MuonResidualsFromTrack muonResidualsFromTrack(globalGeometry, traj, pNavigator(), 1000.);
143 
145 
146  if (muonResidualsFromTrack.trackerNumHits() >= m_minTrackerHits && muonResidualsFromTrack.trackerRedChi2() < m_maxTrackerRedChi2 && (m_allowTIDTEC || !muonResidualsFromTrack.contains_TIDTEC())) {
147  std::vector<DetId> chamberIds = muonResidualsFromTrack.chamberIds();
148 
150 
151  for (std::vector<DetId>::const_iterator chamberId = chamberIds.begin(); chamberId != chamberIds.end(); ++chamberId) {
152 
153  if (chamberId->det() == DetId::Muon && chamberId->subdetId() == MuonSubdetId::DT) {
154  MuonChamberResidual *dt13 = muonResidualsFromTrack.chamberResidual(*chamberId, MuonChamberResidual::kDT13);
155  MuonChamberResidual *dt2 = muonResidualsFromTrack.chamberResidual(*chamberId, MuonChamberResidual::kDT2);
156  DTChamberId id(chamberId->rawId());
157 
158  m_counter_dt++;
159 
160  if (dt13 != NULL && dt13->numHits() >= m_minDT13Hits) {
162 
163  double residual = dt13->global_residual();
164  double resslope = dt13->global_resslope();
165  double chi2 = dt13->chi2();
166  int dof = dt13->ndof();
167 
168  GlobalPoint trackpos = dt13->global_trackpos();
169  double phi = atan2(trackpos.y(), trackpos.x());
170  double z = trackpos.z();
171 
172  assert(1 <= id.sector() && id.sector() <= 14);
173 
174  if (id.station() == 1) m_DTvsz_station1[id.sector()-1]->fill_x(charge, z, residual, chi2, dof);
175  if (id.station() == 2) m_DTvsz_station2[id.sector()-1]->fill_x(charge, z, residual, chi2, dof);
176  if (id.station() == 3) m_DTvsz_station3[id.sector()-1]->fill_x(charge, z, residual, chi2, dof);
177  if (id.station() == 4) m_DTvsz_station4[id.sector()-1]->fill_x(charge, z, residual, chi2, dof);
178 
179  if (id.station() == 1) m_DTvsz_station1[id.sector()-1]->fill_dxdz(charge, z, resslope, chi2, dof);
180  if (id.station() == 2) m_DTvsz_station2[id.sector()-1]->fill_dxdz(charge, z, resslope, chi2, dof);
181  if (id.station() == 3) m_DTvsz_station3[id.sector()-1]->fill_dxdz(charge, z, resslope, chi2, dof);
182  if (id.station() == 4) m_DTvsz_station4[id.sector()-1]->fill_dxdz(charge, z, resslope, chi2, dof);
183 
184  if (id.station() == 1) m_DTvsphi_station1[id.wheel()+2]->fill_x(charge, phi, residual, chi2, dof);
185  if (id.station() == 2) m_DTvsphi_station2[id.wheel()+2]->fill_x(charge, phi, residual, chi2, dof);
186  if (id.station() == 3) m_DTvsphi_station3[id.wheel()+2]->fill_x(charge, phi, residual, chi2, dof);
187  if (id.station() == 4) m_DTvsphi_station4[id.wheel()+2]->fill_x(charge, phi, residual, chi2, dof);
188 
189  if (id.station() == 1) m_DTvsphi_station1[id.wheel()+2]->fill_dxdz(charge, phi, resslope, chi2, dof);
190  if (id.station() == 2) m_DTvsphi_station2[id.wheel()+2]->fill_dxdz(charge, phi, resslope, chi2, dof);
191  if (id.station() == 3) m_DTvsphi_station3[id.wheel()+2]->fill_dxdz(charge, phi, resslope, chi2, dof);
192  if (id.station() == 4) m_DTvsphi_station4[id.wheel()+2]->fill_dxdz(charge, phi, resslope, chi2, dof);
193  }
194 
195  if (dt2 != NULL && dt2->numHits() >= m_minDT2Hits) {
197 
198  double residual = dt2->global_residual();
199  double resslope = dt2->global_resslope();
200  double chi2 = dt2->chi2();
201  int dof = dt2->ndof();
202 
203  GlobalPoint trackpos = dt2->global_trackpos();
204  double phi = atan2(trackpos.y(), trackpos.x());
205  double z = trackpos.z();
206 
207  assert(1 <= id.sector() && id.sector() <= 14);
208 
209  if (id.station() == 1) m_DTvsz_station1[id.sector()-1]->fill_y(charge, z, residual, chi2, dof);
210  if (id.station() == 2) m_DTvsz_station2[id.sector()-1]->fill_y(charge, z, residual, chi2, dof);
211  if (id.station() == 3) m_DTvsz_station3[id.sector()-1]->fill_y(charge, z, residual, chi2, dof);
212  if (id.station() == 4) m_DTvsz_station4[id.sector()-1]->fill_y(charge, z, residual, chi2, dof);
213 
214  if (id.station() == 1) m_DTvsz_station1[id.sector()-1]->fill_dydz(charge, z, resslope, chi2, dof);
215  if (id.station() == 2) m_DTvsz_station2[id.sector()-1]->fill_dydz(charge, z, resslope, chi2, dof);
216  if (id.station() == 3) m_DTvsz_station3[id.sector()-1]->fill_dydz(charge, z, resslope, chi2, dof);
217  if (id.station() == 4) m_DTvsz_station4[id.sector()-1]->fill_dydz(charge, z, resslope, chi2, dof);
218 
219  if (id.station() == 1) m_DTvsphi_station1[id.wheel()+2]->fill_y(charge, phi, residual, chi2, dof);
220  if (id.station() == 2) m_DTvsphi_station2[id.wheel()+2]->fill_y(charge, phi, residual, chi2, dof);
221  if (id.station() == 3) m_DTvsphi_station3[id.wheel()+2]->fill_y(charge, phi, residual, chi2, dof);
222  if (id.station() == 4) m_DTvsphi_station4[id.wheel()+2]->fill_y(charge, phi, residual, chi2, dof);
223 
224  if (id.station() == 1) m_DTvsphi_station1[id.wheel()+2]->fill_dydz(charge, phi, resslope, chi2, dof);
225  if (id.station() == 2) m_DTvsphi_station2[id.wheel()+2]->fill_dydz(charge, phi, resslope, chi2, dof);
226  if (id.station() == 3) m_DTvsphi_station3[id.wheel()+2]->fill_dydz(charge, phi, resslope, chi2, dof);
227  if (id.station() == 4) m_DTvsphi_station4[id.wheel()+2]->fill_dydz(charge, phi, resslope, chi2, dof);
228  }
229  }
230 
231  else if (chamberId->det() == DetId::Muon && chamberId->subdetId() == MuonSubdetId::CSC) {
232  MuonChamberResidual *csc = muonResidualsFromTrack.chamberResidual(*chamberId, MuonChamberResidual::kCSC);
233  CSCDetId id(chamberId->rawId());
234 
235  m_counter_csc++;
236 
237  if (csc != NULL && csc->numHits() >= m_minCSCHits) {
239 
240  double residual = csc->global_residual();
241  double resslope = csc->global_resslope();
242  double chi2 = csc->chi2();
243  int dof = csc->ndof();
244 
245  GlobalPoint trackpos = csc->global_trackpos();
246  double phi = atan2(trackpos.y(), trackpos.x());
247  double R = sqrt(pow(trackpos.x(), 2) + pow(trackpos.y(), 2));
248 
249  int chamber = id.chamber() - 1;
250  if (id.station() > 1 && id.ring() == 1) chamber *= 2;
251 
252  assert(1 <= id.endcap() && id.endcap() <= 2 && 0 <= chamber && chamber <= 35);
253 
254  if (id.station() == 1) m_CSCvsr_me1[id.endcap()-1][chamber]->fill_x(charge, R, residual, chi2, dof);
255  if (id.station() == 2) m_CSCvsr_me2[id.endcap()-1][chamber]->fill_x(charge, R, residual, chi2, dof);
256  if (id.station() == 3) m_CSCvsr_me3[id.endcap()-1][chamber]->fill_x(charge, R, residual, chi2, dof);
257  if (id.station() == 4) m_CSCvsr_me4[id.endcap()-1][chamber]->fill_x(charge, R, residual, chi2, dof);
258 
259  if (id.station() == 1) m_CSCvsr_me1[id.endcap()-1][chamber]->fill_dxdz(charge, R, resslope, chi2, dof);
260  if (id.station() == 2) m_CSCvsr_me2[id.endcap()-1][chamber]->fill_dxdz(charge, R, resslope, chi2, dof);
261  if (id.station() == 3) m_CSCvsr_me3[id.endcap()-1][chamber]->fill_dxdz(charge, R, resslope, chi2, dof);
262  if (id.station() == 4) m_CSCvsr_me4[id.endcap()-1][chamber]->fill_dxdz(charge, R, resslope, chi2, dof);
263 
264  if (id.station() == 1 && id.ring() == 1) m_CSCvsphi_me11[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
265  if (id.station() == 1 && id.ring() == 2) m_CSCvsphi_me12[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
266  if (id.station() == 1 && id.ring() == 3) m_CSCvsphi_me13[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
267  if (id.station() == 1 && id.ring() == 4) m_CSCvsphi_me14[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
268  if (id.station() == 2 && id.ring() == 1) m_CSCvsphi_me21[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
269  if (id.station() == 2 && id.ring() == 2) m_CSCvsphi_me22[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
270  if (id.station() == 3 && id.ring() == 1) m_CSCvsphi_me31[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
271  if (id.station() == 3 && id.ring() == 2) m_CSCvsphi_me32[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
272  if (id.station() == 4 && id.ring() == 1) m_CSCvsphi_me41[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
273  if (id.station() == 4 && id.ring() == 2) m_CSCvsphi_me42[id.endcap()-1]->fill_x(charge, phi, residual, chi2, dof);
274 
275  if (id.station() == 1 && id.ring() == 1) m_CSCvsphi_me11[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
276  if (id.station() == 1 && id.ring() == 2) m_CSCvsphi_me12[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
277  if (id.station() == 1 && id.ring() == 3) m_CSCvsphi_me13[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
278  if (id.station() == 1 && id.ring() == 4) m_CSCvsphi_me14[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
279  if (id.station() == 2 && id.ring() == 1) m_CSCvsphi_me21[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
280  if (id.station() == 2 && id.ring() == 2) m_CSCvsphi_me22[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
281  if (id.station() == 3 && id.ring() == 1) m_CSCvsphi_me31[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
282  if (id.station() == 3 && id.ring() == 2) m_CSCvsphi_me32[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
283  if (id.station() == 4 && id.ring() == 1) m_CSCvsphi_me41[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
284  if (id.station() == 4 && id.ring() == 2) m_CSCvsphi_me42[id.endcap()-1]->fill_dxdz(charge, phi, resslope, chi2, dof);
285  }
286  }
287 
288  else { assert(false); }
289 
290  } // end loop over chambers
291  } // end if track has enough tracker hits
292  } // end if track has acceptable momentum
293  } // end loop over tracks
294 }
295 
297  std::cout << "monitor m_counter_event = " << m_counter_event << std::endl;
298  std::cout << "monitor m_counter_track = " << m_counter_track << std::endl;
299  std::cout << "monitor m_counter_trackpt = " << m_counter_trackpt << std::endl;
300  std::cout << "monitor m_counter_trackokay = " << m_counter_trackokay << std::endl;
301  std::cout << "monitor m_counter_dt = " << m_counter_dt << std::endl;
302  std::cout << "monitor m_counter_13numhits = " << m_counter_13numhits << std::endl;
303  std::cout << "monitor m_counter_2numhits = " << m_counter_2numhits << std::endl;
304  std::cout << "monitor m_counter_csc = " << m_counter_csc << std::endl;
305  std::cout << "monitor m_counter_cscnumhits = " << m_counter_cscnumhits << std::endl;
306 }
307 
308 //
309 // constructors and destructor
310 //
311 
312 // AlignmentMonitorMuonSystemMap1D::AlignmentMonitorMuonSystemMap1D(const AlignmentMonitorMuonSystemMap1D& rhs)
313 // {
314 // // do actual copying here;
315 // }
316 
317 //
318 // assignment operators
319 //
320 // const AlignmentMonitorMuonSystemMap1D& AlignmentMonitorMuonSystemMap1D::operator=(const AlignmentMonitorMuonSystemMap1D& rhs)
321 // {
322 // //An exception safe implementation is
323 // AlignmentMonitorMuonSystemMap1D temp(rhs);
324 // swap(rhs);
325 //
326 // return *this;
327 // }
328 
329 //
330 // const member functions
331 //
332 
333 //
334 // static member functions
335 //
336 
337 //
338 // SEAL definitions
339 //
340 
AlignableNavigator * pNavigator()
T y() const
Definition: PV3DBase.h:57
#define NULL
Definition: scimark2.h:8
void fill_y(char charge, double abscissa, double residy, double chi2, int dof)
double global_residual() const
double charge(const std::vector< uint8_t > &Ampls)
void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks)
Called for each event (by &quot;run()&quot;): may be reimplemented.
GlobalPoint global_trackpos()
const std::vector< DetId > chamberIds() const
void fill_x(char charge, double abscissa, double residx, double chi2, int dof)
int iEvent
Definition: GenABIO.cc:243
static const int CSC
Definition: MuonSubdetId.h:15
Definition: DDAxes.h:10
T sqrt(T t)
Definition: SSEVec.h:28
double pt() const
track transverse momentum
Definition: TrackBase.h:130
void book()
Book or retrieve histograms; MUST be reimplemented.
T z() const
Definition: PV3DBase.h:58
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
void afterAlignment(const edm::EventSetup &iSetup)
#define M_PI
Definition: BFit3D.cc:3
const T & get() const
Definition: EventSetup.h:55
long long int num
Definition: procUtils.cc:71
std::vector< MuonSystemMapPlot1D * > m_plots
void fill_dydz(char charge, double abscissa, double slopey, double chi2, int dof)
double global_resslope() const
tuple cout
Definition: gather_cfg.py:41
int charge() const
track electric charge
Definition: TrackBase.h:112
#define DEFINE_EDM_PLUGIN(factory, type, name)
static const int DT
Definition: MuonSubdetId.h:14
void fill_dxdz(char charge, double abscissa, double slopex, double chi2, int dof)
T x() const
Definition: PV3DBase.h:56
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:40
MuonChamberResidual * chamberResidual(DetId chamberId, int type)
AlignmentMonitorMuonSystemMap1D(const edm::ParameterSet &cfg)
Definition: DDAxes.h:10