CMS 3D CMS Logo

BeamProfile2DBReader.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: CondTools/BeamProfile2DBReader
4 // Class: BeamProfile2DBReader
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 BeamProfile2DBReader : public edm::one::EDAnalyzer<edm::one::SharedResources> {
45 public:
46  explicit BeamProfile2DBReader(const edm::ParameterSet&);
47  ~BeamProfile2DBReader() 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 fX0, fY0, fZ0;
59  double fMeanX, fMeanY, fMeanZ;
62  double fPhi, fAlpha;
63  double fTimeOffset;
64  void init();
65  } theBSfromDB_;
66 
69  TTree* bstree_;
70 
71  // ----------member data ---------------------------
73  std::unique_ptr<std::ofstream> output_;
74 };
75 
76 //
77 // constructors and destructor
78 //
80  : beamSpotToken_(esConsumes()), bstree_(nullptr) {
81  //now do what ever initialization is needed
82  usesResource("TFileService");
83  std::string fileName(iConfig.getUntrackedParameter<std::string>("rawFileName"));
84  if (!fileName.empty()) {
85  output_ = std::make_unique<std::ofstream>(fileName.c_str());
86  if (!output_->good()) {
87  edm::LogError("IOproblem") << "Could not open output file " << fileName << ".";
88  output_.reset();
89  }
90  }
91 }
92 
93 //
94 // member functions
95 //
96 
98  float dummy_double = 0.0;
99  int dummy_int = 0;
100 
101  run = dummy_int;
102  ls = dummy_int;
103  fX0 = dummy_double;
104  fY0 = dummy_double;
105  fZ0 = dummy_double;
106  fMeanX = dummy_double;
107  fMeanY = dummy_double;
108  fMeanZ = dummy_double;
109  fSigmaX = dummy_double;
110  fSigmaY = dummy_double;
111  fSigmaZ = dummy_double;
112  fbetastar = dummy_double;
113  femittance = dummy_double;
114  fPhi = dummy_double;
115  fAlpha = dummy_double;
116  fTimeOffset = dummy_double;
117 }
118 
119 // ------------ method called for each event ------------
121  using namespace edm;
122  std::ostringstream output;
123 
124  // initialize the ntuple
125  theBSfromDB_.init();
126 
127  if (watcher_.check(iSetup)) { // check for new IOV for this run / LS
128 
129  output << " for runs: " << iEvent.id().run() << " - " << iEvent.id().luminosityBlock() << std::endl;
130 
131  // Get SimBeamSpot from EventSetup:
132  const SimBeamSpotObjects* mybeamspot = &iSetup.getData(beamSpotToken_);
133 
134  theBSfromDB_.run = iEvent.id().run();
135  theBSfromDB_.ls = iEvent.id().luminosityBlock();
136  theBSfromDB_.fX0 = mybeamspot->x();
137  theBSfromDB_.fY0 = mybeamspot->y();
138  theBSfromDB_.fZ0 = mybeamspot->z();
139  theBSfromDB_.fMeanX = mybeamspot->meanX();
140  theBSfromDB_.fMeanY = mybeamspot->meanY();
141  theBSfromDB_.fMeanZ = mybeamspot->meanZ();
142  theBSfromDB_.fSigmaX = mybeamspot->sigmaX();
143  theBSfromDB_.fSigmaY = mybeamspot->sigmaY();
144  theBSfromDB_.fSigmaZ = mybeamspot->sigmaZ();
145  theBSfromDB_.fbetastar = mybeamspot->betaStar();
146  theBSfromDB_.femittance = mybeamspot->emittance();
147  theBSfromDB_.fPhi = mybeamspot->phi();
148  theBSfromDB_.fAlpha = mybeamspot->alpha();
149  theBSfromDB_.fTimeOffset = mybeamspot->timeOffset();
150  bstree_->Fill();
151  output << *mybeamspot << std::endl;
152  }
153 
154  // Final output - either message logger or output file:
155  if (output_.get())
156  *output_ << output.str();
157  else
158  edm::LogInfo("BeamProfile2DBReader") << output.str();
159 }
160 
161 // ------------ method called once each job just before starting event loop ------------
163  bstree_ = tFileService->make<TTree>("BSNtuple", "SimBeamSpot analyzer ntuple");
164 
165  //Tree Branches
166  bstree_->Branch("run", &theBSfromDB_.run, "run/I");
167  bstree_->Branch("ls", &theBSfromDB_.ls, "ls/I");
168  bstree_->Branch("BSx0", &theBSfromDB_.fX0, "BSx0/F");
169  bstree_->Branch("BSy0", &theBSfromDB_.fY0, "BSy0/F");
170  bstree_->Branch("BSz0", &theBSfromDB_.fZ0, "BSz0/F");
171  bstree_->Branch("BSmeanX", &theBSfromDB_.fMeanX, "BSmeanX/F");
172  bstree_->Branch("BSmeanY", &theBSfromDB_.fMeanY, "BSmeanY/F");
173  bstree_->Branch("BSmeanZ", &theBSfromDB_.fMeanZ, "BSmeanZ/F");
174  bstree_->Branch("Beamsigmax", &theBSfromDB_.fSigmaX, "Beamsigmax/F");
175  bstree_->Branch("Beamsigmay", &theBSfromDB_.fSigmaY, "Beamsigmay/F");
176  bstree_->Branch("Beamsigmaz", &theBSfromDB_.fSigmaZ, "Beamsigmaz/F");
177  bstree_->Branch("BetaStar", &theBSfromDB_.fbetastar, "BetaStar/F");
178  bstree_->Branch("Emittance", &theBSfromDB_.femittance, "Emittance/F");
179  bstree_->Branch("Phi", &theBSfromDB_.fPhi, "Phi/F");
180  bstree_->Branch("Alpha", &theBSfromDB_.fAlpha, "Alpha/F");
181  bstree_->Branch("TimeOffset", &theBSfromDB_.fTimeOffset, "TimeOffset/F");
182 }
183 
184 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
187  desc.addUntracked<std::string>("rawFileName", {});
188  descriptions.addWithDefaultLabel(desc);
189 }
190 
191 //define this as a plug-in
double sigmaZ() const
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
struct BeamProfile2DBReader::TheBSfromDB theBSfromDB_
double sigmaX() const
get sigmaX, sigmaY, sigmaZ
void analyze(const edm::Event &, const edm::EventSetup &) override
T const & getData(const ESGetToken< T, R > &iToken) const noexcept(false)
Definition: EventSetup.h:119
BeamProfile2DBReader(const edm::ParameterSet &)
edm::Service< TFileService > tFileService
double meanX() const
get meanX, meanY, meanZ position
Log< level::Error, false > LogError
double alpha() const
std::unique_ptr< std::ofstream > output_
T getUntrackedParameter(std::string const &, T const &) const
int iEvent
Definition: GenABIO.cc:224
edm::ESWatcher< SimBeamSpotObjectsRcd > watcher_
const edm::ESGetToken< SimBeamSpotObjects, SimBeamSpotObjectsRcd > beamSpotToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
Log< level::Info, false > LogInfo
double timeOffset() const
double meanY() const
bool check(const edm::EventSetup &iSetup)
Definition: ESWatcher.h:57
double meanZ() const
double emittance() const
HLT enums.
double phi() const
get Phi, Alpha and TimeOffset
double betaStar() const
get BetaStar and Emittance
T * make(const Args &...args) const
make new ROOT object
Definition: TFileService.h:64
Definition: output.py:1
double x() const
get X, Y, Z position
~BeamProfile2DBReader() override=default