Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:09:16

0001 // -*- C++ -*-
0002 //
0003 // Package:    DQMHOAlCaRecoStream
0004 // Class:      DQMHOAlCaRecoStream
0005 //
0006 /**\class DQMHOAlCaRecoStream DQMHOAlCaRecoStream.cc
0007  DQMOffline/DQMHOAlCaRecoStream/src/DQMHOAlCaRecoStream.cc
0008 
0009  Description: <one line class summary>
0010 
0011  Implementation:
0012      <Notes on implementation>
0013 */
0014 //
0015 // Original Author:  Gobinda Majumder
0016 //         Created:  Mon Mar  2 12:33:08 CET 2009
0017 //
0018 //
0019 
0020 // system include files
0021 #include <string>
0022 
0023 // user include files
0024 #include "DataFormats/HcalCalibObjects/interface/HOCalibVariables.h"
0025 
0026 #include "FWCore/Framework/interface/Frameworkfwd.h"
0027 #include "FWCore/Utilities/interface/InputTag.h"
0028 #include "FWCore/Framework/interface/Event.h"
0029 #include "FWCore/Framework/interface/MakerMacros.h"
0030 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0031 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0032 
0033 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0034 #include "DQMServices/Core/interface/DQMStore.h"
0035 #include "DQMServices/Core/interface/DQMStore.h"
0036 #include "FWCore/ServiceRegistry/interface/Service.h"
0037 
0038 //
0039 // class decleration
0040 //
0041 
0042 class DQMHOAlCaRecoStream : public DQMEDAnalyzer {
0043 public:
0044   explicit DQMHOAlCaRecoStream(const edm::ParameterSet &);
0045   ~DQMHOAlCaRecoStream() override;
0046 
0047 private:
0048   void analyze(const edm::Event &, const edm::EventSetup &) override;
0049   void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override;
0050 
0051   // ----------member data ---------------------------
0052 
0053   MonitorElement *hMuonMultipl;
0054   MonitorElement *hMuonMom;
0055   MonitorElement *hMuonEta;
0056   MonitorElement *hMuonPhi;
0057 
0058   MonitorElement *hDirCosine;
0059   MonitorElement *hHOTime;
0060 
0061   MonitorElement *hSigRing[5];
0062   MonitorElement *hPedRing[5];
0063   MonitorElement *hSignal3x3[9];
0064 
0065   int Nevents;
0066   int Nmuons;
0067 
0068   std::string theRootFileName;
0069   std::string folderName_;
0070   double m_sigmaValue;
0071 
0072   double m_lowRadPosInMuch;
0073   double m_highRadPosInMuch;
0074 
0075   int m_nbins;
0076   double m_lowEdge;
0077   double m_highEdge;
0078 
0079   bool saveToFile_;
0080   edm::EDGetTokenT<HOCalibVariableCollection> hoCalibVariableCollectionTag;
0081 };
0082 
0083 //
0084 // constructors and destructor
0085 //
0086 DQMHOAlCaRecoStream::DQMHOAlCaRecoStream(const edm::ParameterSet &iConfig)
0087     : hoCalibVariableCollectionTag(
0088           consumes<HOCalibVariableCollection>(iConfig.getParameter<edm::InputTag>("hoCalibVariableCollectionTag"))) {
0089   // now do what ever initialization is needed
0090 
0091   theRootFileName = iConfig.getUntrackedParameter<std::string>("RootFileName", "tmp.root");
0092   folderName_ = iConfig.getUntrackedParameter<std::string>("folderName");
0093   m_sigmaValue = iConfig.getUntrackedParameter<double>("sigmaval", 0.2);
0094   m_lowRadPosInMuch = iConfig.getUntrackedParameter<double>("lowradposinmuch", 400.0);
0095   m_highRadPosInMuch = iConfig.getUntrackedParameter<double>("highradposinmuch", 480.0);
0096   m_lowEdge = iConfig.getUntrackedParameter<double>("lowedge", -2.0);
0097   m_highEdge = iConfig.getUntrackedParameter<double>("highedge", 6.0);
0098   m_nbins = iConfig.getUntrackedParameter<int>("nbins", 40);
0099   saveToFile_ = iConfig.getUntrackedParameter<bool>("saveToFile", false);
0100 }
0101 
0102 DQMHOAlCaRecoStream::~DQMHOAlCaRecoStream() {
0103   // do anything here that needs to be done at desctruction time
0104   // (e.g. close files, deallocate resources etc.)
0105 }
0106 
0107 //
0108 // member functions
0109 //
0110 
0111 // ------------ method called to for each event  ------------
0112 void DQMHOAlCaRecoStream::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0113   Nevents++;
0114 
0115   edm::Handle<HOCalibVariableCollection> HOCalib;
0116   bool isCosMu = true;
0117 
0118   iEvent.getByToken(hoCalibVariableCollectionTag, HOCalib);
0119 
0120   if (!HOCalib.isValid()) {
0121     LogDebug("") << "DQMHOAlCaRecoStream:: Error! can't get HOCalib product!" << std::endl;
0122     return;
0123   }
0124 
0125   if (isCosMu) {
0126     hMuonMultipl->Fill((*HOCalib).size(), 1.);
0127     if (!(*HOCalib).empty()) {
0128       for (HOCalibVariableCollection::const_iterator hoC = (*HOCalib).begin(); hoC != (*HOCalib).end(); hoC++) {
0129         // OK!!!!
0130         float okt = 2.;
0131         double okx = std::pow((*hoC).trkvx, okt) + std::pow((*hoC).trkvy, okt);
0132         ///////
0133         double dr = std::pow(okx, 0.5);
0134         if (dr < m_lowRadPosInMuch || dr > m_highRadPosInMuch)
0135           continue;
0136 
0137         if ((*hoC).isect < 0)
0138           continue;
0139         if (fabs((*hoC).trkth - acos(-1.) / 2) < 0.000001)
0140           continue;
0141         int ieta = int((std::abs((*hoC).isect) % 10000) / 100.) - 30;
0142 
0143         if (std::abs(ieta) >= 16)
0144           continue;
0145 
0146         Nmuons++;
0147 
0148         hMuonMom->Fill((*hoC).trkmm, 1.0);
0149         hMuonEta->Fill(-log(tan((*hoC).trkth / 2.0)), 1.0);
0150         hMuonPhi->Fill((*hoC).trkph, 1.0);
0151         hDirCosine->Fill((*hoC).hoang, 1.0);
0152         hHOTime->Fill((*hoC).htime, 1.0);
0153 
0154         double energy = (*hoC).hosig[4];
0155         double pedval = (*hoC).hocro;
0156         int iring = 0;
0157         if (ieta >= -15 && ieta <= -11) {
0158           iring = -2;
0159         } else if (ieta >= -10 && ieta <= -5) {
0160           iring = -1;
0161         } else if (ieta >= 5 && ieta <= 10) {
0162           iring = 1;
0163         } else if (ieta >= 11 && ieta <= 15) {
0164           iring = 2;
0165         }
0166 
0167         hSigRing[iring + 2]->Fill(energy, 1.0);
0168         hPedRing[iring + 2]->Fill(pedval, 1.0);
0169 
0170         for (int k = 0; k < 9; k++) {
0171           hSignal3x3[k]->Fill((*hoC).hosig[k]);
0172         }
0173       }  // for (HOCalibVariableCollection::const_iterator hoC=(*HOCalib).begin()
0174     }    // if ((*HOCalib).size() >0 ) {
0175   }      // if (isCosMu) {
0176 }
0177 
0178 // ------------ method called once each job just before starting event loop
0179 // ------------
0180 void DQMHOAlCaRecoStream::bookHistograms(DQMStore::IBooker &ibooker,
0181                                          edm::Run const &irun,
0182                                          edm::EventSetup const &isetup) {
0183   ibooker.setCurrentFolder(folderName_);
0184 
0185   char title[200];
0186   char name[200];
0187 
0188   hMuonMom = ibooker.book1D("hMuonMom", "Muon momentum (GeV)", 50, -100, 100);
0189   hMuonMom->setAxisTitle("Muon momentum (GeV)", 1);
0190 
0191   hMuonEta = ibooker.book1D("hMuonEta", "Pseudo-rapidity of muon", 50, -1.5, 1.5);
0192   hMuonEta->setAxisTitle("Pseudo-rapidity of muon", 1);
0193 
0194   hMuonPhi = ibooker.book1D("hMuonPhi", "Azimuthal angle of muon", 24, -acos(-1), acos(-1));
0195   hMuonPhi->setAxisTitle("Azimuthal angle of muon", 1);
0196 
0197   hMuonMultipl = ibooker.book1D("hMuonMultipl", "Muon Multiplicity", 10, 0.5, 10.5);
0198   hMuonMultipl->setAxisTitle("Muon Multiplicity", 1);
0199 
0200   hDirCosine = ibooker.book1D("hDirCosine", "Direction Cosine of muon at HO tower", 50, -1., 1.);
0201   hDirCosine->setAxisTitle("Direction Cosine of muon at HO tower", 1);
0202 
0203   hHOTime = ibooker.book1D("hHOTime", "HO time distribution", 60, -20, 100.);
0204   hHOTime->setAxisTitle("HO time distribution", 1);
0205 
0206   for (int i = 0; i < 5; i++) {
0207     sprintf(name, "hSigRing_%i", i - 2);
0208     sprintf(title, "HO signal in Ring_%i", i - 2);
0209     hSigRing[i] = ibooker.book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
0210     hSigRing[i]->setAxisTitle(title, 1);
0211 
0212     sprintf(name, "hPedRing_%i", i - 2);
0213     sprintf(title, "HO Pedestal in Ring_%i", i - 2);
0214     hPedRing[i] = ibooker.book1D(name, title, m_nbins, m_lowEdge, m_highEdge);
0215     hPedRing[i]->setAxisTitle(title, 1);
0216   }
0217 
0218   for (int i = -1; i <= 1; i++) {
0219     for (int j = -1; j <= 1; j++) {
0220       int k = 3 * (i + 1) + j + 1;
0221 
0222       sprintf(title, "hSignal3x3_deta%i_dphi%i", i, j);
0223       hSignal3x3[k] = ibooker.book1D(title, title, m_nbins, m_lowEdge, m_highEdge);
0224       hSignal3x3[k]->setAxisTitle(title, 1);
0225     }
0226   }
0227 
0228   Nevents = 0;
0229   Nmuons = 0;
0230 }
0231 
0232 DEFINE_FWK_MODULE(DQMHOAlCaRecoStream);