VMC Examples Version 6.6
Loading...
Searching...
No Matches
VMC::Gflash::SensitiveDetector Class Reference

The calorimeter sensitive detector. More...

#include <SensitiveDetector.h>

Inheritance diagram for VMC::Gflash::SensitiveDetector:

Public Member Functions

 SensitiveDetector (const char *name)
 
 SensitiveDetector (const SensitiveDetector &origin)
 
 SensitiveDetector ()
 
virtual ~SensitiveDetector ()
 
void Initialize ()
 
Bool_t ProcessHits ()
 
void EndOfEvent ()
 
void Register ()
 
virtual void Print (Option_t *option="") const
 
void PrintTotal () const
 
void SetVerboseLevel (Int_t level)
 
HitGetHit (Int_t i) const
 
TClonesArray * GetHitsCollection () const
 Return the hits collection.
 

Private Attributes

TClonesArray * fCaloHitsCollection
 Hits collection.
 
Int_t fCrystalVolId
 The crystal volume Id.
 
Int_t fVerboseLevel
 Verbosity level.
 

Detailed Description

The calorimeter sensitive detector.

Geant4 gflash example adapted to Virtual Monte Carlo.

Date
28/10/2015
Author
I. Hrivnacova; IPN, Orsay

Definition at line 38 of file SensitiveDetector.h.

Constructor & Destructor Documentation

◆ SensitiveDetector() [1/3]

VMC::Gflash::SensitiveDetector::SensitiveDetector ( const char * name)

Standard constructor. Create hits collection.

Parameters
nameThe calorimeter name

Definition at line 39 of file SensitiveDetector.cxx.

41{
42 /// Standard constructor.
43 /// Create hits collection.
44 /// \param name The calorimeter name
45
46 fCaloHitsCollection = new TClonesArray("VMC::Gflash::Hit", 500);
47}
Int_t fVerboseLevel
Verbosity level.
Int_t fCrystalVolId
The crystal volume Id.
TClonesArray * fCaloHitsCollection
Hits collection.

◆ SensitiveDetector() [2/3]

VMC::Gflash::SensitiveDetector::SensitiveDetector ( const SensitiveDetector & origin)

Copy constructor (for cloning on worker thread in MT mode). Create hits collection.

Parameters
originThe source object (on master).

Definition at line 50 of file SensitiveDetector.cxx.

51 : TNamed(origin),
54 fVerboseLevel(origin.fVerboseLevel)
55{
56 /// Copy constructor (for cloning on worker thread in MT mode).
57 /// Create hits collection.
58 /// \param origin The source object (on master).
59
60 fCaloHitsCollection = new TClonesArray("Gflash::Hit");
61}

◆ SensitiveDetector() [3/3]

VMC::Gflash::SensitiveDetector::SensitiveDetector ( )

Default constructor

Definition at line 64 of file SensitiveDetector.cxx.

66{
67 /// Default constructor
68}

◆ ~SensitiveDetector()

VMC::Gflash::SensitiveDetector::~SensitiveDetector ( )
virtual

Destructor

Definition at line 71 of file SensitiveDetector.cxx.

72{
73 /// Destructor
74
77}

Member Function Documentation

◆ Initialize()

void VMC::Gflash::SensitiveDetector::Initialize ( )

Register hits collection in the Root manager; set sensitive volumes.

Definition at line 97 of file SensitiveDetector.cxx.

98{
99 /// Register hits collection in the Root manager;
100 /// set sensitive volumes.
101
102 if (TMCRootManager::Instance()) Register();
103 fCrystalVolId = gMC->VolId("Crystal_log");
104}

◆ ProcessHits()

Bool_t VMC::Gflash::SensitiveDetector::ProcessHits ( )

Account energy deposit and track lengths for each layer in its hit.

Definition at line 107 of file SensitiveDetector.cxx.

