CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/Alignment/CommonAlignmentMonitor/plugins/AlignmentMonitorSegmentDifferences.cc

Go to the documentation of this file.
00001 // -*- C++ -*-
00002 //
00003 // Package:     CommonAlignmentProducer
00004 // Class  :     AlignmentMonitorSegmentDifferences
00005 // 
00006 // Implementation:
00007 //     <Notes on implementation>
00008 //
00009 // Original Author:  Jim Pivarski
00010 //         Created:  Mon Nov 12 13:30:14 CST 2007
00011 //
00012 
00013 // system include files
00014 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorPluginFactory.h"
00015 #include "Alignment/CommonAlignmentMonitor/interface/AlignmentMonitorBase.h"
00016 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00018 #include "FWCore/Framework/interface/EventSetup.h"
00019 #include "Geometry/CommonDetUnit/interface/GlobalTrackingGeometry.h"
00020 #include "Geometry/Records/interface/GlobalTrackingGeometryRecord.h"
00021 
00022 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsFromTrack.h"
00023 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsPositionFitter.h"
00024 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsAngleFitter.h"
00025 #include "Alignment/MuonAlignmentAlgorithms/interface/MuonResidualsTwoBin.h"
00026 
00027 #include <sstream>
00028 
00029 // user include files
00030 
00031 // 
00032 // class definition
00033 // 
00034 
00035 class AlignmentMonitorSegmentDifferences: public AlignmentMonitorBase {
00036 public:
00037   AlignmentMonitorSegmentDifferences(const edm::ParameterSet& cfg);
00038   ~AlignmentMonitorSegmentDifferences() {};
00039 
00040   void book();
00041   void event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& iTrajTracks);
00042   void afterAlignment(const edm::EventSetup &iSetup);
00043 
00044 private:
00045   double m_minTrackPt;
00046   int m_minTrackerHits;
00047   double m_maxTrackerRedChi2;
00048   bool m_allowTIDTEC;
00049   int m_minDT13Hits;
00050   int m_minDT2Hits;
00051   int m_minCSCHits;
00052 
00053   // wheel, sector, stationdiff
00054   TProfile *m_dt13_resid[5][12][3];
00055   TProfile *m_dt13_slope[5][12][3];
00056   TProfile *m_dt2_resid[5][12][2];
00057   TProfile *m_dt2_slope[5][12][2];
00058   TH1F *m_posdt13_resid[5][12][3];
00059   TH1F *m_posdt13_slope[5][12][3];
00060   TH1F *m_posdt2_resid[5][12][2];
00061   TH1F *m_posdt2_slope[5][12][2];
00062   TH1F *m_negdt13_resid[5][12][3];
00063   TH1F *m_negdt13_slope[5][12][3];
00064   TH1F *m_negdt2_resid[5][12][2];
00065   TH1F *m_negdt2_slope[5][12][2];
00066 
00067   // endcap, chamber, stationdiff
00068   TProfile *m_cscouter_resid[2][36][2];
00069   TProfile *m_cscouter_slope[2][36][2];
00070   TProfile *m_cscinner_resid[2][18][3];
00071   TProfile *m_cscinner_slope[2][18][3];
00072   TH1F *m_poscscouter_resid[2][36][2];
00073   TH1F *m_poscscouter_slope[2][36][2];
00074   TH1F *m_poscscinner_resid[2][18][3];
00075   TH1F *m_poscscinner_slope[2][18][3];
00076   TH1F *m_negcscouter_resid[2][36][2];
00077   TH1F *m_negcscouter_slope[2][36][2];
00078   TH1F *m_negcscinner_resid[2][18][3];
00079   TH1F *m_negcscinner_slope[2][18][3];
00080 };
00081 
00082 //
00083 // constants, enums and typedefs
00084 //
00085 
00086 //
00087 // static data member definitions
00088 //
00089 
00090 //
00091 // member functions
00092 //
00093 
00094 AlignmentMonitorSegmentDifferences::AlignmentMonitorSegmentDifferences(const edm::ParameterSet& cfg)
00095    : AlignmentMonitorBase(cfg, "AlignmentMonitorSegmentDifferences")
00096    , m_minTrackPt(cfg.getParameter<double>("minTrackPt"))
00097    , m_minTrackerHits(cfg.getParameter<int>("minTrackerHits"))
00098    , m_maxTrackerRedChi2(cfg.getParameter<double>("maxTrackerRedChi2"))
00099    , m_allowTIDTEC(cfg.getParameter<bool>("allowTIDTEC"))
00100    , m_minDT13Hits(cfg.getParameter<int>("minDT13Hits"))
00101    , m_minDT2Hits(cfg.getParameter<int>("minDT2Hits"))
00102    , m_minCSCHits(cfg.getParameter<int>("minCSCHits"))
00103 {
00104 }
00105 
00106 void AlignmentMonitorSegmentDifferences::book() {
00107    for (int wheel = -2;  wheel <= +2;  wheel++) {
00108       for (int sector = 1;  sector <= 12;  sector++) {
00109          char num[3];
00110          num[0] = ('0' + (sector / 10));
00111          num[1] = ('0' + (sector % 10));
00112          num[2] = 0;
00113          
00114          std::string wheelletter;
00115          if (wheel == -2) wheelletter = "A";
00116          else if (wheel == -1) wheelletter = "B";
00117          else if (wheel ==  0) wheelletter = "C";
00118          else if (wheel == +1) wheelletter = "D";
00119          else if (wheel == +2) wheelletter = "E";
00120 
00121          std::string name, pos, neg;
00122 
00123          name = (std::string("dt13_resid_") + wheelletter + std::string("_") + std::string(num) + std::string("_12"));
00124          m_dt13_resid[wheel+2][sector-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00125          pos = std::string("pos") + name;
00126          m_posdt13_resid[wheel+2][sector-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00127          neg = std::string("neg") + name;
00128          m_negdt13_resid[wheel+2][sector-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00129 
00130          name = (std::string("dt13_resid_") + wheelletter + std::string("_") + std::string(num) + std::string("_23"));
00131          m_dt13_resid[wheel+2][sector-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00132          pos = std::string("pos") + name;
00133          m_posdt13_resid[wheel+2][sector-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00134          neg = std::string("neg") + name;
00135          m_negdt13_resid[wheel+2][sector-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00136 
00137          name = (std::string("dt13_resid_") + wheelletter + std::string("_") + std::string(num) + std::string("_34"));
00138          m_dt13_resid[wheel+2][sector-1][2] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00139          pos = std::string("pos") + name;
00140          m_posdt13_resid[wheel+2][sector-1][2] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00141          neg = std::string("neg") + name;
00142          m_negdt13_resid[wheel+2][sector-1][2] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00143 
00144          name = (std::string("dt2_resid_") + wheelletter + std::string("_") + std::string(num) + std::string("_12"));
00145          m_dt2_resid[wheel+2][sector-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -200., 200., " ");
00146          pos = std::string("pos") + name;
00147          m_posdt2_resid[wheel+2][sector-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -20., 20.);
00148          neg = std::string("neg") + name;
00149          m_negdt2_resid[wheel+2][sector-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -20., 20.);
00150 
00151          name = (std::string("dt2_resid_") + wheelletter + std::string("_") + std::string(num) + std::string("_23"));
00152          m_dt2_resid[wheel+2][sector-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -200., 200., " ");
00153          pos = std::string("pos") + name;
00154          m_posdt2_resid[wheel+2][sector-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -20., 20.);
00155          neg = std::string("neg") + name;
00156          m_negdt2_resid[wheel+2][sector-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -20., 20.);
00157 
00158          name = (std::string("dt13_slope_") + wheelletter + std::string("_") + std::string(num) + std::string("_12"));
00159          m_dt13_slope[wheel+2][sector-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00160          pos = std::string("pos") + name;
00161          m_posdt13_slope[wheel+2][sector-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00162          neg = std::string("neg") + name;
00163          m_negdt13_slope[wheel+2][sector-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00164 
00165          name = (std::string("dt13_slope_") + wheelletter + std::string("_") + std::string(num) + std::string("_23"));
00166          m_dt13_slope[wheel+2][sector-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00167          pos = std::string("pos") + name;
00168          m_posdt13_slope[wheel+2][sector-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00169          neg = std::string("neg") + name;
00170          m_negdt13_slope[wheel+2][sector-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00171 
00172          name = (std::string("dt13_slope_") + wheelletter + std::string("_") + std::string(num) + std::string("_34"));
00173          m_dt13_slope[wheel+2][sector-1][2] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00174          pos = std::string("pos") + name;
00175          m_posdt13_slope[wheel+2][sector-1][2] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00176          neg = std::string("neg") + name;
00177          m_negdt13_slope[wheel+2][sector-1][2] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00178 
00179          name = (std::string("dt2_slope_") + wheelletter + std::string("_") + std::string(num) + std::string("_12"));
00180          m_dt2_slope[wheel+2][sector-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -1000., 1000., " ");
00181          pos = std::string("pos") + name;
00182          m_posdt2_slope[wheel+2][sector-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -100., 100.);
00183          neg = std::string("neg") + name;
00184          m_negdt2_slope[wheel+2][sector-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -100., 100.);
00185 
00186          name = (std::string("dt2_slope_") + wheelletter + std::string("_") + std::string(num) + std::string("_23"));
00187          m_dt2_slope[wheel+2][sector-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -1000., 1000., " ");
00188          pos = std::string("pos") + name;
00189          m_posdt2_slope[wheel+2][sector-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -100., 100.);
00190          neg = std::string("neg") + name;
00191          m_negdt2_slope[wheel+2][sector-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -100., 100.);
00192       }
00193    }
00194 
00195    for (int endcap = 1;  endcap <= 2;  endcap++) {
00196       for (int chamber = 1;  chamber <= 36;  chamber++) {
00197          char num[3];
00198          num[0] = ('0' + (chamber / 10));
00199          num[1] = ('0' + (chamber % 10));
00200          num[2] = 0;
00201 
00202          std::string endcapletter;
00203          if (endcap == 1) endcapletter = "p";
00204          else if (endcap == 2) endcapletter = "m";
00205 
00206          std::string name, pos, neg;
00207 
00208          name = (std::string("cscouter_resid_") + endcapletter + std::string("_") + std::string(num) + std::string("_12"));
00209          m_cscouter_resid[endcap-1][chamber-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00210          pos = std::string("pos") + name;
00211          m_poscscouter_resid[endcap-1][chamber-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00212          neg = std::string("neg") + name;
00213          m_negcscouter_resid[endcap-1][chamber-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00214 
00215          name = (std::string("cscouter_resid_") + endcapletter + std::string("_") + std::string(num) + std::string("_23"));
00216          m_cscouter_resid[endcap-1][chamber-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00217          pos = std::string("pos") + name;
00218          m_poscscouter_resid[endcap-1][chamber-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00219          neg = std::string("neg") + name;
00220          m_negcscouter_resid[endcap-1][chamber-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00221 
00222          name = (std::string("cscouter_slope_") + endcapletter + std::string("_") + std::string(num) + std::string("_12"));
00223          m_cscouter_slope[endcap-1][chamber-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00224          pos = std::string("pos") + name;
00225          m_poscscouter_slope[endcap-1][chamber-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00226          neg = std::string("neg") + name;
00227          m_negcscouter_slope[endcap-1][chamber-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00228 
00229          name = (std::string("cscouter_slope_") + endcapletter + std::string("_") + std::string(num) + std::string("_23"));
00230          m_cscouter_slope[endcap-1][chamber-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00231          pos = std::string("pos") + name;
00232          m_poscscouter_slope[endcap-1][chamber-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00233          neg = std::string("neg") + name;
00234          m_negcscouter_slope[endcap-1][chamber-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00235       }
00236 
00237       for (int chamber = 1;  chamber <= 18;  chamber++) {
00238          char num[3];
00239          num[0] = ('0' + (chamber / 10));
00240          num[1] = ('0' + (chamber % 10));
00241          num[2] = 0;
00242 
00243          std::string endcapletter;
00244          if (endcap == 1) endcapletter = "p";
00245          else if (endcap == 2) endcapletter = "m";
00246 
00247          std::string name, pos, neg;
00248 
00249          name = (std::string("cscinner_resid_") + endcapletter + std::string("_") + std::string(num) + std::string("_12"));
00250          m_cscinner_resid[endcap-1][chamber-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00251          pos = std::string("pos") + name;
00252          m_poscscinner_resid[endcap-1][chamber-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00253          neg = std::string("neg") + name;
00254          m_negcscinner_resid[endcap-1][chamber-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00255 
00256          name = (std::string("cscinner_resid_") + endcapletter + std::string("_") + std::string(num) + std::string("_23"));
00257          m_cscinner_resid[endcap-1][chamber-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00258          pos = std::string("pos") + name;
00259          m_poscscinner_resid[endcap-1][chamber-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00260          neg = std::string("neg") + name;
00261          m_negcscinner_resid[endcap-1][chamber-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00262 
00263          name = (std::string("cscinner_resid_") + endcapletter + std::string("_") + std::string(num) + std::string("_34"));
00264          m_cscinner_resid[endcap-1][chamber-1][2] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00265          pos = std::string("pos") + name;
00266          m_poscscinner_resid[endcap-1][chamber-1][2] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00267          neg = std::string("neg") + name;
00268          m_negcscinner_resid[endcap-1][chamber-1][2] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00269 
00270          name = (std::string("cscinner_slope_") + endcapletter + std::string("_") + std::string(num) + std::string("_12"));
00271          m_cscinner_slope[endcap-1][chamber-1][0] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00272          pos = std::string("pos") + name;
00273          m_poscscinner_slope[endcap-1][chamber-1][0] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00274          neg = std::string("neg") + name;
00275          m_negcscinner_slope[endcap-1][chamber-1][0] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00276 
00277          name = (std::string("cscinner_slope_") + endcapletter + std::string("_") + std::string(num) + std::string("_23"));
00278          m_cscinner_slope[endcap-1][chamber-1][1] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00279          pos = std::string("pos") + name;
00280          m_poscscinner_slope[endcap-1][chamber-1][1] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00281          neg = std::string("neg") + name;
00282          m_negcscinner_slope[endcap-1][chamber-1][1] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00283 
00284          name = (std::string("cscinner_slope_") + endcapletter + std::string("_") + std::string(num) + std::string("_34"));
00285          m_cscinner_slope[endcap-1][chamber-1][2] = bookProfile("/iterN/", name.c_str(), name.c_str(), 20, -1./m_minTrackPt, 1./m_minTrackPt, 1, -100., 100., " ");
00286          pos = std::string("pos") + name;
00287          m_poscscinner_slope[endcap-1][chamber-1][2] = book1D("/iterN/", pos.c_str(), pos.c_str(), 100, -10., 10.);
00288          neg = std::string("neg") + name;
00289          m_negcscinner_slope[endcap-1][chamber-1][2] = book1D("/iterN/", neg.c_str(), neg.c_str(), 100, -10., 10.);
00290       }
00291    }
00292 }
00293 
00294 void AlignmentMonitorSegmentDifferences::event(const edm::Event &iEvent, const edm::EventSetup &iSetup, const ConstTrajTrackPairCollection& trajtracks) {
00295   edm::ESHandle<GlobalTrackingGeometry> globalGeometry;
00296   iSetup.get<GlobalTrackingGeometryRecord>().get(globalGeometry);
00297 
00298   for (ConstTrajTrackPairCollection::const_iterator trajtrack = trajtracks.begin();  trajtrack != trajtracks.end();  ++trajtrack) {
00299     const Trajectory* traj = (*trajtrack).first;
00300     const reco::Track* track = (*trajtrack).second;
00301 
00302     if (track->pt() > m_minTrackPt) {
00303       double qoverpt = (track->charge() > 0 ? 1. : -1.) / track->pt();
00304       MuonResidualsFromTrack muonResidualsFromTrack(globalGeometry, traj, pNavigator(), 1000.);
00305 
00306       if (muonResidualsFromTrack.trackerNumHits() >= m_minTrackerHits  &&  muonResidualsFromTrack.trackerRedChi2() < m_maxTrackerRedChi2  &&  (m_allowTIDTEC  ||  !muonResidualsFromTrack.contains_TIDTEC())) {
00307         std::vector<DetId> chamberIds = muonResidualsFromTrack.chamberIds();
00308 
00309         for (std::vector<DetId>::const_iterator chamberId = chamberIds.begin();  chamberId != chamberIds.end();  ++chamberId) {
00310           if (chamberId->det() == DetId::Muon  &&  chamberId->subdetId() == MuonSubdetId::DT) {
00311             MuonChamberResidual *dt13 = muonResidualsFromTrack.chamberResidual(*chamberId, MuonChamberResidual::kDT13);
00312             MuonChamberResidual *dt2 = muonResidualsFromTrack.chamberResidual(*chamberId, MuonChamberResidual::kDT2);
00313 
00314             if (dt13 != NULL  &&  dt13->numHits() >= m_minDT13Hits) {
00315                DTChamberId thisid(chamberId->rawId());
00316 
00317                for (std::vector<DetId>::const_iterator otherId = chamberIds.begin();  otherId != chamberIds.end();  ++otherId) {
00318                   if (otherId->det() == DetId::Muon  &&  otherId->subdetId() == MuonSubdetId::DT) {
00319                      DTChamberId thatid(otherId->rawId());
00320                      if (thisid.rawId() != thatid.rawId()  &&  thisid.wheel() == thatid.wheel()  &&  thisid.sector() == thatid.sector()) {
00321                         MuonChamberResidual *dt13other = muonResidualsFromTrack.chamberResidual(*otherId, MuonChamberResidual::kDT13);
00322                         if (dt13other != NULL  &&  dt13other->numHits() >= m_minDT13Hits) {
00323                            double slopediff = 1000. * (dt13->global_resslope() - dt13other->global_resslope());
00324 
00325                            double length = dt13->chamberAlignable()->surface().toGlobal(LocalPoint(0,0,0)).perp() - dt13other->chamberAlignable()->surface().toGlobal(LocalPoint(0,0,0)).perp();
00326                            double residdiff = 10. * (dt13->global_residual() + length*dt13->global_resslope() - dt13other->global_residual());
00327 
00328                            if (thisid.station() == 1  &&  thatid.station() == 2) {
00329                               m_dt13_resid[thisid.wheel()+2][thisid.sector()-1][0]->Fill(qoverpt, residdiff);
00330                               m_dt13_slope[thisid.wheel()+2][thisid.sector()-1][0]->Fill(qoverpt, slopediff);
00331                               if (qoverpt > 0) {
00332                                  m_posdt13_resid[thisid.wheel()+2][thisid.sector()-1][0]->Fill(residdiff);
00333                                  m_posdt13_slope[thisid.wheel()+2][thisid.sector()-1][0]->Fill(slopediff);
00334                               }
00335                               else {
00336                                  m_negdt13_resid[thisid.wheel()+2][thisid.sector()-1][0]->Fill(residdiff);
00337                                  m_negdt13_slope[thisid.wheel()+2][thisid.sector()-1][0]->Fill(slopediff);
00338                               }
00339                            }
00340                            else if (thisid.station() == 2  &&  thatid.station() == 3) {
00341                               m_dt13_resid[thisid.wheel()+2][thisid.sector()-1][1]->Fill(qoverpt, residdiff);
00342                               m_dt13_slope[thisid.wheel()+2][thisid.sector()-1][1]->Fill(qoverpt, slopediff);
00343                               if (qoverpt > 0) {
00344                                  m_posdt13_resid[thisid.wheel()+2][thisid.sector()-1][1]->Fill(residdiff);
00345                                  m_posdt13_slope[thisid.wheel()+2][thisid.sector()-1][1]->Fill(slopediff);
00346                               }
00347                               else {
00348                                  m_negdt13_resid[thisid.wheel()+2][thisid.sector()-1][1]->Fill(residdiff);
00349                                  m_negdt13_slope[thisid.wheel()+2][thisid.sector()-1][1]->Fill(slopediff);
00350                               }
00351                            }
00352                            else if (thisid.station() == 3  &&  thatid.station() == 4) {
00353                               m_dt13_resid[thisid.wheel()+2][thisid.sector()-1][2]->Fill(qoverpt, residdiff);
00354                               m_dt13_slope[thisid.wheel()+2][thisid.sector()-1][2]->Fill(qoverpt, slopediff);
00355                               if (qoverpt > 0) {
00356                                  m_posdt13_resid[thisid.wheel()+2][thisid.sector()-1][2]->Fill(residdiff);
00357                                  m_posdt13_slope[thisid.wheel()+2][thisid.sector()-1][2]->Fill(slopediff);
00358                               }
00359                               else {
00360                                  m_negdt13_resid[thisid.wheel()+2][thisid.sector()-1][2]->Fill(residdiff);
00361                                  m_negdt13_slope[thisid.wheel()+2][thisid.sector()-1][2]->Fill(slopediff);
00362                               }
00363                            }
00364 
00365                         } // end other numhits
00366                      } // end this near other
00367                   } // end other is DT
00368                } // end loop over other
00369 
00370             } // end if DT13
00371 
00372             if (dt2 != NULL  &&  dt2->numHits() >= m_minDT2Hits) {
00373                DTChamberId thisid(chamberId->rawId());
00374 
00375                for (std::vector<DetId>::const_iterator otherId = chamberIds.begin();  otherId != chamberIds.end();  ++otherId) {
00376                   if (otherId->det() == DetId::Muon  &&  otherId->subdetId() == MuonSubdetId::DT) {
00377                      DTChamberId thatid(otherId->rawId());
00378                      if (thisid.rawId() != thatid.rawId()  &&  thisid.wheel() == thatid.wheel()  &&  thisid.sector() == thatid.sector()) {
00379                         MuonChamberResidual *dt2other = muonResidualsFromTrack.chamberResidual(*otherId, MuonChamberResidual::kDT2);
00380                         if (dt2other != NULL  &&  dt2other->numHits() >= m_minDT2Hits) {
00381                            double slopediff = 1000. * (dt2->global_resslope() - dt2other->global_resslope());
00382 
00383                            double length = dt2->chamberAlignable()->surface().toGlobal(LocalPoint(0,0,0)).perp() - dt2other->chamberAlignable()->surface().toGlobal(LocalPoint(0,0,0)).perp();
00384                            double residdiff = 10. * (dt2->global_residual() + length*dt2->global_resslope() - dt2other->global_residual());
00385 
00386                            if (thisid.station() == 1  &&  thatid.station() == 2) {
00387                               m_dt2_resid[thisid.wheel()+2][thisid.sector()-1][0]->Fill(qoverpt, residdiff);
00388                               m_dt2_slope[thisid.wheel()+2][thisid.sector()-1][0]->Fill(qoverpt, slopediff);
00389                               if (qoverpt > 0) {
00390                                  m_posdt2_resid[thisid.wheel()+2][thisid.sector()-1][0]->Fill(residdiff);
00391                                  m_posdt2_slope[thisid.wheel()+2][thisid.sector()-1][0]->Fill(slopediff);
00392                               }
00393                               else {
00394                                  m_negdt2_resid[thisid.wheel()+2][thisid.sector()-1][0]->Fill(residdiff);
00395                                  m_negdt2_slope[thisid.wheel()+2][thisid.sector()-1][0]->Fill(slopediff);
00396                               }
00397                            }
00398                            else if (thisid.station() == 2  &&  thatid.station() == 3) {
00399                               m_dt2_resid[thisid.wheel()+2][thisid.sector()-1][1]->Fill(qoverpt, residdiff);
00400                               m_dt2_slope[thisid.wheel()+2][thisid.sector()-1][1]->Fill(qoverpt, slopediff);
00401                               if (qoverpt > 0) {
00402                                  m_posdt2_resid[thisid.wheel()+2][thisid.sector()-1][1]->Fill(residdiff);
00403                                  m_posdt2_slope[thisid.wheel()+2][thisid.sector()-1][1]->Fill(slopediff);
00404                               }
00405                               else {
00406                                  m_negdt2_resid[thisid.wheel()+2][thisid.sector()-1][1]->Fill(residdiff);
00407                                  m_negdt2_slope[thisid.wheel()+2][thisid.sector()-1][1]->Fill(slopediff);
00408                               }
00409                            }
00410                            else if (thisid.station() == 3  &&  thatid.station() == 4) {
00411                               m_dt2_resid[thisid.wheel()+2][thisid.sector()-1][2]->Fill(qoverpt, residdiff);
00412                               m_dt2_slope[thisid.wheel()+2][thisid.sector()-1][2]->Fill(qoverpt, slopediff);
00413                               if (qoverpt > 0) {
00414                                  m_posdt2_resid[thisid.wheel()+2][thisid.sector()-1][2]->Fill(residdiff);
00415                                  m_posdt2_slope[thisid.wheel()+2][thisid.sector()-1][2]->Fill(slopediff);
00416                               }
00417                               else {
00418                                  m_negdt2_resid[thisid.wheel()+2][thisid.sector()-1][2]->Fill(residdiff);
00419                                  m_negdt2_slope[thisid.wheel()+2][thisid.sector()-1][2]->Fill(slopediff);
00420                               }
00421                            }
00422 
00423                         } // end other numhits
00424                      } // end this near other
00425                   } // end other is DT
00426                } // end loop over other
00427 
00428             } // end if DT2
00429           } // end if DT
00430 
00431           else if (chamberId->det() == DetId::Muon  &&  chamberId->subdetId() == MuonSubdetId::CSC) {
00432             MuonChamberResidual *csc = muonResidualsFromTrack.chamberResidual(*chamberId, MuonChamberResidual::kCSC);
00433 
00434             if (csc->numHits() >= m_minCSCHits) {
00435                CSCDetId thisid(chamberId->rawId());
00436 
00437                for (std::vector<DetId>::const_iterator otherId = chamberIds.begin();  otherId != chamberIds.end();  ++otherId) {
00438                   if (otherId->det() == DetId::Muon  &&  otherId->subdetId() == MuonSubdetId::CSC) {
00439                      CSCDetId thatid(otherId->rawId());
00440                      if (thisid.rawId() != thatid.rawId()  &&  thisid.endcap() == thatid.endcap()) {
00441                         MuonChamberResidual *cscother = muonResidualsFromTrack.chamberResidual(*otherId, MuonChamberResidual::kCSC);
00442                         if (cscother != NULL  &&  cscother->numHits() >= m_minCSCHits) {
00443                            double slopediff = 1000. * (csc->global_resslope() - cscother->global_resslope());
00444 
00445                            double length = csc->chamberAlignable()->surface().toGlobal(LocalPoint(0,0,0)).z() - cscother->chamberAlignable()->surface().toGlobal(LocalPoint(0,0,0)).z();
00446                            double residdiff = 10. * (csc->global_residual() + length*csc->global_resslope() - cscother->global_residual());
00447 
00448                            int thischamber = thisid.chamber();
00449                            int thisring = thisid.ring();
00450                            if (thisid.station() == 1  &&  (thisring == 1  ||  thisring == 4)) {
00451                               thischamber = (thischamber - 1) / 2 + 1;
00452                               thisring = 1;
00453                            }
00454 
00455                            if (thisring == thatid.ring()  &&  thischamber == thatid.chamber()) {
00456                               if (thisring == 2  &&  thisid.station() == 1  &&  thatid.station() == 2) {
00457                                  m_cscouter_resid[thisid.endcap()-1][thischamber-1][0]->Fill(qoverpt, residdiff);
00458                                  m_cscouter_slope[thisid.endcap()-1][thischamber-1][0]->Fill(qoverpt, slopediff);
00459                                  if (qoverpt > 0) {
00460                                     m_poscscouter_resid[thisid.endcap()-1][thischamber-1][0]->Fill(residdiff);
00461                                     m_poscscouter_slope[thisid.endcap()-1][thischamber-1][0]->Fill(slopediff);
00462                                  }
00463                                  else {
00464                                     m_negcscouter_resid[thisid.endcap()-1][thischamber-1][0]->Fill(residdiff);
00465                                     m_negcscouter_slope[thisid.endcap()-1][thischamber-1][0]->Fill(slopediff);
00466                                  }
00467                               }
00468                               else if (thisring == 2  &&  thisid.station() == 2  &&  thatid.station() == 3) {
00469                                  m_cscouter_resid[thisid.endcap()-1][thischamber-1][1]->Fill(qoverpt, residdiff);
00470                                  m_cscouter_slope[thisid.endcap()-1][thischamber-1][1]->Fill(qoverpt, slopediff);
00471                                  if (qoverpt > 0) {
00472                                     m_poscscouter_resid[thisid.endcap()-1][thischamber-1][1]->Fill(residdiff);
00473                                     m_poscscouter_slope[thisid.endcap()-1][thischamber-1][1]->Fill(slopediff);
00474                                  }
00475                                  else {
00476                                     m_negcscouter_resid[thisid.endcap()-1][thischamber-1][1]->Fill(residdiff);
00477                                     m_negcscouter_slope[thisid.endcap()-1][thischamber-1][1]->Fill(slopediff);
00478                                  }
00479                               }
00480                               else if (thisring == 1  &&  thisid.station() == 1  &&  thatid.station() == 2) {
00481                                  m_cscinner_resid[thisid.endcap()-1][thischamber-1][0]->Fill(qoverpt, residdiff);
00482                                  m_cscinner_slope[thisid.endcap()-1][thischamber-1][0]->Fill(qoverpt, slopediff);
00483                                  if (qoverpt > 0) {
00484                                     m_poscscinner_resid[thisid.endcap()-1][thischamber-1][0]->Fill(residdiff);
00485                                     m_poscscinner_slope[thisid.endcap()-1][thischamber-1][0]->Fill(slopediff);
00486                                  }
00487                                  else {
00488                                     m_negcscinner_resid[thisid.endcap()-1][thischamber-1][0]->Fill(residdiff);
00489                                     m_negcscinner_slope[thisid.endcap()-1][thischamber-1][0]->Fill(slopediff);
00490                                  }
00491                               }
00492                               else if (thisring == 1  &&  thisid.station() == 2  &&  thatid.station() == 3) {
00493                                  m_cscinner_resid[thisid.endcap()-1][thischamber-1][1]->Fill(qoverpt, residdiff);
00494                                  m_cscinner_slope[thisid.endcap()-1][thischamber-1][1]->Fill(qoverpt, slopediff);
00495                                  if (qoverpt > 0) {
00496                                     m_poscscinner_resid[thisid.endcap()-1][thischamber-1][1]->Fill(residdiff);
00497                                     m_poscscinner_slope[thisid.endcap()-1][thischamber-1][1]->Fill(slopediff);
00498                                  }
00499                                  else {
00500                                     m_negcscinner_resid[thisid.endcap()-1][thischamber-1][1]->Fill(residdiff);
00501                                     m_negcscinner_slope[thisid.endcap()-1][thischamber-1][1]->Fill(slopediff);
00502                                  }
00503                               }
00504                               else if (thisring == 1  &&  thisid.station() == 3  &&  thatid.station() == 4) {
00505                                  m_cscinner_resid[thisid.endcap()-1][thischamber-1][2]->Fill(qoverpt, residdiff);
00506                                  m_cscinner_slope[thisid.endcap()-1][thischamber-1][2]->Fill(qoverpt, slopediff);
00507                                  if (qoverpt > 0) {
00508                                     m_poscscinner_resid[thisid.endcap()-1][thischamber-1][2]->Fill(residdiff);
00509                                     m_poscscinner_slope[thisid.endcap()-1][thischamber-1][2]->Fill(slopediff);
00510                                  }
00511                                  else {
00512                                     m_negcscinner_resid[thisid.endcap()-1][thischamber-1][2]->Fill(residdiff);
00513                                     m_negcscinner_slope[thisid.endcap()-1][thischamber-1][2]->Fill(slopediff);
00514                                  }
00515                               }
00516                            }
00517                         } // end other numhits
00518                      } // end this near other
00519                   } // end other is DT
00520                } // end loop over other
00521 
00522             } // end if csc
00523           } // end if CSC
00524 
00525         } // end loop over chamberIds
00526       } // end if refit is okay
00527     } // end if track pT is within range
00528   } // end loop over tracks
00529 }
00530 
00531 void AlignmentMonitorSegmentDifferences::afterAlignment(const edm::EventSetup &iSetup) {
00532 }
00533 
00534 //
00535 // constructors and destructor
00536 //
00537 
00538 // AlignmentMonitorSegmentDifferences::AlignmentMonitorSegmentDifferences(const AlignmentMonitorSegmentDifferences& rhs)
00539 // {
00540 //    // do actual copying here;
00541 // }
00542 
00543 //
00544 // assignment operators
00545 //
00546 // const AlignmentMonitorSegmentDifferences& AlignmentMonitorSegmentDifferences::operator=(const AlignmentMonitorSegmentDifferences& rhs)
00547 // {
00548 //   //An exception safe implementation is
00549 //   AlignmentMonitorSegmentDifferences temp(rhs);
00550 //   swap(rhs);
00551 //
00552 //   return *this;
00553 // }
00554 
00555 //
00556 // const member functions
00557 //
00558 
00559 //
00560 // static member functions
00561 //
00562 
00563 //
00564 // SEAL definitions
00565 //
00566 
00567 DEFINE_EDM_PLUGIN(AlignmentMonitorPluginFactory, AlignmentMonitorSegmentDifferences, "AlignmentMonitorSegmentDifferences");