Loading...
Searching...
No Matches
RE05CalorimeterHit.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 RE05/src/RE05CalorimeterHit.cc
28/// \brief Implementation of the RE05CalorimeterHit class
29//
30
31#include "RE05CalorimeterHit.hh"
32#include "G4ios.hh"
33#include "G4VVisManager.hh"
34#include "G4Colour.hh"
35#include "G4VisAttributes.hh"
36#include "G4LogicalVolume.hh"
37#include "G4UIcommand.hh"
38#include "G4UnitsTable.hh"
39#include "G4AttValue.hh"
40#include "G4AttDef.hh"
41#include "G4AttCheck.hh"
42
44
45std::map<G4String,G4AttDef> RE05CalorimeterHit::fAttDefs;
46
47//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
48
50: G4VHit(),
51 fZCellID(-1),fPhiCellID(-1),fEdep(0.),fPos(),fRot(),fLogV(0)
52{}
53
54//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
55
57: G4VHit(),
58 fZCellID(z), fPhiCellID(phi),fEdep(0.),fPos(),fRot(),fLogV(logVol)
59{}
60
61//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
62
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
69 : G4VHit()
70{
71 fZCellID = right.fZCellID;
72 fPhiCellID = right.fPhiCellID;
73 fEdep = right.fEdep;
74 fPos = right.fPos;
75 fRot = right.fRot;
76 fLogV = right.fLogV;
77}
78
79//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
80
82{
84 fPhiCellID = right.fPhiCellID;
85 fEdep = right.fEdep;
86 fPos = right.fPos;
87 fRot = right.fRot;
88 fLogV = right.fLogV;
89 return *this;
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
95{
96 return ((fZCellID==right.fZCellID)&&(fPhiCellID==right.fPhiCellID));
97}
98
99//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
100
102{
103 G4VVisManager* pVVisManager = G4VVisManager::GetConcreteInstance();
104 if(pVVisManager)
105 {
106 G4Transform3D trans(fRot,fPos);
107 G4VisAttributes attribs;
108 const G4VisAttributes* pVA = fLogV->GetVisAttributes();
109 if(pVA) attribs = *pVA;
110 G4Colour colour(1.,0.,0.);
111 attribs.SetColour(colour);
112 attribs.SetForceSolid(true);
113 pVVisManager->Draw(*fLogV,attribs,trans);
114 }
115}
116
117//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
118
119const std::map<G4String,G4AttDef>* RE05CalorimeterHit::GetAttDefs() const
120{
121 // G4AttDefs have to have long life. Use static member...
122 if (fAttDefs.empty()) {
123 fAttDefs["HitType"] =
124 G4AttDef("HitType","Type of hit","Physics","","G4String");
125 fAttDefs["ZID"] = G4AttDef("ZID","Z Cell ID","Physics","","G4int");
126 fAttDefs["PhiID"] = G4AttDef("PhiID","Phi Cell ID","Physics","","G4int");
127 fAttDefs["EDep"] =
128 G4AttDef("EDep","Energy defPosited","Physics","G4BestUnit","G4double");
129 }
130 return &fAttDefs;
131}
132
133//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
134
135std::vector<G4AttValue>* RE05CalorimeterHit::CreateAttValues() const
136{
137 // Create expendable G4AttsValues for picking...
138 std::vector<G4AttValue>* attValues = new std::vector<G4AttValue>;
139 attValues->push_back
140 (G4AttValue("HitType","RE05CalorimeterHit",""));
141 attValues->push_back
142 (G4AttValue("ZID",G4UIcommand::ConvertToString(fZCellID),""));
143 attValues->push_back
144 (G4AttValue("PhiID",G4UIcommand::ConvertToString(fPhiCellID),""));
145 attValues->push_back
146 (G4AttValue("EDep",G4BestUnit(fEdep,"Energy"),""));
147 //G4cout << "Checking...\n" << G4AttCheck(attValues, GetAttDefs());
148 return attValues;
149}
150
151//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
152
155
156//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
157
G4ThreadLocal G4Allocator< RE05CalorimeterHit > * RE05CalorimeterHitAllocator
Definition of the RE05CalorimeterHit class.
G4bool operator==(const RE05CalorimeterHit &right) const
const RE05CalorimeterHit & operator=(const RE05CalorimeterHit &right)
const G4LogicalVolume * fLogV
virtual const std::map< G4String, G4AttDef > * GetAttDefs() const
virtual std::vector< G4AttValue > * CreateAttValues() const
static std::map< G4String, G4AttDef > fAttDefs
G4RotationMatrix fRot

Applications | User Support | Publications | Collaboration