File indexing completed on 2024-04-06 12:08:02
0001 #include "DQM/Physics/src/CentralityDQM.h"
0002
0003 #include <memory>
0004
0005
0006 #include "DQMServices/Core/interface/DQMStore.h"
0007 #include "DQMServices/Core/interface/DQMEDAnalyzer.h"
0008
0009
0010 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/EventSetup.h"
0013 #include "FWCore/Framework/interface/LuminosityBlock.h"
0014 #include "FWCore/ServiceRegistry/interface/Service.h"
0015 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0016
0017
0018 #include "DataFormats/VertexReco/interface/Vertex.h"
0019
0020 using namespace edm;
0021 using namespace reco;
0022 using namespace hi;
0023
0024
0025
0026
0027 CentralityDQM::CentralityDQM(const edm::ParameterSet& ps) {
0028 edm::LogInfo("CentralityDQM") << " Starting CentralityDQM "
0029 << "\n";
0030
0031 centralityTag_ = ps.getParameter<InputTag>("centralitycollection");
0032 centralityToken = consumes<reco::Centrality>(centralityTag_);
0033
0034 centralityBinTag_ = (ps.getParameter<edm::InputTag>("centralitybincollection"));
0035 centralityBinToken = consumes<int>(centralityBinTag_);
0036
0037 vertexTag_ = ps.getParameter<InputTag>("vertexcollection");
0038 vertexToken = consumes<std::vector<reco::Vertex> >(vertexTag_);
0039
0040 eventplaneTag_ = ps.getParameter<InputTag>("eventplanecollection");
0041 eventplaneToken = consumes<reco::EvtPlaneCollection>(eventplaneTag_);
0042
0043
0044 }
0045
0046
0047
0048
0049 CentralityDQM::~CentralityDQM() {
0050 edm::LogInfo("CentralityDQM") << " Deleting CentralityDQM "
0051 << "\n";
0052 }
0053
0054
0055
0056
0057 void CentralityDQM::bookHistograms(DQMStore::IBooker& bei, edm::Run const&, edm::EventSetup const&) {
0058
0059
0060 bei.setCurrentFolder("Physics/Centrality/");
0061
0062 h_hiNpix = bei.book1D("h_hiNpix", "h_hiNpix", 750, 0, 75000);
0063 h_hiNpixelTracks = bei.book1D("h_hiNpixelTracks", "hiNpixelTracks", 500, 0, 5000);
0064 h_hiNtracks = bei.book1D("h_hiNtracks", "h_hiNtracks", 500, 0, 5000);
0065 h_hiNtracksPtCut = bei.book1D("h_hiNtracksPtCut", "h_hiNtracksPtCut", 500, 0, 5000);
0066 h_hiNtracksEtaCut = bei.book1D("h_hiNtracksEtaCut", "h_hiNtracksEtaCut", 500, 0, 5000);
0067 h_hiNtracksEtaPtCut = bei.book1D("h_hiNtracksEtaPtCut", "h_hiNtracksEtaPtCut", 500, 0, 5000);
0068
0069 h_hiHF = bei.book1D("h_hiHF", "h_hiHF", 900, 0, 9000);
0070 h_hiHFplus = bei.book1D("h_hiHFplus", "h_hiHFplus", 900, 0, 9000);
0071 h_hiHFminus = bei.book1D("h_hiHFminus", "h_hiHFminus", 900, 0, 9000);
0072 h_hiHFplusEta4 = bei.book1D("h_hiHFplusEta4", "h_hiHFplusEta4", 900, 0, 9000);
0073 h_hiHFminusEta4 = bei.book1D("h_hiHFminusEta4", "h_hiHFminusEta4", 900, 0, 9000);
0074
0075 h_hiHFhit = bei.book1D("h_hiHFhit", "h_hiHFhit", 3000, 0, 300000);
0076 h_hiHFhitPlus = bei.book1D("h_hiHFhitPlus", "h_hiHFhitPlus", 2000, 0, 200000);
0077 h_hiHFhitMinus = bei.book1D("h_hiHFhitMinus", "h_hiHFhitMinus", 2000, 0, 200000);
0078
0079 h_hiEB = bei.book1D("h_hiEB", "h_hiEB", 600, 0, 6000);
0080 h_hiET = bei.book1D("h_hiET", "h_hiET", 600, 0, 6000);
0081 h_hiEE = bei.book1D("h_hiEE", "h_hiEE", 600, 0, 6000);
0082 h_hiEEplus = bei.book1D("h_hiEEplus", "h_hiEEplus", 600, 0, 6000);
0083 h_hiEEminus = bei.book1D("h_hiEEminus", "h_hiEEminus", 600, 0, 6000);
0084 h_hiZDC = bei.book1D("h_hiZDC", "h_hiZDC", 600, 0, 6000);
0085 h_hiZDCplus = bei.book1D("h_hiZDCplus", "h_hiZDCplus", 600, 0, 6000);
0086 h_hiZDCminus = bei.book1D("h_hiZDCminus", "h_hiZDCminus", 600, 0, 6000);
0087 h_hiPF = bei.book1D("h_hiPF", "h_hiPF", 900, 0, 9000);
0088 h_hiPFplus = bei.book1D("h_hiPFplus", "h_hiPFplus", 900, 0, 9000);
0089 h_hiPFminus = bei.book1D("h_hiPFminus", "h_hiPFminus", 900, 0, 9000);
0090
0091 h_vertex_x = bei.book1D("h_vertex_x", "h_vertex_x", 400, -4, 4);
0092 h_vertex_y = bei.book1D("h_vertex_y", "h_vertex_y", 400, -4, 4);
0093 h_vertex_z = bei.book1D("h_vertex_z", "h_vertex_z", 400, -40, 40);
0094
0095 h_cent_bin = bei.book1D("h_cent_bin", "h_cent_bin", 200, 0, 200);
0096
0097 Double_t psirange = 4;
0098 bei.setCurrentFolder("Physics/Centrality/EventPlane/");
0099 h_ep_HFm2 = bei.book1D("h_ep_HFm2", "h_ep_HFm2", 800, -psirange, psirange);
0100 h_ep_HFp2 = bei.book1D("h_ep_HFp2", "h_ep_HFp2", 800, -psirange, psirange);
0101 h_ep_trackmid2 = bei.book1D("h_ep_trackmid2", "h_ep_trackmid2", 800, -psirange, psirange);
0102 h_ep_trackm2 = bei.book1D("h_ep_trackm2", "h_ep_trackm2", 800, -psirange, psirange);
0103 h_ep_trackp2 = bei.book1D("h_ep_trackp2", "h_ep_trackp2", 800, -psirange, psirange);
0104 h_ep_HFm3 = bei.book1D("h_ep_HFm3", "h_ep_HFm3", 800, -psirange, psirange);
0105 h_ep_HFp3 = bei.book1D("h_ep_HFp3", "h_ep_HFp3", 800, -psirange, psirange);
0106 h_ep_trackmid3 = bei.book1D("h_ep_trackmid3", "h_ep_trackmid3", 800, -psirange, psirange);
0107 }
0108
0109
0110
0111
0112 void CentralityDQM::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0113 using namespace edm;
0114 edm::Handle<reco::Centrality> cent;
0115 iEvent.getByToken(centralityToken, cent);
0116
0117 edm::Handle<int> cbin;
0118 iEvent.getByToken(centralityBinToken, cbin);
0119
0120
0121
0122 edm::Handle<reco::EvtPlaneCollection> ep;
0123 iEvent.getByToken(eventplaneToken, ep);
0124
0125 if (cent.isValid()) {
0126 int hibin = -999;
0127 if (cbin.isValid())
0128 hibin = *cbin;
0129
0130
0131
0132 h_cent_bin->Fill(hibin);
0133 h_hiNpix->Fill(cent->multiplicityPixel());
0134 h_hiNpixelTracks->Fill(cent->NpixelTracks());
0135 h_hiNtracks->Fill(cent->Ntracks());
0136
0137 h_hiNtracksPtCut->Fill(cent->NtracksPtCut());
0138 h_hiNtracksEtaCut->Fill(cent->NtracksEtaCut());
0139 h_hiNtracksEtaPtCut->Fill(cent->NtracksEtaPtCut());
0140
0141 h_hiHF->Fill(cent->EtHFtowerSum());
0142 h_hiHFplus->Fill(cent->EtHFtowerSumPlus());
0143 h_hiHFminus->Fill(cent->EtHFtowerSumMinus());
0144 h_hiHFplusEta4->Fill(cent->EtHFtruncatedPlus());
0145 h_hiHFminusEta4->Fill(cent->EtHFtruncatedMinus());
0146
0147 h_hiHFhit->Fill(cent->EtHFhitSum());
0148 h_hiHFhitPlus->Fill(cent->EtHFhitSumPlus());
0149 h_hiHFhitMinus->Fill(cent->EtHFhitSumMinus());
0150
0151 h_hiZDC->Fill(cent->zdcSum());
0152 h_hiZDCplus->Fill(cent->zdcSumPlus());
0153 h_hiZDCminus->Fill(cent->zdcSumMinus());
0154
0155 h_hiEEplus->Fill(cent->EtEESumPlus());
0156 h_hiEEminus->Fill(cent->EtEESumMinus());
0157 h_hiEE->Fill(cent->EtEESum());
0158 h_hiEB->Fill(cent->EtEBSum());
0159 h_hiET->Fill(cent->EtMidRapiditySum());
0160
0161 h_hiPF->Fill(cent->EtPFhfSum());
0162 h_hiPFplus->Fill(cent->EtPFhfSumPlus());
0163 h_hiPFminus->Fill(cent->EtPFhfSumMinus());
0164
0165 edm::Handle<std::vector<reco::Vertex> > vertex;
0166 iEvent.getByToken(vertexToken, vertex);
0167 h_vertex_x->Fill(vertex->begin()->x());
0168 h_vertex_y->Fill(vertex->begin()->y());
0169 h_vertex_z->Fill(vertex->begin()->z());
0170 }
0171
0172 if (ep.isValid()) {
0173 EvtPlaneCollection::const_iterator rp = ep->begin();
0174
0175 h_ep_HFm2->Fill(rp[HFm2].angle(0));
0176 h_ep_HFp2->Fill(rp[HFp2].angle(0));
0177 h_ep_trackmid2->Fill(rp[trackmid2].angle(0));
0178 h_ep_trackm2->Fill(rp[trackm2].angle(0));
0179 h_ep_trackp2->Fill(rp[trackp2].angle(0));
0180
0181 h_ep_HFm3->Fill(rp[HFm3].angle(0));
0182 h_ep_HFp3->Fill(rp[HFp3].angle(0));
0183 h_ep_trackmid3->Fill(rp[trackmid3].angle(0));
0184 }
0185 }