Belle II Software development
SensitiveDetectorDebugHelper.h
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#pragma once
10
11#include <vxd/simulation/SensorTraversal.h>
12#include <vxd/dataobjects/VxdID.h>
13#include <vxd/dataobjects/VXDTrueHit.h>
14#include <vxd/dataobjects/VXDSimHit.h>
15
16#include <TFile.h>
17#include <TTree.h>
18
19namespace Belle2 {
24 namespace VXD {
25
40 public:
44 void startTraversal(VxdID sensorID, const SensorTraversal& traversal);
46 void addTrueHit(const VXDTrueHit* truehit);
48 void addSimHit(const VXDSimHit* simhit, int startPoint, int endPoint);
51 {
52 m_tree->Fill();
53 }
54 private:
63 {
64 m_tree->Write();
65 m_file->Close();
66 }
67
69 TFile* m_file {0};
71 TTree* m_tree {0};
72
73 enum {
75 MAX_STEPS = 5000,
77 MAX_HITS = 500,
79 MAX_EDEP = 10000,
80
89 };
90
124 };
125 }
127} //Belle2 namespace
Class to keep track of the traversal of the sensitive volume for one track.
Class VXDSimHit - Geant4 simulated hit for the VXD.
Definition VXDSimHit.h:28
Class VXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition VXDTrueHit.h:34
void finishTraversal()
Finish the entry by calling Fill()
static SensitiveDetectorDebugHelper & getInstance()
Singleton class: get instance.
TTree * m_tree
ROOT Tree to write information to.
@ MAX_EDEP
assume maximum number of energy deposition points per sensor traversal
@ MAX_HITS
assume maximum number of simhits per sensor traversal
@ SIZE_EDEP
number of values per energy deposition point
@ MAX_STEPS
assume maximum number of steps per sensor traversal
void addSimHit(const VXDSimHit *simhit, int startPoint, int endPoint)
Write the information normally used when creating SimHits to the entry.
void startTraversal(VxdID sensorID, const SensorTraversal &traversal)
Start writing a new sensor traversal: Geant4 steps will be added to the ROOT entry.
struct Belle2::VXD::SensitiveDetectorDebugHelper::info m_info
object to store all variables
SensitiveDetectorDebugHelper()
Singleton class: private constructor.
void addTrueHit(const VXDTrueHit *truehit)
Write the information normally used when creating TrueHits to the entry.
SensitiveDetectorDebugHelper(const SensitiveDetectorDebugHelper &)=delete
Singleton class: no copy constructor.
TFile * m_file
ROOT File to write information to.
SensitiveDetectorDebugHelper & operator=(const SensitiveDetectorDebugHelper &)=delete
Singleton class: no assignment operator.
~SensitiveDetectorDebugHelper()
Destructor: write root file.
Class to uniquely identify a any structure of the PXD and SVD.
Definition VxdID.h:32
Namespace to provide code needed by both Vertex Detectors, PXD and SVD, and also testbeam telescopes.
Definition GeoCache.h:33
Abstract base class for different kinds of events.
double stepInfo[SIZE_STEP *MAX_STEPS]
values for the step points
int sensorID
id of the sensor for the current track
int primary
whether or not the track is from a primary particle
double simhitInfo[SIZE_HITS *MAX_HITS]
values for the simhits
int contained
whether or not the track was contained in the senitive volume
double edepInfo[SIZE_EDEP *MAX_EDEP]
values for the energy deposition