Belle II Software  release-08-01-10
BkgSensitiveDetector.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/background/BkgSensitiveDetector.h>
10 #include <simulation/background/BkgNeutronWeight.h>
11 #include <simulation/dataobjects/BeamBackHit.h>
12 #include <framework/gearbox/Unit.h>
13 #include <framework/gearbox/Const.h>
14 
15 #include <G4Track.hh>
16 #include <G4Step.hh>
17 
18 #include <framework/datastore/StoreArray.h>
19 #include <framework/datastore/RelationArray.h>
20 #include <TVector3.h>
21 
22 using namespace std;
23 
24 namespace Belle2 {
31  BkgSensitiveDetector::BkgSensitiveDetector(const char* subDett, int iden):
32  Simulation::SensitiveDetectorBase("BKG", Const::invalidDetector), m_trackID(0), m_startPos(0., 0., 0.), m_startMom(0., 0., 0.),
33  m_startTime(0), m_startEnergy(0), m_energyDeposit(0), m_trackLength(0.)
34  {
35  // registration of store arrays and relations
36 
37  StoreArray<MCParticle> mcParticles;
38  StoreArray<BeamBackHit> beamBackHits;
39 
40  beamBackHits.registerInDataStore();
41  mcParticles.registerRelationTo(beamBackHits);
42 
43  // additional registration of MCParticle relation
44 
45  RelationArray relation(mcParticles, beamBackHits);
47 
48  m_eclrepscale = 0;
49  std::string subDet = subDett;
50  if (subDet == "IR") m_subDet = 0;
51  else if (subDet == "PXD") m_subDet = 1;
52  else if (subDet == "SVD") m_subDet = 2;
53  else if (subDet == "CDC") m_subDet = 3;
54  else if (subDet == "ARICH") m_subDet = 4;
55  else if (subDet == "TOP") m_subDet = 5;
56  else if (subDet == "ECL") m_subDet = 6;
57  else if (subDet == "EKLM") m_subDet = 7;
58  else if (subDet == "BKLM") m_subDet = 8;
59  else m_subDet = 99;
60 
61  if (m_subDet == 6) {
62  const int nc[] = {3, 4, 6, 9, 2, 9, 6, 4}; // the number of crystals in a sector per ring
63  const int indx[] = {96, 288, 864, 1151, 7776, 8064, 8544, 8736}; // regions with the same number of crystals in a ring
64  std::vector<int> indxv(indx, indx + sizeof(indx) / sizeof(int));
65  m_eclrepscale = nc[upper_bound(indxv.begin(), indxv.end(), iden) - indxv.begin()];
66  }
67 
68  m_identifier = iden;
69 
70  }
71 
72 
73  G4bool BkgSensitiveDetector::step(G4Step* aStep, G4TouchableHistory*)
74  {
75 
76  const G4StepPoint& preStep = *aStep->GetPreStepPoint();
77  G4Track& track = *aStep->GetTrack();
78 
79  if (m_trackID != track.GetTrackID()) {
80  //TrackID changed, store track informations
81  m_trackID = track.GetTrackID();
82  //Get world position
83  const G4ThreeVector& worldPosition = preStep.GetPosition();
84  m_startPos.SetXYZ(worldPosition.x() * Unit::mm / Unit::cm, worldPosition.y() * Unit::mm / Unit::cm,
85  worldPosition.z() * Unit::mm / Unit::cm);
86  //Get momentum
87  const G4ThreeVector& momentum = preStep.GetMomentum() ;
88  m_startMom.SetXYZ(momentum.x() * Unit::MeV, momentum.y() * Unit::MeV,
89  momentum.z() * Unit::MeV);
90  //Get time
91  m_startTime = preStep.GetGlobalTime();
92  //Get energy
93  m_startEnergy = preStep.GetKineticEnergy() * Unit::MeV;
94  //Reset energy deposit;
95  m_energyDeposit = 0;
96  //Reset track lenght;
97  m_trackLength = 0;
98  }
99  //Update energy deposit
100  m_energyDeposit += aStep->GetTotalEnergyDeposit() * Unit::MeV;
101 
102  m_trackLength += aStep->GetStepLength() * Unit::mm;
103 
104  //Save Hit if track leaves volume or is killed
105  if (track.GetNextVolume() != track.GetVolume() || track.GetTrackStatus() >= fStopAndKill) {
106  int pdgCode = track.GetDefinition()->GetPDGEncoding();
107  double endEnergy = track.GetKineticEnergy() * Unit::MeV;
108  double neutWeight = 0;
109  if (pdgCode == Const::neutron.getPDGCode()) {
111  neutWeight = wt.getWeight(m_startEnergy / Unit::MeV);
112  }
113  StoreArray<BeamBackHit> beamBackHits;
114  int nentr = beamBackHits.getEntries();
115  int id = m_identifier;
116  if (m_subDet == 6) { // ECL
117  int sector = preStep.GetTouchable()->GetCopyNumber(1);
118  if ((sector == 0) && (m_eclrepscale == 2) && ((id & 1) != 0)) id += 144; // barrel odd-numbered diodes in sector 0 -> sector 144
119  id += sector * m_eclrepscale;
120  }
121  beamBackHits.appendNew(m_subDet, id, pdgCode, m_trackID, m_startPos, m_startMom, m_startTime, m_startEnergy, endEnergy,
122  m_energyDeposit, m_trackLength, neutWeight);
123 
124 
125  // create relation to MCParticle
126  StoreArray<MCParticle> mcParticles;
127  RelationArray relBeamBackHitToMCParticle(mcParticles, beamBackHits);
128  relBeamBackHitToMCParticle.add(m_trackID, nentr);
129 
130  //Reset TrackID
131  m_trackID = 0;
132  }
133 
134  return true;
135  }
136 
137 
139 } // end of namespace Belle2
The class to get the weighting factor for a 1-MeV-equivalent neutron flux on silicon.
double m_energyDeposit
Energy deposited in volume.
ROOT::Math::XYZVector m_startPos
Particle position at the entrance in volume.
double m_startEnergy
Particle energy at the entrance in volume.
double m_trackLength
Length of the track in the volume.
int m_eclrepscale
Replica (=sector) scale in ECL.
int m_subDet
Subdetector id number.
ROOT::Math::XYZVector m_startMom
Particle momentum at the entrance in volume.
int m_identifier
Identifier of subdetector component.
This class provides a set of constants for the framework.
Definition: Const.h:34
static const ParticleType neutron
neutron particle
Definition: Const.h:666
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
static const double mm
[millimeters]
Definition: Unit.h:70
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
double getWeight(double ke)
Get weighting factor to convert a neutron to its 1-MeV equivalent.
static BkgNeutronWeight & getInstance()
Return a reference to the singleton BkgNeutronWeight instance.
bool step(G4Step *aStep, G4TouchableHistory *) override
Process each step and calculate variables defined in PXDSimHit.
Abstract base class for different kinds of events.