Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 11:56:17

0001 // \file AlignmentRcdScan.cpp
0002 //
0003 // \author    : Andreas Mussgiller
0004 // Revision   : $Revision: 1.1 $
0005 // last update: $Date: 2010/06/01 07:45:46 $
0006 // by         : $Author: mussgill $
0007 
0008 #include <string>
0009 #include <map>
0010 #include <vector>
0011 
0012 #include <TMath.h>
0013 
0014 #include "FWCore/Framework/interface/EventSetup.h"
0015 #include "FWCore/Framework/interface/Event.h"
0016 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0017 #include "FWCore/Framework/interface/MakerMacros.h"
0018 #include "FWCore/Framework/interface/ESHandle.h"
0019 #include "CondFormats/Alignment/interface/Alignments.h"
0020 #include "CondFormats/Alignment/interface/AlignTransform.h"
0021 #include "DataFormats/DetId/interface/DetId.h"
0022 
0023 #include "CondFormats/AlignmentRecord/interface/TrackerAlignmentRcd.h"
0024 #include "CondFormats/AlignmentRecord/interface/DTAlignmentRcd.h"
0025 #include "CondFormats/AlignmentRecord/interface/CSCAlignmentRcd.h"
0026 
0027 #include "FWCore/Framework/interface/ESWatcher.h"
0028 
0029 class AlignmentRcdScan : public edm::one::EDAnalyzer<> {
0030 public:
0031   enum Mode { Unknown = 0, Tk = 1, DT = 2, CSC = 3 };
0032 
0033   explicit AlignmentRcdScan(const edm::ParameterSet& iConfig);
0034   ~AlignmentRcdScan();
0035 
0036   virtual void analyze(const edm::Event& evt, const edm::EventSetup& evtSetup);
0037 
0038 private:
0039   void inspectRecord(const std::string& rcdname, const edm::Event& evt, const Alignments* alignments);
0040 
0041   int mode_;
0042   bool verbose_;
0043 
0044   edm::ESWatcher<TrackerAlignmentRcd> watchTk_;
0045   edm::ESWatcher<DTAlignmentRcd> watchDT_;
0046   edm::ESWatcher<CSCAlignmentRcd> watchCSC_;
0047 
0048   const edm::ESGetToken<Alignments, TrackerAlignmentRcd> tkAliToken_;
0049   const edm::ESGetToken<Alignments, DTAlignmentRcd> dtAliToken_;
0050   const edm::ESGetToken<Alignments, CSCAlignmentRcd> cscAliToken_;
0051 
0052   Alignments* refAlignments_;
0053 };
0054 
0055 AlignmentRcdScan::AlignmentRcdScan(const edm::ParameterSet& iConfig)
0056     : verbose_(iConfig.getUntrackedParameter<bool>("verbose")),
0057       tkAliToken_(esConsumes()),
0058       dtAliToken_(esConsumes()),
0059       cscAliToken_(esConsumes()),
0060       refAlignments_(0) {
0061   std::string modestring = iConfig.getUntrackedParameter<std::string>("mode");
0062   if (modestring == "Tk") {
0063     mode_ = Tk;
0064   } else if (modestring == "DT") {
0065     mode_ = DT;
0066   } else if (modestring == "CSC") {
0067     mode_ = CSC;
0068   } else {
0069     mode_ = Unknown;
0070   }
0071 
0072   if (mode_ == Unknown) {
0073     throw cms::Exception("BadConfig") << "Mode " << modestring << " not known";
0074   }
0075 }
0076 
0077 AlignmentRcdScan::~AlignmentRcdScan() {
0078   if (refAlignments_)
0079     delete refAlignments_;
0080 }
0081 
0082 void AlignmentRcdScan::analyze(const edm::Event& evt, const edm::EventSetup& evtSetup) {
0083   if (mode_ == Tk && watchTk_.check(evtSetup)) {
0084     const Alignments* alignments = &evtSetup.getData(tkAliToken_);
0085     inspectRecord("TrackerAlignmentRcd", evt, alignments);
0086   }
0087   if (mode_ == DT && watchDT_.check(evtSetup)) {
0088     const Alignments* alignments = &evtSetup.getData(dtAliToken_);
0089     inspectRecord("DTAlignmentRcd", evt, alignments);
0090   }
0091   if (mode_ == CSC && watchCSC_.check(evtSetup)) {
0092     const Alignments* alignments = &evtSetup.getData(cscAliToken_);
0093     inspectRecord("CSCAlignmentRcd", evt, alignments);
0094   }
0095 }
0096 
0097 void AlignmentRcdScan::inspectRecord(const std::string& rcdname, const edm::Event& evt, const Alignments* alignments) {
0098   edm::LogPrint("inspectRecord") << rcdname << " content starting from run " << evt.run();
0099 
0100   if (verbose_ == false) {
0101     edm::LogPrint("inspectRecord") << std::endl;
0102     return;
0103   }
0104 
0105   edm::LogPrint("inspectRecord") << " with " << alignments->m_align.size() << " entries" << std::endl;
0106 
0107   if (refAlignments_) {
0108     edm::LogPrint("inspectRecord") << "  Compared to previous record:" << std::endl;
0109 
0110     double meanX = 0;
0111     double rmsX = 0;
0112     double meanY = 0;
0113     double rmsY = 0;
0114     double meanZ = 0;
0115     double rmsZ = 0;
0116     double meanR = 0;
0117     double rmsR = 0;
0118     double dPhi;
0119     double meanPhi = 0;
0120     double rmsPhi = 0;
0121 
0122     std::vector<AlignTransform>::const_iterator iref = refAlignments_->m_align.begin();
0123     for (std::vector<AlignTransform>::const_iterator i = alignments->m_align.begin(); i != alignments->m_align.end();
0124          ++i, ++iref) {
0125       meanX += i->translation().x() - iref->translation().x();
0126       rmsX += pow(i->translation().x() - iref->translation().x(), 2);
0127 
0128       meanY += i->translation().y() - iref->translation().y();
0129       rmsY += pow(i->translation().y() - iref->translation().y(), 2);
0130 
0131       meanZ += i->translation().z() - iref->translation().z();
0132       rmsZ += pow(i->translation().z() - iref->translation().z(), 2);
0133 
0134       meanR += i->translation().perp() - iref->translation().perp();
0135       rmsR += pow(i->translation().perp() - iref->translation().perp(), 2);
0136 
0137       dPhi = i->translation().phi() - iref->translation().phi();
0138       if (dPhi > M_PI)
0139         dPhi -= 2.0 * M_PI;
0140       if (dPhi < -M_PI)
0141         dPhi += 2.0 * M_PI;
0142 
0143       meanPhi += dPhi;
0144       rmsPhi += dPhi * dPhi;
0145     }
0146 
0147     meanX /= alignments->m_align.size();
0148     rmsX /= alignments->m_align.size();
0149     meanY /= alignments->m_align.size();
0150     rmsY /= alignments->m_align.size();
0151     meanZ /= alignments->m_align.size();
0152     rmsZ /= alignments->m_align.size();
0153     meanR /= alignments->m_align.size();
0154     rmsR /= alignments->m_align.size();
0155     meanPhi /= alignments->m_align.size();
0156     rmsPhi /= alignments->m_align.size();
0157 
0158     edm::LogPrint("inspectRecord") << "    mean X shift:   " << std::setw(12) << std::scientific << std::setprecision(3)
0159                                    << meanX << " (RMS = " << sqrt(rmsX) << ")" << std::endl;
0160     edm::LogPrint("inspectRecord") << "    mean Y shift:   " << std::setw(12) << std::scientific << std::setprecision(3)
0161                                    << meanY << " (RMS = " << sqrt(rmsY) << ")" << std::endl;
0162     edm::LogPrint("inspectRecord") << "    mean Z shift:   " << std::setw(12) << std::scientific << std::setprecision(3)
0163                                    << meanZ << " (RMS = " << sqrt(rmsZ) << ")" << std::endl;
0164     edm::LogPrint("inspectRecord") << "    mean R shift:   " << std::setw(12) << std::scientific << std::setprecision(3)
0165                                    << meanR << " (RMS = " << sqrt(rmsR) << ")" << std::endl;
0166     edm::LogPrint("inspectRecord") << "    mean Phi shift: " << std::setw(12) << std::scientific << std::setprecision(3)
0167                                    << meanPhi << " (RMS = " << sqrt(rmsPhi) << ")" << std::endl;
0168 
0169     delete refAlignments_;
0170   }
0171 
0172   refAlignments_ = new Alignments(*alignments);
0173 
0174   edm::LogPrint("inspectRecord") << std::endl;
0175 }
0176 
0177 //define this as a plug-in
0178 DEFINE_FWK_MODULE(AlignmentRcdScan);