File indexing completed on 2024-07-05 03:36:31
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019 #include <fstream>
0020 #include <memory>
0021 #include <sstream>
0022
0023
0024 #include "CondFormats/BeamSpotObjects/interface/SimBeamSpotHLLHCObjects.h"
0025 #include "CondFormats/DataRecord/interface/SimBeamSpotHLLHCObjectsRcd.h"
0026 #include "FWCore/Framework/interface/ESHandle.h"
0027 #include "FWCore/Framework/interface/ESWatcher.h"
0028 #include "FWCore/Framework/interface/Event.h"
0029 #include "FWCore/Framework/interface/EventSetup.h"
0030 #include "FWCore/Framework/interface/Frameworkfwd.h"
0031 #include "FWCore/Framework/interface/MakerMacros.h"
0032 #include "FWCore/Framework/interface/one/EDAnalyzer.h"
0033 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0034
0035
0036 #include "CommonTools/UtilAlgos/interface/TFileService.h"
0037 #include "FWCore/ServiceRegistry/interface/Service.h"
0038 #include <TTree.h>
0039
0040
0041
0042
0043
0044 class BeamProfileHLLHC2DBReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
0045 public:
0046 explicit BeamProfileHLLHC2DBReader(const edm::ParameterSet&);
0047 ~BeamProfileHLLHC2DBReader() override = default;
0048
0049 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0050
0051 private:
0052 void beginJob() override;
0053 void analyze(const edm::Event&, const edm::EventSetup&) override;
0054
0055 struct TheBSfromDB {
0056 int run;
0057 int ls;
0058 double fMeanX, fMeanY, fMeanZ;
0059 double fEProton, fCrabFrequency, fRF800;
0060 double fCrossingAngle, fCrabbingAngleCrossing, fCrabbingAngleSeparation;
0061 double fBetaCrossingPlane, fBetaSeparationPlane;
0062 double fHorizontalEmittance, fVerticalEmittance;
0063 double fBunchLength, fTimeOffset;
0064 void init();
0065 } theBSfromDB_;
0066
0067 const edm::ESGetToken<SimBeamSpotHLLHCObjects, SimBeamSpotHLLHCObjectsRcd> beamSpotToken_;
0068 edm::Service<TFileService> tFileService;
0069 TTree* bstree_;
0070
0071
0072 edm::ESWatcher<SimBeamSpotHLLHCObjectsRcd> watcher_;
0073 std::unique_ptr<std::ofstream> output_;
0074 };
0075
0076
0077 BeamProfileHLLHC2DBReader::BeamProfileHLLHC2DBReader(const edm::ParameterSet& iConfig)
0078 : beamSpotToken_(esConsumes()), bstree_(nullptr) {
0079
0080 usesResource("TFileService");
0081 std::string fileName(iConfig.getUntrackedParameter<std::string>("rawFileName"));
0082 if (!fileName.empty()) {
0083 output_ = std::make_unique<std::ofstream>(fileName.c_str());
0084 if (!output_->good()) {
0085 edm::LogError("IOproblem") << "Could not open output file " << fileName << ".";
0086 output_.reset();
0087 }
0088 }
0089 }
0090
0091
0092 void BeamProfileHLLHC2DBReader::TheBSfromDB::init() {
0093 float dummy_double = 0.0;
0094 int dummy_int = 0;
0095
0096 run = dummy_int;
0097 ls = dummy_int;
0098 fMeanX = dummy_double;
0099 fMeanY = dummy_double;
0100 fMeanZ = dummy_double;
0101 fEProton = dummy_double;
0102 fCrabFrequency = dummy_double;
0103 fRF800 = dummy_double;
0104 fCrossingAngle = dummy_double;
0105 fCrabbingAngleCrossing = dummy_double;
0106 fCrabbingAngleSeparation = dummy_double;
0107 fBetaCrossingPlane = dummy_double;
0108 fBetaSeparationPlane = dummy_double;
0109 fHorizontalEmittance = dummy_double;
0110 fVerticalEmittance = dummy_double;
0111 fBunchLength = dummy_double;
0112 fTimeOffset = dummy_double;
0113 }
0114
0115
0116 void BeamProfileHLLHC2DBReader::analyze(const edm::Event& iEvent, const edm::EventSetup& iSetup) {
0117 using namespace edm;
0118 std::ostringstream output;
0119
0120
0121 theBSfromDB_.init();
0122
0123 if (watcher_.check(iSetup)) {
0124
0125 output << " for runs: " << iEvent.id().run() << " - " << iEvent.id().luminosityBlock() << std::endl;
0126
0127
0128 const SimBeamSpotHLLHCObjects* mybeamspot = &iSetup.getData(beamSpotToken_);
0129
0130 theBSfromDB_.run = iEvent.id().run();
0131 theBSfromDB_.ls = iEvent.id().luminosityBlock();
0132 theBSfromDB_.fMeanX = mybeamspot->meanX();
0133 theBSfromDB_.fMeanY = mybeamspot->meanY();
0134 theBSfromDB_.fMeanZ = mybeamspot->meanZ();
0135 theBSfromDB_.fEProton = mybeamspot->eProton();
0136 theBSfromDB_.fCrabFrequency = mybeamspot->crabFrequency();
0137 theBSfromDB_.fRF800 = mybeamspot->rf800();
0138 theBSfromDB_.fCrossingAngle = mybeamspot->crossingAngle();
0139 theBSfromDB_.fCrabbingAngleCrossing = mybeamspot->crabbingAngleCrossing();
0140 theBSfromDB_.fCrabbingAngleSeparation = mybeamspot->crabbingAngleSeparation();
0141 theBSfromDB_.fBetaCrossingPlane = mybeamspot->betaCrossingPlane();
0142 theBSfromDB_.fBetaSeparationPlane = mybeamspot->betaSeparationPlane();
0143 theBSfromDB_.fHorizontalEmittance = mybeamspot->horizontalEmittance();
0144 theBSfromDB_.fVerticalEmittance = mybeamspot->verticalEmittance();
0145 theBSfromDB_.fBunchLength = mybeamspot->bunchLenght();
0146 theBSfromDB_.fTimeOffset = mybeamspot->timeOffset();
0147 bstree_->Fill();
0148 output << *mybeamspot << std::endl;
0149 }
0150
0151
0152 if (output_.get())
0153 *output_ << output.str();
0154 else
0155 edm::LogInfo("BeamProfileHLLHC2DBReader") << output.str();
0156 }
0157
0158
0159 void BeamProfileHLLHC2DBReader::beginJob() {
0160 bstree_ = tFileService->make<TTree>("BSNtuple", "SimBeamSpotHLLHC analyzer ntuple");
0161
0162
0163 bstree_->Branch("run", &theBSfromDB_.run, "run/I");
0164 bstree_->Branch("ls", &theBSfromDB_.ls, "ls/I");
0165 bstree_->Branch("MeanX", &theBSfromDB_.fMeanX, "MeanX/F");
0166 bstree_->Branch("MeanY", &theBSfromDB_.fMeanY, "MeanY/F");
0167 bstree_->Branch("MeanZ", &theBSfromDB_.fMeanZ, "MeanZ/F");
0168 bstree_->Branch("EProton", &theBSfromDB_.fEProton, "EProton/F");
0169 bstree_->Branch("CrabFrequency", &theBSfromDB_.fCrabFrequency, "CrabFrequency/F");
0170 bstree_->Branch("RF800", &theBSfromDB_.fRF800, "RF800/O");
0171 bstree_->Branch("CrossingAngle", &theBSfromDB_.fCrossingAngle, "CrossingAngle/F");
0172 bstree_->Branch("CrabbingAngleCrossing", &theBSfromDB_.fCrabbingAngleCrossing, "CrabbingAngleCrossing/F");
0173 bstree_->Branch("CrabbingAngleSeparation", &theBSfromDB_.fCrabbingAngleSeparation, "CrabbingAngleSeparation/F");
0174 bstree_->Branch("BetaCrossingPlane", &theBSfromDB_.fBetaCrossingPlane, "BetaCrossingPlane/F");
0175 bstree_->Branch("BetaSeparationPlane", &theBSfromDB_.fBetaSeparationPlane, "BetaSeparationPlane/F");
0176 bstree_->Branch("HorizontalEmittance", &theBSfromDB_.fHorizontalEmittance, "HorizontalEmittance/F");
0177 bstree_->Branch("VerticalEmittance", &theBSfromDB_.fVerticalEmittance, "VerticalEmittance/F");
0178 bstree_->Branch("BunchLength", &theBSfromDB_.fBunchLength, "BunchLength/F");
0179 bstree_->Branch("TimeOffset", &theBSfromDB_.fTimeOffset, "TimeOffset/F");
0180 }
0181
0182
0183 void BeamProfileHLLHC2DBReader::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0184 edm::ParameterSetDescription desc;
0185 desc.addUntracked<std::string>("rawFileName", {});
0186 descriptions.addWithDefaultLabel(desc);
0187 }
0188
0189
0190 DEFINE_FWK_MODULE(BeamProfileHLLHC2DBReader);