VMC Examples Version 6.6
Loading...
Searching...
No Matches
A01DriftChamberSD.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 A01DriftChamberSD.cxx
11/// \brief Implementation of the A01DriftChamberSD 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 "A01DriftChamberHit.h"
25#include "A01DriftChamberSD.h"
26
27using namespace std;
28
29/// \cond CLASSIMP
30ClassImp(A01DriftChamberSD)
31 /// \endcond
32
33 using namespace std;
34
35//_____________________________________________________________________________
36A01DriftChamberSD::A01DriftChamberSD(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 /// Creates hits collection.
46 /// \param name The calorimeter hits collection name
47 /// \param volName The sensitive volume name
48
49 fHitsCollection = new TClonesArray("A01DriftChamberHit", 500);
50}
51
52//_____________________________________________________________________________
54 : TNamed(origin),
55 fHitsCollection(0),
56 fVolName(origin.fVolName),
57 fVolId(origin.fVolId),
58 fWriteHits(origin.fWriteHits),
59 fVerboseLevel(origin.fVerboseLevel)
60{
61 /// Copy constructor (for clonig on worker thread in MT mode).
62 /// Creates hits collection.
63 /// \param origin The source object (on master).
64
65 fHitsCollection = new TClonesArray("A01DriftChamberHit", 500);
66}
67
68//_____________________________________________________________________________
70 : TNamed(),
71 fHitsCollection(0),
72 fVolName(),
73 fVolId(0),
74 fWriteHits(true),
75 fVerboseLevel(1)
76{
77 /// Default constructor
78}
79
80//_____________________________________________________________________________
82{
83 /// Destructor
84
85 if (fHitsCollection) fHitsCollection->Delete();
86 delete fHitsCollection;
87}
88
89//
90// private methods
91//
92
93//_____________________________________________________________________________
95{
96 /// \return The hit for the specified layer.
97 /// \param i The layer number
98
99 return (A01DriftChamberHit*)fHitsCollection->At(i);
100}
101
102//
103// public methods
104//
105
106//_____________________________________________________________________________
108{
109 /// Register hits collection in the Root manager;
110 /// set sensitive volumes.
111
112 if (TMCRootManager::Instance()) Register();
113
114 fVolId = gMC->VolId(fVolName.Data());
115}
116
117//_____________________________________________________________________________
119{
120 /// Account the hit time, local and global position for each layer in its hit.
121
122 Int_t copyNo;
123 Int_t id = gMC->CurrentVolID(copyNo);
124 if (id != fVolId) return false;
125
126 Double_t charge = gMC->TrackCharge();
127 if (charge == 0.) return false;
128
129 if (!gMC->IsTrackEntering()) return false;
130
131 // get copyNo in mother
132 gMC->CurrentVolOffID(1, copyNo);
133
134 // global and local hit position
135 Double_t globalPos[3];
136 Double_t localPos[3];
137 gMC->TrackPosition(globalPos[0], globalPos[1], globalPos[2]);
138 gMC->Gmtod(globalPos, localPos, 1);
139
140 // Debug printing
141 // cout << "** Drift Chamber: Create hit in DriftChamber copyNo, worldPos[cm]
142 // "
143 // << copyNo << ", "
144 // << globalPos[0] << ", " << globalPos[1] << ", " << globalPos[2] <<
145 // endl;
146 // cout << "gMC->CurrentVolOffName(1): " << gMC->CurrentVolOffName(1) << endl;
147 // cout << "localPos[cm] "
148 // << localPos[0] << ", " << localPos[1] << ", " << localPos[2] << endl;
149
150 Int_t nofHits = fHitsCollection->GetEntriesFast();
151 A01DriftChamberHit* hit =
152 new ((*fHitsCollection)[nofHits]) A01DriftChamberHit(copyNo);
153
154 hit->SetWorldPos(TVector3(globalPos[0], globalPos[1], globalPos[2]));
155 hit->SetLocalPos(TVector3(localPos[0], localPos[1], localPos[2]));
156 hit->SetTime(gMC->TrackTime());
157
158 return true;
159}
160
161//_____________________________________________________________________________
163{
164 /// Print hits collection (if verbose) and reset hits afterwards.
165
166 if (fVerboseLevel > 0) Print();
167
168 // Reset hits collection
169 fHitsCollection->Clear();
170}
171
172//_____________________________________________________________________________
174{
175 /// Register the hits collection in Root manager.
176
177 if (fWriteHits) {
178 TMCRootManager::Instance()->Register(
179 GetName(), "TClonesArray", &fHitsCollection);
180 }
181}
182
183//_____________________________________________________________________________
184void A01DriftChamberSD::Print(Option_t* /*option*/) const
185{
186 /// Print the hits collection.
187
188 Int_t nofHits = fHitsCollection->GetEntriesFast();
189 cout << GetName() << " has " << nofHits << " hits." << endl;
190 if (fVerboseLevel > 1) {
191 for (Int_t i2 = 0; i2 < 5; ++i2) {
192 for (Int_t i = 0; i < nofHits; ++i) {
193 A01DriftChamberHit* hit = GetHit(i);
194 if (hit->GetLayerID() == i2) hit->Print();
195 }
196 }
197 }
198}
Definition of the A01DriftChamberHit class.
Definition of the A01DriftChamberSD class.
The drift chamber hit.
void SetLocalPos(const TVector3 &pos)
virtual void Print(Option_t *option="") const
void SetTime(Double_t t)
void SetWorldPos(const TVector3 &pos)
The calorimeter sensitive detector.
Bool_t fWriteHits
Option to write hits.
Int_t fVolId
The calorimeter volume Id.
TClonesArray * fHitsCollection
Hits collection.
A01DriftChamberHit * GetHit(Int_t i) const
Int_t fVerboseLevel
Verbosity level.
TString fVolName
The sensitive volume name.
virtual void Print(Option_t *option="") const