Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
B5DriftChamberSD.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: B5DriftChamberSD.cc 76474 2013-11-11 10:36:34Z gcosmo $
27 //
28 /// \file B5DriftChamberSD.cc
29 /// \brief Implementation of the B5DriftChamber class
30 
31 #include "B5DriftChamberSD.hh"
32 #include "B5DriftChamberHit.hh"
33 
34 #include "G4HCofThisEvent.hh"
35 #include "G4TouchableHistory.hh"
36 #include "G4Track.hh"
37 #include "G4Step.hh"
38 #include "G4SDManager.hh"
39 #include "G4ios.hh"
40 
41 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
42 
44 : G4VSensitiveDetector(name), fHitsCollection(0), fHCID(-1)
45 {
46  collectionName.insert("driftChamberColl");
47 }
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
50 
52 {}
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55 
57 {
58  fHitsCollection
60  if (fHCID<0)
61  { fHCID = G4SDManager::GetSDMpointer()->GetCollectionID(fHitsCollection); }
62  hce->AddHitsCollection(fHCID,fHitsCollection);
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
68 {
69  G4double charge = step->GetTrack()->GetDefinition()->GetPDGCharge();
70  if (charge==0.) return true;
71 
72  G4StepPoint* preStepPoint = step->GetPreStepPoint();
73 
74  G4TouchableHistory* touchable
76  G4VPhysicalVolume* motherPhysical = touchable->GetVolume(1); // mother
77  G4int copyNo = motherPhysical->GetCopyNo();
78 
79  G4ThreeVector worldPos = preStepPoint->GetPosition();
80  G4ThreeVector localPos
81  = touchable->GetHistory()->GetTopTransform().TransformPoint(worldPos);
82 
83  B5DriftChamberHit* hit = new B5DriftChamberHit(copyNo);
84  hit->SetWorldPos(worldPos);
85  hit->SetLocalPos(localPos);
86  hit->SetTime(preStepPoint->GetGlobalTime());
87 
88  fHitsCollection->insert(hit);
89 
90  return true;
91 }
92 
93 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
virtual G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist)
G4ParticleDefinition * GetDefinition() const
G4THitsCollection< B5DriftChamberHit > B5DriftChamberHitsCollection
G4int GetCollectionID(G4String colName)
Definition: G4SDManager.cc:131
Definition of the B5DriftChamberSD class.
G4VPhysicalVolume * GetVolume(G4int depth=0) const
virtual ~B5DriftChamberSD()
const XML_Char * name
const G4VTouchable * GetTouchable() const
void SetTime(G4double t)
void SetLocalPos(G4ThreeVector xyz)
int G4int
Definition: G4Types.hh:78
void SetWorldPos(G4ThreeVector xyz)
virtual void Initialize(G4HCofThisEvent *HCE)
B5DriftChamberSD(G4String name)
G4StepPoint * GetPreStepPoint() const
const G4ThreeVector & GetPosition() const
bool G4bool
Definition: G4Types.hh:79
Definition: G4Step.hh:76
void AddHitsCollection(G4int HCID, G4VHitsCollection *aHC)
G4ThreeVector TransformPoint(const G4ThreeVector &vec) const
static G4SDManager * GetSDMpointer()
Definition: G4SDManager.cc:40
virtual G4int GetCopyNo() const =0
G4double GetGlobalTime() const
const G4AffineTransform & GetTopTransform() const
G4CollectionNameVector collectionName
const G4NavigationHistory * GetHistory() const
double G4double
Definition: G4Types.hh:76
G4Track * GetTrack() const
G4double GetPDGCharge() const
Definition of the B5DriftChamberHit class.