VMC Examples Version 6.6
Loading...
Searching...
No Matches
A01HodoscopeSD.cxx
Go to the documentation of this file.
1//------------------------------------------------
2// The Virtual Monte Carlo examples
3// Copyright (C) 2007 - 2014 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 A01HodoscopeSD.cxx
11/// \brief Implementation of the A01HodoscopeSD class
12///
13/// Geant4 example A01 adapted to Virtual Monte Carlo \n
14///
15/// \date 12/05/2012
16/// \author I. Hrivnacova; IPN, Orsay
17
18#include <Riostream.h>
19#include <TLorentzVector.h>
20#include <TMCRootManager.h>
21#include <TTree.h>
22#include <TVirtualMC.h>
23
24#include "A01HodoscopeHit.h"
25#include "A01HodoscopeSD.h"
26
27using namespace std;
28
29/// \cond CLASSIMP
30ClassImp(A01HodoscopeSD)
31 /// \endcond
32
33 using namespace std;
34
35//_____________________________________________________________________________
36A01HodoscopeSD::A01HodoscopeSD(const char* name, const char* volName)
37 : TNamed(name, ""),
38 fHitsCollection(0),
39 fVolName(volName),
40 fVolId(0),
41 fWriteHits(true),
42 fVerboseLevel(1)
43{
44 /// Standard constructor.
45 /// Create hits collection and an empty hit for each layer.
46 /// \param name The calorimeter hits collection name
47 /// \param volName The sensitive volume name
48
49 // Create hits collection and an empty hit for each layer
50 fHitsCollection = new TClonesArray("A01HodoscopeHit", 500);
51 // cout << "Hodoscope nofHits: " << fHitsCollection->GetEntriesFast() << endl;
52}
53
54//_____________________________________________________________________________
56 : TNamed(origin),
57 fHitsCollection(0),
58 fVolName(origin.fVolName),
59 fVolId(origin.fVolId),
60 fWriteHits(origin.fWriteHits),
61 fVerboseLevel(origin.fVerboseLevel)
62{
63 /// Copy constructor (for clonig on worker thread in MT mode).
64 /// Create hits collection and an empty hit for each layer.
65 /// \param origin The source object (on master).
66
67 fHitsCollection = new TClonesArray("A01HodoscopeHit", 500);
68 // cout << "Hodoscope nofHits: " << fHitsCollection->GetEntriesFast() << endl;
69}
70
71//_____________________________________________________________________________
73 : TNamed(),
74 fHitsCollection(0),
75 fVolName(),
76 fVolId(0),
77 fWriteHits(true),
78 fVerboseLevel(1)
79{
80 /// Default constructor
81}
82
83//_____________________________________________________________________________
85{
86 /// Destructor
87
88 if (fHitsCollection) fHitsCollection->Delete();
89 delete fHitsCollection;
90}
91
92//
93// private methods
94//
95
96//_____________________________________________________________________________
98{
99 /// \return The hit for the specified layer.
100 /// \param i The layer number
101
102 return (A01HodoscopeHit*)fHitsCollection->At(i);
103}
104
105//
106// public methods
107//
108
109//_____________________________________________________________________________
111{
112 /// Register hits collection in the Root manager;
113 /// set sensitive volumes.
114
115 if (TMCRootManager::Instance()) Register();
116
117 fVolId = gMC->VolId(fVolName.Data());
118}
119
120//_____________________________________________________________________________
122{
123 /// Account hit time; create a new hit per detector cell if it does not yet
124 /// exist
125
126 Int_t copyNo;
127 Int_t id = gMC->CurrentVolID(copyNo);
128 if (id != fVolId) return false;
129
130 Double_t edep = gMC->Edep();
131 if (edep == 0.) return false;
132
133 Double_t hitTime = gMC->TrackTime();
134
135 // check if this finger already has a hit
136 Int_t ix = -1;
137 Int_t nofHits = fHitsCollection->GetEntriesFast();
138 for (Int_t i = 0; i < nofHits; i++) {
139 A01HodoscopeHit* hit = GetHit(i);
140 if (hit->GetID() == copyNo) {
141 ix = i;
142 break;
143 }
144 }
145
146 if (ix >= 0) {
147 // if it has, then take the earlier time
148 A01HodoscopeHit* hit = GetHit(ix);
149 if (hit->GetTime() > hitTime) {
150 hit->SetTime(hitTime);
151 }
152 }
153 else {
154 // Debug printing
155 // cout << "** Hodoscope: Create hit in nofHits, copyNo, hitTime[s] "
156 // << nofHits << ", " << copyNo << ", " << hitTime << endl;
157 // cout << "gMC->CurrentVolName(): " << gMC->CurrentVolName() << endl;
158
159 // if not, create a new hit and set it to the collection
160 A01HodoscopeHit* hit =
161 new ((*fHitsCollection)[nofHits]) A01HodoscopeHit(copyNo, hitTime);
162 hit->SetVolId(id);
163
164 // get transformation
165 // add later
166 }
167
168 return true;
169}
170
171//_____________________________________________________________________________
173{
174 /// Print hits collection (if verbose) and reset hits afterwards.
175
176 if (fVerboseLevel > 0) Print();
177
178 // Reset hits collection
179 fHitsCollection->Clear();
180}
181
182//_____________________________________________________________________________
184{
185 /// Register the hits collection in Root manager.
186
187 if (fWriteHits) {
188 TMCRootManager::Instance()->Register(
189 GetName(), "TClonesArray", &fHitsCollection);
190 }
191}
192
193//_____________________________________________________________________________
194void A01HodoscopeSD::Print(Option_t* /*option*/) const
195{
196 /// Print the hits collection.
197
198 Int_t nofHits = fHitsCollection->GetEntriesFast();
199 cout << GetName() << " has " << nofHits << " hits." << endl;
200
201 if (fVerboseLevel > 1) {
202 for (Int_t i = 0; i < nofHits; i++) (*fHitsCollection)[i]->Print();
203 }
204}
205
206/*
207//_____________________________________________________________________________
208void A01HodoscopeSD::PrintTotal() const
209{
210/// Print the total values for all layers.
211
212 Double_t totEAbs=0.;
213 Double_t totLAbs=0.;
214 Double_t totEGap=0.;
215 Double_t totLGap=0.;
216
217 Int_t nofHits = fHitsCollection->GetEntriesFast();
218 for (Int_t i=0; i<nofHits; i++) {
219 totEAbs += GetHit(i)->GetEdepAbs();
220 totLAbs += GetHit(i)->GetTrakAbs();
221 totEGap += GetHit(i)->GetEdepGap();
222 totLGap += GetHit(i)->GetTrakGap();
223 }
224
225 cout << " Absorber: total energy (MeV): "
226 << setw(7) << totEAbs * 1.0e03
227 << " total track length (cm): "
228 << setw(7) << totLAbs
229 << endl
230 << " Gap: total energy (MeV): "
231 << setw(7) << totEGap * 1.0e03
232 << " total track length (cm): "
233 << setw(7) << totLGap
234 << endl;
235}
236*/
Definition of the A01HodoscopeHit class.
Definition of the A01HodoscopeSD class.
The hodoscope hit.
Double_t GetTime() const
void SetVolId(Int_t volId)
void SetTime(Double_t t)
Int_t GetID() const
The calorimeter sensitive detector.
TString fVolName
The sensitive volume Name.
Bool_t fWriteHits
Option to write hits.
virtual void Print(Option_t *option="") const
Int_t fVolId
The calorimeter volume Id.
A01HodoscopeHit * GetHit(Int_t i) const
virtual ~A01HodoscopeSD()
TClonesArray * fHitsCollection
Hits collection.
Int_t fVerboseLevel
Verbosity level.