00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorPluginFactory.h"
00012 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorBase.h"
00013 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsFromTrack.h"
00014 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsPositionFitter.h"
00015 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsAngleFitter.h"
00016 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsTwoBin.h"
00017
00018 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00019 #include "FWCore/Utilities/interface/InputTag.h"
00020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00021 #include "FWCore/Framework/interface/EventSetup.h"
00022 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00023 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00024
00025 #include <sstream>
00026
00027
00028 class AlignmentMonitorSegmentDifferences: public AlignmentMonitorBase
00029 {
00030 public:
00031 AlignmentMonitorSegmentDifferences(const edm::ParameterSet& cfg);
00032 ~AlignmentMonitorSegmentDifferences() {}
00033
00034 void book();
00035
00036 void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& iTrajTracks);
00037 void processMuonResidualsFromTrack(MuonResidualsFromTrack &mrft);
00038
00039 void afterAlignment(const edm::EventSetup &iSetup) {}
00040
00041 private:
00042
00043 edm::InputTag m_muonCollectionTag;
00044 double m_minTrackPt;
00045 double m_minTrackP;
00046 double m_maxDxy;
00047 int m_minTrackerHits;
00048 double m_maxTrackerRedChi2;
00049 bool m_allowTIDTEC;
00050 bool m_minNCrossedChambers;
00051 int m_minDT13Hits;
00052 int m_minDT2Hits;
00053 int m_minCSCHits;
00054 bool m_doDT;
00055 bool m_doCSC;
00056
00057
00058 TProfile *m_dt13_resid[5][12][3];
00059 TProfile *m_dt13_slope[5][12][3];
00060 TProfile *m_dt2_resid[5][12][2];
00061 TProfile *m_dt2_slope[5][12][2];
00062 TH1F *m_posdt13_resid[5][12][3];
00063 TH1F *m_posdt13_slope[5][12][3];
00064 TH1F *m_posdt2_resid[5][12][2];
00065 TH1F *m_posdt2_slope[5][12][2];
00066 TH1F *m_negdt13_resid[5][12][3];
00067 TH1F *m_negdt13_slope[5][12][3];
00068 TH1F *m_negdt2_resid[5][12][2];
00069 TH1F *m_negdt2_slope[5][12][2];
00070
00071
00072 TProfile *m_cscouter_resid[2][36][2];
00073 TProfile *m_cscouter_slope[2][36][2];
00074 TProfile *m_cscinner_resid[2][18][3];
00075 TProfile *m_cscinner_slope[2][18][3];
00076 TH1F *m_poscscouter_resid[2][36][2];
00077 TH1F *m_poscscouter_slope[2][36][2];
00078 TH1F *m_poscscinner_resid[2][18][3];
00079 TH1F *m_poscscinner_slope[2][18][3];
00080 TH1F *m_negcscouter_resid[2][36][2];
00081 TH1F *m_negcscouter_slope[2][36][2];
00082 TH1F *m_negcscinner_resid[2][18][3];
00083 TH1F *m_negcscinner_slope[2][18][3];
00084
00085
00086 TH1F *m_x_pos_dt1_csc1_resid[2][12];
00087 TH1F *m_x_pos_dt1_csc2_resid[2][12];
00088 TH1F *m_x_pos_dt2_csc1_resid[2][12];
00089 TH1F *m_x_neg_dt1_csc1_resid[2][12];
00090 TH1F *m_x_neg_dt1_csc2_resid[2][12];
00091 TH1F *m_x_neg_dt2_csc1_resid[2][12];
00092 };
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106 AlignmentMonitorSegmentDifferences::AlignmentMonitorSegmentDifferences(const edm::ParameterSet& cfg)
00107 : AlignmentMonitorBase(cfg, "AlignmentMonitorSegmentDifferences")
00108 , m_muonCollectionTag(cfg.getParameter<edm::InputTag>("muonCollectionTag"))
00109 , m_minTrackPt(cfg.getParameter<double>("minTrackPt"))
00110 , m_minTrackP(cfg.getParameter<double>("minTrackP"))
00111 , m_maxDxy(cfg.getParameter<double>("maxDxy"))
00112 , m_minTrackerHits(cfg.getParameter<int>("minTrackerHits"))
00113 , m_maxTrackerRedChi2(cfg.getParameter<double>("maxTrackerRedChi2"))
00114 , m_allowTIDTEC(cfg.getParameter<bool>("allowTIDTEC"))
00115 , m_minNCrossedChambers(cfg.getParameter<int>("minNCrossedChambers"))
00116 , m_minDT13Hits(cfg.getParameter<int>("minDT13Hits"))
00117 , m_minDT2Hits(cfg.getParameter<int>("minDT2Hits"))
00118 , m_minCSCHits(cfg.getParameter<int>("minCSCHits"))
00119 , m_doDT(cfg.getParameter<bool>("doDT"))
00120 , m_doCSC(cfg.getParameter<bool>("doCSC"))
00121 {
00122 }
00123
00124 void AlignmentMonitorSegmentDifferences::book()
00125 {
00126 char name[222], pos[222], neg[222];
00127
00128 double max_curv = 1./m_minTrackPt;
00129
00130 if (m_doDT) for (int wheel = -2; wheel <= +2; wheel++)
00131 {
00132 char wheel_label[][2]={"A","B","C","D","E"};
00133 for (int sector = 1; sector <= 12; sector++)
00134 {
00135 char wheel_sector[50];
00136 sprintf(wheel_sector,"%s_%02d", wheel_label[wheel+2], sector );
00137
00138 int nb = 100;
00139 double wnd = 25.;
00140
00141 sprintf(name, "dt13_resid_%s_12", wheel_sector);
00142 sprintf(pos,"pos%s", name);
00143 sprintf(neg,"neg%s", name);
00144 m_dt13_resid[wheel+2][sector-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00145 m_posdt13_resid[wheel+2][sector-1][0] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00146 m_negdt13_resid[wheel+2][sector-1][0] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00147
00148 sprintf(name, "dt13_resid_%s_23", wheel_sector);
00149 sprintf(pos,"pos%s", name);
00150 sprintf(neg,"neg%s", name);
00151 m_dt13_resid[wheel+2][sector-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00152 m_posdt13_resid[wheel+2][sector-1][1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00153 m_negdt13_resid[wheel+2][sector-1][1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00154
00155 sprintf(name, "dt13_resid_%s_34", wheel_sector);
00156 sprintf(pos,"pos%s", name);
00157 sprintf(neg,"neg%s", name);
00158 m_dt13_resid[wheel+2][sector-1][2] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00159 m_posdt13_resid[wheel+2][sector-1][2] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00160 m_negdt13_resid[wheel+2][sector-1][2] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00161
00162 sprintf(name, "dt2_resid_%s_12", wheel_sector);
00163 sprintf(pos,"pos%s", name);
00164 sprintf(neg,"neg%s", name);
00165 m_dt2_resid[wheel+2][sector-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -200., 200., " ");
00166 m_posdt2_resid[wheel+2][sector-1][0] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00167 m_negdt2_resid[wheel+2][sector-1][0] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00168
00169 sprintf(name, "dt2_resid_%s_23", wheel_sector);
00170 sprintf(pos,"pos%s", name);
00171 sprintf(neg,"neg%s", name);
00172 m_dt2_resid[wheel+2][sector-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -200., 200., " ");
00173 m_posdt2_resid[wheel+2][sector-1][1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00174 m_negdt2_resid[wheel+2][sector-1][1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00175
00176 sprintf(name, "dt13_slope_%s_12", wheel_sector);
00177 sprintf(pos,"pos%s", name);
00178 sprintf(neg,"neg%s", name);
00179 m_dt13_slope[wheel+2][sector-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00180 m_posdt13_slope[wheel+2][sector-1][0] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00181 m_negdt13_slope[wheel+2][sector-1][0] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00182
00183 sprintf(name, "dt13_slope_%s_23", wheel_sector);
00184 sprintf(pos,"pos%s", name);
00185 sprintf(neg,"neg%s", name);
00186 m_dt13_slope[wheel+2][sector-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00187 m_posdt13_slope[wheel+2][sector-1][1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00188 m_negdt13_slope[wheel+2][sector-1][1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00189
00190 sprintf(name, "dt13_slope_%s_34", wheel_sector);
00191 sprintf(pos,"pos%s", name);
00192 sprintf(neg,"neg%s", name);
00193 m_dt13_slope[wheel+2][sector-1][2] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00194 m_posdt13_slope[wheel+2][sector-1][2] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00195 m_negdt13_slope[wheel+2][sector-1][2] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00196
00197 sprintf(name, "dt2_slope_%s_12", wheel_sector);
00198 sprintf(pos,"pos%s", name);
00199 sprintf(neg,"neg%s", name);
00200 m_dt2_slope[wheel+2][sector-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -1000., 1000., " ");
00201 m_posdt2_slope[wheel+2][sector-1][0] = book1D("/iterN/", pos, pos, nb, -100., 100.);
00202 m_negdt2_slope[wheel+2][sector-1][0] = book1D("/iterN/", neg, neg, nb, -100., 100.);
00203
00204 sprintf(name, "dt2_slope_%s_23", wheel_sector);
00205 sprintf(pos,"pos%s", name);
00206 sprintf(neg,"neg%s", name);
00207 m_dt2_slope[wheel+2][sector-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -1000., 1000., " ");
00208 m_posdt2_slope[wheel+2][sector-1][1] = book1D("/iterN/", pos, pos, nb, -100., 100.);
00209 m_negdt2_slope[wheel+2][sector-1][1] = book1D("/iterN/", neg, neg, nb, -100., 100.);
00210 }
00211 }
00212
00213 if (m_doCSC) for (int endcap = 1; endcap <= 2; endcap++)
00214 {
00215 std::string endcapletter;
00216 if (endcap == 1) endcapletter = "p";
00217 else if (endcap == 2) endcapletter = "m";
00218
00219 for (int chamber = 1; chamber <= 36; chamber++)
00220 {
00221 char ec_chamber[50];
00222 sprintf(ec_chamber,"%s_%02d", endcapletter.c_str(), chamber );
00223
00224 int nb = 100;
00225 double wnd = 60.;
00226
00227 sprintf(name, "cscouter_resid_%s_12",ec_chamber);
00228 sprintf(pos,"pos%s", name);
00229 sprintf(neg,"neg%s", name);
00230 m_cscouter_resid[endcap-1][chamber-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00231 m_poscscouter_resid[endcap-1][chamber-1][0] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00232 m_negcscouter_resid[endcap-1][chamber-1][0] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00233
00234 sprintf(name, "cscouter_resid_%s_23",ec_chamber);
00235 sprintf(pos,"pos%s", name);
00236 sprintf(neg,"neg%s", name);
00237 m_cscouter_resid[endcap-1][chamber-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00238 m_poscscouter_resid[endcap-1][chamber-1][1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00239 m_negcscouter_resid[endcap-1][chamber-1][1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00240
00241 sprintf(name, "cscouter_slope_%s_12",ec_chamber);
00242 sprintf(pos,"pos%s", name);
00243 sprintf(neg,"neg%s", name);
00244 m_cscouter_slope[endcap-1][chamber-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00245 m_poscscouter_slope[endcap-1][chamber-1][0] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00246 m_negcscouter_slope[endcap-1][chamber-1][0] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00247
00248 sprintf(name, "cscouter_slope_%s_23",ec_chamber);
00249 sprintf(pos,"pos%s", name);
00250 sprintf(neg,"neg%s", name);
00251 m_cscouter_slope[endcap-1][chamber-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00252 m_poscscouter_slope[endcap-1][chamber-1][1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00253 m_negcscouter_slope[endcap-1][chamber-1][1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00254 }
00255
00256 for (int chamber = 1; chamber <= 18; chamber++)
00257 {
00258 char ec_chamber[50];
00259 sprintf(ec_chamber,"%s_%02d", endcapletter.c_str(), chamber );
00260
00261 int nb = 100;
00262 double wnd = 40.;
00263
00264 sprintf(name, "cscinner_resid_%s_12",ec_chamber);
00265 sprintf(pos,"pos%s", name);
00266 sprintf(neg,"neg%s", name);
00267 m_cscinner_resid[endcap-1][chamber-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00268 m_poscscinner_resid[endcap-1][chamber-1][0] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00269 m_negcscinner_resid[endcap-1][chamber-1][0] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00270
00271 sprintf(name, "cscinner_resid_%s_23",ec_chamber);
00272 sprintf(pos,"pos%s", name);
00273 sprintf(neg,"neg%s", name);
00274 m_cscinner_resid[endcap-1][chamber-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00275 m_poscscinner_resid[endcap-1][chamber-1][1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00276 m_negcscinner_resid[endcap-1][chamber-1][1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00277
00278 sprintf(name, "cscinner_resid_%s_34",ec_chamber);
00279 sprintf(pos,"pos%s", name);
00280 sprintf(neg,"neg%s", name);
00281 m_cscinner_resid[endcap-1][chamber-1][2] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00282 m_poscscinner_resid[endcap-1][chamber-1][2] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00283 m_negcscinner_resid[endcap-1][chamber-1][2] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00284
00285 sprintf(name, "cscinner_slope_%s_12",ec_chamber);
00286 sprintf(pos,"pos%s", name);
00287 sprintf(neg,"neg%s", name);
00288 m_cscinner_slope[endcap-1][chamber-1][0] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00289 m_poscscinner_slope[endcap-1][chamber-1][0] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00290 m_negcscinner_slope[endcap-1][chamber-1][0] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00291
00292 sprintf(name, "cscinner_slope_%s_23",ec_chamber);
00293 sprintf(pos,"pos%s", name);
00294 sprintf(neg,"neg%s", name);
00295 m_cscinner_slope[endcap-1][chamber-1][1] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00296 m_poscscinner_slope[endcap-1][chamber-1][1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00297 m_negcscinner_slope[endcap-1][chamber-1][1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00298
00299 sprintf(name, "cscinner_slope_%s_34",ec_chamber);
00300 sprintf(pos,"pos%s", name);
00301 sprintf(neg,"neg%s", name);
00302 m_cscinner_slope[endcap-1][chamber-1][2] = bookProfile("/iterN/", name, name, 20, -max_curv, max_curv, 1, -100., 100., " ");
00303 m_poscscinner_slope[endcap-1][chamber-1][2] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00304 m_negcscinner_slope[endcap-1][chamber-1][2] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00305 }
00306 }
00307
00308
00309 for (int e = 1; e<=2; e++)
00310 for (int s = 1; s <= 12; s++)
00311 {
00312 char endcap_sector[50];
00313 if (e == 1) sprintf(endcap_sector,"Wp2S%02d", s);
00314 if (e == 2) sprintf(endcap_sector,"Wm2S%02d", s);
00315
00316 int nb = 200;
00317 double wnd = 100.;
00318
00319 sprintf(pos,"pos_x_dt1_csc1_%s", endcap_sector);
00320 sprintf(neg,"neg_x_dt1_csc1_%s", endcap_sector);
00321 m_x_pos_dt1_csc1_resid[e-1][s-1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00322 m_x_neg_dt1_csc1_resid[e-1][s-1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00323
00324 sprintf(pos,"pos_x_dt1_csc2_%s", endcap_sector);
00325 sprintf(neg,"neg_x_dt1_csc2_%s", endcap_sector);
00326 m_x_pos_dt1_csc2_resid[e-1][s-1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00327 m_x_neg_dt1_csc2_resid[e-1][s-1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00328
00329 sprintf(pos,"pos_x_dt2_csc1_%s", endcap_sector);
00330 sprintf(neg,"neg_x_dt2_csc1_%s", endcap_sector);
00331 m_x_pos_dt2_csc1_resid[e-1][s-1] = book1D("/iterN/", pos, pos, nb, -wnd, wnd);
00332 m_x_neg_dt2_csc1_resid[e-1][s-1] = book1D("/iterN/", neg, neg, nb, -wnd, wnd);
00333 }
00334
00335 }
00336
00337
00338 void AlignmentMonitorSegmentDifferences::event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& trajtracks)
00339 {
00340 edm::ESHandle<GlobalTrackingGeometry> globalGeometry;
00341 iSetup.get<GlobalTrackingGeometryRecord>().get(globalGeometry);
00342
00343 edm::Handle<reco::BeamSpot> beamSpot;
00344 iEvent.getByLabel(m_beamSpotTag, beamSpot);
00345
00346 if (m_muonCollectionTag.label().empty())
00347 {
00348 for (ConstTrajTrackPairCollection::const_iterator trajtrack = trajtracks.begin(); trajtrack != trajtracks.end(); ++trajtrack)
00349 {
00350 const Trajectory* traj = (*trajtrack).first;
00351 const reco::Track* track = (*trajtrack).second;
00352
00353 if (track->pt() > m_minTrackPt && track->p() > m_minTrackP && fabs(track->dxy(beamSpot->position())) < m_maxDxy )
00354 {
00355 MuonResidualsFromTrack muonResidualsFromTrack(globalGeometry, traj, track, pNavigator(), 1000.);
00356 processMuonResidualsFromTrack(muonResidualsFromTrack);
00357 }
00358 }
00359 }
00360 else
00361 {
00362 edm::Handle<reco::MuonCollection> muons;
00363 iEvent.getByLabel(m_muonCollectionTag, muons);
00364
00365 for (reco::MuonCollection::const_iterator muon = muons->begin(); muon != muons->end(); ++muon)
00366 {
00367 if ( !(muon->isTrackerMuon() && muon->innerTrack().isNonnull() ) ) continue;
00368
00369 if (m_minTrackPt < muon->pt() && m_minTrackP < muon->p() && fabs(muon->innerTrack()->dxy(beamSpot->position())) < m_maxDxy)
00370 {
00371 MuonResidualsFromTrack muonResidualsFromTrack(globalGeometry, &(*muon), pNavigator(), 100.);
00372 processMuonResidualsFromTrack(muonResidualsFromTrack);
00373 }
00374 }
00375 }
00376 }
00377
00378
00379 void AlignmentMonitorSegmentDifferences::processMuonResidualsFromTrack(MuonResidualsFromTrack &mrft)
00380 {
00381 if (mrft.trackerNumHits() < m_minTrackerHits) return;
00382 if (!m_allowTIDTEC && mrft.contains_TIDTEC()) return;
00383 if (mrft.normalizedChi2() > m_maxTrackerRedChi2) return;
00384
00385 int nMuChambers = 0;
00386 std::vector<DetId> chamberIds = mrft.chamberIds();
00387 for (unsigned ch=0; ch < chamberIds.size(); ch++) if (chamberIds[ch].det() == DetId::Muon) nMuChambers++;
00388 if (nMuChambers < m_minNCrossedChambers ) return;
00389
00390 double qoverpt = (mrft.getTrack()->charge() > 0 ? 1. : -1.) / mrft.getTrack()->pt();
00391 double qoverpz = 0.;
00392 if (fabs(mrft.getTrack()->pz()) > 0.01) qoverpz = mrft.getTrack()->charge() / fabs(mrft.getTrack()->pz());
00393
00394 for (std::vector<DetId>::const_iterator chamberId = chamberIds.begin(); chamberId != chamberIds.end(); ++chamberId)
00395 {
00396 if (chamberId->det() != DetId::Muon ) continue;
00397
00398
00399 if (m_doDT && chamberId->subdetId() == MuonSubdetId::DT)
00400 {
00401 MuonChamberResidual *dt13 = mrft.chamberResidual(*chamberId, MuonChamberResidual::kDT13);
00402 MuonChamberResidual *dt2 = mrft.chamberResidual(*chamberId, MuonChamberResidual::kDT2);
00403
00404 if (dt13 != NULL && dt13->numHits() >= m_minDT13Hits)
00405 {
00406 DTChamberId thisid(chamberId->rawId());
00407 for (std::vector<DetId>::const_iterator otherId = chamberIds.begin(); otherId != chamberIds.end(); ++otherId)
00408 {
00409 if (otherId->det() == DetId::Muon && otherId->subdetId() == MuonSubdetId::DT)
00410 {
00411 DTChamberId thatid(otherId->rawId());
00412 if (thisid.rawId() != thatid.rawId() && thisid.wheel() == thatid.wheel() && thisid.sector() == thatid.sector())
00413 {
00414 MuonChamberResidual *dt13other = mrft.chamberResidual(*otherId, MuonChamberResidual::kDT13);
00415 if (dt13other != NULL && dt13other->numHits() >= m_minDT13Hits)
00416 {
00417 double slopediff = 1000. * (dt13->global_resslope() - dt13other->global_resslope());
00418
00419
00420
00421 double residdiff = 10. * (dt13->global_residual() - dt13other->global_residual());
00422
00423 int st = 0;
00424 if (thatid.station() - thisid.station() == 1) st = thisid.station();
00425 if (st>0)
00426 {
00427 m_dt13_resid[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(qoverpt, residdiff);
00428 m_dt13_slope[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(qoverpt, slopediff);
00429 if (qoverpt > 0)
00430 {
00431 m_posdt13_resid[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(residdiff);
00432 m_posdt13_slope[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(slopediff);
00433 }
00434 else
00435 {
00436 m_negdt13_resid[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(residdiff);
00437 m_negdt13_slope[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(slopediff);
00438 }
00439 }
00440 }
00441 }
00442 }
00443
00444
00445
00446 if ( !(abs(thisid.wheel()) == 2 && (thisid.station() == 1 || thisid.station() == 2)) ) continue;
00447 if (otherId->det() == DetId::Muon && otherId->subdetId() == MuonSubdetId::CSC)
00448 {
00449 CSCDetId thatid(otherId->rawId());
00450
00451 if ( !( (thatid.station()==1 && thatid.ring()==3) || (thatid.station()==2 && thatid.ring()==2) ) ) continue;
00452
00453 MuonChamberResidual *cscother = mrft.chamberResidual(*otherId, MuonChamberResidual::kCSC);
00454 if (cscother != NULL && cscother->numHits() >= m_minCSCHits)
00455 {
00456
00457 double csc_scale = dt13->chamberAlignable()->surface().toGlobal(align::LocalPoint(dt13->trackx(), dt13->tracky(),0)).perp() /
00458 cscother->chamberAlignable()->surface().toGlobal(align::LocalPoint(cscother->trackx(), cscother->tracky(),0)).perp();
00459 double residdiff = 10. * (dt13->global_residual() - cscother->global_residual() * csc_scale);
00460 if (thisid.station() == 1 && thatid.station()==1)
00461 {
00462 if (qoverpt > 0) m_x_pos_dt1_csc1_resid[thatid.endcap()-1][thisid.sector()-1]->Fill(residdiff);
00463 else m_x_neg_dt1_csc1_resid[thatid.endcap()-1][thisid.sector()-1]->Fill(residdiff);
00464 }
00465 else if (thisid.station() == 1 && thatid.station()==2)
00466 {
00467 if (qoverpt > 0) m_x_pos_dt1_csc2_resid[thatid.endcap()-1][thisid.sector()-1]->Fill(residdiff);
00468 else m_x_neg_dt1_csc2_resid[thatid.endcap()-1][thisid.sector()-1]->Fill(residdiff);
00469 }
00470 else if (thisid.station() == 2 && thatid.station()==1)
00471 {
00472 if (qoverpt > 0) m_x_pos_dt2_csc1_resid[thatid.endcap()-1][thisid.sector()-1]->Fill(residdiff);
00473 else m_x_neg_dt2_csc1_resid[thatid.endcap()-1][thisid.sector()-1]->Fill(residdiff);
00474 }
00475 }
00476 }
00477 }
00478 }
00479
00480
00481 if (dt2 != NULL && dt2->numHits() >= m_minDT2Hits)
00482 {
00483 DTChamberId thisid(chamberId->rawId());
00484 for (std::vector<DetId>::const_iterator otherId = chamberIds.begin(); otherId != chamberIds.end(); ++otherId)
00485 {
00486 if (otherId->det() == DetId::Muon && otherId->subdetId() == MuonSubdetId::DT)
00487 {
00488 DTChamberId thatid(otherId->rawId());
00489 if (thisid.rawId() != thatid.rawId() && thisid.wheel() == thatid.wheel() && thisid.sector() == thatid.sector())
00490 {
00491 MuonChamberResidual *dt2other = mrft.chamberResidual(*otherId, MuonChamberResidual::kDT2);
00492 if (dt2other != NULL && dt2other->numHits() >= m_minDT2Hits)
00493 {
00494 double slopediff = 1000. * (dt2->global_resslope() - dt2other->global_resslope());
00495
00496
00497
00498 double residdiff = 10. * (dt2->global_residual() - dt2other->global_residual());
00499
00500 int st = 0;
00501 if (thatid.station() - thisid.station() == 1) st = thisid.station();
00502 if (st>0)
00503 {
00504 m_dt2_resid[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(qoverpt, residdiff);
00505 m_dt2_slope[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(qoverpt, slopediff);
00506 if (qoverpt > 0)
00507 {
00508 m_posdt2_resid[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(residdiff);
00509 m_posdt2_slope[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(slopediff);
00510 }
00511 else
00512 {
00513 m_negdt2_resid[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(residdiff);
00514 m_negdt2_slope[thisid.wheel()+2][thisid.sector()-1][st-1]->Fill(slopediff);
00515 }
00516 }
00517 }
00518 }
00519 }
00520 }
00521 }
00522 }
00523
00524
00525 else if (m_doCSC && chamberId->subdetId() == MuonSubdetId::CSC)
00526 {
00527 MuonChamberResidual *csc = mrft.chamberResidual(*chamberId, MuonChamberResidual::kCSC);
00528 if (csc->numHits() >= m_minCSCHits)
00529 {
00530 CSCDetId thisid(chamberId->rawId());
00531 for (std::vector<DetId>::const_iterator otherId = chamberIds.begin(); otherId != chamberIds.end(); ++otherId)
00532 {
00533 if (otherId->det() == DetId::Muon && otherId->subdetId() == MuonSubdetId::CSC)
00534 {
00535 CSCDetId thatid(otherId->rawId());
00536 if (thisid.rawId() != thatid.rawId() && thisid.endcap() == thatid.endcap())
00537 {
00538 MuonChamberResidual *cscother = mrft.chamberResidual(*otherId, MuonChamberResidual::kCSC);
00539 if (cscother != NULL && cscother->numHits() >= m_minCSCHits)
00540 {
00541 double slopediff = 1000. * (csc->global_resslope() - cscother->global_resslope());
00542
00543
00544
00545 double residdiff = 10. * (csc->global_residual() - cscother->global_residual());
00546
00547 int thischamber = thisid.chamber();
00548 int thisring = thisid.ring();
00549 if (thisid.station() == 1 && (thisring == 1 || thisring == 4))
00550 {
00551 thischamber = (thischamber - 1) / 2 + 1;
00552 thisring = 1;
00553 }
00554
00555 if (thisring == thatid.ring() && thischamber == thatid.chamber())
00556 {
00557 bool inner = (thisring == 1);
00558 bool outer = (thisring == 2);
00559 int st = 0;
00560 if (thatid.station() - thisid.station() == 1 && (inner || thisid.station()<3) ) st = thisid.station();
00561
00562 if (outer && st>0)
00563 {
00564 m_cscouter_resid[thisid.endcap()-1][thischamber-1][st-1]->Fill(qoverpz, residdiff);
00565 m_cscouter_slope[thisid.endcap()-1][thischamber-1][st-1]->Fill(qoverpz, slopediff);
00566 if (qoverpz > 0)
00567 {
00568 m_poscscouter_resid[thisid.endcap()-1][thischamber-1][st-1]->Fill(residdiff);
00569 m_poscscouter_slope[thisid.endcap()-1][thischamber-1][st-1]->Fill(slopediff);
00570 }
00571 else
00572 {
00573 m_negcscouter_resid[thisid.endcap()-1][thischamber-1][st-1]->Fill(residdiff);
00574 m_negcscouter_slope[thisid.endcap()-1][thischamber-1][st-1]->Fill(slopediff);
00575 }
00576 }
00577 if (inner && st>0)
00578 {
00579 m_cscinner_resid[thisid.endcap()-1][thischamber-1][st-1]->Fill(qoverpz, residdiff);
00580 m_cscinner_slope[thisid.endcap()-1][thischamber-1][st-1]->Fill(qoverpz, slopediff);
00581 if (qoverpz > 0)
00582 {
00583 m_poscscinner_resid[thisid.endcap()-1][thischamber-1][st-1]->Fill(residdiff);
00584 m_poscscinner_slope[thisid.endcap()-1][thischamber-1][st-1]->Fill(slopediff);
00585 }
00586 else
00587 {
00588 m_negcscinner_resid[thisid.endcap()-1][thischamber-1][st-1]->Fill(residdiff);
00589 m_negcscinner_slope[thisid.endcap()-1][thischamber-1][st-1]->Fill(slopediff);
00590 }
00591 }
00592 }
00593 }
00594 }
00595 }
00596 }
00597
00598 }
00599 }
00600
00601 }
00602 }
00603
00604
00605 DEFINE_EDM_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorSegmentDifferences, "AlignmentMonitorSegmentDifferences");