File indexing completed on 2024-04-06 12:31:58
0001
0002
0003
0004
0005
0006 #include "FWCore/Framework/interface/Frameworkfwd.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/ESHandle.h"
0009 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0010
0011 #include "FWCore/Framework/interface/Event.h"
0012 #include "FWCore/Framework/interface/MakerMacros.h"
0013
0014 #include "CondTools/RunInfo/interface/LHCInfoCombined.h"
0015
0016 #include "TFile.h"
0017 #include "TH1D.h"
0018 #include "TH2D.h"
0019
0020
0021
0022 class CTPPSLHCInfoPlotter : public edm::one::EDAnalyzer<> {
0023 public:
0024 explicit CTPPSLHCInfoPlotter(const edm::ParameterSet &);
0025
0026 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0027
0028 private:
0029 void analyze(const edm::Event &, const edm::EventSetup &) override;
0030 void endJob() override;
0031
0032 const edm::ESGetToken<LHCInfo, LHCInfoRcd> lhcInfoToken_;
0033 const edm::ESGetToken<LHCInfoPerLS, LHCInfoPerLSRcd> lhcInfoPerLSToken_;
0034 const edm::ESGetToken<LHCInfoPerFill, LHCInfoPerFillRcd> lhcInfoPerFillToken_;
0035 const bool useNewLHCInfo_;
0036
0037 std::string outputFile_;
0038
0039 TH1D *h_beamEnergy_;
0040 TH1D *h_xangle_;
0041 TH1D *h_betaStar_;
0042 TH2D *h2_betaStar_vs_xangle_;
0043
0044 TH1D *h_fill_;
0045 TH1D *h_run_;
0046 };
0047
0048
0049
0050 using namespace std;
0051 using namespace edm;
0052
0053
0054
0055 CTPPSLHCInfoPlotter::CTPPSLHCInfoPlotter(const edm::ParameterSet &iConfig)
0056 : lhcInfoToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoLabel")))),
0057 lhcInfoPerLSToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoPerLSLabel")))),
0058 lhcInfoPerFillToken_(esConsumes(ESInputTag("", iConfig.getParameter<std::string>("lhcInfoPerFillLabel")))),
0059 useNewLHCInfo_(iConfig.getParameter<bool>("useNewLHCInfo")),
0060 outputFile_(iConfig.getParameter<string>("outputFile")),
0061
0062 h_beamEnergy_(new TH1D("h_beamEnergy", ";beam energy (GeV)", 81, -50., 8050.)),
0063 h_xangle_(new TH1D("h_xangle", ";(half) crossing angle (#murad)", 201, -0.5, 200.5)),
0064 h_betaStar_(new TH1D("h_betaStar", ";#beta^{*} (m)", 101, -0.005, 1.005)),
0065 h2_betaStar_vs_xangle_(new TH2D("h2_betaStar_vs_xangle",
0066 ";(half) crossing angle (#murad);#beta^{*} (m)",
0067 201,
0068 -0.5,
0069 200.5,
0070 101,
0071 -0.005,
0072 1.005)),
0073
0074 h_fill_(new TH1D("h_fill", ";fill", 6001, 3999.5, 10000.5)),
0075 h_run_(new TH1D("h_run", ";run", 6000, 270E3, 430E3)) {}
0076
0077
0078
0079 void CTPPSLHCInfoPlotter::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0080 edm::ParameterSetDescription desc;
0081
0082 desc.add<std::string>("lhcInfoLabel", "")->setComment("label of the LHCInfo record");
0083 desc.add<std::string>("lhcInfoPerLSLabel", "")->setComment("label of the LHCInfoPerLS record");
0084 desc.add<std::string>("lhcInfoPerFillLabel", "")->setComment("label of the LHCInfoPerFill record");
0085 desc.add<bool>("useNewLHCInfo", false)->setComment("flag whether to use new LHCInfoPer* records or old LHCInfo");
0086
0087 desc.add<std::string>("outputFile", "")->setComment("output file");
0088
0089 descriptions.add("ctppsLHCInfoPlotterDefault", desc);
0090 }
0091
0092
0093
0094 void CTPPSLHCInfoPlotter::analyze(const edm::Event &iEvent, const edm::EventSetup &iSetup) {
0095 LHCInfoCombined lhcInfoCombined(iSetup, lhcInfoPerLSToken_, lhcInfoPerFillToken_, lhcInfoToken_, useNewLHCInfo_);
0096
0097 h_beamEnergy_->Fill(lhcInfoCombined.energy);
0098 h_xangle_->Fill(lhcInfoCombined.crossingAngle());
0099 h_betaStar_->Fill(lhcInfoCombined.betaStarX);
0100
0101 h2_betaStar_vs_xangle_->Fill(lhcInfoCombined.crossingAngle(),
0102 lhcInfoCombined.betaStarX);
0103
0104
0105 h_fill_->Fill(lhcInfoCombined.fillNumber);
0106 h_run_->Fill(iEvent.id().run());
0107 }
0108
0109
0110
0111 void CTPPSLHCInfoPlotter::endJob() {
0112 auto f_out = std::make_unique<TFile>(outputFile_.c_str(), "recreate");
0113
0114 h_beamEnergy_->Write();
0115 h_xangle_->Write();
0116 h_betaStar_->Write();
0117 h2_betaStar_vs_xangle_->Write();
0118
0119 h_fill_->Write();
0120 h_run_->Write();
0121 }
0122
0123
0124
0125 DEFINE_FWK_MODULE(CTPPSLHCInfoPlotter);