Belle II Software development
SensitiveDetectorDebugHelper Class Reference

Small helper class to facilitate debugging of VXD::SensitiveDetector implementation. More...

#include <SensitiveDetectorDebugHelper.h>

Classes

struct  info
 struct with all the branches needed More...
 

Public Member Functions

void startTraversal (VxdID sensorID, const SensorTraversal &traversal)
 Start writing a new sensor traversal: Geant4 steps will be added to the ROOT entry.
 
void addTrueHit (const VXDTrueHit *truehit)
 Write the information normally used when creating TrueHits to the entry.
 
void addSimHit (const VXDSimHit *simhit, int startPoint, int endPoint)
 Write the information normally used when creating SimHits to the entry.
 
void finishTraversal ()
 Finish the entry by calling Fill()
 

Static Public Member Functions

static SensitiveDetectorDebugHelpergetInstance ()
 Singleton class: get instance.
 

Private Types

enum  {
  MAX_STEPS = 5000 ,
  MAX_HITS = 500 ,
  MAX_EDEP = 10000 ,
  SIZE_STEP = 9 ,
  SIZE_TRUE = 20 ,
  SIZE_HITS = 13 ,
  SIZE_EDEP = 2
}
 

Private Member Functions

 SensitiveDetectorDebugHelper ()
 Singleton class: private constructor.
 
 SensitiveDetectorDebugHelper (const SensitiveDetectorDebugHelper &)=delete
 Singleton class: no copy constructor.
 
SensitiveDetectorDebugHelperoperator= (const SensitiveDetectorDebugHelper &)=delete
 Singleton class: no assignment operator.
 
 ~SensitiveDetectorDebugHelper ()
 Destructor: write root file.
 

Private Attributes

TFile * m_file {0}
 ROOT File to write information to.
 
TTree * m_tree {0}
 ROOT Tree to write information to.
 
struct Belle2::VXD::SensitiveDetectorDebugHelper::info m_info
 object to store all variables
 

Detailed Description

Small helper class to facilitate debugging of VXD::SensitiveDetector implementation.

This file will write a root file containing the Geant4 Step information as well as the information used for TrueHit and SimHit creation to check whether these are consistent. Normally this class should not be in use, a define can be set in VXD::SensitiveDetectorBase to enable its use.

This class is designed to make access very easy in python using numpy, e.g. the SimHits can be obtained as numpy array from the tree using

‍tree.GetEntry(i) simhits = np.array(tree.simhitInfo).reshape((-1, 13)) # HIT_SIZE=13

Definition at line 41 of file SensitiveDetectorDebugHelper.h.

Member Enumeration Documentation

◆ anonymous enum

anonymous enum
private
Enumerator
MAX_STEPS 

assume maximum number of steps per sensor traversal

MAX_HITS 

assume maximum number of simhits per sensor traversal

MAX_EDEP 

assume maximum number of energy deposition points per sensor traversal

SIZE_STEP 

number of values per step

SIZE_TRUE 

number of values per truehit

SIZE_HITS 

number of values per simhit

SIZE_EDEP 

number of values per energy deposition point

Definition at line 75 of file SensitiveDetectorDebugHelper.h.

75 {
77 MAX_STEPS = 5000,
79 MAX_HITS = 500,
81 MAX_EDEP = 10000,
82
84 SIZE_STEP = 9,
86 SIZE_TRUE = 20,
88 SIZE_HITS = 13,
90 SIZE_EDEP = 2
91 };
@ 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

Constructor & Destructor Documentation

◆ SensitiveDetectorDebugHelper()

Singleton class: private constructor.

Definition at line 26 of file SensitiveDetectorDebugHelper.cc.

