Belle II Software  release-06-02-00
B2BIIConvertBeamParamsModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include <b2bii/modules/B2BIIMdstInput/B2BIIConvertBeamParamsModule.h>
10 #include <belle_legacy/benergy/BeamEnergy.h>
11 #include <belle_legacy/ip/IpProfile.h>
12 #include <belle_legacy/panther/panther.h>
13 #include <belle_legacy/tables/belletdf.h>
14 
15 #include <framework/database/Database.h>
16 #include <framework/database/EventDependency.h>
17 #include <framework/database/IntervalOfValidity.h>
18 #include <framework/dbobjects/BeamParameters.h>
19 #include <framework/utilities/FileSystem.h>
20 #include <mdst/dbobjects/BeamSpot.h>
21 #include <mdst/dbobjects/CollisionBoostVector.h>
22 #include <mdst/dbobjects/CollisionInvariantMass.h>
23 
24 #include <list>
25 
26 namespace Belle {
30  int check_beginrun(bool) { return 0; }
31 }
32 
33 namespace {
35  TVector3 CLHEPtoROOT(const HepPoint3D& point)
36  {
37  return TVector3(point.x(), point.y(), point.z());
38  }
40  TMatrixDSym CLHEPtoROOT(const HepSymMatrix& matrix)
41  {
42  TMatrixDSym result(matrix.num_col());
43  for (int i = 0; i < matrix.num_row(); ++i) {
44  for (int j = 0; j < matrix.num_col(); ++j) {
45  result(i, j) = matrix(i + 1, j + 1);
46  }
47  }
48  return result;
49  }
50 }
51 
52 namespace Belle2 {
58  REG_MODULE(B2BIIConvertBeamParams);
59 
61  {
62  addParam("mcFlag", m_mcFlag, "which mc flag to use", m_mcFlag);
63  addParam("missingBenergy", m_missingBenergy, "where to store information about runs with missing database info", m_missingBenergy);
64  addParam("missingIp", m_missingIp, "where to store information about runs with missing IP profile info", m_missingIp);
65  addParam("smearEnergy", m_SmearEnergy, "Allow beam energy smearing for MC generation", true);
66  addParam("smearDirection", m_SmearDirection, "Allow beam direction smearing for MC generation", true);
67  addParam("smearVertex", m_SmearVertex, "Allow IP position smearing for MC generation", true);
68  addParam("generateCMS", m_GenerateCMS, "Generate events in CMS, not lab system.", false);
69  addParam("storeBeamParameters", m_storeBeamParameters, "Store the BeamParameters payloads in the localDB", true);
70  addParam("storeCollisionInvariantMass", m_storeCollisionInvariantMass, "Store the CollisionInvariantMass payloads in the localDB",
71  true);
72  addParam("storeCollisionBoostVector", m_storeCollisionBoostVector, "Store the CollisionBoostVector payloads in the localDB", true);
73  addParam("storeBeamSpot", m_storeBeamSpot, "Store the BeamSpot payloads in the localDB", true);
74  }
75 
77  {
78  BsInit(0);
79  m_event.isRequired();
80  }
81 
83  {
84  // Fudge a runhead record with all the info we need to be able to get the
85  // Beamenergy/IpProfile from the database
86  BsClrTab(BBS_CLEAR);
87  Belle::Belle_runhead_Manager& runMgr = Belle::Belle_runhead_Manager::get_manager();
88  runMgr.add();
89  Belle::Belle_runhead& runhead = runMgr[0];
90  runhead.ExpMC(m_mcFlag);
91  runhead.ExpNo(m_event->getExperiment());
92  runhead.RunNo(m_event->getRun());
93  B2INFO("Obtaining values for exp " << m_event->getExperiment() << ", run " << m_event->getRun());
94  // and then get the values from the database
95  Belle::BeamEnergy::flag_database(0);
96  Belle::IpProfile::begin_run();
97  Belle::Ip_profile_general_Manager& ipMgr = Belle::Ip_profile_general_Manager::get_manager();
98 
99  if (!Belle::BeamEnergy::usable()) {
101  if (!lock.lock()) {
102  B2ERROR("No BeamEnergy for exp " << m_event->getExperiment() << ", run " << m_event->getRun());
103  } else {
104  std::ofstream file(m_missingBenergy.c_str(), std::ios::app);
105  file << m_event->getExperiment() << "," << m_event->getRun() << std::endl;
106  }
107  return;
108  }
109 
110  const double Eher = Belle::BeamEnergy::E_HER();
111  const double Eler = Belle::BeamEnergy::E_LER();
112  const double crossingAngle = Belle::BeamEnergy::Cross_angle();
113  const double angleLer = M_PI; //parallel to negative z axis (different from Belle II!)
114  const double angleHer = crossingAngle; //in positive z and x direction, verified to be consistent with Upsilon(4S) momentum
115 
116  // Beam energy spread taken from the hard-coded values in the Belle evtgen:
117  // /sw/belle/belle/b20090127_0910/src/util/evtgenutil/basf_if/evtgen.cc
118  // The beam energy spread seems not have been logged in the belle DB at all.
119  // If you are reading this and you know how to access the run-dependent beam energy spread,
120  // please fix this (April 2021)!
121  std::vector<double> covarianceHer = {0.00513 * 0.00513}; // Energy spread only. No idea about the direction spread
122  std::vector<double> covarianceLer = {0.002375 * 0.002375}; // Energy spread only. No idea about the direction spread
123 
124  IntervalOfValidity iov(m_event->getExperiment(), m_event->getRun(), m_event->getExperiment(), m_event->getRun());
125 
126  BeamParameters beamParams;
127  beamParams.setLER(Eler, angleLer, 0, covarianceLer);
128  beamParams.setHER(Eher, angleHer, 0, covarianceHer);
129 
130  // set the flags according to the module settings.
131  int flags = 0;
132  if (m_GenerateCMS)
134  if (m_SmearEnergy)
136  if (m_SmearDirection)
138  if (m_SmearVertex)
140  beamParams.setGenerationFlags(flags);
141 
142  CollisionBoostVector collisionBoostVector;
143  CollisionInvariantMass collisionInvM;
144  TLorentzVector cms = beamParams.getLER() + beamParams.getHER();
145  collisionBoostVector.setBoost(cms.BoostVector(), TMatrixTSym<double>(3));
146  //note: maybe we could use Belle::BeamEnergy::E_beam_corr(), Belle::BeamEnergy::E_beam_err()
147  collisionInvM.setMass(cms.M(), 0.0 , 0.0);
148 
149  // Boost vector and invariant mass are not intra-run dependent, store now
151  Database::Instance().storeData("CollisionBoostVector", &collisionBoostVector, iov);
153  Database::Instance().storeData("CollisionInvariantMass", &collisionInvM, iov);
154 
155  // and now we continue with the vertex
156  if (!Belle::IpProfile::usable()) {
158  if (!lock.lock()) {
159  B2ERROR("No IpProfile for exp " << m_event->getExperiment() << ", run " << m_event->getRun());
160  } else {
161  std::ofstream file(m_missingIp.c_str(), std::ios::app);
162  file << m_event->getExperiment() << "," << m_event->getRun() << std::endl;
163  }
164  std::vector<double> covariance; // 0 = no error
165  beamParams.setVertex(TVector3(
166  std::numeric_limits<double>::quiet_NaN(),
167  std::numeric_limits<double>::quiet_NaN(),
168  std::numeric_limits<double>::quiet_NaN()), covariance);
170  Database::Instance().storeData("BeamParameters", &beamParams, iov);
171 
172  BeamSpot beamSpot;
173  beamSpot.setIP(
174  TVector3(std::numeric_limits<double>::quiet_NaN(),
175  std::numeric_limits<double>::quiet_NaN(),
176  std::numeric_limits<double>::quiet_NaN()
177  ), TMatrixTSym<double>(3)
178  );
179  if (m_storeBeamSpot) {
180  Database::Instance().storeData("BeamSpot", &beamSpot, iov);
181  B2INFO("No IpProfile, created BeamSpot Payload");
182  }
183  return;
184  }
185 
186  int nbins = std::max(ipMgr.count() - 1, 1);
187  B2INFO("exp " << m_event->getExperiment() << ", run " << m_event->getRun() << ": " << nbins << " IpProfile bins");
188 
189  // need to keep the objects alive until we created the payloads
190  std::list<BeamParameters> beamparamsList;
191  std::list<BeamSpot> beamspotList;
192  // and we want them to be intra-run dependent
193  std::unique_ptr<EventDependency> beamparams;
194  std::unique_ptr<EventDependency> beamspots;
195 
196  Belle::Belle_event_Manager& evtMgr = Belle::Belle_event_Manager::get_manager();
197  if (!evtMgr.count()) evtMgr.add();
198  for (int i = 0; i < nbins; ++i) {
199  int evtNr = i * BIN_EVENTS;
200  // fudge belle event record to get the correct run dependent ip profile
201  Belle::Belle_event& evt = evtMgr[0];
202  evt.ExpMC(m_mcFlag);
203  evt.ExpNo(m_event->getExperiment());
204  evt.RunNo(m_event->getRun());
205  evt.EvtNo(evtNr);
206  Belle::IpProfile::set_evtbin_number();
207  B2ASSERT("something wrong: " << Belle::IpProfile::EvtBinNo() << "!=" << i, Belle::IpProfile::EvtBinNo() == i);
208 
209  // and convert
210  HepPoint3D ip = Belle::IpProfile::e_position();
211  HepSymMatrix ipErr = Belle::IpProfile::e_position_err();
212 
213  beamParams.setVertex(CLHEPtoROOT(ip));
214  beamParams.setCovVertex(CLHEPtoROOT(ipErr));
215  beamparamsList.emplace_back(beamParams);
216 
217  BeamSpot beamSpot;
218  beamSpot.setIP(CLHEPtoROOT(ip), CLHEPtoROOT(ipErr));
219  beamspotList.emplace_back(beamSpot);
220 
221  if (!beamparams) {
222  beamparams = std::make_unique<EventDependency>(&beamparamsList.back(), false);
223  beamspots = std::make_unique<EventDependency>(&beamspotList.back(), false);
224  } else {
225  beamparams->add(evtNr, &beamparamsList.back());
226  beamspots->add(evtNr, &beamspotList.back());
227  }
228  }
229  if (beamparamsList.size() < 1) {
230  B2ERROR("Something is wrong with exp " << m_event->getExperiment() << ", run " << m_event->getRun() << ": no bins found");
231  return;
232  }
233  if (beamparamsList.size() == 1) {
234  // just one bin? no need to store event dependency
236  Database::Instance().storeData("BeamParameters", &beamparamsList.back(), iov);
237  if (m_storeBeamSpot)
238  Database::Instance().storeData("BeamSpot", &beamspotList.back(), iov);
239  B2INFO("Created event independent payload");
240  } else {
241  // otherwise store full information
243  Database::Instance().storeData("BeamParameters", beamparams.get(), iov);
244  if (m_storeBeamSpot)
245  Database::Instance().storeData("BeamSpot", beamspots.get(), iov);
246  B2INFO("Created event dependent payload with " << beamparamsList.size() << " entries");
247  }
248  }
250 } //Belle2 namespace
bool m_SmearDirection
Smear beam direction when generating initial events.
std::string m_missingBenergy
Where to store information about runs without beam energy information.
std::string m_missingIp
Where to store information about runs without IP profile information.
bool m_storeCollisionInvariantMass
Store the CollisionInvariantMass payloads in the localDB.
StoreObjPtr< EventMetaData > m_event
Event metadata.
bool m_SmearVertex
Smear vertex position when generating initial events.
bool m_storeBeamParameters
Store the BeamParameters payloads in the localDB.
bool m_GenerateCMS
Generate events in CMS, not lab system.
bool m_storeBeamSpot
Store the BeamSpot payloads in the localDB.
bool m_SmearEnergy
Smear energy when generating initial events.
bool m_storeCollisionBoostVector
Store the CollisionBoostVector payloads in the localDB.
int m_mcFlag
mcFlag to use when getting belle database content
This class contains the nominal beam parameters and the parameters used for smearing of the primary v...
void setHER(double energy, double angleX, double angleY, const std::vector< double > &cov)
Set the HER FourVector and error matrix from beam energy, angle and covariance entries.
void setLER(double energy, double angleX, double angleY, const std::vector< double > &cov)
Set the LER FourVector and error matrix from beam energy, angle and covariance entries.
void setVertex(const TVector3 &vertex, const std::vector< double > &cov)
Set the vertex position and error matrix.
void setCovVertex(const TMatrixDSym &cov)
Set the covariance matrix of the vertex position.
This class contains the beam spot position and size modeled as a gaussian distribution in space.
Definition: BeamSpot.h:22
void setIP(const TVector3 &ipPosition, const TMatrixDSym &covariance)
Set the IP position and its error matrix.
Definition: BeamSpot.h:59
This class contains the measured average boost vector vec(beta) = (beta_x, beta_y,...
void setBoost(const TVector3 &boost, const TMatrixDSym &covariance)
Set the boost vector and its error matrix.
This class contains the measured average center-of-mass energy, which is equal to the invariant mass ...
void setMass(double mass, double error, double spread)
Set the CMS energy and its uncertainty.
Helper class for locking a file.
Definition: FileSystem.h:97
bool lock(int timeout=300, bool ignoreErrors=false)
Try to lock the file.
Definition: FileSystem.cc:183
A class that describes the interval of experiments/runs for which an object in the database is valid.
@ c_generateCMS
generate initial event in CMS instead of lab
@ c_smearBeamEnergy
smear energy of HER and LER (but not direction)
@ c_smearBeamDirection
smear direction of HER and LER (but not energy)
const TLorentzVector & getLER() const
Get 4vector of the low energy beam.
const TLorentzVector & getHER() const
Get 4vector of the high energy beam.
void setGenerationFlags(int flags)
Set the generation flags to be used for event generation (ORed combination of EGenerationFlags)
void initialize() override
Initialize phanter banks.
void beginRun() override
Set run info in panther and load IPProfile/Benergy and convert the values to payloads.
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
static Database & Instance()
Instance of a singleton Database.
Definition: Database.cc:41
bool storeData(const std::string &name, TObject *object, const IntervalOfValidity &iov)
Store an object in the database.
Definition: Database.cc:140
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.