CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
AlignmentMonitorMuonVsCurvature.cc
Go to the documentation of this file.
1 /*
2  * Package: CommonAlignmentProducer
3  * Class : AlignmentMonitorMuonVsCurvature
4  *
5  * Original Author: Jim Pivarski
6  * Created: Fri Feb 19 21:45:02 CET 2010
7  *
8  * $Id:$
9  */
10 
14 
21 
22 #include <sstream>
23 #include "TProfile.h"
24 #include "TH2F.h"
25 #include "TH1F.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;
44  void processMuonResidualsFromTrack(MuonResidualsFromTrack &mrft, const Trajectory *traj = nullptr);
45 
46 private:
47  // es token
53 
54  // parameters
56  double m_minTrackPt;
57  double m_minTrackP;
62  double m_maxDxy;
66  int m_layer;
68  bool m_doDT;
69  bool m_doCSC;
70 
72 
73  // DT [wheel][station][sector]
76 
77  // CSC [endcap*station][ring][chamber]
80 
83 };
84 
87  : AlignmentMonitorBase(cfg, iC, "AlignmentMonitorMuonVsCurvature"),
88  m_esTokenGBTGeom(iC.esConsumes()),
89  m_esTokenDetId(iC.esConsumes(edm::ESInputTag("", "MuonDetIdAssociator"))),
90  m_esTokenProp(iC.esConsumes(edm::ESInputTag("", "SteppingHelixPropagatorAny"))),
91  m_esTokenMF(iC.esConsumes()),
92  m_esTokenBuilder(iC.esConsumes(MuonResidualsFromTrack::builderESInputTag())),
93  m_muonCollectionTag(cfg.getParameter<edm::InputTag>("muonCollectionTag")),
94  m_minTrackPt(cfg.getParameter<double>("minTrackPt")),
95  m_minTrackP(cfg.getParameter<double>("minTrackP")),
96  m_minTrackerHits(cfg.getParameter<int>("minTrackerHits")),
97  m_maxTrackerRedChi2(cfg.getParameter<double>("maxTrackerRedChi2")),
98  m_allowTIDTEC(cfg.getParameter<bool>("allowTIDTEC")),
99  m_minNCrossedChambers(cfg.getParameter<int>("minNCrossedChambers")),
100  m_maxDxy(cfg.getParameter<double>("maxDxy")),
101  m_minDT13Hits(cfg.getParameter<int>("minDT13Hits")),
102  m_minDT2Hits(cfg.getParameter<int>("minDT2Hits")),
103  m_minCSCHits(cfg.getParameter<int>("minCSCHits")),
104  m_layer(cfg.getParameter<int>("layer")),
105  m_propagator(cfg.getParameter<std::string>("propagator")),
106  m_doDT(cfg.getParameter<bool>("doDT")),
107  m_doCSC(cfg.getParameter<bool>("doCSC")) {}
108 
110  // DT
111  std::string wheelname[5] = {"wheelm2_", "wheelm1_", "wheelz_", "wheelp1_", "wheelp2_"};
112  if (m_doDT)
113  for (int wheel = -2; wheel <= 2; wheel++)
114  for (int station = 1; station <= 4; station++)
115  for (int sector = 1; sector <= 14; sector++) {
116  if (station != 4 && sector > 12)
117  continue;
118 
119  char stationname[20];
120  sprintf(stationname, "st%d_", station);
121 
122  char sectorname[20];
123  sprintf(sectorname, "sector%02d_", sector);
124 
125  for (int component = 0; component < kNumComponents; component++) {
126  std::stringstream th2f_name, tprofile_name;
127  th2f_name << "th2f_" << wheelname[wheel + 2] << stationname << sectorname;
128  tprofile_name << "tprofile_" << wheelname[wheel + 2] << stationname << sectorname;
129 
130  double yminmax = 50., xminmax = 0.05;
131  if (m_minTrackPt > 0.)
132  xminmax = 1. / m_minTrackPt;
133  int ynbins = 50;
134  if (component == kDeltaX) {
135  th2f_name << "deltax";
136  tprofile_name << "deltax";
137  } else if (component == kDeltaDxDz) {
138  th2f_name << "deltadxdz";
139  tprofile_name << "deltadxdz";
140  }
141 
142  th2f_wheel_st_sector[wheel + 2][station - 1][sector - 1][component] =
143  book2D("/iterN/", th2f_name.str(), "", 30, -xminmax, xminmax, ynbins, -yminmax, yminmax);
144  tprofile_wheel_st_sector[wheel + 2][station - 1][sector - 1][component] =
145  bookProfile("/iterN/", tprofile_name.str(), "", 30, -xminmax, xminmax);
146  }
147  }
148 
149  // CSC
150  std::string stname[8] = {"Ep_S1_", "Ep_S2_", "Ep_S3_", "Ep_S4_", "Em_S1_", "Em_S2_", "Em_S3_", "Em_S4_"};
151  if (m_doCSC)
152  for (int station = 0; station < 8; station++)
153  for (int ring = 1; ring <= 3; ring++)
154  for (int chamber = 1; chamber <= 36; chamber++) {
155  int st = station % 4 + 1;
156  if (st > 1 && ring > 2)
157  continue; // only station 1 has more then 2 rings
158  if (st > 1 && ring == 1 && chamber > 18)
159  continue; // ring 1 stations 1,2,3 have 18 chambers
160 
161  char ringname[20];
162  sprintf(ringname, "R%d_", ring);
163 
164  char chname[20];
165  sprintf(chname, "C%02d_", chamber);
166 
167  for (int component = 0; component < kNumComponents; component++) {
168  std::stringstream componentname;
169  double yminmax = 50., xminmax = 0.05;
170  if (m_minTrackPt > 0.)
171  xminmax = 1. / m_minTrackPt;
172  if (ring == 1)
173  xminmax *= 0.5;
174  if (component == kDeltaX) {
175  componentname << "deltax";
176  } else if (component == kDeltaDxDz) {
177  componentname << "deltadxdz";
178  }
179 
180  std::stringstream th2f_name, tprofile_name;
181  th2f_name << "th2f_" << stname[station] << ringname << chname << componentname.str();
182  tprofile_name << "tprofile_" << stname[station] << ringname << chname << componentname.str();
183 
184  th2f_st_ring_chamber[station][ring - 1][chamber - 1][component] =
185  book2D("/iterN/", th2f_name.str(), "", 30, -xminmax, xminmax, 100, -yminmax, yminmax);
186  tprofile_st_ring_chamber[station][ring - 1][chamber - 1][component] =
187  bookProfile("/iterN/", tprofile_name.str(), "", 30, -xminmax, xminmax);
188  }
189  }
190 
191  th1f_trackerRedChi2 = book1D("/iterN/", "trackerRedChi2", "Refit tracker reduced chi^2", 100, 0., 30.);
193  book1D("/iterN/", "trackerRedChi2Diff", "Fit-minus-refit tracker reduced chi^2", 100, -5., 5.);
194 }
195 
197  const edm::EventSetup &iSetup,
198  const ConstTrajTrackPairCollection &trajtracks) {
200  iEvent.getByLabel(m_beamSpotTag, beamSpot);
201 
202  const GlobalTrackingGeometry *globalGeometry = &iSetup.getData(m_esTokenGBTGeom);
203  const DetIdAssociator *muonDetIdAssociator_ = &iSetup.getData(m_esTokenDetId);
204  const Propagator *prop = &iSetup.getData(m_esTokenProp);
206  auto builder = iSetup.getHandle(m_esTokenBuilder);
207 
208  if (m_muonCollectionTag.label().empty()) // use trajectories
209  {
210  for (ConstTrajTrackPairCollection::const_iterator trajtrack = trajtracks.begin(); trajtrack != trajtracks.end();
211  ++trajtrack) {
212  const Trajectory *traj = (*trajtrack).first;
213  const reco::Track *track = (*trajtrack).second;
214 
215  if (track->pt() > m_minTrackPt && track->p() > m_minTrackP && fabs(track->dxy(beamSpot->position())) < m_maxDxy) {
216  MuonResidualsFromTrack muonResidualsFromTrack(
217  builder, magneticField, globalGeometry, muonDetIdAssociator_, prop, traj, track, pNavigator(), 1000.);
218  processMuonResidualsFromTrack(muonResidualsFromTrack, traj);
219  } // end if track pT is within range
220  } // end loop over tracks
221  } else {
223  iEvent.getByLabel(m_muonCollectionTag, muons);
224 
225  for (reco::MuonCollection::const_iterator muon = muons->begin(); muon != muons->end(); ++muon) {
226  if (!(muon->isTrackerMuon() && muon->innerTrack().isNonnull()))
227  continue;
228 
229  if (m_minTrackPt < muon->pt() && m_minTrackP < muon->p() &&
230  fabs(muon->innerTrack()->dxy(beamSpot->position())) < m_maxDxy) {
231  MuonResidualsFromTrack muonResidualsFromTrack(globalGeometry, &(*muon), pNavigator(), 100.);
232  processMuonResidualsFromTrack(muonResidualsFromTrack);
233  }
234  }
235  }
236 }
237 
239  const Trajectory *traj) {
240  if (mrft.trackerNumHits() < m_minTrackerHits)
241  return;
242  if (!m_allowTIDTEC && mrft.contains_TIDTEC())
243  return;
244 
245  int nMuChambers = 0;
246  std::vector<DetId> chamberIds = mrft.chamberIds();
247  for (unsigned ch = 0; ch < chamberIds.size(); ch++)
248  if (chamberIds[ch].det() == DetId::Muon)
249  nMuChambers++;
250  if (nMuChambers < m_minNCrossedChambers)
251  return;
252 
253  th1f_trackerRedChi2->Fill(mrft.trackerRedChi2());
255 
256  if (mrft.normalizedChi2() > m_maxTrackerRedChi2)
257  return;
258 
259  double qoverpt = mrft.getTrack()->charge() / mrft.getTrack()->pt();
260  double qoverpz = 0.;
261  if (fabs(mrft.getTrack()->pz()) > 0.01)
262  qoverpz = mrft.getTrack()->charge() / fabs(mrft.getTrack()->pz());
263 
264  for (std::vector<DetId>::const_iterator chamberId = chamberIds.begin(); chamberId != chamberIds.end(); ++chamberId) {
265  if (chamberId->det() != DetId::Muon)
266  continue;
267 
268  if (m_doDT && chamberId->subdetId() == MuonSubdetId::DT) {
269  DTChamberId dtid(chamberId->rawId());
271 
272  if (dt13 != nullptr && dt13->numHits() >= m_minDT13Hits) {
273  int wheel = dtid.wheel() + 2;
274  int station = dtid.station() - 1;
275  int sector = dtid.sector() - 1;
276 
277  double resid_x = 10. * dt13->global_residual();
278  double resid_dxdz = 1000. * dt13->global_resslope();
279 
280  if (fabs(resid_x) < 100. && fabs(resid_dxdz) < 100.) {
281  th2f_wheel_st_sector[wheel][station][sector][kDeltaX]->Fill(qoverpt, resid_x);
282  tprofile_wheel_st_sector[wheel][station][sector][kDeltaX]->Fill(qoverpt, resid_x);
283  th2f_wheel_st_sector[wheel][station][sector][kDeltaDxDz]->Fill(qoverpt, resid_dxdz);
284  tprofile_wheel_st_sector[wheel][station][sector][kDeltaDxDz]->Fill(qoverpt, resid_dxdz);
285  }
286  } // if it's a good segment
287  } // if DT
288 
289  if (m_doCSC && chamberId->subdetId() == MuonSubdetId::CSC) {
290  CSCDetId cscid(chamberId->rawId());
292 
293  if (csc != nullptr && csc->numHits() >= m_minCSCHits) {
294  int station = 4 * cscid.endcap() + cscid.station() - 5;
295  int ring = cscid.ring() - 1;
296  if (cscid.station() == 1 && cscid.ring() == 4)
297  ring = 0; // join ME1/a to ME1/b
298  int chamber = cscid.chamber() - 1;
299 
300  double resid_x = 10. * csc->global_residual();
301  double resid_dxdz = 1000. * csc->global_resslope();
302 
303  if (fabs(resid_x) < 100. && fabs(resid_dxdz) < 100.) {
304  th2f_st_ring_chamber[station][ring][chamber][kDeltaX]->Fill(qoverpz, resid_x);
305  tprofile_st_ring_chamber[station][ring][chamber][kDeltaX]->Fill(qoverpz, resid_x);
306  th2f_st_ring_chamber[station][ring][chamber][kDeltaDxDz]->Fill(qoverpz, resid_dxdz);
307  tprofile_st_ring_chamber[station][ring][chamber][kDeltaDxDz]->Fill(qoverpz, resid_dxdz);
308  }
309  } // if it's a good segment
310  } // if CSC
311 
312  } // end loop over chamberIds
313 }
314 
double p() const
momentum vector magnitude
Definition: TrackBase.h:631
AlignableNavigator * pNavigator()
tuple cfg
Definition: looper.py:296
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
const MuonResidualsFromTrack::BuilderToken m_esTokenBuilder
TH2F * book2D(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX, int nchY, double lowY, double highY)
const edm::ESGetToken< Propagator, TrackingComponentsRecord > m_esTokenProp
TH2F * th2f_st_ring_chamber[8][3][36][kNumComponents]
tuple magneticField
TProfile * tprofile_st_ring_chamber[8][3][36][kNumComponents]
AlignmentMonitorMuonVsCurvature(const edm::ParameterSet &cfg, edm::ConsumesCollector iC)
const edm::ESGetToken< MagneticField, IdealMagneticFieldRecord > m_esTokenMF
TProfile * bookProfile(std::string dir, std::string name, std::string title, int nchX, double lowX, double highX, int nchY=1, double lowY=0., double highY=0., const char *option="s")
const std::vector< DetId > chamberIds() const
void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection &iTrajTracks) override
Called for each event (by &quot;run()&quot;): may be reimplemented.
bool getData(T &iHolder) const
Definition: EventSetup.h:128
void processMuonResidualsFromTrack(MuonResidualsFromTrack &mrft, const Trajectory *traj=nullptr)
int iEvent
Definition: GenABIO.cc:224
const edm::ESGetToken< GlobalTrackingGeometry, GlobalTrackingGeometryRecord > m_esTokenGBTGeom
double pt() const
track transverse momentum
Definition: TrackBase.h:637
const reco::Track * getTrack()
void book() override
Book or retrieve histograms; MUST be reimplemented.
std::vector< ConstTrajTrackPair > ConstTrajTrackPairCollection
const edm::InputTag m_beamSpotTag
TProfile * tprofile_wheel_st_sector[5][4][14][kNumComponents]
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
double pz() const
z coordinate of momentum vector
Definition: TrackBase.h:646
const edm::ESGetToken< DetIdAssociator, DetIdAssociatorRecord > m_esTokenDetId
std::string const & label() const
Definition: InputTag.h:36
tuple muons
Definition: patZpeak.py:39
int charge() const
track electric charge
Definition: TrackBase.h:596
#define DEFINE_EDM_PLUGIN(factory, type, name)
ESHandle< T > getHandle(const ESGetToken< T, R > &iToken) const
Definition: EventSetup.h:157
static constexpr int DT
Definition: MuonSubdetId.h:11
TH2F * th2f_wheel_st_sector[5][4][14][kNumComponents]
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608
static constexpr int CSC
Definition: MuonSubdetId.h:12
ESGetTokenH3DDVariant esConsumes(std::string const &Reccord, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
MuonChamberResidual * chamberResidual(DetId chamberId, int type)