Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-05-10 02:20:58

0001 #include "FWCore/Framework/interface/Frameworkfwd.h"
0002 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0003 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0004 #include "FWCore/Utilities/interface/EDGetToken.h"
0005 #include "DataFormats/Common/interface/ValidHandle.h"
0006 
0007 #include "FWCore/Framework/interface/Event.h"
0008 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
0009 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0010 #include "FWCore/Utilities/interface/InputTag.h"
0011 
0012 #include "FWCore/ServiceRegistry/interface/Service.h"
0013 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0014 
0015 #include "TFile.h"
0016 #include "TH1.h"
0017 #include "TH2.h"
0018 
0019 #include <CLHEP/Units/PhysicalConstants.h>
0020 #include <CLHEP/Units/SystemOfUnits.h>
0021 
0022 #include "FWCore/Framework/interface/MakerMacros.h"
0023 
0024 class VtxTester : public edm::one::EDAnalyzer<> {
0025 public:
0026   //
0027   explicit VtxTester(const edm::ParameterSet&);
0028   virtual ~VtxTester() {}
0029 
0030   void beginJob() override {}
0031   void analyze(const edm::Event&, const edm::EventSetup&) override;
0032   void endJob() override {}
0033 
0034   static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0035 
0036 private:
0037   edm::EDGetTokenT<edm::HepMCProduct> srcHepMCToken_;
0038 
0039   TH1D* fVtxHistTime_;
0040   TH1D* fVtxHistz_;
0041   TH2D* fVtxHistzTime_;
0042   TH1D* fVtxHistx_;
0043   TH1D* fVtxHisty_;
0044   TH2D* fVtxHistxy_;
0045   TH1D* fPhiHist_;
0046   TH1D* fEtaHist_;
0047 };
0048 
0049 VtxTester::VtxTester(const edm::ParameterSet& cfg)
0050     : srcHepMCToken_(consumes<edm::HepMCProduct>(cfg.getParameter<edm::InputTag>("src"))) {
0051   edm::Service<TFileService> fs;
0052 
0053   fVtxHistTime_ = fs->make<TH1D>("VtxHistTime", "#vtx, t [ns]", 60, -0.9, 0.9);
0054   fVtxHistz_ = fs->make<TH1D>("VtxHistz", "#vtx, z [mm]", 400, -200., 200.);
0055   fVtxHistzTime_ = fs->make<TH2D>("VtxHistzTime", "#vtx time [ns] vs z [mm]", 400, -200., 200., 60, -0.9, 0.9);
0056   fVtxHistx_ = fs->make<TH1D>("VtxHistx", "#vtx, x [mm]", 200, -1., 1.);
0057   fVtxHisty_ = fs->make<TH1D>("VtxHisty", "#vtx, y [mm]", 200, -1., 1.);
0058   fVtxHistxy_ = fs->make<TH2D>("VtxHistxy", "#vtx y vs x [mm]", 200, -1., 1., 200, -1., 1.);
0059 
0060   fPhiHist_ = fs->make<TH1D>("PhiHist", "#vtx phi", 70, -3.5, 3.5);
0061   fEtaHist_ = fs->make<TH1D>("EtaHist", "#vtx eta", 120, -6., 6.);
0062 }
0063 
0064 void VtxTester::analyze(const edm::Event& e, const edm::EventSetup&) {
0065   auto theGenEvent = makeValid(e.getHandle(srcHepMCToken_));
0066 
0067   const HepMC::GenEvent* Evt = theGenEvent->GetEvent();
0068 
0069   // take only 1st vertex for now - it's been tested only of PGuns...
0070 
0071   HepMC::GenEvent::vertex_const_iterator Vtx = Evt->vertices_begin();
0072 
0073   fVtxHistTime_->Fill((*Vtx)->position().t() * CLHEP::mm / CLHEP::c_light);
0074   fVtxHistz_->Fill((*Vtx)->position().z() * CLHEP::mm);
0075   fVtxHistzTime_->Fill((*Vtx)->position().z() * CLHEP::mm, (*Vtx)->position().t() * CLHEP::mm / CLHEP::c_light);
0076   fVtxHistx_->Fill((*Vtx)->position().x() * CLHEP::mm);
0077   fVtxHisty_->Fill((*Vtx)->position().y() * CLHEP::mm);
0078   fVtxHistxy_->Fill((*Vtx)->position().x() * CLHEP::mm, (*Vtx)->position().y() * CLHEP::mm);
0079 
0080   for (HepMC::GenEvent::particle_const_iterator Part = Evt->particles_begin(); Part != Evt->particles_end(); Part++) {
0081     HepMC::FourVector Mom = (*Part)->momentum();
0082     double Phi = Mom.phi();
0083     double Eta = -log(tan(Mom.theta() / 2.));
0084 
0085     fPhiHist_->Fill(Phi);
0086     fEtaHist_->Fill(Eta);
0087   }
0088 
0089   return;
0090 }
0091 
0092 void VtxTester::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0093   edm::ParameterSetDescription desc;
0094   desc.add<edm::InputTag>("src", edm::InputTag("generatorSmeared"))
0095       ->setComment("Input generated HepMC event after vtx smearing");
0096   descriptions.add("vtxTester", desc);
0097 }
0098 
0099 DEFINE_FWK_MODULE(VtxTester);