108{
109 /// Account energy deposit and track lengths for each layer in its hit.
110
111 Int_t copyNo;
112 Int_t id = gMC->CurrentVolID(copyNo);
113
114 if (id != fCrystalVolId) return false;
115
116 // cout << "In crystal " << copyNo << endl;
117
118 Double_t edep = gMC->Edep();
119 // cout << " edep [keV]" << edep*1.e06 << endl;
120 Double_t xpos;
121 Double_t ypos;
122 Double_t zpos;
123 gMC->TrackPosition(xpos, ypos, zpos);
124
125 Int_t next = fCaloHitsCollection->GetEntriesFast();
126 Hit* caloHit = new ((*fCaloHitsCollection)[next]) Hit();
127 caloHit->SetEdep(edep);
128 caloHit->SetPos(TVector3(xpos, ypos, zpos));
129 caloHit->SetCrystalNum(copyNo);
130
131 return true;
132}

◆ EndOfEvent()

void VMC::Gflash::SensitiveDetector::EndOfEvent ( )

Print hits collection (if verbose) and reset hits afterwards.

Definition at line 135 of file SensitiveDetector.cxx.

136{
137 /// Print hits collection (if verbose) and reset hits afterwards.
138
139 if (fVerboseLevel > 1) Print();
140
141 // Reset hits collection
142 fCaloHitsCollection->Clear();
143}
virtual void Print(Option_t *option="") const

◆ Register()

void VMC::Gflash::SensitiveDetector::Register ( )

Register the hits collection in Root manager.

Definition at line 146 of file SensitiveDetector.cxx.

147{
148 /// Register the hits collection in Root manager.
149
150 TMCRootManager::Instance()->Register(
151 "hits", "TClonesArray", &fCaloHitsCollection);
152}

◆ Print()

void VMC::Gflash::SensitiveDetector::Print ( Option_t * option = "") const
virtual

Print the hits collection.

Definition at line 155 of file SensitiveDetector.cxx.

156{
157 /// Print the hits collection.
158
159 Int_t nofHits = fCaloHitsCollection->GetEntriesFast();
160
161 cout << "\n-------->Hits Collection: in this event: " << endl;
162
163 for (Int_t i = 0; i < nofHits; i++) (*fCaloHitsCollection)[i]->Print();
164}

◆ PrintTotal()

void VMC::Gflash::SensitiveDetector::PrintTotal ( ) const

Print the total values for all layers.

Definition at line 167 of file SensitiveDetector.cxx.

168{
169 /// Print the total values for all layers.
170
171 Double_t totEdep = 0.;
172
173 Int_t nofHits = fCaloHitsCollection->GetEntriesFast();
174 for (Int_t i = 0; i < nofHits; i++) {
175 totEdep += GetHit(i)->GetEdep();
176 }
177
178 cout << " Calorimeter: total energy (MeV): " << setw(7) << totEdep * 1.0e03
179 << endl;
180}
Double_t GetEdep() const
Definition Hit.h:52

◆ SetVerboseLevel()

void VMC::Gflash::SensitiveDetector::SetVerboseLevel ( Int_t level)
inline

Set verbose level

Parameters
levelThe new verbose level value

Definition at line 72 of file SensitiveDetector.h.

73{
74 fVerboseLevel = level;
75}

◆ GetHit()

Hit * VMC::Gflash::SensitiveDetector::GetHit ( Int_t i) const
Returns
The hit for the specified layer.
Parameters
iThe layer number

Definition at line 84 of file SensitiveDetector.cxx.

85{
86 /// \return The hit for the specified layer.
87 /// \param i The layer number
88
89 return (Hit*)fCaloHitsCollection->At(i);
90}

◆ GetHitsCollection()

TClonesArray * VMC::Gflash::SensitiveDetector::GetHitsCollection ( ) const
inline

Return the hits collection.

Definition at line 78 of file SensitiveDetector.h.

79{
81}

Member Data Documentation

◆ fCaloHitsCollection

TClonesArray* VMC::Gflash::SensitiveDetector::fCaloHitsCollection
private

Hits collection.

Definition at line 63 of file SensitiveDetector.h.

◆ fCrystalVolId

Int_t VMC::Gflash::SensitiveDetector::fCrystalVolId
private

The crystal volume Id.

Definition at line 64 of file SensitiveDetector.h.

◆ fVerboseLevel

Int_t VMC::Gflash::SensitiveDetector::fVerboseLevel
private

Verbosity level.

Definition at line 65 of file SensitiveDetector.h.


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