27 {
28 m_file = new TFile("debug-vxdsensitivedetector.root", "RECREATE");
29 m_tree = new TTree("sensortraversal", "Debug info for VXD Sensitive Detector implementation");
30//Due to laziness, define some macros to create branches with less typing
31#define ADDBRANCH__(name,var,type) m_tree->Branch(#name,&m_info.var,#type)
32#define ADDBARRAY__(name,var,n,type) ADDBRANCH__(name,var,name[n]/type)
33#define ADDBRANCH(x,type) ADDBRANCH__(x,x,x/type)
34#define ADDBARRAY(x,n,type) ADDBARRAY__(x,x,n,type)
35 //and create all branches
36 ADDBRANCH(sensorID, I);
37 ADDBRANCH(pdg, I);
38 ADDBRANCH(contained, I);
39 ADDBRANCH(primary, I);
40 ADDBRANCH(length, D);
41 ADDBRANCH(stepN, I);
42 ADDBARRAY(stepInfo, stepN, D);
43 ADDBRANCH(trueN, I);
44 ADDBARRAY(trueInfo, trueN, D);
45 ADDBRANCH(edepN, I);
46 ADDBARRAY(edepInfo, edepN, D);
47 ADDBRANCH(simhitN, I);
48 ADDBARRAY(simhitInfo, simhitN, D);
49 }
TTree * m_tree
ROOT Tree to write information to.
TFile * m_file
ROOT File to write information to.

◆ ~SensitiveDetectorDebugHelper()

~SensitiveDetectorDebugHelper ( )
inlineprivate

Destructor: write root file.

Definition at line 64 of file SensitiveDetectorDebugHelper.h.

65 {
66 m_tree->Write();
67 m_file->Close();
68 }

Member Function Documentation

◆ addSimHit()

void addSimHit ( const VXDSimHit simhit,
int  startPoint,
int  endPoint 
)

Write the information normally used when creating SimHits to the entry.

Definition at line 103 of file SensitiveDetectorDebugHelper.cc.

104 {
106 m_info.simhitInfo[m_info.simhitN + 0] = simhit->getPosIn().X();
107 m_info.simhitInfo[m_info.simhitN + 1] = simhit->getPosIn().Y();
108 m_info.simhitInfo[m_info.simhitN + 2] = simhit->getPosIn().Z();
109 m_info.simhitInfo[m_info.simhitN + 3] = simhit->getPosOut().X();
110 m_info.simhitInfo[m_info.simhitN + 4] = simhit->getPosOut().Y();
111 m_info.simhitInfo[m_info.simhitN + 5] = simhit->getPosOut().Z();
112 m_info.simhitInfo[m_info.simhitN + 6] = simhit->getGlobalTime();
113 m_info.simhitInfo[m_info.simhitN + 7] = startPoint;
114 m_info.simhitInfo[m_info.simhitN + 8] = endPoint;
115 //Add optimized energy deposition profile
117 for (auto step : simhit->getElectronProfile()) {
118 m_info.edepInfo[m_info.edepN] = step.first;
119 ++m_info.edepN;
120 m_info.edepInfo[m_info.edepN] = step.second;
121 ++m_info.edepN;
122 }
124 //And also add electron profile sampled by distance
125 for (auto step : simhit->getElectronsConstantDistance(50 * Unit::um)) {
126 m_info.edepInfo[m_info.edepN] = step.first;
127 ++m_info.edepN;
128 m_info.edepInfo[m_info.edepN] = step.second;
129 ++m_info.edepN;
130 }
132 //And sampled by electrons
133 for (auto step : simhit->getElectronsConstantNumber(1000)) {
134 m_info.edepInfo[m_info.edepN] = step.first;
135 ++m_info.edepN;
136 m_info.edepInfo[m_info.edepN] = step.second;
137 ++m_info.edepN;
138 }
141 }
static const double um
[micrometers]
Definition: Unit.h:71
struct Belle2::VXD::SensitiveDetectorDebugHelper::info m_info
object to store all variables
double simhitInfo[SIZE_HITS *MAX_HITS]
values for the simhits
double edepInfo[SIZE_EDEP *MAX_EDEP]
values for the energy deposition

◆ addTrueHit()

void addTrueHit ( const VXDTrueHit truehit)

Write the information normally used when creating TrueHits to the entry.

