Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:32

0001 // -*- C++ -*-
0002 //
0003 // Package:    SiStripBaselineAnalyzer
0004 // Class:      SiStripBaselineAnalyzer
0005 //
0006 /**\class SiStripBaselineAnalyzer SiStripBaselineAnalyzer.cc Validation/SiStripAnalyzer/src/SiStripBaselineAnalyzer.cc
0007 
0008  Description: <one line class summary>
0009 
0010  Implementation:
0011      <Notes on implementation>
0012 */
0013 //
0014 // Original Author:  Ivan Amos Cali
0015 //         Created:  Mon Jul 28 14:10:52 CEST 2008
0016 //
0017 //
0018 
0019 // system include files
0020 #include <memory>
0021 #include <iostream>
0022 
0023 // user include files
0024 #include "FWCore/Framework/interface/Frameworkfwd.h"
0025 #include "FWCore/Framework/interface/MakerMacros.h"
0026 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0027 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0028 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0029 
0030 #include "FWCore/Framework/interface/Event.h"
0031 #include "FWCore/Framework/interface/MakerMacros.h"
0032 #include "FWCore/Framework/interface/ESHandle.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034 #include "DataFormats/Common/interface/DetSet.h"
0035 #include "DataFormats/Common/interface/DetSetVector.h"
0036 #include "DataFormats/Common/interface/DetSetVectorNew.h"
0037 
0038 #include "DataFormats/SiStripDigi/interface/SiStripProcessedRawDigi.h"
0039 #include "DataFormats/SiStripDigi/interface/SiStripRawDigi.h"
0040 #include "DataFormats/SiStripCluster/interface/SiStripCluster.h"
0041 #include "DataFormats/SiStripCluster/interface/SiStripClusterCollection.h"
0042 
0043 #include "FWCore/ServiceRegistry/interface/Service.h"
0044 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0045 #include "CommonTools/Utils/interface/TFileDirectory.h"
0046 
0047 //ROOT inclusion
0048 #include "TH1F.h"
0049 #include "TProfile.h"
0050 
0051 //
0052 // class decleration
0053 //
0054 
0055 class SiStripBaselineComparator : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0056 public:
0057   explicit SiStripBaselineComparator(const edm::ParameterSet&);
0058   ~SiStripBaselineComparator() override;
0059   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0060 
0061 private:
0062   void beginJob() override;
0063   void analyze(const edm::Event&, const edm::EventSetup&) override;
0064   void endJob() override;
0065 
0066   edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster> > srcClusters_;
0067   edm::EDGetTokenT<edmNew::DetSetVector<SiStripCluster> > srcClusters2_;
0068 
0069   edm::Service<TFileService> fs_;
0070 
0071   TH1F* h1_nOldClusters_;
0072   TH1F* h1_nMatchedClusters_;
0073   TH1F* h1_nSplitClusters_;
0074   TProfile* h1_matchingMult_;
0075   TProfile* h1_matchedWidth_;
0076   TProfile* h1_matchedCharges_;
0077 
0078   const uint16_t nbins_clusterSize_ = 128;
0079   const uint16_t min_clusterSize_ = 0;
0080   const uint16_t max_clusterSize_ = 128;
0081 };
0082 
0083 SiStripBaselineComparator::SiStripBaselineComparator(const edm::ParameterSet& conf) {
0084   usesResource(TFileService::kSharedResource);
0085 
0086   srcClusters_ = consumes<edmNew::DetSetVector<SiStripCluster> >(conf.getParameter<edm::InputTag>("srcClusters"));
0087   srcClusters2_ = consumes<edmNew::DetSetVector<SiStripCluster> >(conf.getParameter<edm::InputTag>("srcClusters2"));
0088 
0089   h1_nOldClusters_ = fs_->make<TH1F>(
0090       "nOldClusters", "nOldClusters;ClusterSize", nbins_clusterSize_, min_clusterSize_, max_clusterSize_);
0091   h1_nMatchedClusters_ = fs_->make<TH1F>(
0092       "nMatchedClusters", "nMatchedClusters;ClusterSize", nbins_clusterSize_, min_clusterSize_, max_clusterSize_);
0093   h1_nSplitClusters_ = fs_->make<TH1F>("nSplitClusters",
0094                                        "nMatchedClusters;ClusterSize; n Split Clusters",
0095                                        nbins_clusterSize_,
0096                                        min_clusterSize_,
0097                                        max_clusterSize_);
0098   h1_matchingMult_ = fs_->make<TProfile>("matchingMult",
0099                                          "matchingMult;ClusterSize; average number of clusters if split",
0100                                          nbins_clusterSize_,
0101                                          min_clusterSize_,
0102                                          max_clusterSize_);
0103   h1_matchedWidth_ = fs_->make<TProfile>("matchedWidth",
0104                                          "matchedWidth;ClusterSize;average #Delta Width",
0105                                          nbins_clusterSize_,
0106                                          min_clusterSize_,
0107                                          max_clusterSize_);
0108   h1_matchedCharges_ = fs_->make<TProfile>("matchedCharges",
0109                                            "matchedCharges;ClusterSize;cluster Charge ratio",
0110                                            nbins_clusterSize_,
0111                                            min_clusterSize_,
0112                                            max_clusterSize_);
0113 }
0114 
0115 SiStripBaselineComparator::~SiStripBaselineComparator() = default;
0116 
0117 void SiStripBaselineComparator::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0118   edm::ParameterSetDescription desc;
0119   desc.add<edm::InputTag>("srcClusters", edm::InputTag("siStripClusters"));
0120   desc.add<edm::InputTag>("srcClusters2", edm::InputTag("moddedsiStripClusters"));
0121   descriptions.add("siStripBaselineComparator", desc);
0122 }
0123 
0124 void SiStripBaselineComparator::analyze(const edm::Event& e, const edm::EventSetup& es) {
0125   edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters;
0126   e.getByToken(srcClusters_, clusters);
0127   edm::Handle<edmNew::DetSetVector<SiStripCluster> > clusters2;
0128   e.getByToken(srcClusters2_, clusters2);
0129 
0130   for (auto const& set : *clusters) {
0131     for (auto const& clus : set) {
0132       h1_nOldClusters_->Fill(clus.amplitudes().size(), 1);
0133       int nMatched = 0;
0134       const int charge1 = std::accumulate(clus.amplitudes().begin(), clus.amplitudes().end(), 0);
0135       std::vector<int> matchedWidths;
0136       std::vector<int> matchedCharges;
0137 
0138       //scan other set of clusters
0139       for (auto const& set2 : *clusters2) {
0140         if (set.id() != set2.id())
0141           continue;
0142         for (auto const& clus2 : set2) {
0143           const int charge2 = std::accumulate(clus2.amplitudes().begin(), clus2.amplitudes().end(), 0);
0144           if ((clus.firstStrip() <= clus2.firstStrip()) &&
0145               (clus2.firstStrip() < clus.firstStrip() + clus.amplitudes().size())) {
0146             matchedWidths.push_back(clus2.amplitudes().size());
0147             matchedCharges.push_back(charge2);
0148             if (nMatched == 0) {
0149               h1_nMatchedClusters_->Fill(clus.amplitudes().size(), 1);
0150             } else if (nMatched == 1) {
0151               h1_nSplitClusters_->Fill(clus.amplitudes().size(), 1);
0152             }
0153             ++nMatched;
0154           }
0155         }
0156       }
0157       for (int i = 0; i < nMatched; i++) {
0158         if (matchedWidths.at(i) - clus.amplitudes().size() < 1000)
0159           h1_matchedWidth_->Fill(clus.amplitudes().size(), matchedWidths.at(i) - clus.amplitudes().size());
0160         if (charge1 != 0 && matchedCharges.at(i) / (float)charge1 < 10000000)
0161           h1_matchedCharges_->Fill(clus.amplitudes().size(), matchedCharges.at(i) / (float)charge1);
0162       }
0163       if (nMatched > 1)
0164         h1_matchingMult_->Fill(clus.amplitudes().size(), nMatched);
0165     }
0166   }
0167 }
0168 
0169 // ------------ method called once each job just before starting event loop  ------------
0170 void SiStripBaselineComparator::beginJob() {}
0171 
0172 // ------------ method called once each job just after ending the event loop  ------------
0173 void SiStripBaselineComparator::endJob() {}
0174 
0175 //define this as a plug-in
0176 DEFINE_FWK_MODULE(SiStripBaselineComparator);