Loading...
Searching...
No Matches
DriftChamberHit.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//
27/// \file B5/src/DriftChamberHit.cc
28/// \brief Implementation of the B5::DriftChamberHit class
29
30#include "DriftChamberHit.hh"
31
32#include "G4VVisManager.hh"
33#include "G4VisAttributes.hh"
34#include "G4Circle.hh"
35#include "G4Colour.hh"
36#include "G4AttDefStore.hh"
37#include "G4AttDef.hh"
38#include "G4AttValue.hh"
39#include "G4UIcommand.hh"
40#include "G4UnitsTable.hh"
41#include "G4SystemOfUnits.hh"
42#include "G4ios.hh"
43
44namespace B5
45{
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
54: fLayerID(layerID)
55{}
56
57//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
58
59G4bool DriftChamberHit::operator==(const DriftChamberHit &/*right*/) const
60{
61 return false;
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
65
67{
68 auto visManager = G4VVisManager::GetConcreteInstance();
69 if (! visManager) return;
70
71 G4Circle circle(fWorldPos);
72 circle.SetScreenSize(2);
73 circle.SetFillStyle(G4Circle::filled);
74 G4VisAttributes attribs(G4Colour::Yellow());
75 circle.SetVisAttributes(attribs);
76 visManager->Draw(circle);
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
81const std::map<G4String,G4AttDef>* DriftChamberHit::GetAttDefs() const
82{
83 G4bool isNew;
84 auto store = G4AttDefStore::GetInstance("DriftChamberHit",isNew);
85
86 if (isNew) {
87 (*store)["HitType"]
88 = G4AttDef("HitType","Hit Type","Physics","","G4String");
89
90 (*store)["ID"]
91 = G4AttDef("ID","ID","Physics","","G4int");
92
93 (*store)["Time"]
94 = G4AttDef("Time","Time","Physics","G4BestUnit","G4double");
95
96 (*store)["Pos"]
97 = G4AttDef("Pos", "Position", "Physics","G4BestUnit","G4ThreeVector");
98 }
99
100 return store;
101}
102
103//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
104
105std::vector<G4AttValue>* DriftChamberHit::CreateAttValues() const
106{
107 auto values = new std::vector<G4AttValue>;
108
109 values
110 ->push_back(G4AttValue("HitType","DriftChamberHit",""));
111 values
112 ->push_back(G4AttValue("ID",G4UIcommand::ConvertToString(fLayerID),""));
113 values
114 ->push_back(G4AttValue("Time",G4BestUnit(fTime,"Time"),""));
115 values
116 ->push_back(G4AttValue("Pos",G4BestUnit(fWorldPos,"Length"),""));
117
118 return values;
119}
120
121//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
122
124{
125 G4cout << " Layer[" << fLayerID << "] : time " << fTime/ns
126 << " (nsec) --- local (x,y) " << fLocalPos.x()
127 << ", " << fLocalPos.y() << G4endl;
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131
132}
Definition of the B5::DriftChamberHit class.
Drift chamber hit.
G4ThreeVector fWorldPos
G4ThreeVector fLocalPos
void Print() override
std::vector< G4AttValue > * CreateAttValues() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4bool operator==(const DriftChamberHit &right) const
void Draw() override
DriftChamberHit()=default
G4ThreadLocal G4Allocator< DriftChamberHit > * DriftChamberHitAllocator

Applications | User Support | Publications | Collaboration