Belle II Software development
ARICHBackgroundModule.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// Own header.
9
10#include <arich/modules/arichBackground/ARICHBackgroundModule.h>
11
12#include <time.h>
13
14// Hit classes
15
16#include <simulation/dataobjects/BeamBackHit.h>
17#include <arich/dataobjects/ARICHSimHit.h>
18#include <arich/dataobjects/ARICHDigit.h>
19
20#include <mdst/dataobjects/MCParticle.h>
21
22// framework - DataStore
23#include <framework/datastore/RelationIndex.h>
24#include <framework/datastore/RelationArray.h>
25
26// framework aux
27#include <framework/logging/Logger.h>
28#include <framework/gearbox/Const.h>
29#include <framework/geometry/B2Vector3.h>
30
31using namespace boost;
32
33namespace Belle2 {
38 namespace arich {
39 //-----------------------------------------------------------------
41 //-----------------------------------------------------------------
42
43 REG_MODULE(ARICHBackground);
44
45
46 //-----------------------------------------------------------------
47 // Implementation
48 //-----------------------------------------------------------------
49
50 ARICHBackgroundModule::ARICHBackgroundModule() : Module(), m_phpos(TVector3()), m_phmom(TVector3()),
51 m_phVtx(TVector3()), m_phMmom(TVector3()), m_phMvtx(TVector3()), m_phPvtx(TVector3()), m_phPmom(TVector3()),
52 m_phGMvtx(TVector3()), m_phGMmom(TVector3()), m_modOrig(TVector3()), m_source(0), m_phPDG(0), m_phMPDG(0), m_phPPDG(0),
53 m_phGMPDG(0), m_type(0), m_edep(0.0), m_ttime(0.0), m_moduleID(0), m_phnw(0.0), m_trackLength(0.0), m_energy(0.0),
54 m_outputFile(NULL),
55 m_outputTree(NULL)
56 {
57 // Set description()
58 setDescription("ARICHBackground module. Used to extract information relevant for ARICH background from background files");
59
60 // Add parameters
61 addParam("FileName", m_filename, "output file name", std::string("mytree.root"));
62 addParam("BkgTag", m_bkgTag, "background type tag (appended to hits in the output tree", 0);
63
64 }
65
67 {
68
69 }
70
72 {
73 m_outputFile = new TFile(m_filename.c_str(), "RECREATE");
74
75 m_outputTree = new TTree("m_outputTree", "tree of arich background hits");
76 m_outputTree->Branch("source", &m_source, "m_source/I");
77 m_outputTree->Branch("phpos", "TVector3", &m_phpos);
78 m_outputTree->Branch("phmom", "TVector3", &m_phmom);
79 m_outputTree->Branch("phVtx", "TVector3", &m_phVtx);
80 m_outputTree->Branch("phMmom", "TVector3", &m_phMmom);
81 m_outputTree->Branch("phMvtx", "TVector3", &m_phMvtx);
82 m_outputTree->Branch("phPvtx", "TVector3", &m_phPvtx);
83 m_outputTree->Branch("phPmom", "TVector3", &m_phPmom);
84 m_outputTree->Branch("phGMvtx", "TVector3", &m_phGMvtx);
85 m_outputTree->Branch("phGMmom", "TVector3", &m_phGMmom);
86 m_outputTree->Branch("modOrig", "TVector3", &m_modOrig);
87 m_outputTree->Branch("phPDG", &m_phPDG, "m_phPDG/I");
88 m_outputTree->Branch("phMPDG", &m_phMPDG, "m_phMPDG/I");
89 m_outputTree->Branch("phPPDG", &m_phPPDG, "m_phPPDG/I");
90 m_outputTree->Branch("phGMPDG", &m_phGMPDG, "m_phGMPDG/I");
91 m_outputTree->Branch("type", &m_type, "m_type/I");
92 m_outputTree->Branch("edep", &m_edep, "m_edep/D");
93 m_outputTree->Branch("ttime", &m_ttime, "m_ttime/D");
94 m_outputTree->Branch("moduleID", &m_moduleID, "m_moduleID/I");
95 m_outputTree->Branch("phnw", &m_phnw, "m_phnw/D");
96 m_outputTree->Branch("trackLength", &m_trackLength, "m_trackLength/D");
97 m_outputTree->Branch("energy", &m_energy, "m_energy/D");
98
100 }
101
103 {
104 // Print run number
105 B2INFO("ARICHBackground: Processing. ");
106
107 }
108
110 {
111 for (const BeamBackHit& arichBeamBkgHit : m_BeamBackHits) {
112
113 int subdet = arichBeamBkgHit.getSubDet();
114 if (subdet != 4) continue;
115 m_type = 0;
116 if (arichBeamBkgHit.getIdentifier() == 1) m_type = 1;
117 m_edep = arichBeamBkgHit.getEnergyDeposit();
118 m_ttime = arichBeamBkgHit.getTime();
119 m_phPDG = arichBeamBkgHit.getPDG();
120 m_phpos = arichBeamBkgHit.getPosition();
121 m_phmom = arichBeamBkgHit.getMomentum();
122 m_moduleID = m_arichgp->getDetectorPlane().pointSlotID(m_phpos.X(), m_phpos.Y());
123 double r = m_arichgp->getDetectorPlane().getSlotR(m_moduleID);
124 double phi = m_arichgp->getDetectorPlane().getSlotPhi(m_moduleID);
125 m_modOrig = TVector3(r * cos(phi), r * sin(phi), 0);
126 m_energy = arichBeamBkgHit.getEnergy();
127
128 if (m_phPDG == Const::neutron.getPDGCode()) {
129 m_phnw = arichBeamBkgHit.getNeutronWeight();
130 m_trackLength = arichBeamBkgHit.getTrackLength();
131 }
132 m_phVtx = TVector3(0, 0, 0); m_phMvtx = TVector3(0, 0, 0); m_phMmom = TVector3(0, 0, 0);
133 m_phMPDG = -1; m_phPPDG = -1; m_phGMPDG = -1;
134 m_phPvtx = TVector3(0, 0, 0); m_phPmom = TVector3(0, 0, 0); m_phGMvtx = TVector3(0, 0, 0); m_phGMmom = TVector3(0, 0, 0);
135
137 if (relBeamBackHitToMCParticle.getFirstElementTo(arichBeamBkgHit)) {
138 const MCParticle* currParticle = relBeamBackHitToMCParticle.getFirstElementTo(arichBeamBkgHit)->from;
139 m_phVtx = B2Vector3D(currParticle->getVertex());
140 const MCParticle* mother = currParticle->getMother();
141 int mm = 0;
142 while (mother) {
143 if (mm == 0) {
144 m_phMPDG = mother->getPDG();
145 m_phMvtx = B2Vector3D(mother->getVertex());
146 m_phMmom = B2Vector3D(mother->getMomentum());
147 }
148 if (mm == 1) {
149 m_phGMPDG = mother->getPDG();
150 m_phGMvtx = B2Vector3D(mother->getVertex());
151 m_phGMmom = B2Vector3D(mother->getMomentum());
152 }
153 const MCParticle* pommother = mother->getMother();
154 if (!pommother) {
155 m_phPPDG = mother->getPDG();
156 m_phPvtx = B2Vector3D(mother->getVertex());
157 m_phPmom = B2Vector3D(mother->getMomentum());
158 }
159 mother = pommother;
160 mm++;
161 }
162
163 }
164 m_outputTree->Fill();
165 }
166
167 // these assignments are the same in every iteration of the loop
168 // so just assign them once outside the loop
169 m_type = 2;
170 m_phPDG = 0;
171 m_edep = 0;
172 for (const ARICHSimHit& arichsimhit : m_ARICHSimHits) {
173 m_moduleID = arichsimhit.getModuleID();
174 m_ttime = arichsimhit.getGlobalTime();
175
176 m_phVtx = TVector3(0, 0, 0); m_phMvtx = TVector3(0, 0, 0); m_phMmom = TVector3(0, 0, 0);
177 m_phMPDG = -1; m_phPPDG = -1; m_phGMPDG = -1; m_phmom = TVector3(0, 0, 0);
178 m_phPvtx = TVector3(0, 0, 0); m_phPmom = TVector3(0, 0, 0); m_phGMvtx = TVector3(0, 0, 0); m_phGMmom = TVector3(0, 0, 0);
180 if (relSimHitToMCParticle.getFirstElementTo(arichsimhit)) {
181 const MCParticle* currParticle = relSimHitToMCParticle.getFirstElementTo(arichsimhit)->from;
182 m_phVtx = B2Vector3D(currParticle->getVertex());
183 const MCParticle* mother = currParticle->getMother();
184 int mm = 0;
185 while (mother) {
186 if (mm == 0) {
187 m_phMPDG = mother->getPDG();
188 m_phMvtx = B2Vector3D(mother->getVertex());
189 m_phMmom = B2Vector3D(mother->getMomentum());
190 }
191 if (mm == 1) {
192 m_phGMPDG = mother->getPDG();
193 m_phGMvtx = B2Vector3D(mother->getVertex());
194 m_phGMmom = B2Vector3D(mother->getMomentum());
195 }
196 const MCParticle* pommother = mother->getMother();
197 if (!pommother) {
198 m_phPPDG = mother->getPDG();
199 m_phPvtx = B2Vector3D(mother->getVertex());
200 m_phPmom = B2Vector3D(mother->getMomentum());
201 }
202 mother = pommother;
203 mm++;
204 }
205
206 }
207 m_outputTree->Fill();
208 }
209 }
210
212 {
213 // Announce
214 B2INFO("ARICHBackground finished.");
215
216 //m_outputTree->Write();
217 m_outputFile->Write();
218 m_outputFile->Close();
219 }
220
221 } // end arich namespace
223} // end Belle2 namespace
Class ARICHSimHit - Geant4 simulated hit for ARICH.
Definition: ARICHSimHit.h:29
Class BeamBackHit - Stores hits from beam backgound simulation.
Definition: BeamBackHit.h:28
static const ParticleType neutron
neutron particle
Definition: Const.h:675
A Class to store the Monte Carlo particle information.
Definition: MCParticle.h:32
ROOT::Math::XYZVector getVertex() const
Return production vertex position, shorthand for getProductionVertex().
Definition: MCParticle.h:183
int getPDG() const
Return PDG code of particle.
Definition: MCParticle.h:112
ROOT::Math::XYZVector getMomentum() const
Return momentum.
Definition: MCParticle.h:198
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
Provides access to fast ( O(log n) ) bi-directional lookups on a specified relation.
Definition: RelationIndex.h:76
const Element * getFirstElementTo(const TO &to) const
Return a pointer to the first relation Element of the given object.
virtual void event()
Event processor.
ROOT::Math::XYZVector m_phPvtx
primary particle vertex
StoreArray< BeamBackHit > m_BeamBackHits
StoreArray for BeamBackHits.
virtual void initialize()
Initialize the Module.
virtual void beginRun()
Called when entering a new run.
int m_source
hit source (RBB_HER, ...)
DBObjPtr< ARICHGeometryConfig > m_arichgp
Geometry parametrization.
int m_phMPDG
hit particle mother PDG code
ROOT::Math::XYZVector m_phVtx
hit particle vertex position
ROOT::Math::XYZVector m_phPmom
primary particle momentum
StoreArray< ARICHSimHit > m_ARICHSimHits
StoreArray for ARICHSimHits.
int m_phGMPDG
hit particle grand mother PDG code
ROOT::Math::XYZVector m_phmom
hit momentum
virtual void terminate()
Termination action.
ROOT::Math::XYZVector m_phGMvtx
hit particle grand mother vertex
ROOT::Math::XYZVector m_phpos
hit position
int m_type
hit particle type; 0 hit in board, 1 hit in HAPD bottom, 2 photon hit
int m_phPPDG
hit particle primary PDG code
double m_trackLength
particle track lenght in hit volume
ROOT::Math::XYZVector m_modOrig
HAPD module position.
std::string m_filename
Output file name.
ROOT::Math::XYZVector m_phGMmom
hit particle grand mother momentum
StoreArray< MCParticle > m_MCParticles
StoreArray for MCParticles.
ROOT::Math::XYZVector m_phMvtx
hit particle mother vertex
ROOT::Math::XYZVector m_phMmom
hit particle mother momentum
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
B2Vector3< double > B2Vector3D
typedef for common usage with double
Definition: B2Vector3.h:516
MCParticle * getMother() const
Returns a pointer to the mother particle.
Definition: MCParticle.h:600
Abstract base class for different kinds of events.
const FROM * from
pointer of the element from which the relation points.