Belle II Software development
BeamBkgNeutronModule.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 <background/modules/BeamBkgNeutron/BeamBkgNeutronModule.h>
10
11#include <time.h>
12
13// Hit classes
14
15#include <simulation/dataobjects/BeamBackHit.h>
16#include <mdst/dataobjects/MCParticle.h>
17#include <pxd/dataobjects/PXDSimHit.h>
18#include <svd/dataobjects/SVDSimHit.h>
19#include <cdc/dataobjects/CDCSimHit.h>
20#include <arich/dataobjects/ARICHSimHit.h>
21#include <top/dataobjects/TOPSimHit.h>
22#include <ecl/dataobjects/ECLSimHit.h>
23#include <klm/dataobjects/KLMSimHit.h>
24#include <simulation/dataobjects/MCParticleTrajectory.h>
25
26// framework - DataStore
27#include <framework/datastore/StoreArray.h>
28#include <framework/datastore/RelationIndex.h>
29#include <framework/datastore/RelationArray.h>
30
31// framework aux
32#include <framework/logging/Logger.h>
33#include <framework/gearbox/Const.h>
34#include <framework/geometry/B2Vector3.h>
35
36using namespace std;
37using namespace boost;
38
39namespace Belle2 {
44 //-----------------------------------------------------------------
46 //-----------------------------------------------------------------
47
48 REG_MODULE(BeamBkgNeutron);
49
50
51 //-----------------------------------------------------------------
52 // Implementation
53 //-----------------------------------------------------------------
54
56 {
57 // Set description()
58 setDescription("BeamBkgNeutronModule module. Used to extract information relevant for the neutron background from background files");
59
60 // Add parameters
61 addParam("FileName", m_filename, "output file name", string("mytree.root"));
62 }
63
65 {
66 // create/open the file
67 m_ff = new TFile(m_filename.c_str(), "RECREATE");
68
69 // create trees
70 m_tree1 = new TTree("tree1", "BeamBackHit data");
71 m_tree2 = new TTree("tree2", "SimHits data");
72
73 // create branches
74 m_tree1->Branch("event", &m_iEvent, "event/I");
75 m_tree1->Branch("subDet", &m_subDet, "subDet/I");
76 m_tree1->Branch("iden", &m_iden, "iden/I");
77 m_tree1->Branch("PDG", &m_PDG, "PDG/I");
78 m_tree1->Branch("trackID", &m_trackID, "trackID/I");
79 m_tree1->Branch("momentumX", &m_momentumX, "momentumX/F");
80 m_tree1->Branch("momentumY", &m_momentumY, "momentumY/F");
81 m_tree1->Branch("momentumZ", &m_momentumZ, "momentumZ/F");
82 m_tree1->Branch("positionX", &m_positionX, "positionX/F");
83 m_tree1->Branch("positionY", &m_positionY, "positionY/F");
84 m_tree1->Branch("positionZ", &m_positionZ, "positionZ/F");
85 m_tree1->Branch("E_start", &m_E_start, "E_start/F");
86 m_tree1->Branch("E_end", &m_E_end, "E_end/F");
87 m_tree1->Branch("eDep", &m_eDep, "eDep/F");
88 m_tree1->Branch("trackLength", &m_trackLength, "trackLength/F");
89 m_tree1->Branch("nWeight", &m_nWeight, "nWeight/F");
90 m_tree1->Branch("E_init", &m_E_init, "E_init/F");
91 m_tree1->Branch("mass", &m_mass, "mass/F");
92 m_tree1->Branch("lifeTime", &m_lifeTime, "lifeTime/F");
93 m_tree1->Branch("vtxProdX", &m_vtxProdX, "vtxProdX/F");
94 m_tree1->Branch("vtxProdY", &m_vtxProdY, "vtxProdY/F");
95 m_tree1->Branch("vtxProdZ", &m_vtxProdZ, "vtxProdZ/F");
96 m_tree1->Branch("trj_x", &m_trj_x);
97 m_tree1->Branch("trj_y", &m_trj_y);
98 m_tree1->Branch("trj_z", &m_trj_z);
99 m_tree1->Branch("trj_px", &m_trj_px);
100 m_tree1->Branch("trj_py", &m_trj_py);
101 m_tree1->Branch("trj_pz", &m_trj_pz);
102
103 m_tree2->Branch("nSimHits", m_nSimHits, "nSimHits[13]/I");
104 m_tree2->Branch("hitPDG", m_hitPDG, "hitPDG[13]/I");
105 m_tree2->Branch("momPDG", m_momPDG, "momPDG[13]/I");
106// std::string xyzvector_typedef = "std::vector<ROOT::Math::DisplacementVector3D<ROOT::Math::Cartesian3D<double>,ROOT::Math::DefaultCoordinateSystemTag>>";
107// tt->Branch("vtxProd", "xyzvector_typedef.c_str()", &m_vtxProd);
108 }
109
111 {
112 // Print run number
113 B2INFO("BeamBkgNeutron: Processing. ");
114 }
115
117 {
118 // MCParticles
119 StoreArray<MCParticle> McParticles;
120
121 // SimHits
122 StoreArray<PXDSimHit> PXDSimHits;
123 StoreArray<SVDSimHit> SVDSimHits;
124 StoreArray<CDCSimHit> CDCSimHits;
125 StoreArray<ARICHSimHit> ARICHSimHits;
126 StoreArray<TOPSimHit> TOPSimHits;
127 StoreArray<ECLSimHit> ECLSimHits;
128 StoreArray<KLMSimHit> KLMSimHits;
129
130 for (Int_t i = 0; i < 13; i++) {
131 m_nSimHits[i] = 0;
132 m_hitPDG[i] = 0;
133 m_momPDG[i] = 0;
134 }
135
136 m_nSimHits[1] = PXDSimHits.getEntries();
137 m_nSimHits[2] = SVDSimHits.getEntries();
138 m_nSimHits[3] = CDCSimHits.getEntries();
139 m_nSimHits[4] = ARICHSimHits.getEntries();
140 m_nSimHits[5] = TOPSimHits.getEntries();
141 m_nSimHits[6] = ECLSimHits.getEntries();
142 m_nSimHits[7] = KLMSimHits.getEntries();
143
144 // loop over KLM simHits
145 for (Int_t hit = 0; hit < m_nSimHits[7]; hit++) {
146 // get KLMSimHit
147// EKLMSimHit* simHit = EKLMSimHits[hit];
148 const KLMSimHit* simHit = KLMSimHits[hit];
149 Float_t posZ = simHit->getPositionZ();
150 if (280.0 < posZ && posZ < 288.0) m_nSimHits[9]++; // FWD EKLM innermost layer
151 else if (400.0 < posZ && posZ < 406.0) m_nSimHits[10]++; // FWD EKLM outermost layer
152 else if (-194.0 < posZ && posZ < -186.0) m_nSimHits[11]++; // BWD EKLM innermost layer
153 else if (-294.0 < posZ && posZ < -286.0) m_nSimHits[12]++; // BWD EKLM outermost layer
154 }
155
156 RelationIndex<MCParticle, PXDSimHit> relPXDSimHitToMCParticle(McParticles, PXDSimHits);
157 RelationIndex<MCParticle, SVDSimHit> relSVDSimHitToMCParticle(McParticles, SVDSimHits);
158 RelationIndex<MCParticle, CDCSimHit> relCDCSimHitToMCParticle(McParticles, CDCSimHits);
159 RelationIndex<MCParticle, ARICHSimHit> relARICHSimHitToMCParticle(McParticles, ARICHSimHits);
160 RelationIndex<MCParticle, TOPSimHit> relTOPSimHitToMCParticle(McParticles, TOPSimHits);
161 RelationIndex<MCParticle, ECLSimHit> relECLSimHitToMCParticle(McParticles, ECLSimHits);
162 RelationIndex<MCParticle, KLMSimHit> relKLMSimHitToMCParticle(McParticles, KLMSimHits);
163
164 Int_t detID;
165
166 //--- PXD
167 detID = 1;
168 // loop over simhits
169 for (Int_t iHit = 0; iHit < m_nSimHits[detID]; iHit++) {
170 const PXDSimHit* simHit = PXDSimHits[iHit];
171 // get related MCparticle
172 if (relPXDSimHitToMCParticle.getFirstElementTo(simHit)) {
173 const MCParticle* currParticle = relPXDSimHitToMCParticle.getFirstElementTo(simHit)->from;
174 m_hitPDG[detID] = currParticle->getPDG();
175 if (!currParticle->isPrimaryParticle()) {
176 const MCParticle* momParticle = currParticle->getMother();
177 m_momPDG[detID] = momParticle->getPDG();
178 }
179 }
180 }
181
182 //--- SVD
183 detID = 2;
184 // loop over simhits
185 for (Int_t iHit = 0; iHit < m_nSimHits[detID]; iHit++) {
186 const SVDSimHit* simHit = SVDSimHits[iHit];
187 // get related MCparticle
188 if (relSVDSimHitToMCParticle.getFirstElementTo(simHit)) {
189 const MCParticle* currParticle = relSVDSimHitToMCParticle.getFirstElementTo(simHit)->from;
190 m_hitPDG[detID] = currParticle->getPDG();
191 if (!currParticle->isPrimaryParticle()) {
192 const MCParticle* momParticle = currParticle->getMother();
193 m_momPDG[detID] = momParticle->getPDG();
194 }
195 }
196 }
197
198 //--- CDC
199 detID = 3;
200 // loop over simhits
201 for (Int_t iHit = 0; iHit < m_nSimHits[detID]; iHit++) {
202 const CDCSimHit* simHit = CDCSimHits[iHit];
203 // get related MCparticle
204 if (relCDCSimHitToMCParticle.getFirstElementTo(simHit)) {
205 const MCParticle* currParticle = relCDCSimHitToMCParticle.getFirstElementTo(simHit)->from;
206 m_hitPDG[detID] = currParticle->getPDG();
207 if (!currParticle->isPrimaryParticle()) {
208 const MCParticle* momParticle = currParticle->getMother();
209 m_momPDG[detID] = momParticle->getPDG();
210 }
211 }
212 }
213
214 //--- ARICH
215 detID = 4;
216 // loop over simhits
217 for (Int_t iHit = 0; iHit < m_nSimHits[detID]; iHit++) {
218 const ARICHSimHit* simHit = ARICHSimHits[iHit];
219 // get related MCparticle
220 if (relARICHSimHitToMCParticle.getFirstElementTo(simHit)) {
221 const MCParticle* currParticle = relARICHSimHitToMCParticle.getFirstElementTo(simHit)->from;
222 m_hitPDG[detID] = currParticle->getPDG();
223 if (!currParticle->isPrimaryParticle()) {
224 const MCParticle* momParticle = currParticle->getMother();
225 m_momPDG[detID] = momParticle->getPDG();
226 }
227 }
228 }
229
230 //--- TOP
231 detID = 5;
232 // loop over simhits
233 for (Int_t iHit = 0; iHit < m_nSimHits[detID]; iHit++) {
234 const TOPSimHit* simHit = TOPSimHits[iHit];
235 // get related MCparticle
236 if (relTOPSimHitToMCParticle.getFirstElementTo(simHit)) {
237 const MCParticle* currParticle = relTOPSimHitToMCParticle.getFirstElementTo(simHit)->from;
238 m_hitPDG[detID] = currParticle->getPDG();
239 if (!currParticle->isPrimaryParticle()) {
240 const MCParticle* momParticle = currParticle->getMother();
241 m_momPDG[detID] = momParticle->getPDG();
242 }
243 }
244 }
245
246 //--- ECL
247 detID = 6;
248 // loop over simhits
249 for (Int_t iHit = 0; iHit < m_nSimHits[detID]; iHit++) {
250 const ECLSimHit* simHit = ECLSimHits[iHit];
251 // get related MCparticle
252 if (relECLSimHitToMCParticle.getFirstElementTo(simHit)) {
253 const MCParticle* currParticle = relECLSimHitToMCParticle.getFirstElementTo(simHit)->from;
254 m_hitPDG[detID] = currParticle->getPDG();
255 if (!currParticle->isPrimaryParticle()) {
256 const MCParticle* momParticle = currParticle->getMother();
257 m_momPDG[detID] = momParticle->getPDG();
258 }
259 }
260 }
261 //--- KLM
262 detID = 7;
263 // loop over simhits
264 for (Int_t iHit = 0; iHit < m_nSimHits[detID]; iHit++) {
265 const KLMSimHit* simHit = KLMSimHits[iHit];
266 // get related MCparticle
267 if (relKLMSimHitToMCParticle.getFirstElementTo(simHit)) {
268 const MCParticle* currParticle = relKLMSimHitToMCParticle.getFirstElementTo(simHit)->from;
269 m_hitPDG[detID] = currParticle->getPDG();
270 if (!currParticle->isPrimaryParticle()) {
271 const MCParticle* momParticle = currParticle->getMother();
272 m_momPDG[detID] = momParticle->getPDG();
273 }
274 }
275 }
276
277 // fill the tree
278 m_tree2->Fill();
279
280 // BeamBkgHits
281 StoreArray<BeamBackHit> BeamBackHits;
282 RelationIndex<MCParticle, BeamBackHit> relBeamBackHitToMCParticle(McParticles, BeamBackHits);
283
284 Int_t nHits = BeamBackHits.getEntries();
285
286 // loop over bkgHits
287 for (Int_t iHit = 0; iHit < nHits; iHit++) {
288 // get one bkgHit
289 const BeamBackHit* bkgHit = BeamBackHits[iHit];
290 if (bkgHit->getPDG() == 2112) { /*neutron*/
292 m_iden = bkgHit->getIdentifier();
293 m_subDet = bkgHit->getSubDet();
294 m_PDG = bkgHit->getPDG();
295 m_trackID = bkgHit->getTrackID();
296 auto position = bkgHit->getPosition();
297 m_positionX = position.X();
298 m_positionY = position.Y();
299 m_positionZ = position.Z();
300 auto momentum = bkgHit->getMomentum();
301 m_momentumX = momentum.X();
302 m_momentumY = momentum.Y();
303 m_momentumZ = momentum.Z();
304 m_E_start = bkgHit->getEnergy();
305 m_E_end = bkgHit->getEnergyAtExit();
306 m_eDep = bkgHit->getEnergyDeposit();
307 m_trackLength = bkgHit->getTrackLength();
308 m_nWeight = bkgHit->getNeutronWeight();
309
310 m_trj_x.clear();
311 m_trj_y.clear();
312 m_trj_z.clear();
313 m_trj_px.clear();
314 m_trj_py.clear();
315 m_trj_pz.clear();
316
317 // get related MCparticle
318 if (relBeamBackHitToMCParticle.getFirstElementTo(bkgHit)) {
319 const MCParticle* currParticle = relBeamBackHitToMCParticle.getFirstElementTo(bkgHit)->from;
320 auto m_vtxProd = currParticle->getVertex();
321 m_E_init = currParticle->getEnergy();
322 m_mass = currParticle->getMass();
323 m_lifeTime = currParticle->getLifetime();
324 m_vtxProdX = m_vtxProd.X();
325 m_vtxProdY = m_vtxProd.Y();
326 m_vtxProdZ = m_vtxProd.Z();
327
328 // get trajectory stored for this particle
329 const auto mcTrajectories = currParticle->getRelationsTo<MCParticleTrajectory>();
330 for (auto rel : mcTrajectories.relations()) {
331 // the trajectories with negative weight are from secondary daughters which were ignored so we don't use them
332 if (rel.weight <= 0) continue;
333 // get the trajectory
334 const MCParticleTrajectory& trajectory = dynamic_cast<const MCParticleTrajectory&>(*rel.object);
335 // loop over each point
336 for (const MCTrajectoryPoint& pt : trajectory) {
337 m_trj_x.push_back(pt.x);
338 m_trj_y.push_back(pt.y);
339 m_trj_z.push_back(pt.z);
340 m_trj_px.push_back(pt.px);
341 m_trj_py.push_back(pt.py);
342 m_trj_pz.push_back(pt.pz);
343 }
344 }
345 }
346
347 // fill the tree
348 m_tree1->Fill();
349 }
350 }
351 // increase the entry counter
352 m_iEntry++;
353 }
354
356 {
357 // CPU time end
358
359 // Announce
360 B2INFO("BeamBkgNeutron finished.");
361
362 // write
363 m_ff->cd();
364 m_tree1->Write();
365 m_tree2->Write();
366 // close the tree
367 m_ff->Close();
368 }
369
371} // end Belle2 namespace
Class ARICHSimHit - Geant4 simulated hit for ARICH.
Definition ARICHSimHit.h:29
Class BeamBackHit - Stores hits from beam background simulation.
Definition BeamBackHit.h:28
double getNeutronWeight() const
get the effective neutron weight
double getEnergy() const
Get energy of the particle.
int getTrackID() const
the traci ID of the particle
Definition BeamBackHit.h:92
double getEnergyDeposit() const
Get particle energy deposit in sensitive volume.
double getEnergyAtExit() const
Get energy of the particle.
int getPDG() const
Get the lund code of the particle that hit the sensitive area.
Definition BeamBackHit.h:89
double getTrackLength() const
the length of the track in the volume
ROOT::Math::XYZVector getMomentum() const
Get momentum of the particle hit.
ROOT::Math::XYZVector getPosition() const
Get global position of the particle hit.
Definition BeamBackHit.h:95
int getIdentifier() const
Get the identifier of subdetector component in which hit occurred.
Definition BeamBackHit.h:83
int getSubDet() const
Det the index of subdetector in which hit occurred.
Definition BeamBackHit.h:86
Float_t m_mass
McParticle mass [GeV].
std::vector< Float_t > m_trj_pz
Z component of trajectory momentum.
Float_t m_lifeTime
McParticle lifetime [ns].
Float_t m_E_init
McParticle energy [GeV].
Float_t m_vtxProdX
McParticle prod.
Int_t m_iEntry
Entry identifier.
std::vector< Float_t > m_trj_py
Y component of trajectory momentum.
Float_t m_momentumY
Y component of momentum.
std::vector< Float_t > m_trj_z
Z component of trajectory position.
std::vector< Float_t > m_trj_y
Y component of trajectory position.
Float_t m_positionY
Y component of position.
Float_t m_positionX
X component of position.
Float_t m_vtxProdZ
McParticle prod.
Float_t m_vtxProdY
McParticle prod.
Float_t m_momentumZ
Z component of momentum.
std::vector< Float_t > m_trj_x
X component of trajectory position.
Float_t m_E_start
Starting energy.
Int_t m_nSimHits[13]
Array with number of SimHits.
std::string m_filename
Output file name.
Int_t m_momPDG[13]
Array with PDG momentum.
Int_t m_hitPDG[13]
Array with PDG hits.
Float_t m_eDep
Deposited energy.
Int_t m_iEvent
Event identifier.
TFile * m_ff
Output root file.
Float_t m_positionZ
Z component of position.
std::vector< Float_t > m_trj_px
X component of trajectory momentum.
Float_t m_momentumX
X component of momentum.
Example Detector.
Definition CDCSimHit.h:21
ClassECLSimHit - Geant4 simulated hit for the ECL.
Definition ECLSimHit.h:29
KLM simulation hit.
Definition KLMSimHit.h:31
float getPositionZ() const
Get hit global position z coordinate.
Definition KLMSimHit.h:405
Class to save the full simulated trajectory of a particle.
A Class to store the Monte Carlo particle information.
Definition MCParticle.h:32
float getEnergy() const
Return particle energy in GeV.
Definition MCParticle.h:136
float getLifetime() const
Return the lifetime in ns.
Definition MCParticle.h:166
float getMass() const
Return the particle mass in GeV.
Definition MCParticle.h:124
ROOT::Math::XYZVector getVertex() const
Return production vertex position, shorthand for getProductionVertex().
Definition MCParticle.h:172
int getPDG() const
Return PDG code of particle.
Definition MCParticle.h:101
void setDescription(const std::string &description)
Sets the description of the module.
Definition Module.cc:214
Module()
Constructor.
Definition Module.cc:30
Class PXDSimHit - Geant4 simulated hit for the PXD.
Definition PXDSimHit.h:24
Provides access to fast ( O(log n) ) bi-directional lookups on a specified relation.
const Element * getFirstElementTo(const TO &to) const
Return a pointer to the first relation Element of the given object.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
Class SVDSimHit - Geant4 simulated hit for the SVD.
Definition SVDSimHit.h:26
Accessor to arrays stored in the data store.
Definition StoreArray.h:113
int getEntries() const
Get the number of objects in the array.
Definition StoreArray.h:216
Class to store simulated hits of Cherenkov photons on PMT's input for digitization module (TOPDigitiz...
Definition TOPSimHit.h:27
virtual void initialize() override
Initialize the Module.
virtual void event() override
Event processor.
virtual void terminate() override
Termination action.
virtual void beginRun() override
Called when entering a new run.
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:559
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
MCParticle * getMother() const
Returns a pointer to the mother particle.
Definition MCParticle.h:590
bool isPrimaryParticle() const
Check if particle is a primary particle which was created by the generator (and not,...
Definition MCParticle.h:585
Abstract base class for different kinds of events.
STL namespace.
Small struct to encode a position/momentum without additional overhead.
const FROM * from
pointer of the element from which the relation points.