Belle II Software development
InitialParticleGeneration.h
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#pragma once
10
11#include <framework/datastore/StoreObjPtr.h>
12#include <framework/database/DBObjPtr.h>
13#include <framework/utilities/MultivariateNormalGenerator.h>
14#include <framework/dbobjects/BeamParameters.h>
15#include <framework/dataobjects/MCInitialParticles.h>
16
17#include <framework/utilities/ConditionalGaussGenerator.h>
18
19#include <TMatrixDSym.h>
20
21namespace Belle2 {
26
36 public:
41 explicit InitialParticleGeneration(int allowedFlags = 0): m_allowedFlags(allowedFlags | MCInitialParticles::c_generateCMS) {};
42
45
47 ROOT::Math::XYZVector getVertexConditional();
48
66 ROOT::Math::XYZVector updateVertex(bool force = false);
67
69 const BeamParameters& getBeamParameters() const { return *m_beamParams; }
70
72 void initialize();
73
75 void setAllowedFlags(int allowedFlags)
76 {
78 }
79
81 ConditionalGaussGenerator initConditionalGenerator(const ROOT::Math::PxPyPzEVector& pHER, const ROOT::Math::PxPyPzEVector& pLER,
82 const TMatrixDSym& covHER, const TMatrixDSym& covLER);
83
85 double getNominalEcms() { return m_beamParams->getMass(); }
86
88 double getNominalEcmsSpread() { return m_generateLorentzTransformation.getX0spread(); }
89
92
93 private:
94
99 MCInitialParticles& generate(int allowedFlags);
100
102 TMatrixDSym adjustCovMatrix(TMatrixDSym cov) const;
103
104
110 ROOT::Math::XYZVector generateVertex(const ROOT::Math::XYZVector& initial, const TMatrixDSym& cov,
111 MultivariateNormalGenerator& gen) const;
117 ROOT::Math::PxPyPzEVector generateBeam(const ROOT::Math::PxPyPzEVector& initial, const TMatrixDSym& cov,
118 MultivariateNormalGenerator& gen) const;
129
132
135 };
136
137
139} //Belle2 namespace
This class contains the nominal beam parameters and the parameters used for smearing of the primary v...
Class implementing n-dimensional random number generator from Gaussian distribution where the first c...
Class for accessing objects in the database.
Definition DBObjPtr.h:21
InitialParticleGeneration(int allowedFlags=0)
constructor
const BeamParameters & getBeamParameters() const
Return reference to nominal beam parameters.
DBObjPtr< BeamParameters > m_beamParams
Datastore object containing the nominal beam parameters.
double getNominalEcmsSpread()
Get spread of CMS collision energy calculated from beam parameters.
int m_allowedFlags
Allowed generation flags.
MultivariateNormalGenerator m_generateHER
Generator for HER.
StoreObjPtr< MCInitialParticles > m_event
Datastore object containing the generated event.
MultivariateNormalGenerator m_generateLER
Generator for LER.
MultivariateNormalGenerator m_generateVertex
Generator for Vertex.
ConditionalGaussGenerator m_generateLorentzTransformation
Generator of the Lorentz transformation.
double getNominalEcms()
Get the CMS energy of collisions.
void setAllowedFlags(int allowedFlags)
Set allowed flags.
const ConditionalGaussGenerator & getLorentzGenerator()
Get the generator for the Lorentz transformation.
This class contains the initial state for the given event.
@ c_generateCMS
generate initial event in CMS instead of lab
Class to generate normal distributed, correlated random numbers given the mean values and the covaria...
Type-safe access to single objects in the data store.
Definition StoreObjPtr.h:96
ROOT::Math::XYZVector getVertexConditional()
Generate vertex position and possibly update the generator of Lorentz transformation.
void initialize()
function to be executed on initialize()
ROOT::Math::XYZVector updateVertex(bool force=false)
Update the vertex position:
MCInitialParticles & generate()
Generate a new event.
ROOT::Math::XYZVector generateVertex(const ROOT::Math::XYZVector &initial, const TMatrixDSym &cov, MultivariateNormalGenerator &gen) const
generate the vertex
ConditionalGaussGenerator initConditionalGenerator(const ROOT::Math::PxPyPzEVector &pHER, const ROOT::Math::PxPyPzEVector &pLER, const TMatrixDSym &covHER, const TMatrixDSym &covLER)
Initialize the conditional generator using HER & LER 4-vectors and HER & LER covariance matrices desc...
TMatrixDSym adjustCovMatrix(TMatrixDSym cov) const
adjust smearing covariance matrix based on the generation flags
ROOT::Math::PxPyPzEVector generateBeam(const ROOT::Math::PxPyPzEVector &initial, const TMatrixDSym &cov, MultivariateNormalGenerator &gen) const
generate 4 vector for one beam
Abstract base class for different kinds of events.