Belle II Software  release-06-02-00
MCParticleGenerator.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 <simulation/kernel/MCParticleGenerator.h>
10 
11 #include <mdst/dataobjects/MCParticle.h>
12 #include <framework/gearbox/Unit.h>
13 #include <framework/logging/Logger.h>
14 #include <framework/datastore/StoreArray.h>
15 #include <simulation/kernel/UserInfo.h>
16 #include <geometry/GeometryManager.h>
17 
18 #include <G4ParticleTable.hh>
19 #include <G4VUserPrimaryParticleInformation.hh>
20 #include <G4ParticleDefinition.hh>
21 #include <G4ParticleTypes.hh>
22 #include <G4VPhysicalVolume.hh>
23 #include <G4LogicalVolume.hh>
24 #include <G4VSolid.hh>
25 #include <CLHEP/Vector/LorentzVector.h>
26 
27 #include <TLorentzVector.h>
28 
29 using namespace std;
30 using namespace Belle2;
31 using namespace Belle2::Simulation;
32 
33 
34 MCParticleGenerator::MCParticleGenerator(const string& mcCollectionName, MCParticleGraph& mcParticleGraph) :
35  G4VPrimaryGenerator(), m_mcCollectionName(mcCollectionName), m_mcParticleGraph(mcParticleGraph)
36 {
37  G4VPhysicalVolume* topVolume = geometry::GeometryManager::getInstance().getTopVolume();
38  if (topVolume) {
39  m_topSolid = topVolume->GetLogicalVolume()->GetSolid();
40  }
41 }
42 
43 
45 {
46 
47 }
48 
49 
51 {
52  //Get the MCParticle collection from the DataStore
54 
55  //Prepare the MCParticle graph
57 
58  //Loop over the primary particles. The MCParticle collection has to be
59  //sorted breadth first: primary particles come first and then the daughters
60  // ... at least that would be the case if there is only one generator, but if
61  // more than one generator is used that is not necessarily true. So ignore
62  // that last statement.
63  //
64  // Let's add all top level particles, i.e. without mother, and addParticle
65  // will be called recursively for all daughters
66  int nPart = mcParticles.getEntries();
67  for (int iPart = 0; iPart < nPart; iPart++) {
68  MCParticle* currParticle = mcParticles[iPart];
69  if (currParticle->getMother() != NULL) continue;
70 
71  //Add primary particle (+ daughters) and the vertex
72  addParticle(*currParticle, event, NULL, 0, true);
73  }
74 }
75 
76 G4PrimaryVertex* MCParticleGenerator::determineVertex(const MCParticleGraph::GraphParticle& p, double& productionTimeShift)
77 {
78  // We want to determine the simulation vertex for the given particle.
79  // So lets see where it is
80  TVector3 mcProdVtx = p.getProductionVertex();
81  G4ThreeVector pos(mcProdVtx.X() / Unit::mm * CLHEP::mm,
82  mcProdVtx.Y() / Unit::mm * CLHEP::mm,
83  mcProdVtx.Z() / Unit::mm * CLHEP::mm);
84  productionTimeShift = 0;
85  // and make sure its inside the simulation volume
86  if (m_topSolid && m_topSolid->Inside(pos) == kOutside) {
87  // Ok, we're outside the simulation volume. Check if we can actually
88  // get to the simulation volume in time, no fields outside the volume
89  TLorentzVector mcPartMom4 = p.get4Vector();
90  CLHEP::HepLorentzVector mom(mcPartMom4.X() / Unit::MeV * CLHEP::MeV,
91  mcPartMom4.Y() / Unit::MeV * CLHEP::MeV,
92  mcPartMom4.Z() / Unit::MeV * CLHEP::MeV,
93  mcPartMom4.E() / Unit::MeV * CLHEP::MeV);
94  const auto dir = mom.vect().unit();
95  double distance = m_topSolid->DistanceToIn(pos, dir);
96 
97  if (distance == kInfinity) {
98  B2DEBUG(10, "Particle starts outside simulation volume and isn't intersecting, skipping");
99  return nullptr;
100  }
101 
102  const double speed = (p.getMass() == 0 ? mom.beta() : 1) * CLHEP::c_light;
103  const double flightLength = speed * p.getLifetime();
104  if (distance > flightLength) {
105  B2DEBUG(10, "Particle starts outside simulation volume and doesn't reach it before decaying, skipping");
106  return nullptr;
107  }
108  pos += dir * distance;
109  productionTimeShift = (distance / speed);
110  }
111  return new G4PrimaryVertex(pos.x(), pos.y(), pos.z(), p.getProductionTime() + productionTimeShift);
112 }
113 
115  G4Event* event,
116  G4PrimaryParticle* lastG4Mother,
117  int motherIndex,
118  bool useTime)
119 {
120  G4PrimaryParticle* g4Mother = lastG4Mother;
121 
122  //Check if the particle should be added to Geant4
123  //Only add particle if its PDG code is known to Geant4 and it is not flagged as virtual or initial.
124  bool addToG4 = (!mcParticle.isVirtual()) && (!mcParticle.isInitial());
125 
126  //MS: distinguish optical photon from the rest of particles
127  bool opticalPhoton = mcParticle.getPDG() == 0 && mcParticle.getEnergy() < 10.0 * Unit::eV;
128  //SD: add neutral long-lived particle
129  bool neutral_llp = mcParticle.getCharge() == 0.0 && mcParticle.getLifetime() > 0.0;
130 
131  G4ParticleDefinition* pdef = NULL;
132  if (opticalPhoton) {
133  pdef = G4OpticalPhoton::OpticalPhotonDefinition();
134  } else {
135  pdef = G4ParticleTable::GetParticleTable()->FindParticle(mcParticle.getPDG());
136  if (pdef == NULL && neutral_llp) {
137  pdef = G4ParticleTable::GetParticleTable()->FindParticle("LongLivedNeutralParticle");
138  }
139  }
140 
141  if (pdef == NULL) {
142  // RG, DK when the skipped particle has a long lived and it is not flagged as virtual
143  if (addToG4 && mcParticle.getLifetime() != 0.0)
144  B2WARNING("PDG code " << mcParticle.getPDG() << " unknown to Geant4 - particle skipped.");
145  addToG4 = false;
146  }
147 
148  //Add particle to MCParticle graph
150  graphParticle.setPDG(mcParticle.getPDG());
151  graphParticle.setStatus(mcParticle.getStatus());
152  graphParticle.setMass(mcParticle.getMass());
153  graphParticle.setEnergy(mcParticle.getEnergy());
154  graphParticle.setValidVertex(mcParticle.hasValidVertex());
155  graphParticle.setProductionTime(mcParticle.getProductionTime());
156  graphParticle.setDecayTime(mcParticle.getDecayTime());
157  graphParticle.setProductionVertex(mcParticle.getProductionVertex());
158  graphParticle.setMomentum(mcParticle.getMomentum());
159  graphParticle.setDecayVertex(mcParticle.getDecayVertex());
160  graphParticle.setFirstDaughter(mcParticle.getFirstDaughter());
161  graphParticle.setLastDaughter(mcParticle.getLastDaughter());
162  //Add decay to MCParticle graph
163  if (motherIndex > 0) graphParticle.comesFrom(m_mcParticleGraph[motherIndex - 1]);
164 
165  //Create a new Geant4 Primary particle and
166  //store the link to the GraphMCParticle object as user info.
167  G4PrimaryParticle* newPart = nullptr;
168  // And maybe we need to generate a new primary Geant4 vertex for this particle.
169  G4PrimaryVertex* vertex = nullptr;
170  // And in the rare case that the particle started outside the simulation volume
171  // we need an adjusted production time
172  double productionTimeShift{0};
173 
174  // Ok, we want to add the particle and it doesn't have a mother so lets see
175  // if we can find a proper vertex
176  if (addToG4 and lastG4Mother == nullptr) {
177  vertex = determineVertex(graphParticle, productionTimeShift);
178  // we don't have a mother and no vertex so skip this particle
179  if (not vertex) addToG4 = false;
180  }
181  // So now we have a particle to add and either a mother or a vertex.
182  // Let's create the particle in Geant4
183  if (addToG4) {
184  const TLorentzVector mcPartMom4 = graphParticle.get4Vector();
185  newPart = new G4PrimaryParticle(pdef,
186  mcPartMom4.X() / Unit::MeV * CLHEP::MeV,
187  mcPartMom4.Y() / Unit::MeV * CLHEP::MeV,
188  mcPartMom4.Z() / Unit::MeV * CLHEP::MeV,
189  mcPartMom4.E() / Unit::MeV * CLHEP::MeV);
190  newPart->SetMass(graphParticle.getMass() / Unit::MeV * CLHEP::MeV);
191  if (opticalPhoton) {
192  TVector3 polarization = graphParticle.getDecayVertex(); // temporary stored here
193  newPart->SetPolarization(polarization.X(), polarization.Y(), polarization.Z());
194  }
195  newPart->SetUserInformation(new ParticleInfo(graphParticle));
196 
197  //Set propagation time only if useTime is true, the MCparticle has a valid vertex and has children.
198  useTime &= graphParticle.hasValidVertex() && mcParticle.getFirstDaughter() > 0;
199  if (useTime) {
200  //ProperTime is in particle eigentime, so convert lab lifetime to eigentime
201  //correct for any possible shift between simulation vertex and generation vertex
202  double propertime = (graphParticle.getLifetime() - productionTimeShift) / mcPartMom4.Gamma();
203  newPart->SetProperTime(propertime);
204  }
205 
206  if (lastG4Mother) {
207  lastG4Mother->SetDaughter(newPart);
208  // This particle will be produced at the decay vertex of its mother so we **want** Geant4 to
209  // overwrite production time/vertex, so as convention having NaN here means we'll
210  // take over the values from Geant4 instead of the ones from the generator
211  graphParticle.setProductionTime(std::numeric_limits<float>::quiet_NaN());
212  } else {
213  vertex->SetPrimary(newPart);
214  event->AddPrimaryVertex(vertex);
215  }
216 
217  g4Mother = newPart;
218 
219  //Do not store the generator info in the MCParticles block unless Geant4 creates a track in the detector.
220  graphParticle.setIgnore();
221 
222  } else {
223  B2DEBUG(100, "The particle " << mcParticle.getIndex() << " (PDG " << mcParticle.getPDG() << ") was not added to Geant4");
224  }
225 
226  //Add all children
227  int currMotherIndex = m_mcParticleGraph.size();
228  for (MCParticle* daughter : mcParticle.getDaughters()) {
229  addParticle(*daughter, event, g4Mother, currMotherIndex, useTime);
230  }
231 }
232 
233 
234 //===================================================================
235 // Protected methods
236 //===================================================================
Class to represent Particle data in graph.
void comesFrom(GraphParticle &mother)
Tells the graph that this particle is a decay product of mother.
void setFirstDaughter(int daughter)
Set the 1-based index of the first daughter, 0 means no daughters.
void setLastDaughter(int daughter)
Set the 1-based index of the last daughter, 0 means no daughters.
void setIgnore(bool ignore=true)
Set or remove the ignore flag.
Class to build, validate and sort a particle decay chain.
size_t size() const
Return the number of particles in the graph.
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
void setMomentum(const TVector3 &momentum)
Set particle momentum.
Definition: MCParticle.h:414
float getEnergy() const
Return particle energy in GeV.
Definition: MCParticle.h:147
void setDecayTime(float time)
Set decay time.
Definition: MCParticle.h:387
void setMass(float mass)
Set particle mass.
Definition: MCParticle.h:363
int getIndex() const
Get 1-based index of the particle in the corresponding MCParticle list.
Definition: MCParticle.h:230
void setProductionVertex(const TVector3 &vertex)
Set production vertex position.
Definition: MCParticle.h:393
float getLifetime() const
Return the lifetime in ns.
Definition: MCParticle.h:177
void setEnergy(float energy)
Set energy.
Definition: MCParticle.h:369
void setDecayVertex(const TVector3 &vertex)
Set decay vertex.
Definition: MCParticle.h:441
std::vector< Belle2::MCParticle * > getDaughters() const
Get vector of all daughter particles, empty vector if none.
Definition: MCParticle.cc:50
float getMass() const
Return the particle mass in GeV.
Definition: MCParticle.h:135
float getDecayTime() const
Return the decay time in ns.
Definition: MCParticle.h:168
void setValidVertex(bool valid)
Set indication wether vertex and time information is valid or just default.
Definition: MCParticle.h:375
TVector3 getMomentum() const
Return momentum.
Definition: MCParticle.h:198
bool hasValidVertex() const
Indication whether vertex and time information is useful or just default.
Definition: MCParticle.h:153
unsigned int getStatus(unsigned short int bitmask=USHRT_MAX) const
Return status code of particle.
Definition: MCParticle.h:122
TVector3 getDecayVertex() const
Return decay vertex.
Definition: MCParticle.h:219
int getLastDaughter() const
Get 1-based index of last daughter, 0 if no daughters.
Definition: MCParticle.h:259
float getCharge() const
Return the particle charge defined in TDatabasePDG.
Definition: MCParticle.cc:34
void setPDG(int pdg)
Set PDG code of the particle.
Definition: MCParticle.h:332
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
float getProductionTime() const
Return production time in ns.
Definition: MCParticle.h:159
TLorentzVector get4Vector() const
Return 4Vector of particle.
Definition: MCParticle.h:207
TVector3 getProductionVertex() const
Return production vertex position.
Definition: MCParticle.h:189
void setStatus(unsigned short int status)
Set Status code for the particle.
Definition: MCParticle.h:343
int getFirstDaughter() const
Get 1-based index of first daughter, 0 if no daughters.
Definition: MCParticle.h:251
void setProductionTime(float time)
Set production time.
Definition: MCParticle.h:381
std::string m_mcCollectionName
Name of the MCParticle collection.
virtual ~MCParticleGenerator()
The destructor of the MCParticleGenerator class.
MCParticleGraph & m_mcParticleGraph
Reference to the MCParticle graph.
void addParticle(const MCParticle &mcParticle, G4Event *event, G4PrimaryParticle *lastG4Mother, int motherIndex, bool useTime)
Takes a MCParticle and creates a primary particle for Geant4.
virtual void GeneratePrimaryVertex(G4Event *event)
Create G4 primary particles from MCParticle list.
G4VSolid * m_topSolid
Pointer to the top level solid to check if particles are inside the simulation volume.
G4PrimaryVertex * determineVertex(const MCParticleGraph::GraphParticle &p, double &productionTimeShift)
Create a simulation vertex for the given particle.
UserInfo class which is used to attach additional information to Geant4 particles and tracks.
Definition: UserInfo.h:36
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
static const double mm
[millimeters]
Definition: Unit.h:70
static const double eV
[electronvolt]
Definition: Unit.h:112
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
static GeometryManager & getInstance()
Return a reference to the instance.
G4VPhysicalVolume * getTopVolume()
Return a pointer to the top volume.
bool isInitial() const
Check if particle is an initial particle such as ISR.
Definition: MCParticle.h:572
MCParticle * getMother() const
Returns a pointer to the mother particle.
Definition: MCParticle.h:582
void clear()
Reset particles and decay information to make the class reusable.
GraphParticle & addParticle()
Add new particle to the graph.
bool isVirtual() const
Check if particle is virtual.
Definition: MCParticle.h:557
Abstract base class for different kinds of events.