Definition at line 78 of file SensitiveDetectorDebugHelper.cc.

79 {
80 m_info.trueInfo[0] = truehit->getEntryU();
81 m_info.trueInfo[1] = truehit->getEntryV();
82 m_info.trueInfo[2] = truehit->getEntryW();
83 m_info.trueInfo[3] = truehit->getU();
84 m_info.trueInfo[4] = truehit->getV();
85 m_info.trueInfo[5] = truehit->getW();
86 m_info.trueInfo[6] = truehit->getExitU();
87 m_info.trueInfo[7] = truehit->getExitV();
88 m_info.trueInfo[8] = truehit->getExitW();
89 m_info.trueInfo[9] = truehit->getEntryMomentum().X();
90 m_info.trueInfo[10] = truehit->getEntryMomentum().Y();
91 m_info.trueInfo[11] = truehit->getEntryMomentum().Z();
92 m_info.trueInfo[12] = truehit->getMomentum().X();
93 m_info.trueInfo[13] = truehit->getMomentum().Y();
94 m_info.trueInfo[14] = truehit->getMomentum().Z();
95 m_info.trueInfo[15] = truehit->getExitMomentum().X();
96 m_info.trueInfo[16] = truehit->getExitMomentum().Y();
97 m_info.trueInfo[17] = truehit->getExitMomentum().Z();
98 m_info.trueInfo[18] = truehit->getEnergyDep();
99 m_info.trueInfo[19] = truehit->getGlobalTime();
101 }

◆ finishTraversal()

void finishTraversal ( )
inline

Finish the entry by calling Fill()

Definition at line 52 of file SensitiveDetectorDebugHelper.h.

53 {
54 m_tree->Fill();
55 }

◆ getInstance()

SensitiveDetectorDebugHelper & getInstance ( )
static

Singleton class: get instance.

Definition at line 20 of file SensitiveDetectorDebugHelper.cc.

21 {
22 static SensitiveDetectorDebugHelper instance;
23 return instance;
24 }
SensitiveDetectorDebugHelper()
Singleton class: private constructor.

◆ startTraversal()

void startTraversal ( VxdID  sensorID,
const SensorTraversal traversal 
)

Start writing a new sensor traversal: Geant4 steps will be added to the ROOT entry.

Definition at line 50 of file SensitiveDetectorDebugHelper.cc.

51 {
52 m_info.sensorID = sensorID;
53 m_info.pdg = traversal.getPDGCode();
54 m_info.contained = traversal.isContained();
55 m_info.primary = traversal.isPrimary();
56 m_info.length = traversal.getLength();
57 m_info.stepN = 0;
58 m_info.trueN = 0;
59 m_info.edepN = 0;
60 m_info.simhitN = 0;
61
62 //add all steps if there is enough space
63 assert(traversal.size() <= MAX_STEPS);
64 for (const StepInformation& step : traversal) {
65 m_info.stepInfo[m_info.stepN + 0] = step.position.x();
66 m_info.stepInfo[m_info.stepN + 1] = step.position.y();
67 m_info.stepInfo[m_info.stepN + 2] = step.position.z();
68 m_info.stepInfo[m_info.stepN + 3] = step.momentum.x();
69 m_info.stepInfo[m_info.stepN + 4] = step.momentum.y();
70 m_info.stepInfo[m_info.stepN + 5] = step.momentum.z();
71 m_info.stepInfo[m_info.stepN + 6] = step.electrons;
72 m_info.stepInfo[m_info.stepN + 7] = step.time;
73 m_info.stepInfo[m_info.stepN + 8] = step.length;
75 }
76 }
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
int contained
whether or not the track was contained in the senitive volume

Member Data Documentation

◆ m_file

TFile* m_file {0}
private

ROOT File to write information to.

Definition at line 71 of file SensitiveDetectorDebugHelper.h.

◆ m_tree

TTree* m_tree {0}
private

ROOT Tree to write information to.

Definition at line 73 of file SensitiveDetectorDebugHelper.h.


The documentation for this class was generated from the following files: