Belle II Software  release-06-02-00
GeneratedVertexDisplacerModule.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 <generators/modules/GeneratedVertexDisplacerModule.h>
10 #include <framework/core/ModuleParam.templateDetails.h>
11 
12 
13 
14 using namespace Belle2;
15 using namespace std;
16 
17 //-----------------------------------------------------------------
18 // Register the Module
19 //-----------------------------------------------------------------
20 REG_MODULE(GeneratedVertexDisplacer)
21 
22 //-----------------------------------------------------------------
23 // Implementation
24 //-----------------------------------------------------------------
25 
27 {
28  // Set module properties
29  setDescription(
30  R"DOC(""Takes a list of PDG values and lifetime paramters to displaces the vertex of MCParticles with matching PDG value corresponding to the given lifetime parameter. Can be used betwenerator and the detector simulation.)DOC");
31 
32  // Parameter definitions
33  addParam("lifetimeOption", m_lifetimeOption,
34  "Set the lifetime option, either 'fixed', 'flat' or 'exponential'. Option is set globally.", string("fixed"));
35  addParam("lifetime", m_lifetime,
36  "Set the numerical value of the lifetime c*tau in units of [cm]. Example: (lifetime_1, lifetime_2, ... ,lifetime_N). Lifetime_X is set for PDG_X.",
37  std::vector<float> {0});
38  addParam("pdgVal", m_pdgVals,
39  "PDG values of MCParticles that should be displaced. Subsequent daughters of these will be displaced accordingly. Example: (PDG_1, PDG_2, ... , PDG_N). PDG_X corresponds to lifetime_X.",
40  std::vector<int> { -999});
41  addParam("maxDecayTime", m_maxDecayTime,
42  "Set the maximal allowed decay time in c*tau [cm] for the option 'flat'. The 'exponential' option is not bound. Default 300cm.",
43  static_cast<float>(300));
44  addParam("ctau", m_ctau,
45  "Input unit option. True (default): module expects lifetime as ctau [cm]. False: lifetime as tau [ns].",
46  true);
47 
48 }
49 
50 
52 {
53 
54  B2DEBUG(0, "Initialize GeneratedVertexDisplacerModule");
55  m_mcparticles.isRequired(m_particleList);
56 
57  if (m_pdgVals.size() != m_lifetime.size()) {
58  B2FATAL("List of PDG values and lifetime parameters have different sizes. Specify a lifetime for each PDG value.");
59  }
60 
61  if ((m_lifetimeOption.compare("fixed") != 0) && (m_lifetimeOption.compare("flat") != 0)
62  && (m_lifetimeOption.compare("exponential") != 0)) {
63  B2FATAL("Lifetime option must be 0 (fixed), 1 (flat) or 2 (exponential).");
64  }
65 }
66 
67 
69 {
70 
71  if (m_mcparticles.getEntries() < 1) {
72  B2WARNING("MC particle container empty. Calling next event.");
73  return;
74  }
75 
76  for (int mcp_index = 0; mcp_index < m_mcparticles.getEntries(); mcp_index++) {
77 
78  MCParticle& mcp = *m_mcparticles[mcp_index];
79 
81 
82  int mcp_pdg = mcp.getPDG();
83 
84  for (unsigned int param_index = 0; param_index < m_pdgVals.size(); param_index++) {
85  if (m_pdgVals.at(param_index) == mcp_pdg) {
86  B2DEBUG(0, "Displacing Vertex of particle with pdg_id: " << mcp_pdg << " with lifetime: " << m_lifetime.at(
87  param_index) << " and option: " << m_lifetimeOption);
88  displace(mcp, m_lifetime.at(param_index));
89  B2DEBUG(0, "new vertex (x,y,z) at:" << "(" << mcp.getDecayVertex().X() << "," << mcp.getDecayVertex().Y() << "," <<
90  mcp.getDecayVertex().Z() << ")" << " with decaylength=" << std::sqrt(std::pow(mcp.getDecayVertex().X(),
91  2) + std::pow(mcp.getDecayVertex().Y(), 2) + std::pow(mcp.getDecayVertex().Z(), 2))) ;
92  }
93  }
94  }
95 }
96 
97 
99 {
100  TLorentzVector* displacementVector = new TLorentzVector();
101  getDisplacement(particle, lifetime, *displacementVector);
102 
103  TVector3 newDecayVertex = particle.getProductionVertex();
104  newDecayVertex.SetX(newDecayVertex.X() + displacementVector->X());
105  newDecayVertex.SetY(newDecayVertex.Y() + displacementVector->Y());
106  newDecayVertex.SetZ(newDecayVertex.Z() + displacementVector->Z());
107 
108  particle.setDecayVertex(newDecayVertex);
109  particle.setDecayTime(particle.getProductionTime() + displacementVector->T());
110  particle.setValidVertex(true);
111 
112  // displace direct daughters
113  if (particle.getNDaughters()) {
114  displaceDaughter(*displacementVector, particle.getDaughters());
115  }
116 
117  delete displacementVector;
118 }
119 
120 
121 void GeneratedVertexDisplacerModule::displaceDaughter(TLorentzVector& motherDisplacementVector, std::vector<MCParticle*> daughters)
122 {
123  for (unsigned int daughter_index = 0; daughter_index < daughters.size(); daughter_index++) {
124 
125  MCParticle* daughter_mcp = daughters.at(daughter_index);
126  int daughter_mcpArrayIndex = daughter_mcp->getArrayIndex();
127 
128  // getDaughters returns a copied list, need to change parameters of the original particle in the MCParticle StoreArray.
129  MCParticle& mcp = *m_mcparticles[daughter_mcpArrayIndex];
130  mcp.setProductionVertex(mcp.getProductionVertex() + motherDisplacementVector.Vect());
131  mcp.setProductionTime(mcp.getProductionTime() + motherDisplacementVector.T());
132  mcp.setValidVertex(true);
133 
134  // Displace subsequent daughters
135  if (daughter_mcp->getNDaughters()) {
136  displaceDaughter(motherDisplacementVector, daughter_mcp->getDaughters());
137  }
138  }
139 }
140 
141 
143  const MCParticle& particle, float lifetime, TLorentzVector& displacement)
144 {
145  TLorentzVector fourVector_mcp = particle.get4Vector();
146  float decayLength_mcp = 0;
147  if (!m_ctau) lifetime *= Const::speedOfLight;
148 
149  if (m_lifetimeOption.compare("fixed") == 0) decayLength_mcp = lifetime;
150  else if (m_lifetimeOption.compare("flat") == 0) decayLength_mcp = gRandom->Uniform(0, m_maxDecayTime);
151  else {
152  decayLength_mcp = -1 * std::log(gRandom->Uniform(0, 1.)) * lifetime * fourVector_mcp.Gamma() * fourVector_mcp.Beta();
153 
154  if (!particle.getMass()) {
155  B2WARNING("Displacing a particle with zero mass. Forcing Gamma=1 for the decay time.");
156  decayLength_mcp = -1 * std::log(gRandom->Uniform(0, 1.)) * lifetime;
157  }
158  }
159 
160  float pMag = fourVector_mcp.P();
161  displacement.SetX(decayLength_mcp * fourVector_mcp.X() / pMag);
162  displacement.SetY(decayLength_mcp * fourVector_mcp.Y() / pMag);
163  displacement.SetZ(decayLength_mcp * fourVector_mcp.Z() / pMag);
164  displacement.SetT(decayLength_mcp / (Const::speedOfLight * fourVector_mcp.Beta()));
165 }
166 
167 
169 {
170 
171 }
172 
173 
static const double speedOfLight
[cm/ns]
Definition: Const.h:575
"Takes a list of PDG values and lifetime paramters to displaces the vertex of MCParticles with matchi...
void getDisplacement(const MCParticle &particle, float lifetime, TLorentzVector &displacement)
Helper function to calculate the numerical value of the vertex displacement (x,y,z,...
void displaceDaughter(TLorentzVector &motherDisplacementVector, std::vector< MCParticle * > daughters)
Helper function to loop over subsequent daughters and displaces their vertices corresponding to their...
virtual void initialize() override
Register input and output data, initialises the module.
virtual void event() override
Method is called for each event.
virtual void terminate() override
Terminates the module.
void displace(MCParticle &particle, float lifetime)
Helper function to displace the mother particles (corresponding to given pdf values).
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
@ c_Initial
bit 5: Particle is initial such as e+ or e- and not going to Geant4
Definition: MCParticle.h:57
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
void setProductionVertex(const TVector3 &vertex)
Set production vertex position.
Definition: MCParticle.h:393
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
Definition: MCParticle.cc:50
bool hasStatus(unsigned short int bitmask) const
Return if specific status bit is set.
Definition: MCParticle.h:129
int getArrayIndex() const
Get 0-based index of the particle in the corresponding MCParticle list.
Definition: MCParticle.h:244
void setValidVertex(bool valid)
Set indication wether vertex and time information is valid or just default.
Definition: MCParticle.h:375
int getNDaughters() const
Return number of daughter MCParticles.
Definition: MCParticle.cc:65
TVector3 getDecayVertex() const
Return decay vertex.
Definition: MCParticle.h:219
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
float getProductionTime() const
Return production time in ns.
Definition: MCParticle.h:159
TVector3 getProductionVertex() const
Return production vertex position.
Definition: MCParticle.h:189
void setProductionTime(float time)
Set production time.
Definition: MCParticle.h:381
Base class for Modules.
Definition: Module.h:72
#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.