Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:33:14

0001 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0002 #include "DQMServices/Core/interface/DQMStore.h"
0003 #include "DataFormats/CaloRecHit/interface/CaloCluster.h"
0004 #include "DataFormats/CaloRecHit/interface/CaloClusterFwd.h"
0005 #include "DataFormats/CaloTowers/interface/CaloTowerCollection.h"
0006 #include "DataFormats/CaloTowers/interface/CaloTowerDetId.h"
0007 #include "DataFormats/Common/interface/Handle.h"
0008 #include "DataFormats/DetId/interface/DetId.h"
0009 #include "DataFormats/EcalDetId/interface/EcalSubdetector.h"
0010 #include "DataFormats/HcalDetId/interface/HcalDetId.h"
0011 #include "DataFormats/HcalDetId/interface/HcalSubdetector.h"
0012 #include "DataFormats/Math/interface/Vector3D.h"
0013 #include "DataFormats/Math/interface/deltaR.h"
0014 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0015 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
0016 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementCluster.h"
0017 #include "DataFormats/ParticleFlowReco/interface/PFBlockElementTrack.h"
0018 #include "DataFormats/ParticleFlowReco/interface/PFCluster.h"
0019 #include "DataFormats/ParticleFlowReco/interface/PFClusterFwd.h"
0020 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
0021 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0022 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFraction.h"
0023 #include "FWCore/Framework/interface/Event.h"
0024 #include "FWCore/Framework/interface/EventSetup.h"
0025 #include "FWCore/Framework/interface/Frameworkfwd.h"
0026 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0027 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0028 #include "FWCore/PluginManager/interface/ModuleDef.h"
0029 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0030 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0031 
0032 #include <algorithm>
0033 #include <cmath>
0034 #include <iostream>
0035 #include <string>
0036 #include <utility>
0037 #include <vector>
0038 
0039 #ifdef PFLOW_DEBUG
0040 #define LOGVERB(x) edm::LogVerbatim(x)
0041 #else
0042 #define LOGVERB(x) LogTrace(x)
0043 #endif
0044 
0045 class PFCaloGPUComparisonTask : public DQMEDAnalyzer {
0046 public:
0047   PFCaloGPUComparisonTask(edm::ParameterSet const& conf);
0048   ~PFCaloGPUComparisonTask() override = default;
0049   void analyze(edm::Event const& e, edm::EventSetup const& c) override;
0050   void bookHistograms(DQMStore::IBooker&, edm::Run const&, edm::EventSetup const&) override;
0051 
0052 private:
0053   edm::EDGetTokenT<reco::PFClusterCollection> pfClusterTok_ref_;
0054   edm::EDGetTokenT<reco::PFClusterCollection> pfClusterTok_target_;
0055 
0056   MonitorElement* pfCluster_Multiplicity_GPUvsCPU_;
0057   MonitorElement* pfCluster_Energy_GPUvsCPU_;
0058   MonitorElement* pfCluster_RecHitMultiplicity_GPUvsCPU_;
0059   MonitorElement* pfCluster_Layer_GPUvsCPU_;
0060   MonitorElement* pfCluster_Depth_GPUvsCPU_;
0061   MonitorElement* pfCluster_Eta_GPUvsCPU_;
0062   MonitorElement* pfCluster_Phi_GPUvsCPU_;
0063   MonitorElement* pfCluster_DuplicateMatches_GPUvsCPU_;
0064 
0065   std::string pfCaloGPUCompDir_;
0066 };
0067 
0068 PFCaloGPUComparisonTask::PFCaloGPUComparisonTask(const edm::ParameterSet& conf)
0069     : pfClusterTok_ref_{consumes<reco::PFClusterCollection>(
0070           conf.getUntrackedParameter<edm::InputTag>("pfClusterToken_ref"))},
0071       pfClusterTok_target_{
0072           consumes<reco::PFClusterCollection>(conf.getUntrackedParameter<edm::InputTag>("pfClusterToken_target"))},
0073       pfCaloGPUCompDir_{conf.getUntrackedParameter<std::string>("pfCaloGPUCompDir")} {}
0074 
0075 void PFCaloGPUComparisonTask::bookHistograms(DQMStore::IBooker& ibooker,
0076                                              edm::Run const& irun,
0077                                              edm::EventSetup const& isetup) {
0078   const char* histo;
0079 
0080   ibooker.setCurrentFolder("ParticleFlow/" + pfCaloGPUCompDir_);
0081 
0082   histo = "pfCluster_Multiplicity_GPUvsCPU";
0083   pfCluster_Multiplicity_GPUvsCPU_ = ibooker.book2D(histo, histo, 100, 0, 2000, 100, 0, 2000);
0084 
0085   histo = "pfCluster_Energy_GPUvsCPU";
0086   pfCluster_Energy_GPUvsCPU_ = ibooker.book2D(histo, histo, 100, 0, 500, 100, 0, 500);
0087 
0088   histo = "pfCluster_RecHitMultiplicity_GPUvsCPU";
0089   pfCluster_RecHitMultiplicity_GPUvsCPU_ = ibooker.book2D(histo, histo, 100, 0, 100, 100, 0, 100);
0090 
0091   histo = "pfCluster_Layer_GPUvsCPU";
0092   pfCluster_Layer_GPUvsCPU_ = ibooker.book2D(histo, histo, 100, 0, 100, 100, 0, 100);
0093 
0094   histo = "pfCluster_Depth_GPUvsCPU";
0095   pfCluster_Depth_GPUvsCPU_ = ibooker.book2D(histo, histo, 100, 0, 100, 100, 0, 100);
0096 
0097   histo = "pfCluster_Eta_GPUvsCPU";
0098   pfCluster_Eta_GPUvsCPU_ = ibooker.book2D(histo, histo, 100, 0, 100, 100, 0, 100);
0099 
0100   histo = "pfCluster_Phi_GPUvsCPU";
0101   pfCluster_Phi_GPUvsCPU_ = ibooker.book2D(histo, histo, 100, 0, 100, 100, 0, 100);
0102 
0103   histo = "pfCluster_DuplicateMatches_GPUvsCPU";
0104   pfCluster_DuplicateMatches_GPUvsCPU_ = ibooker.book1D(histo, histo, 100, 0., 1000);
0105 }
0106 void PFCaloGPUComparisonTask::analyze(edm::Event const& event, edm::EventSetup const& c) {
0107   edm::Handle<reco::PFClusterCollection> pfClusters_ref;
0108   event.getByToken(pfClusterTok_ref_, pfClusters_ref);
0109 
0110   edm::Handle<reco::PFClusterCollection> pfClusters_target;
0111   event.getByToken(pfClusterTok_target_, pfClusters_target);
0112 
0113   //
0114   // Compare per-event PF cluster multiplicity
0115 
0116   if (pfClusters_ref->size() != pfClusters_target->size())
0117     LOGVERB("PFCaloGPUComparisonTask") << " PFCluster multiplicity " << pfClusters_ref->size() << " "
0118                                        << pfClusters_target->size();
0119   pfCluster_Multiplicity_GPUvsCPU_->Fill((float)pfClusters_ref->size(), (float)pfClusters_target->size());
0120 
0121   //
0122   // Find matching PF cluster pairs
0123   std::vector<int> matched_idx;
0124   matched_idx.reserve(pfClusters_ref->size());
0125   for (unsigned i = 0; i < pfClusters_ref->size(); ++i) {
0126     bool matched = false;
0127     for (unsigned j = 0; j < pfClusters_target->size(); ++j) {
0128       if (pfClusters_ref->at(i).seed() == pfClusters_target->at(j).seed()) {
0129         if (!matched) {
0130           matched = true;
0131           matched_idx.push_back((int)j);
0132         } else {
0133           edm::LogWarning("PFCaloGPUComparisonTask") << "Found duplicate match";
0134           pfCluster_DuplicateMatches_GPUvsCPU_->Fill((int)j);
0135         }
0136       }
0137     }
0138     if (!matched)
0139       matched_idx.push_back(-1);  // if you don't find a match, put a dummy number
0140   }
0141 
0142   //
0143   // Plot matching PF cluster variables
0144   for (unsigned i = 0; i < pfClusters_ref->size(); ++i) {
0145     if (matched_idx[i] >= 0) {
0146       unsigned int j = matched_idx[i];
0147       int ref_energy_bin = pfCluster_Energy_GPUvsCPU_->getTH2F()->GetXaxis()->FindBin(pfClusters_ref->at(i).energy());
0148       int target_energy_bin =
0149           pfCluster_Energy_GPUvsCPU_->getTH2F()->GetXaxis()->FindBin(pfClusters_target->at(j).energy());
0150       if (ref_energy_bin != target_energy_bin)
0151         edm::LogPrint("PFCaloGPUComparisonTask")
0152             << "Off-diagonal energy bin entries: " << pfClusters_ref->at(i).energy() << " "
0153             << pfClusters_ref->at(i).eta() << " " << pfClusters_ref->at(i).phi() << " "
0154             << pfClusters_target->at(j).energy() << " " << pfClusters_target->at(j).eta() << " "
0155             << pfClusters_target->at(j).phi() << std::endl;
0156       pfCluster_Energy_GPUvsCPU_->Fill(pfClusters_ref->at(i).energy(), pfClusters_target->at(j).energy());
0157       pfCluster_Layer_GPUvsCPU_->Fill(pfClusters_ref->at(i).layer(), pfClusters_target->at(j).layer());
0158       pfCluster_Eta_GPUvsCPU_->Fill(pfClusters_ref->at(i).eta(), pfClusters_target->at(j).eta());
0159       pfCluster_Phi_GPUvsCPU_->Fill(pfClusters_ref->at(i).phi(), pfClusters_target->at(j).phi());
0160       pfCluster_Depth_GPUvsCPU_->Fill(pfClusters_ref->at(i).depth(), pfClusters_target->at(j).depth());
0161       pfCluster_RecHitMultiplicity_GPUvsCPU_->Fill((float)pfClusters_ref->at(i).recHitFractions().size(),
0162                                                    (float)pfClusters_target->at(j).recHitFractions().size());
0163     }
0164   }
0165 }
0166 
0167 #include "FWCore/Framework/interface/MakerMacros.h"
0168 DEFINE_FWK_MODULE(PFCaloGPUComparisonTask);