Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:08:27

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiPixelPhase1Summary
0004 // Class:      SiPixelBarycenter
0005 //
0006 /**\class 
0007 
0008  Description: Create the pixel barycenter plots
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Danilo Meuser
0015 //         Created:  26th May 2021
0016 //
0017 //
0018 #include "DQM/SiPixelPhase1Summary/interface/SiPixelBarycenter.h"
0019 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0020 #include "DQM/SiPixelPhase1Summary/interface/SiPixelBarycenterHelper.h"
0021 
0022 #include <string>
0023 #include <iostream>
0024 
0025 using namespace std;
0026 using namespace edm;
0027 
0028 SiPixelBarycenter::SiPixelBarycenter(const edm::ParameterSet& iConfig)
0029     : DQMEDHarvester(iConfig),
0030       alignmentToken_(esConsumes<edm::Transition::EndRun>()),
0031       gprToken_(esConsumes<edm::Transition::EndRun>()),
0032       trackerTopologyToken_(esConsumes<edm::Transition::EndRun>()) {
0033   LogInfo("PixelDQM") << "SiPixelBarycenter::SiPixelBarycenter: Got DQM BackEnd interface" << endl;
0034 }
0035 
0036 void SiPixelBarycenter::beginRun(edm::Run const& run, edm::EventSetup const& eSetup) {}
0037 
0038 void SiPixelBarycenter::dqmEndJob(DQMStore::IBooker& iBooker, DQMStore::IGetter& iGetter) {}
0039 
0040 void SiPixelBarycenter::dqmEndRun(DQMStore::IBooker& iBooker,
0041                                   DQMStore::IGetter& iGetter,
0042                                   const edm::Run& iRun,
0043                                   edm::EventSetup const& c) {
0044   bookBarycenterHistograms(iBooker);
0045 
0046   const Alignments* alignments = &c.getData(alignmentToken_);
0047   const Alignments* gpr = &c.getData(gprToken_);
0048   const TrackerTopology* tTopo = &c.getData(trackerTopologyToken_);
0049 
0050   fillBarycenterHistograms(iBooker, iGetter, alignments->m_align, gpr->m_align, *tTopo);
0051 }
0052 //------------------------------------------------------------------
0053 // Used to book the barycenter histograms
0054 //------------------------------------------------------------------
0055 void SiPixelBarycenter::bookBarycenterHistograms(DQMStore::IBooker& iBooker) {
0056   iBooker.cd();
0057 
0058   iBooker.setCurrentFolder("PixelPhase1/Barycenter");
0059   //Book one histogram for each subdetector
0060   for (const std::string& subdetector : subdetectors_) {
0061     barycenters_[subdetector] =
0062         iBooker.book1D("barycenters_" + subdetector,
0063                        "Position of the barycenter for " + subdetector + ";Coordinate;Position [mm]",
0064                        3,
0065                        0.5,
0066                        3.5);
0067     barycenters_[subdetector]->setBinLabel(1, "x");
0068     barycenters_[subdetector]->setBinLabel(2, "y");
0069     barycenters_[subdetector]->setBinLabel(3, "z");
0070   }
0071 
0072   //Reset the iBooker
0073   iBooker.setCurrentFolder("PixelPhase1/");
0074 }
0075 
0076 //------------------------------------------------------------------
0077 // Fill the Barycenter histograms
0078 //------------------------------------------------------------------
0079 void SiPixelBarycenter::fillBarycenterHistograms(DQMStore::IBooker& iBooker,
0080                                                  DQMStore::IGetter& iGetter,
0081                                                  const std::vector<AlignTransform>& input,
0082                                                  const std::vector<AlignTransform>& GPR,
0083                                                  const TrackerTopology& tTopo) {
0084   const auto GPR_translation_pixel = GPR[0].translation();
0085   const std::map<DQMBarycenter::coordinate, float> GPR_pixel = {{DQMBarycenter::t_x, GPR_translation_pixel.x()},
0086                                                                 {DQMBarycenter::t_y, GPR_translation_pixel.y()},
0087                                                                 {DQMBarycenter::t_z, GPR_translation_pixel.z()}};
0088 
0089   DQMBarycenter::TkAlBarycenters barycenters;
0090   barycenters.computeBarycenters(input, tTopo, GPR_pixel);
0091 
0092   auto Xbarycenters = barycenters.getX();
0093   auto Ybarycenters = barycenters.getY();
0094   auto Zbarycenters = barycenters.getZ();
0095 
0096   //Fill histogram for each subdetector
0097   for (std::size_t i = 0; i < subdetectors_.size(); ++i) {
0098     barycenters_[subdetectors_[i]]->setBinContent(1, Xbarycenters[i]);
0099     barycenters_[subdetectors_[i]]->setBinContent(2, Ybarycenters[i]);
0100     barycenters_[subdetectors_[i]]->setBinContent(3, Zbarycenters[i]);
0101 
0102     // zero the errors for better comparison display
0103     barycenters_[subdetectors_[i]]->setBinError(1, 0.);
0104     barycenters_[subdetectors_[i]]->setBinError(2, 0.);
0105     barycenters_[subdetectors_[i]]->setBinError(3, 0.);
0106   }
0107 }
0108 
0109 //define this as a plug-in
0110 DEFINE_FWK_MODULE(SiPixelBarycenter);