Loading...
Searching...
No Matches
DriftChamberHit.hh
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//
27/// \file B5/include/DriftChamberHit.hh
28/// \brief Definition of the B5::DriftChamberHit class
29
30#ifndef B5DriftChamberHit_h
31#define B5DriftChamberHit_h 1
32
33#include "G4VHit.hh"
34#include "G4THitsCollection.hh"
35#include "G4Allocator.hh"
36#include "G4ThreeVector.hh"
37#include "G4LogicalVolume.hh"
38#include "G4Transform3D.hh"
39#include "G4RotationMatrix.hh"
40
41class G4AttDef;
42class G4AttValue;
43
44namespace B5
45{
46
47/// Drift chamber hit
48///
49/// It records:
50/// - the layer ID
51/// - the particle time
52/// - the particle local and global positions
53
54class DriftChamberHit : public G4VHit
55{
56 public:
57 DriftChamberHit() = default;
58 DriftChamberHit(G4int layerID);
60 ~DriftChamberHit() override = default;
61
63 G4bool operator==(const DriftChamberHit &right) const;
64
65 inline void *operator new(size_t);
66 inline void operator delete(void *aHit);
67
68 void Draw() override;
69 const std::map<G4String,G4AttDef>* GetAttDefs() const override;
70 std::vector<G4AttValue>* CreateAttValues() const override;
71 void Print() override;
72
73 void SetLayerID(G4int z) { fLayerID = z; }
74 G4int GetLayerID() const { return fLayerID; }
75
76 void SetTime(G4double t) { fTime = t; }
77 G4double GetTime() const { return fTime; }
78
79 void SetLocalPos(G4ThreeVector xyz) { fLocalPos = xyz; }
80 G4ThreeVector GetLocalPos() const { return fLocalPos; }
81
82 void SetWorldPos(G4ThreeVector xyz) { fWorldPos = xyz; }
83 G4ThreeVector GetWorldPos() const { return fWorldPos; }
84
85 private:
86 G4int fLayerID = -1;
87 G4double fTime = 0.;
88 G4ThreeVector fLocalPos;
89 G4ThreeVector fWorldPos;
90};
91
93
95
96inline void* DriftChamberHit::operator new(size_t)
97{
100 }
101 return (void*)DriftChamberHitAllocator->MallocSingle();
102}
103
104inline void DriftChamberHit::operator delete(void* aHit)
105{
106 DriftChamberHitAllocator->FreeSingle((DriftChamberHit*) aHit);
107}
108
109}
110
111//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
112
113#endif
Drift chamber hit.
~DriftChamberHit() override=default
G4double GetTime() const
G4ThreeVector fWorldPos
void SetLayerID(G4int z)
G4ThreeVector GetLocalPos() const
G4ThreeVector fLocalPos
DriftChamberHit & operator=(const DriftChamberHit &right)=default
void Print() override
void SetLocalPos(G4ThreeVector xyz)
std::vector< G4AttValue > * CreateAttValues() const override
DriftChamberHit(const DriftChamberHit &right)=default
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4ThreeVector GetWorldPos() const
G4bool operator==(const DriftChamberHit &right) const
G4int GetLayerID() const
void SetTime(G4double t)
void SetWorldPos(G4ThreeVector xyz)
void Draw() override
DriftChamberHit()=default
G4ThreadLocal G4Allocator< DriftChamberHit > * DriftChamberHitAllocator
G4THitsCollection< DriftChamberHit > DriftChamberHitsCollection

Applications | User Support | Publications | Collaboration