VMC Examples Version 6.6
Loading...
Searching...
No Matches
SensitiveDetector.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Virtual Monte Carlo examples
3// Copyright (C) 2007 - 2015 Ivana Hrivnacova
4// All rights reserved.
5//
6// For the licensing terms see geant4_vmc/LICENSE.
7// Contact: root-vmc@cern.ch
8//-------------------------------------------------
9
10/// \file Gflash/src/SensitiveDetector.cxx
11/// \brief Implementation of the Gflash::SensitiveDetector class
12///
13/// Geant4 gflash example adapted to Virtual Monte Carlo.
14///
15/// \date 28/10/2015
16/// \author I. Hrivnacova; IPN, Orsay
17
18#include "SensitiveDetector.h"
19#include "Hit.h"
20
21#include <Riostream.h>
22#include <TLorentzVector.h>
23#include <TMCRootManager.h>
24#include <TTree.h>
25#include <TVirtualMC.h>
26
27/// \cond CLASSIMP
29 /// \endcond
30
31 using namespace std;
32
33namespace VMC
34{
35namespace Gflash
36{
37
38//_____________________________________________________________________________
40 : TNamed(name, ""), fCaloHitsCollection(0), fCrystalVolId(), fVerboseLevel(1)
41{
42 /// Standard constructor.
43 /// Create hits collection.
44 /// \param name The calorimeter name
45
46 fCaloHitsCollection = new TClonesArray("VMC::Gflash::Hit", 500);
47}
48
49//_____________________________________________________________________________
51 : TNamed(origin),
52 fCaloHitsCollection(0),
53 fCrystalVolId(),
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}
62
63//_____________________________________________________________________________
65 : TNamed(), fCaloHitsCollection(0), fCrystalVolId(), fVerboseLevel(0)
66{
67 /// Default constructor
68}
69
70//_____________________________________________________________________________
72{
73 /// Destructor
74
77}
78
79//
80// private methods
81//
82
83//_____________________________________________________________________________
85{
86 /// \return The hit for the specified layer.
87 /// \param i The layer number
88
89 return (Hit*)fCaloHitsCollection->At(i);
90}
91
92//
93// public methods
94//
95
96//_____________________________________________________________________________
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}
105
106//_____________________________________________________________________________
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}
133
134//_____________________________________________________________________________
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}
144
145//_____________________________________________________________________________
147{
148 /// Register the hits collection in Root manager.
149
150 TMCRootManager::Instance()->Register(
151 "hits", "TClonesArray", &fCaloHitsCollection);
152}
153
154//_____________________________________________________________________________
155void SensitiveDetector::Print(Option_t* /*option*/) const
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}
165
166//_____________________________________________________________________________
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}
181
182} // namespace Gflash
183} // namespace VMC
The calorimeter hit.
Definition Hit.h:37
Double_t GetEdep() const
Definition Hit.h:52
void SetPos(TVector3 xyz)
Definition Hit.h:48
void SetEdep(Double_t de)
Definition Hit.h:47
void SetCrystalNum(Int_t num)
Definition Hit.h:49
The calorimeter sensitive detector.
Int_t fVerboseLevel
Verbosity level.
Int_t fCrystalVolId
The crystal volume Id.
virtual void Print(Option_t *option="") const
TClonesArray * fCaloHitsCollection
Hits collection.