Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2023-03-17 11:28:41

0001 // author: Mike Schmitt, University of Florida
0002 // first version 11/7/2007
0003 
0004 #include "DQMServices/Core/interface/DQMStore.h"
0005 #include "DataFormats/Common/interface/Handle.h"
0006 #include "DataFormats/Common/interface/RefToBase.h"
0007 #include "DataFormats/HepMCCandidate/interface/GenParticle.h"
0008 #include "DataFormats/METReco/interface/CaloMET.h"
0009 #include "DataFormats/METReco/interface/CaloMETCollection.h"
0010 #include "DataFormats/METReco/interface/GenMET.h"
0011 #include "DataFormats/METReco/interface/GenMETCollection.h"
0012 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidate.h"
0013 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateFwd.h"
0014 #include "DataFormats/ParticleFlowReco/interface/PFBlock.h"
0015 #include "DataFormats/ParticleFlowReco/interface/PFBlockElement.h"
0016 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
0017 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0018 #include "FWCore/Framework/interface/ESHandle.h"
0019 #include "FWCore/Framework/interface/Event.h"
0020 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0021 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0022 #include "FWCore/ServiceRegistry/interface/Service.h"
0023 
0024 #include <algorithm>
0025 #include <cmath>
0026 #include <fstream>
0027 #include <iostream>
0028 #include <map>
0029 #include <memory>
0030 #include <ostream>
0031 #include <string>
0032 #include <vector>
0033 
0034 class PFTester : public edm::one::EDAnalyzer<> {
0035 public:
0036   typedef dqm::legacy::DQMStore DQMStore;
0037   typedef dqm::legacy::MonitorElement MonitorElement;
0038 
0039   explicit PFTester(const edm::ParameterSet &);
0040   ~PFTester() override;
0041 
0042   void analyze(const edm::Event &, const edm::EventSetup &) override;
0043   void beginJob() override;
0044   void endJob() override;
0045 
0046 private:
0047   // DAQ Tools
0048   DQMStore *dbe_;
0049   std::map<std::string, MonitorElement *> me;
0050 
0051   // Inputs from Configuration File
0052   std::string outputFile_;
0053   edm::EDGetTokenT<reco::PFCandidateCollection> inputPFlowLabel_tok_;
0054 };
0055 
0056 #include "FWCore/Framework/interface/MakerMacros.h"
0057 DEFINE_FWK_MODULE(PFTester);
0058 
0059 using namespace edm;
0060 using namespace std;
0061 using namespace reco;
0062 
0063 PFTester::PFTester(const edm::ParameterSet &iConfig) {
0064   inputPFlowLabel_tok_ = consumes<reco::PFCandidateCollection>(iConfig.getParameter<std::string>("InputPFlowLabel"));
0065   outputFile_ = iConfig.getUntrackedParameter<std::string>("OutputFile");
0066 
0067   if (!outputFile_.empty())
0068     edm::LogInfo("OutputInfo") << " ParticleFLow Task histograms will be saved to '" << outputFile_.c_str() << "'";
0069   else
0070     edm::LogInfo("OutputInfo") << " ParticleFlow Task histograms will NOT be saved";
0071 }
0072 
0073 PFTester::~PFTester() {}
0074 
0075 void PFTester::beginJob() {
0076   // get ahold of back-end interface
0077   dbe_ = edm::Service<DQMStore>().operator->();
0078 
0079   if (dbe_) {
0080     dbe_->setCurrentFolder("PFTask/PFCandidates");
0081 
0082     me["CandidateEt"] = dbe_->book1D("CandidateEt", "CandidateEt", 1000, 0, 1000);
0083     me["CandidateEta"] = dbe_->book1D("CandidateEta", "CandidateEta", 200, -5, 5);
0084     me["CandidatePhi"] = dbe_->book1D("CandidatePhi", "CandidatePhi", 200, -M_PI, M_PI);
0085     me["CandidateCharge"] = dbe_->book1D("CandidateCharge", "CandidateCharge", 5, -2, 2);
0086     me["PFCandidateType"] = dbe_->book1D("PFCandidateType", "PFCandidateType", 10, 0, 10);
0087 
0088     dbe_->setCurrentFolder("PFTask/PFBlocks");
0089 
0090     me["NumElements"] = dbe_->book1D("NumElements", "NumElements", 25, 0, 25);
0091     me["NumTrackElements"] = dbe_->book1D("NumTrackElements", "NumTrackElements", 5, 0, 5);
0092     me["NumPS1Elements"] = dbe_->book1D("NumPS1Elements", "NumPS1Elements", 5, 0, 5);
0093     me["NumPS2Elements"] = dbe_->book1D("NumPS2Elements", "NumPS2Elements", 5, 0, 5);
0094     me["NumECALElements"] = dbe_->book1D("NumECALElements", "NumECALElements", 5, 0, 5);
0095     me["NumHCALElements"] = dbe_->book1D("NumHCALElements", "NumHCALElements", 5, 0, 5);
0096     me["NumMuonElements"] = dbe_->book1D("NumMuonElements", "NumMuonElements", 5, 0, 5);
0097 
0098     dbe_->setCurrentFolder("PFTask/PFTracks");
0099 
0100     me["TrackCharge"] = dbe_->book1D("TrackCharge", "TrackCharge", 5, -2, 2);
0101     me["TrackNumPoints"] = dbe_->book1D("TrackNumPoints", "TrackNumPoints", 100, 0, 100);
0102     me["TrackNumMeasurements"] = dbe_->book1D("TrackNumMeasurements", "TrackNumMeasurements", 100, 0, 100);
0103     me["TrackImpactParameter"] = dbe_->book1D("TrackImpactParameter", "TrackImpactParameter", 1000, 0, 1);
0104 
0105     dbe_->setCurrentFolder("PFTask/PFClusters");
0106   }
0107 }
0108 
0109 void PFTester::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0110   // Data to Retrieve from the Event
0111   const PFCandidateCollection *pflow_candidates;
0112 
0113   // ==========================================================
0114   // Retrieve!
0115   // ==========================================================
0116 
0117   {
0118     // Get Particle Flow Candidates
0119     Handle<PFCandidateCollection> pflow_hnd;
0120     iEvent.getByToken(inputPFlowLabel_tok_, pflow_hnd);
0121     pflow_candidates = pflow_hnd.product();
0122   }
0123 
0124   if (!pflow_candidates) {
0125     edm::LogInfo("OutputInfo") << " failed to retrieve data required by ParticleFlow Task";
0126     edm::LogInfo("OutputInfo") << " ParticleFlow Task cannot continue...!";
0127     return;
0128   }
0129 
0130   // ==========================================================
0131   // Analyze!
0132   // ==========================================================
0133 
0134   // Loop Over Particle Flow Candidates
0135   PFCandidateCollection::const_iterator pf;
0136   for (pf = pflow_candidates->begin(); pf != pflow_candidates->end(); pf++) {
0137     const PFCandidate *particle = &(*pf);
0138 
0139     // Fill Histograms for Candidate Methods
0140     me["CandidateEt"]->Fill(particle->et());
0141     me["CandidateEta"]->Fill(particle->eta());
0142     me["CandidatePhi"]->Fill(particle->phi());
0143     me["CandidateCharge"]->Fill(particle->charge());
0144     me["CandidatePdgId"]->Fill(particle->pdgId());
0145 
0146     // Fill Histograms for PFCandidate Specific Methods
0147     me["PFCandidateType"]->Fill(particle->particleId());
0148     // particle->elementsInBlocks();
0149 
0150     // Get the PFBlock and Elements
0151     // JW: Returns vector of blocks now ,TO BE FIXED ----
0152     /*PFBlock block = *(particle->block());
0153     OwnVector<PFBlockElement> elements = block.elements();
0154     int numElements = elements.size();
0155     int numTrackElements = 0;
0156     int numPS1Elements = 0;
0157     int numPS2Elements = 0;
0158     int numECALElements = 0;
0159     int numHCALElements = 0;
0160     int numMuonElements = 0;
0161 
0162     // Loop over Elements in Block
0163     OwnVector<PFBlockElement>::const_iterator element;
0164     for (element = elements.begin(); element != elements.end(); element++) {
0165 
0166       int element_type = element->type();
0167       // Element is a Tracker Track
0168       if (element_type == PFBlockElement::TRACK) {
0169 
0170         // Get General Information about the Track
0171         PFRecTrack track = *(element->trackRefPF());
0172         me["TrackCharge"]->Fill(track.charge());
0173         me["TrackNumPoints"]->Fill(track.nTrajectoryPoints());
0174         me["TrackNumMeasurements"]->Fill(track.nTrajectoryMeasurements());
0175 
0176         // Loop Over Points in the Track
0177         vector<PFTrajectoryPoint> points = track.trajectoryPoints();
0178         vector<PFTrajectoryPoint>::iterator point;
0179         for (point = points.begin(); point != points.end(); point++) {
0180           int point_layer = point->layer();
0181           double x = point->positionXYZ().x();
0182           double y = point->positionXYZ().y();
0183           double z = point->positionXYZ().z();
0184           //switch (point_layer) {
0185           //case PFTrajectoryPoint::ClosestApproach:
0186           // Fill the Track's D0
0187           if (point_layer == PFTrajectoryPoint::ClosestApproach) {
0188             me["TrackImpactParameter"]->Fill(sqrt(x*x + y*y + z*z));
0189           }
0190         }
0191         numTrackElements++;
0192       }
0193 
0194       // Element is an ECAL Cluster
0195       else if (element_type == PFBlockElement::ECAL) {
0196         numECALElements++;
0197       }
0198       // Element is a HCAL Cluster
0199       else if (element_type == PFBlockElement::HCAL) {
0200         numHCALElements++;
0201       }
0202       // Element is a Muon Track
0203       else if (element_type == PFBlockElement::MUON) {
0204         numMuonElements++;
0205       }
0206       // Fill the Respective Elements Sizes
0207       me["NumElements"]->Fill(numElements);
0208       me["NumTrackElements"]->Fill(numTrackElements);
0209       me["NumPS1Elements"]->Fill(numPS1Elements);
0210       me["NumPS2Elements"]->Fill(numPS2Elements);
0211       me["NumECALElements"]->Fill(numECALElements);
0212       me["NumHCALElements"]->Fill(numHCALElements);
0213       me["NumMuonElements"]->Fill(numMuonElements);
0214     } ----------------------------------------------  */
0215   }
0216 }
0217 
0218 void PFTester::endJob() {
0219   // Store the DAQ Histograms
0220   if (!outputFile_.empty() && dbe_)
0221     dbe_->save(outputFile_);
0222 }