Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Private Attributes | List of all members
RE01CalorimeterSD Class Reference

#include <Doxymodules_runAndEvent.h>

Inheritance diagram for RE01CalorimeterSD:
G4VSensitiveDetector

Public Member Functions

 RE01CalorimeterSD (G4String name)
 
virtual ~RE01CalorimeterSD ()
 
virtual void Initialize (G4HCofThisEvent *HCE)
 

Protected Member Functions

virtual G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist)
 

Private Attributes

RE01CalorimeterHitsCollectionfCalCollection
 
int fCellID [20][48]
 
const int fNumberOfCellsInZ
 
const int fNumberOfCellsInPhi
 

Detailed Description

Definition at line 22 of file Doxymodules_runAndEvent.h.

Constructor & Destructor Documentation

◆ RE01CalorimeterSD()

RE01CalorimeterSD::RE01CalorimeterSD ( G4String  name)

Definition at line 43 of file RE01CalorimeterSD.cc.

46{
47 // Initialize data member.
48 for(G4int j=0;j<fNumberOfCellsInZ;j++)
49 for(G4int k=0;k<fNumberOfCellsInPhi;k++)
50 {
51 fCellID[j][k] = -1;
52 }
53 //
54 G4String HCname;
55 collectionName.insert(HCname="calCollection");
56}
RE01CalorimeterHitsCollection * fCalCollection

◆ ~RE01CalorimeterSD()

RE01CalorimeterSD::~RE01CalorimeterSD ( )
virtual

Definition at line 59 of file RE01CalorimeterSD.cc.

60{;}

Member Function Documentation

◆ Initialize()

void RE01CalorimeterSD::Initialize ( G4HCofThisEvent HCE)
virtual

Definition at line 63 of file RE01CalorimeterSD.cc.

64{
66 (SensitiveDetectorName,collectionName[0]);
67 for(G4int j=0;j<fNumberOfCellsInZ;j++)
68 for(G4int k=0;k<fNumberOfCellsInPhi;k++)
69 {
70 fCellID[j][k] = -1;
71 }
72
73 static G4int HCID = -1;
74 if(HCID<0)
75 { HCID = GetCollectionID(0); }
76 HCE->AddHitsCollection( HCID, fCalCollection );
77}
G4THitsCollection< RE01CalorimeterHit > RE01CalorimeterHitsCollection

◆ ProcessHits()

G4bool RE01CalorimeterSD::ProcessHits ( G4Step aStep,
G4TouchableHistory ROhist 
)
protectedvirtual

Definition at line 80 of file RE01CalorimeterSD.cc.

81{
82//***** RE05CalorimeterSD has been migrated to Geant4 version 10 that does not
83//***** support Readout Geometry in multi-threaded mode. Now RE05CalorimeterSD
84//***** is assigned to a dedicaed parallel world. The pointer "aStep" points to
85//***** a G4Step object for the parallel world.
86
87 G4double edep = aStep->GetTotalEnergyDeposit();
88 if(edep==0.) return false;
89
90 const G4VTouchable* ROhist = aStep->GetPreStepPoint()->GetTouchable();
91 G4int copyIDinZ = ROhist->GetReplicaNumber();
92 G4int copyIDinPhi = ROhist->GetReplicaNumber(1);
93
94 if(fCellID[copyIDinZ][copyIDinPhi]==-1)
95 {
97 (ROhist->GetVolume()->GetLogicalVolume(),copyIDinZ,copyIDinPhi);
98 calHit->SetEdep( edep );
99 G4AffineTransform aTrans = ROhist->GetHistory()->GetTopTransform();
100 aTrans.Invert();
101 calHit->SetPos(aTrans.NetTranslation());
102 calHit->SetRot(aTrans.NetRotation());
103 calHit->SetTrackInformation(aStep->GetTrack());
104 G4int icell = fCalCollection->insert( calHit );
105 fCellID[copyIDinZ][copyIDinPhi] = icell - 1;
106 if(verboseLevel>0)
107 { G4cout << " New Calorimeter Hit on CellID "
108 << copyIDinZ << " " << copyIDinPhi << G4endl; }
109 }
110 else
111 {
112 (*fCalCollection)[fCellID[copyIDinZ][copyIDinPhi]]->AddEdep(edep);
113 (*fCalCollection)[fCellID[copyIDinZ][copyIDinPhi]]
114 ->SetTrackInformation(aStep->GetTrack());
115 if(verboseLevel>0)
116 { G4cout << " Energy added to CellID "
117 << copyIDinZ << " " << copyIDinPhi << G4endl; }
118 }
119
120 return true;
121}
void SetEdep(G4double de)
void SetPos(G4ThreeVector xyz)
void SetTrackInformation(const G4Track *aTrack)
void SetRot(G4RotationMatrix rmat)

Member Data Documentation

◆ fCalCollection

RE01CalorimeterHitsCollection* RE01CalorimeterSD::fCalCollection
private

Definition at line 54 of file RE01CalorimeterSD.hh.

◆ fCellID

int RE01CalorimeterSD::fCellID[20][48]
private

Definition at line 55 of file RE01CalorimeterSD.hh.

◆ fNumberOfCellsInZ

const int RE01CalorimeterSD::fNumberOfCellsInZ
private

Definition at line 56 of file RE01CalorimeterSD.hh.

◆ fNumberOfCellsInPhi

const int RE01CalorimeterSD::fNumberOfCellsInPhi
private

Definition at line 57 of file RE01CalorimeterSD.hh.


The documentation for this class was generated from the following files:

Applications | User Support | Publications | Collaboration