VMC Examples Version 6.6
Loading...
Searching...
No Matches
A01HadCalorimeterSD.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 A01HadCalorimeterSD.cxx
11/// \brief Implementation of the A01HadCalorimeterSD 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 "A01HadCalorHit.h"
25#include "A01HadCalorimeterSD.h"
26
27using namespace std;
28
29/// \cond CLASSIMP
30ClassImp(A01HadCalorimeterSD)
31 /// \endcond
32
33 using namespace std;
34
37
38//_____________________________________________________________________________
40 : TNamed(name, ""),
41 fCalCollection(0),
42 fVolId(0),
43 fWriteHits(true),
44 fVerboseLevel(1)
45{
46 /// Standard constructor.
47 /// Create hits collection and an empty hit for each layer.
48 /// \param name The calorimeter hits collection name
49
51 new TClonesArray("A01HadCalorHit", fgkNofColumns * fgkNofRows);
52 Int_t counter = 0;
53 for (Int_t iColumn = 0; iColumn < fgkNofColumns; ++iColumn) {
54 for (Int_t iRow = 0; iRow < fgkNofRows; ++iRow) {
55 new ((*fCalCollection)[counter++]) A01HadCalorHit();
56 }
57 }
58}
59
60//_____________________________________________________________________________
62 : TNamed(origin),
63 fCalCollection(0),
64 fVolId(origin.fVolId),
65 fWriteHits(origin.fWriteHits),
66 fVerboseLevel(origin.fVerboseLevel)
67{
68 /// Copy constructor (for clonig on worker thread in MT mode).
69 /// Create hits collection and an empty hit for each layer.
70 /// \param origin The source object (on master).
71
73 new TClonesArray("A01HadCalorHit", fgkNofColumns * fgkNofRows);
74 Int_t counter = 0;
75 for (Int_t iColumn = 0; iColumn < fgkNofColumns; ++iColumn) {
76 for (Int_t iRow = 0; iRow < fgkNofRows; ++iRow) {
77 new ((*fCalCollection)[counter++]) A01HadCalorHit();
78 }
79 }
80}
81
82//_____________________________________________________________________________
84 : TNamed(), fCalCollection(0), fVolId(0), fWriteHits(true), fVerboseLevel(1)
85{
86 /// Default constructor
87}
88
89//_____________________________________________________________________________
91{
92 /// Destructor
93
94 if (fCalCollection) fCalCollection->Delete();
95 delete fCalCollection;
96}
97
98//
99// private methods
100//
101
102//_____________________________________________________________________________
104{
105 /// \return The hit for the specified layer.
106 /// \param i The layer number
107
108 return (A01HadCalorHit*)fCalCollection->At(i);
109}
110
111//_____________________________________________________________________________
113{
114 /// Reset all hits in the hits collection.
115
116 for (Int_t i = 0; i < fCalCollection->GetEntriesFast(); i++)
117 GetHit(i)->Reset();
118}
119
120//
121// public methods
122//
123
124//_____________________________________________________________________________
126{
127 /// Register hits collection in the Root manager;
128 /// set sensitive volumes.
129
130 if (TMCRootManager::Instance()) Register();
131
132 fVolId = gMC->VolId("HadCalScintiLogical");
133}
134
135//_____________________________________________________________________________
137{
138 /// Account energy deposit for each layer in its hit.
139
140 Int_t copyNo;
141 Int_t id = gMC->CurrentVolID(copyNo);
142 if (id != fVolId) return false;
143
144 Double_t edep = gMC->Edep();
145 if (edep == 0.) return false;
146
147 Int_t rowNo;
148 gMC->CurrentVolOffID(2, rowNo);
149 Int_t columnNo;
150 gMC->CurrentVolOffID(3, columnNo);
151 // VMC adopts Root numbering of divisions starting from 1
152 rowNo--;
153 columnNo--;
154 Int_t hitID = fgkNofRows * columnNo + rowNo;
155
156 A01HadCalorHit* hit = GetHit(hitID);
157 if (!hit) {
158 std::cerr << "No hit found for layer with "
159 << "rowNo = " << rowNo << " columnNo = " << columnNo << endl;
160 return false;
161 }
162
163 // check if it is the first touch
164 if (hit->GetColumnID() < 0) {
165 // Debug printing (to check hits indexing)
166 // cout << "HadCalorimeter: First Add in hit in (column, row): "
167 // << columnNo << ", " << rowNo << endl;
168 // cout << "gMC->CurrentVolOffName(2), gMC->CurrentVolOffName(3): "
169 // << gMC->CurrentVolOffName(2) << ", " << gMC->CurrentVolOffName(3) <<
170 // endl;
171
172 // fill volume information
173 hit->SetRowID(rowNo);
174 hit->SetColumnID(columnNo);
175
176 // get transformation
177 // add later
178 }
179
180 // add energy deposition
181 hit->AddEdep(edep);
182
183 return true;
184}
185
186//_____________________________________________________________________________
188{
189 /// Print hits collection (if verbose) and reset hits afterwards.
190
191 if (fVerboseLevel > 0) PrintTotal();
192
193 // Reset hits collection
194 ResetHits();
195}
196
197//_____________________________________________________________________________
199{
200 /// Register the hits collection in Root manager.
201
202 if (fWriteHits) {
203 TMCRootManager::Instance()->Register(
204 GetName(), "TClonesArray", &fCalCollection);
205 }
206}
207
208//_____________________________________________________________________________
209void A01HadCalorimeterSD::Print(Option_t* /*option*/) const
210{
211 /// Print the hits collection.
212
213 Int_t nofHits = fCalCollection->GetEntriesFast();
214
215 cout << "\n-------->Hits Collection: in this event: " << endl;
216
217 if (fVerboseLevel > 1) {
218 for (Int_t i = 0; i < nofHits; i++) (*fCalCollection)[i]->Print();
219 }
220}
221
222//_____________________________________________________________________________
224{
225 /// Print the total values for all layers.
226
227 Int_t nofHits = 0;
228 Double_t totalEdep = 0.;
229 for (Int_t i = 0; i < fgkNofColumns * fgkNofRows; ++i) {
230 Double_t edep = GetHit(i)->GetEdep();
231 if (edep > 0.) {
232 nofHits++;
233 totalEdep += edep;
234 }
235 }
236 cout << GetName() << " has " << nofHits << " hits. Total Edep is "
237 << totalEdep * 1e03 << " (MeV)" << endl;
238}
Definition of the A01HadCalorHit class.
Definition of the A01HadCalorimeterSD class.
The hadron calorimeter hit.
Int_t GetColumnID() const
Double_t GetEdep() const
void SetRowID(Int_t volId)
void SetColumnID(Int_t z)
void AddEdep(Double_t de)
The hadron calorimeter sensitive detector.
A01HadCalorHit * GetHit(Int_t i) const
static const Int_t fgkNofRows
Int_t fVerboseLevel
Verbosity level.
virtual void Print(Option_t *option="") const
Int_t fVolId
The calorimeter volume Id.
static const Int_t fgkNofColumns
Bool_t fWriteHits
Option to write hits.
TClonesArray * fCalCollection
Hits collection.