Belle II Software  release-08-01-10
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/geometry/B2Vector3.h>
20 #include <framework/utilities/FileSystem.h>
21 #include <mdst/dbobjects/BeamSpot.h>
22 #include <mdst/dbobjects/CollisionBoostVector.h>
23 #include <mdst/dbobjects/CollisionInvariantMass.h>
24 
25 #include <list>
26 
27 namespace Belle {
31  int check_beginrun(bool) { return 0; }
32 }
33 
34 namespace {
36  TVector3 CLHEPtoROOT(const HepPoint3D& point)
37  {
38  return TVector3(point.x(), point.y(), point.z());
39  }
41  TMatrixDSym CLHEPtoROOT(const HepSymMatrix& matrix)
42  {
43  TMatrixDSym result(matrix.num_col());
44  for (int i = 0; i < matrix.num_row(); ++i) {
45  for (int j = 0; j < matrix.num_col(); ++j) {
46  result(i, j) = matrix(i + 1, j + 1);
47  }
48  }
49  return result;
50  }
51 }
52 
53 namespace Belle2 {
59  REG_MODULE(B2BIIConvertBeamParams);
60 
62  {
63  addParam("mcFlag", m_mcFlag, "which mc flag to use", m_mcFlag);
64  addParam("missingBenergy", m_missingBenergy, "where to store information about runs with missing database info", m_missingBenergy);
65  addParam("missingIp", m_missingIp, "where to store information about runs with missing IP profile info", m_missingIp);
66  addParam("smearEnergy", m_SmearEnergy, "Allow beam energy smearing for MC generation", true);
67  addParam("smearDirection", m_SmearDirection, "Allow beam direction smearing for MC generation", true);
68  addParam("smearVertex", m_SmearVertex, "Allow IP position smearing for MC generation", true);
69  addParam("generateCMS", m_GenerateCMS, "Generate events in CMS, not lab system.", false);
70  addParam("storeBeamParameters", m_storeBeamParameters, "Store the BeamParameters payloads in the localDB", true);
71  addParam("storeCollisionInvariantMass", m_storeCollisionInvariantMass, "Store the CollisionInvariantMass payloads in the localDB",
72  true);
73  addParam("storeCollisionBoostVector", m_storeCollisionBoostVector, "Store the CollisionBoostVector payloads in the localDB", true);
74  addParam("storeBeamSpot", m_storeBeamSpot, "Store the BeamSpot payloads in the localDB", true);
75  }
76 
78  {
79  BsInit(0);
80  m_event.isRequired();
81  }
82 
84  {
85  // Fudge a runhead record with all the info we need to be able to get the
86  // Beamenergy/IpProfile from the database
87  BsClrTab(BBS_CLEAR);
88  Belle::Belle_runhead_Manager& runMgr = Belle::Belle_runhead_Manager::get_manager();
89  runMgr.add();
90  Belle::Belle_runhead& runhead = runMgr[0];
91  runhead.ExpMC(m_mcFlag);
92  runhead.ExpNo(m_event->getExperiment());
93  runhead.RunNo(m_event->getRun());
94  B2INFO("Obtaining values for exp " << m_event->getExperiment() << ", run " << m_event->getRun());
95  // and then get the values from the database
96  Belle::BeamEnergy::flag_database(0);
97  Belle::IpProfile::begin_run();
98  Belle::Ip_profile_general_Manager& ipMgr = Belle::Ip_profile_general_Manager::get_manager();
99 
100  if (!Belle::BeamEnergy::usable()) {
102  if (!lock.lock()) {
103  B2ERROR("No BeamEnergy for exp " << m_event->getExperiment() << ", run " << m_event->getRun());
104  } else {
105  std::ofstream file(m_missingBenergy.c_str(), std::ios::app);
106  file << m_event->getExperiment() << "," << m_event->getRun() << std::endl;
107  }
108  return;
109  }
110 
111  const double Eher = Belle::BeamEnergy::E_HER();
112  const double Eler = Belle::BeamEnergy::E_LER();
113  const double crossingAngle = Belle::BeamEnergy::Cross_angle();
114  const double angleLer = M_PI; //parallel to negative z axis (different from Belle II!)
115  const double angleHer = crossingAngle; //in positive z and x direction, verified to be consistent with Upsilon(4S) momentum
116 
117  // Beam energy spread taken from the hard-coded values in the Belle evtgen:
118  // /sw/belle/belle/b20090127_0910/src/util/evtgenutil/basf_if/evtgen.cc
119  // The beam energy spread seems not have been logged in the belle DB at all.
120  // If you are reading this and you know how to access the run-dependent beam energy spread,
121  // please fix this (April 2021)!
122  std::vector<double> covarianceHer = {0.00513 * 0.00513}; // Energy spread only. No idea about the direction spread
123  std::vector<double> covarianceLer = {0.002375 * 0.002375}; // Energy spread only. No idea about the direction spread
124 
125  IntervalOfValidity iov(m_event->getExperiment(), m_event->getRun(), m_event->getExperiment(), m_event->getRun());
126 
127  BeamParameters beamParams;
128  beamParams.setLER(Eler, angleLer, 0, covarianceLer);
129  beamParams.setHER(Eher, angleHer, 0, covarianceHer);
130 
131  // set the flags according to the module settings.
132  int flags = 0;
133  if (m_GenerateCMS)
135  if (m_SmearEnergy)
137  if (m_SmearDirection)
139  if (m_SmearVertex)
141  beamParams.setGenerationFlags(flags);
142 
143  CollisionBoostVector collisionBoostVector;
144  CollisionInvariantMass collisionInvM;
145  ROOT::Math::PxPyPzEVector cms = beamParams.getLER() + beamParams.getHER();
146  collisionBoostVector.setBoost(B2Vector3D(-cms.BoostToCM()), TMatrixTSym<double>(3));
147  //note: maybe we could use Belle::BeamEnergy::E_beam_corr(), Belle::BeamEnergy::E_beam_err()
148  collisionInvM.setMass(cms.M(), 0.0, 0.0);
149 
150  // Boost vector and invariant mass are not intra-run dependent, store now
152  Database::Instance().storeData("CollisionBoostVector", &collisionBoostVector, iov);
154  Database::Instance().storeData("CollisionInvariantMass", &collisionInvM, iov);
155 
156  // and now we continue with the vertex
157  if (!Belle::IpProfile::usable()) {
159  if (!lock.lock()) {
160  B2ERROR("No IpProfile for exp " << m_event->getExperiment() << ", run " << m_event->getRun());
161  } else {
162  std::ofstream file(m_missingIp.c_str(), std::ios::app);
163  file << m_event->getExperiment() << "," << m_event->getRun() << std::endl;
164  }
165  std::vector<double> covariance; // 0 = no error
166  beamParams.setVertex(ROOT::Math::XYZVector(
167  std::numeric_limits<double>::quiet_NaN(),
168  std::numeric_limits<double>::quiet_NaN(),
169  std::numeric_limits<double>::quiet_NaN()), covariance);
171  Database::Instance().storeData("BeamParameters", &beamParams, iov);
172 
173  BeamSpot beamSpot;
174  beamSpot.setIP(
175  TVector3(std::numeric_limits<double>::quiet_NaN(),
176  std::numeric_limits<double>::quiet_NaN(),
177  std::numeric_limits<double>::quiet_NaN()
178  ), TMatrixTSym<double>(3)
179  );
180  if (m_storeBeamSpot) {
181  Database::Instance().storeData("BeamSpot", &beamSpot, iov);
182  B2INFO("No IpProfile, created BeamSpot Payload");
183  }
184  return;
185  }
186 
187  int nbins = std::max(ipMgr.count() - 1, 1);
188  B2INFO("exp " << m_event->getExperiment() << ", run " << m_event->getRun() << ": " << nbins << " IpProfile bins");
189 
190  // need to keep the objects alive until we created the payloads
191  std::list<BeamParameters> beamparamsList;
192  std::list<BeamSpot> beamspotList;
193  // and we want them to be intra-run dependent
194  std::unique_ptr<EventDependency> beamparams;
195  std::unique_ptr<EventDependency> beamspots;
196 
197  Belle::Belle_event_Manager& evtMgr = Belle::Belle_event_Manager::get_manager();
198  if (!evtMgr.count()) evtMgr.add();
199  for (int i = 0; i < nbins; ++i) {
200  int evtNr = i * BIN_EVENTS;
201  // fudge belle event record to get the correct run dependent ip profile
202  Belle::Belle_event& evt = evtMgr[0];
203  evt.ExpMC(m_mcFlag);
204  evt.ExpNo(m_event->getExperiment());
205  evt.RunNo(m_event->getRun());
206  evt.EvtNo(evtNr);
207  Belle::IpProfile::set_evtbin_number();
208  B2ASSERT("something wrong: " << Belle::IpProfile::EvtBinNo() << "!=" << i, Belle::IpProfile::EvtBinNo() == i);
209 
210  // and convert
211  HepPoint3D ip = Belle::IpProfile::e_position();
212  HepSymMatrix ipErr = Belle::IpProfile::e_position_err();
213 
214  beamParams.setVertex(ROOT::Math::XYZVector(ip.x(), ip.y(), ip.z()));
215  beamParams.setCovVertex(CLHEPtoROOT(ipErr));
216  beamparamsList.emplace_back(beamParams);
217 
218  BeamSpot beamSpot;
219  beamSpot.setIP(CLHEPtoROOT(ip), CLHEPtoROOT(ipErr));
220  beamspotList.emplace_back(beamSpot);
221 
222  if (!beamparams) {
223  beamparams = std::make_unique<EventDependency>(&beamparamsList.back(), false);
224  beamspots = std::make_unique<EventDependency>(&beamspotList.back(), false);
225  } else {
226  beamparams->add(evtNr, &beamparamsList.back());
227  beamspots->add(evtNr, &beamspotList.back());
228  }
229  }
230  if (beamparamsList.size() < 1) {
231  B2ERROR("Something is wrong with exp " << m_event->getExperiment() << ", run " << m_event->getRun() << ": no bins found");
232  return;
233  }
234  if (beamparamsList.size() == 1) {
235  // just one bin? no need to store event dependency
237  Database::Instance().storeData("BeamParameters", &beamparamsList.back(), iov);
238  if (m_storeBeamSpot)
239  Database::Instance().storeData("BeamSpot", &beamspotList.back(), iov);
240  B2INFO("Created event independent payload");
241  } else {
242  // otherwise store full information
244  Database::Instance().storeData("BeamParameters", beamparams.get(), iov);
245  if (m_storeBeamSpot)
246  Database::Instance().storeData("BeamSpot", beamspots.get(), iov);
247  B2INFO("Created event dependent payload with " << beamparamsList.size() << " entries");
248  }
249  }
251 } //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 ROOT::Math::XYZVector &vertex, const std::vector< double > &cov)
Set the vertex position and error matrix.
void setGenerationFlags(int flags) override
Set the generation flags to be used for event generation (ORed combination of EGenerationFlags).
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:186
A class that describes the interval of experiments/runs for which an object in the database is valid.
const ROOT::Math::PxPyPzEVector & getLER() const
Get 4vector of the low energy beam.
const ROOT::Math::PxPyPzEVector & getHER() const
Get 4vector of the high energy beam.
@ 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)
REG_MODULE(arichBtest)
Register the Module.
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:42
bool storeData(const std::string &name, TObject *object, const IntervalOfValidity &iov)
Store an object in the database.
Definition: Database.cc:141
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
Abstract base class for different kinds of events.