Belle II Software  release-08-01-10
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 #ifndef VXD_SIMULATION_SENSITIVEDETECTORDEBUGHELPER_H
11 #define VXD_SIMULATION_SENSITIVEDETECTORDEBUGHELPER_H
12 
13 #include <vxd/simulation/SensorTraversal.h>
14 #include <vxd/dataobjects/VxdID.h>
15 #include <vxd/dataobjects/VXDTrueHit.h>
16 #include <vxd/dataobjects/VXDSimHit.h>
17 
18 #include <TFile.h>
19 #include <TTree.h>
20 
21 namespace Belle2 {
26  namespace VXD {
27 
42  public:
46  void startTraversal(VxdID sensorID, const SensorTraversal& traversal);
48  void addTrueHit(const VXDTrueHit* truehit);
50  void addSimHit(const VXDSimHit* simhit, int startPoint, int endPoint);
53  {
54  m_tree->Fill();
55  }
56  private:
65  {
66  m_tree->Write();
67  m_file->Close();
68  }
69 
71  TFile* m_file {0};
73  TTree* m_tree {0};
74 
75  enum {
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  };
92 
94  struct info {
96  int sensorID;
98  int pdg;
102  int primary;
104  double length;
105 
107  int stepN;
110 
112  int trueN;
115 
117  int edepN;
120 
122  int simhitN;
125  } m_info;
126  };
127  }
129 } //Belle2 namespace
130 #endif // VXD_SIMULATION_DEBUGFILE_H
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:30
Class VXDTrueHit - Records of tracks that either enter or leave the sensitive volume.
Definition: VXDTrueHit.h:36
Small helper class to facilitate debugging of VXD::SensitiveDetector implementation.
void finishTraversal()
Finish the entry by calling Fill()
static SensitiveDetectorDebugHelper & getInstance()
Singleton class: get instance.
TTree * m_tree
ROOT Tree to write information to.
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.
@ 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
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: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