Loading...
Searching...
No Matches
HadCalorimeterHit.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/HadCalorimeterHit.cc
28/// \brief Implementation of the B5::HadCalorimeterHit class
29
30#include "HadCalorimeterHit.hh"
31#include "DetectorConstruction.hh"
32
33#include "G4VVisManager.hh"
34#include "G4VisAttributes.hh"
35#include "G4Box.hh"
36#include "G4Colour.hh"
37#include "G4AttDefStore.hh"
38#include "G4AttDef.hh"
39#include "G4AttValue.hh"
40#include "G4UIcommand.hh"
41#include "G4UnitsTable.hh"
42#include "G4SystemOfUnits.hh"
43#include "G4ios.hh"
44
45namespace B5
46{
47
48//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
49
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
53
54HadCalorimeterHit::HadCalorimeterHit(G4int columnID,G4int rowID)
55: fColumnID(columnID), fRowID(rowID)
56{}
57
58//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
59
61{
62 return ( fColumnID==right.fColumnID && fRowID==right.fRowID );
63}
64
65//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66
68{
69 auto visManager = G4VVisManager::GetConcreteInstance();
70 if (! visManager || (fEdep==0.)) return;
71
72 // Draw a calorimeter cell with depth propotional to the energy deposition
73 G4Transform3D trans(fRot.inverse(),fPos);
74 G4VisAttributes attribs;
75 attribs.SetColour(G4Colour::Red());
76 attribs.SetForceSolid(true);
77 G4Box box("dummy",15.*cm,15.*cm,1.*m*fEdep/(0.1*GeV));
78 visManager->Draw(box,attribs,trans);
79}
80
81//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
82
83const std::map<G4String,G4AttDef>* HadCalorimeterHit::GetAttDefs() const
84{
85 G4bool isNew;
86 auto store = G4AttDefStore::GetInstance("HadCalorimeterHit",isNew);
87
88 if (isNew) {
89 (*store)["HitType"]
90 = G4AttDef("HitType","Hit Type","Physics","","G4String");
91
92 (*store)["Column"]
93 = G4AttDef("Column","Column ID","Physics","","G4int");
94
95 (*store)["Row"]
96 = G4AttDef("Row","Row ID","Physics","","G4int");
97
98 (*store)["Energy"]
99 = G4AttDef("Energy","Energy Deposited","Physics","G4BestUnit",
100 "G4double");
101
102 (*store)["Pos"]
103 = G4AttDef("Pos", "Position", "Physics","G4BestUnit",
104 "G4ThreeVector");
105 }
106 return store;
107}
108
109//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
110
111std::vector<G4AttValue>* HadCalorimeterHit::CreateAttValues() const
112{
113 auto values = new std::vector<G4AttValue>;
114
115 values
116 ->push_back(G4AttValue("HitType","HadCalorimeterHit",""));
117 values
118 ->push_back(G4AttValue("Column",G4UIcommand::ConvertToString(fColumnID),
119 ""));
120 values
121 ->push_back(G4AttValue("Row",G4UIcommand::ConvertToString(fRowID),""));
122 values
123 ->push_back(G4AttValue("Energy",G4BestUnit(fEdep,"Energy"),""));
124 values
125 ->push_back(G4AttValue("Pos",G4BestUnit(fPos,"Length"),""));
126
127 return values;
128}
129
130//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
131
133{
134 G4cout << " Cell[" << fRowID << ", " << fColumnID << "] "
135 << fEdep/MeV << " (MeV) " << fPos << G4endl;
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
139
140}
Definition of the B5::HadCalorimeterHit class.
Hadron Calorimeter hit.
HadCalorimeterHit()=default
G4bool operator==(const HadCalorimeterHit &right) const
std::vector< G4AttValue > * CreateAttValues() const override
const std::map< G4String, G4AttDef > * GetAttDefs() const override
G4ThreadLocal G4Allocator< HadCalorimeterHit > * HadCalorimeterHitAllocator

Applications | User Support | Publications | Collaboration