CMS 3D CMS Logo

BeamProfileHLLHC2DBReader.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CondTools/BeamProfileHLLHC2DBReader
4 // Class: BeamProfileHLLHC2DBReader
5 //
13 //
14 // Original Author: Francesco Brivio
15 // Created: 11 June 2023
16 //
17 
18 // system include files
19 #include <fstream>
20 #include <memory>
21 #include <sstream>
22 
23 // user include files
34 
35 // For ROOT
38 #include <TTree.h>
39 
40 //
41 // class declaration
42 //
43 
44 class BeamProfileHLLHC2DBReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
45 public:
47  ~BeamProfileHLLHC2DBReader() override = default;
48 
49  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
50 
51 private:
52  void beginJob() override;
53  void analyze(const edm::Event&, const edm::EventSetup&) override;
54 
55  struct TheBSfromDB {
56  int run;
57  int ls;
58  double fMeanX, fMeanY, fMeanZ;
64  void init();
65  } theBSfromDB_;
66 
69  TTree* bstree_;
70 
71  // ----------member data ---------------------------
73  std::unique_ptr<std::ofstream> output_;
74 };
75 
76 // ------------ constructor ------------
78  : beamSpotToken_(esConsumes()), bstree_(nullptr) {
79  //now do what ever initialization is needed
80  usesResource("TFileService");
81  std::string fileName(iConfig.getUntrackedParameter<std::string>("rawFileName"));
82  if (!fileName.empty()) {
83  output_ = std::make_unique<std::ofstream>(fileName.c_str());
84  if (!output_->good()) {
85  edm::LogError("IOproblem") << "Could not open output file " << fileName << ".";
86  output_.reset();
87  }
88  }
89 }
90 
91 // ------------ SimBeamSpotHLLHCObjects initialization ------------
93  float dummy_double = 0.0;
94  int dummy_int = 0;
95 
96  run = dummy_int;
97  ls = dummy_int;
98  fMeanX = dummy_double;
99  fMeanY = dummy_double;
100  fMeanZ = dummy_double;
101  fEProton = dummy_double;
102  fCrabFrequency = dummy_double;
103  fRF800 = dummy_double;
104  fCrossingAngle = dummy_double;
105  fCrabbingAngleCrossing = dummy_double;
106  fCrabbingAngleSeparation = dummy_double;
107  fBetaCrossingPlane = dummy_double;
108  fBetaSeparationPlane = dummy_double;
109  fHorizontalEmittance = dummy_double;
110  fVerticalEmittance = dummy_double;
111  fBunchLength = dummy_double;
112  fTimeOffset = dummy_double;
113 }
114 
115 // ------------ method called for each event ------------
117  using namespace edm;
118  std::ostringstream output;
119 
120  // initialize the ntuple
121  theBSfromDB_.init();
122 
123  if (watcher_.check(iSetup)) { // check for new IOV for this run / LS
124 
125  output << " for runs: " << iEvent.id().run() << " - " << iEvent.id().luminosityBlock() << std::endl;
126 
127  // Get SimBeamSpotHLLHCObjects from EventSetup:
128  const SimBeamSpotHLLHCObjects* mybeamspot = &iSetup.getData(beamSpotToken_);
129 
130  theBSfromDB_.run = iEvent.id().run();
131  theBSfromDB_.ls = iEvent.id().luminosityBlock();
132  theBSfromDB_.fMeanX = mybeamspot->meanX();
133  theBSfromDB_.fMeanY = mybeamspot->meanY();
134  theBSfromDB_.fMeanZ = mybeamspot->meanZ();
135  theBSfromDB_.fEProton = mybeamspot->eProton();
136  theBSfromDB_.fCrabFrequency = mybeamspot->crabFrequency();
137  theBSfromDB_.fRF800 = mybeamspot->rf800();
138  theBSfromDB_.fCrossingAngle = mybeamspot->crossingAngle();
145  theBSfromDB_.fBunchLength = mybeamspot->bunchLenght();
146  theBSfromDB_.fTimeOffset = mybeamspot->timeOffset();
147  bstree_->Fill();
148  output << *mybeamspot << std::endl;
149  }
150 
151  // Final output - either message logger or output file:
152  if (output_.get())
153  *output_ << output.str();
154  else
155  edm::LogInfo("BeamProfileHLLHC2DBReader") << output.str();
156 }
157 
158 // ------------ method called once each job just before starting event loop ------------
160  bstree_ = tFileService->make<TTree>("BSNtuple", "SimBeamSpotHLLHC analyzer ntuple");
161 
162  //Tree Branches
163  bstree_->Branch("run", &theBSfromDB_.run, "run/I");
164  bstree_->Branch("ls", &theBSfromDB_.ls, "ls/I");
165  bstree_->Branch("MeanX", &theBSfromDB_.fMeanX, "MeanX/F");
166  bstree_->Branch("MeanY", &theBSfromDB_.fMeanY, "MeanY/F");
167  bstree_->Branch("MeanZ", &theBSfromDB_.fMeanZ, "MeanZ/F");
168  bstree_->Branch("EProton", &theBSfromDB_.fEProton, "EProton/F");
169  bstree_->Branch("CrabFrequency", &theBSfromDB_.fCrabFrequency, "CrabFrequency/F");
170  bstree_->Branch("RF800", &theBSfromDB_.fRF800, "RF800/O");
171  bstree_->Branch("CrossingAngle", &theBSfromDB_.fCrossingAngle, "CrossingAngle/F");
172  bstree_->Branch("CrabbingAngleCrossing", &theBSfromDB_.fCrabbingAngleCrossing, "CrabbingAngleCrossing/F");
173  bstree_->Branch("CrabbingAngleSeparation", &theBSfromDB_.fCrabbingAngleSeparation, "CrabbingAngleSeparation/F");
174  bstree_->Branch("BetaCrossingPlane", &theBSfromDB_.fBetaCrossingPlane, "BetaCrossingPlane/F");
175  bstree_->Branch("BetaSeparationPlane", &theBSfromDB_.fBetaSeparationPlane, "BetaSeparationPlane/F");
176  bstree_->Branch("HorizontalEmittance", &theBSfromDB_.fHorizontalEmittance, "HorizontalEmittance/F");
177  bstree_->Branch("VerticalEmittance", &theBSfromDB_.fVerticalEmittance, "VerticalEmittance/F");
178  bstree_->Branch("BunchLength", &theBSfromDB_.fBunchLength, "BunchLength/F");
179  bstree_->Branch("TimeOffset", &theBSfromDB_.fTimeOffset, "TimeOffset/F");
180 }
181 
182 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
185  desc.addUntracked<std::string>("rawFileName", {});
186  descriptions.addWithDefaultLabel(desc);
187 }
188 
189 //define this as a plug-in
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
~BeamProfileHLLHC2DBReader() override=default
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
Log< level::Error, false > LogError
edm::Service< TFileService > tFileService
double eProton() const
get EProton, fCrabFrequency, RF800
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< std::ofstream > output_
int iEvent
Definition: GenABIO.cc:224
double meanX() const
get meanX, meanY, meanZ position
void analyze(const edm::Event &, const edm::EventSetup &) override
double crossingAngle() const
set Crossing and Crabbing angles
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
BeamProfileHLLHC2DBReader(const edm::ParameterSet &)
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Log< level::Info, false > LogInfo
double betaCrossingPlane() const
get BetaStar and Emittance
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
struct BeamProfileHLLHC2DBReader::TheBSfromDB theBSfromDB_
HLT enums.
double bunchLenght() const
get BunchLength and TimeOffset
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
Definition: output.py:1
edm::ESWatcher< SimBeamSpotHLLHCObjectsRcd > watcher_
const edm::ESGetToken< SimBeamSpotHLLHCObjects, SimBeamSpotHLLHCObjectsRcd > beamSpotToken_