Belle II Software  release-05-02-19
InitialParticleGeneration.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Ritter *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #include <generators/utilities/InitialParticleGeneration.h>
12 #include <framework/gearbox/Const.h>
13 #include <framework/logging/Logger.h>
14 
15 namespace Belle2 {
22  {
23  m_event.registerInDataStore();
24  }
25 
26  TVector3 InitialParticleGeneration::generateVertex(const TVector3& initial, const TMatrixDSym& cov,
27  MultivariateNormalGenerator& gen)
28  {
29  if (m_event->hasGenerationFlags(BeamParameters::c_smearVertex)) {
30  if (!gen.size()) gen.setMeanCov(initial, cov);
31  return gen.generateVec3();
32  }
33  return initial;
34  }
35 
36  TLorentzVector InitialParticleGeneration::generateBeam(const TLorentzVector& initial, const TMatrixDSym& cov,
38  {
39  if (!(m_event->getGenerationFlags() & BeamParameters::c_smearBeam)) {
40  //no smearing, nothing to do
41  return initial;
42  }
43  if (gen.size() != 3) gen.setMeanCov(TVector3(initial.E(), 0, 0), cov);
44  // generate the beam parameters: energy and horizontal angle and vertical
45  // angle with respect to nominal direction
46  // p[0] = energy
47  // p[1] = horizontal angle
48  // p[2] = vertical angle
49  Eigen::VectorXd p = gen.generate();
50  // if we don't smear beam energy set it to nominal values
51  if (!m_event->hasGenerationFlags(BeamParameters::c_smearBeamEnergy)) {
52  p[0] = initial.E();
53  }
54  // if we don't smear beam direction set angles to 0
55  if (!m_event->hasGenerationFlags(BeamParameters::c_smearBeamDirection)) {
56  p[1] = 0;
57  p[2] = 0;
58  }
59  // calculate TLorentzvector from p
60  const double pz = std::sqrt(p[0] * p[0] - Const::electronMass * Const::electronMass);
61  // Rotate around theta_x and theta_y. Same as result.RotateX(p[1]);
62  // result.RotateY(p[2]) but as we know it's a 0,0,pz vector this can be
63  // optimized.
64  const double sx = sin(p[1]);
65  const double cx = cos(p[1]);
66  const double sy = sin(p[2]);
67  const double cy = cos(p[2]);
68  const double px = sy * cx * pz;
69  const double py = -sx * pz;
70  TLorentzVector result(px, py, cx * cy * pz, p[0]);
71  // And rotate into the direction of the nominal beam
72  result.RotateY(initial.Theta());
73  return result;
74  }
75 
76  MCInitialParticles& InitialParticleGeneration::generate()
77  {
78  return generate(m_allowedFlags);
79  }
80 
81  MCInitialParticles& InitialParticleGeneration::generate(int allowedFlags)
82  {
83  if (!m_event) {
84  m_event.create();
85  }
86  if (!m_beamParams.isValid()) {
87  B2FATAL("Cannot generate beam without valid BeamParameters");
88  }
89  if (m_beamParams.hasChanged()) {
93  }
94  m_event->setGenerationFlags(m_beamParams->getGenerationFlags() & allowedFlags);
95  TLorentzVector her = generateBeam(m_beamParams->getHER(), m_beamParams->getCovHER(), m_generateHER);
96  TLorentzVector ler = generateBeam(m_beamParams->getLER(), m_beamParams->getCovLER(), m_generateLER);
97  TVector3 vtx = generateVertex(m_beamParams->getVertex(), m_beamParams->getCovVertex(), m_generateVertex);
98  m_event->set(her, ler, vtx);
99  //Check if we want to go to CMS, if so boost both
100  if (m_beamParams->hasGenerationFlags(BeamParameters::c_generateCMS)) {
101  m_event->setGenerationFlags(0);
102  her = m_event->getLabToCMS() * her;
103  ler = m_event->getLabToCMS() * ler;
104  m_event->set(her, ler, vtx);
105  m_event->setGenerationFlags(m_beamParams->getGenerationFlags() & allowedFlags);
106  }
107  return *m_event;
108  }
109 
110  TVector3 InitialParticleGeneration::updateVertex(bool force)
111  {
112  if (!m_beamParams.isValid()) {
113  B2FATAL("Cannot generate beam without valid BeamParameters");
114  }
115  if (!m_event) {
116  // generate a new mc initial particle without smearing except for the vertex
118  return m_event->getVertex();
119  }
120  if (!m_beamParams->hasGenerationFlags(BeamParameters::c_smearVertex) or
121  (m_event->hasGenerationFlags(BeamParameters::c_smearVertex) and not force)) {
122  // already has a smeared vertex or smearing disallowed. nothing to do
123  return {0, 0, 0};
124  }
125  // Ok, now we need to just update the vertex, make sure the parameters are up to date ...
126  if (m_beamParams.hasChanged()) {
130  }
131  // Add the smear vertex flag
132  m_event->setGenerationFlags(m_event->getGenerationFlags() | BeamParameters::c_smearVertex);
133  // And generate a vertex ...
134  auto previous = m_event->getVertex();
135  auto vtx = generateVertex(m_beamParams->getVertex(), m_beamParams->getCovVertex(), m_generateVertex);
136  m_event->setVertex(vtx);
137  // Everything done
138  return vtx - previous;
139  }
141 } //Belle2 namespace
Belle2::InitialParticleGeneration::generateVertex
TVector3 generateVertex(const TVector3 &initial, const TMatrixDSym &cov, MultivariateNormalGenerator &gen)
generate the vertex
Definition: InitialParticleGeneration.cc:34
Belle2::MultivariateNormalGenerator::generate
Eigen::VectorXd generate() const
Generate a set of correlated random numbers with the previouly set mean and covariance.
Definition: MultivariateNormalGenerator.h:89
Belle2::MCInitialParticles::c_smearBeamEnergy
@ c_smearBeamEnergy
smear energy of HER and LER (but not direction)
Definition: MCInitialParticles.h:44
Belle2::MCInitialParticles::c_smearBeam
@ c_smearBeam
smear the full beam momentum (energy and direction)
Definition: MCInitialParticles.h:48
Belle2::InitialParticleGeneration::updateVertex
TVector3 updateVertex(bool force=false)
Update the vertex position:
Definition: InitialParticleGeneration.cc:118
Belle2::InitialParticleGeneration::m_event
StoreObjPtr< MCInitialParticles > m_event
Datastore object containing the generated event.
Definition: InitialParticleGeneration.h:100
Belle2::MCInitialParticles::c_smearBeamDirection
@ c_smearBeamDirection
smear direction of HER and LER (but not energy)
Definition: MCInitialParticles.h:46
Belle2::InitialParticleGeneration::m_generateHER
MultivariateNormalGenerator m_generateHER
Generator for HER.
Definition: InitialParticleGeneration.h:102
Belle2::InitialParticleGeneration::m_generateVertex
MultivariateNormalGenerator m_generateVertex
Generator for Vertex.
Definition: InitialParticleGeneration.h:106
Belle2::InitialParticleGeneration::generateBeam
TLorentzVector generateBeam(const TLorentzVector &initial, const TMatrixDSym &cov, MultivariateNormalGenerator &gen)
generate 4 vector for one beam
Definition: InitialParticleGeneration.cc:44
Belle2::Const::electronMass
static const double electronMass
electron mass
Definition: Const.h:558
Belle2::MCInitialParticles::c_smearVertex
@ c_smearVertex
smear vertex
Definition: MCInitialParticles.h:50
Belle2::MultivariateNormalGenerator
Class to generate normal distributed, correlated random numbers given the mean values and the covaria...
Definition: MultivariateNormalGenerator.h:62
Belle2::MCInitialParticles::c_generateCMS
@ c_generateCMS
generate initial event in CMS instead of lab
Definition: MCInitialParticles.h:42
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::InitialParticleGeneration::m_generateLER
MultivariateNormalGenerator m_generateLER
Generator for LER.
Definition: InitialParticleGeneration.h:104
Belle2::MultivariateNormalGenerator::size
size_t size() const
Return the number of elements to be generated on generate()
Definition: MultivariateNormalGenerator.h:104
Belle2::InitialParticleGeneration::m_allowedFlags
int m_allowedFlags
Allowed generation flags.
Definition: InitialParticleGeneration.h:108
Belle2::InitialParticleGeneration::m_beamParams
DBObjPtr< BeamParameters > m_beamParams
Datastore object containing the nominal beam parameters.
Definition: InitialParticleGeneration.h:98
Belle2::InitialParticleGeneration::initialize
void initialize()
function to be executed on initialize()
Definition: InitialParticleGeneration.cc:29
Belle2::InitialParticleGeneration::generate
MCInitialParticles & generate()
Generate a new event.
Definition: InitialParticleGeneration.cc:84
Belle2::MultivariateNormalGenerator::setMeanCov
bool setMeanCov(int n, const double *mean, const double *cov)
set the mean and covariance for the distribution with array interface: mean and covariance are passed...
Definition: MultivariateNormalGenerator.cc:22
Belle2::MultivariateNormalGenerator::reset
void reset()
reset the generator setting the size to 0.
Definition: MultivariateNormalGenerator.cc:16