ATLAS Offline Software
Public Member Functions | Private Attributes | Static Private Attributes | List of all members
MuonWallSD Class Reference

#include <MuonWallSD.h>

Inheritance diagram for MuonWallSD:
Collaboration diagram for MuonWallSD:

Public Member Functions

 MuonWallSD (const std::string &name, const std::string &hitCollectionName, int verbose)
 
 ~MuonWallSD ()
 
void StartOfAthenaEvent ()
 
void Initialize (G4HCofThisEvent *) override final
 
G4bool ProcessHits (G4Step *, G4TouchableHistory *) override final
 
void EndOfAthenaEvent ()
 

Private Attributes

const TileTBIDm_tileTBID
 
int m_nhits [s_nCell]
 
TileSimHitm_hit [s_nCell]
 
Identifier m_id [s_nCell]
 
SG::WriteHandle< TileHitVectorm_HitColl
 

Static Private Attributes

static const int s_nCellMu = 14
 
static const int s_nCellS = 4
 
static const int s_nCell = s_nCellMu+s_nCellS
 

Detailed Description

Definition at line 32 of file MuonWallSD.h.

Constructor & Destructor Documentation

◆ MuonWallSD()

MuonWallSD::MuonWallSD ( const std::string &  name,
const std::string &  hitCollectionName,
int  verbose 
)

Definition at line 27 of file MuonWallSD.cxx.

28  : G4VSensitiveDetector(name),
29  m_nhits(),
30  m_hit(),
31  m_HitColl(hitCollectionName)
32 {
33  verboseLevel = std::max(verboseLevel, verbose);
34 
35  ISvcLocator* svcLocator = Gaudi::svcLocator();
36 
37  StoreGateSvc* detStore(nullptr);
38  if (svcLocator->service("DetectorStore", detStore).isFailure()) {
39  G4ExceptionDescription description;
40  description << "Constructor: DetectorStoreSvc not found!";
41  G4Exception("MuonWallSD", "NoDetStore", FatalException, description);
42  abort();
43  } else if (verboseLevel >= 5) {
44  G4cout << "DetectorStoreSvc initialized" << G4endl;
45  }
46 
47  if (detStore->retrieve(m_tileTBID).isFailure()) {
48  G4ExceptionDescription description;
49  description << "Constructor: No TileTBID helper!";
50  G4Exception("MuonWallSD", "NoTileTBIDHelper", FatalException, description);
51  abort();
52  } else if (verboseLevel >= 5) {
53  G4cout << "TileTBID helper retrieved" << G4endl;
54  }
55 
57  for (int channel=0; channel<s_nCellMu; ++channel) {
59  }
60 
62  for (int channel=0; channel<s_nCellS; ++channel) {
64  }
65 }

◆ ~MuonWallSD()

MuonWallSD::~MuonWallSD ( )

Definition at line 67 of file MuonWallSD.cxx.

67  {
68 }

Member Function Documentation

◆ EndOfAthenaEvent()

void MuonWallSD::EndOfAthenaEvent ( )

Definition at line 147 of file MuonWallSD.cxx.

147  {
148  for (int ind = 0; ind < s_nCell; ++ind) {
149  int nhit = m_nhits[ind];
150  if (nhit > 0) {
151  if (verboseLevel >= 5) {
152  G4cout << "Cell id=" << m_tileTBID->to_string(m_id[ind])
153  << " nhit=" << nhit
154  << " ene=" << m_hit[ind]->energy() << G4endl;
155  }
156  m_HitColl->Insert(TileHit(m_hit[ind]));
157  delete m_hit[ind];
158  } else if (verboseLevel >= 10) {
159  G4cout << "Cell id=" << m_tileTBID->to_string(m_id[ind])
160  << " nhit=0" << G4endl;
161  }
162  }
163 
164  if (verboseLevel >= 5) {
165  G4cout << "Total number of hits is " << m_HitColl->size() << G4endl;
166  }
167  return ;
168 }

◆ Initialize()

void MuonWallSD::Initialize ( G4HCofThisEvent *  )
finaloverride

Definition at line 78 of file MuonWallSD.cxx.

78  {
79  if (verboseLevel >= 5) {
80  G4cout << "MuonWallSD::Initialize()" << G4endl;
81  }
82 
83  if (!m_HitColl.isValid()) {
84  m_HitColl = std::make_unique<TileHitVector>(m_HitColl.name());
85  }
86 }

◆ ProcessHits()

G4bool MuonWallSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
finaloverride

Definition at line 88 of file MuonWallSD.cxx.

88  {
89  if (verboseLevel >= 10) {
90  G4cout << "MuonWallSD::ProcessHits" << G4endl;
91  }
92 
93  const G4TouchableHistory* theTouchable = (G4TouchableHistory*)(aStep->GetPreStepPoint()->GetTouchable());
94  const G4VPhysicalVolume* physVol = theTouchable->GetVolume();
95  const G4LogicalVolume* logiVol = physVol->GetLogicalVolume();
96  const G4String nameLogiVol = logiVol->GetName();
97  const G4int nScinti = physVol->GetCopyNo();
98 
99  const G4double edep = aStep->GetTotalEnergyDeposit() * aStep->GetTrack()->GetWeight();
100  G4double stepl = 0.;
101 
102  if (aStep->GetTrack()->GetDefinition()->GetPDGCharge() != 0.){ // FIXME not-equal check on double
103 
104  stepl = aStep->GetStepLength();
105  }
106 
107  if ((edep == 0.) && (stepl == 0.)) { //FIXME equality check on double
108 
109  return false;
110  }
111 
112  int ind;
113 
114  if(nameLogiVol.find("MuScintillatorLayer") !=G4String::npos) {
115  // for muon wall, nScinti-1 is the correct indice.
116  ind = nScinti-1;
117  } else if(nameLogiVol.find("S1") !=G4String::npos) {
118  ind = s_nCellMu + 0;
119  } else if(nameLogiVol.find("S2") !=G4String::npos) {
120  ind = s_nCellMu + 1;
121  } else if(nameLogiVol.find("S3") !=G4String::npos) {
122  ind = s_nCellMu + 2;
123  } else {
124  ind = s_nCellMu + 3;
125  }
126 
127  if (verboseLevel >= 10) {
128  G4cout << ((m_nhits[ind] > 0)?"Additional hit in ":"First hit in ")
129  << ((ind<s_nCellMu)?"MuonWall ":"beam counter S")
130  << ((ind<s_nCellMu)?nScinti:(ind-s_nCellMu+1))
131  << " time=" << aStep->GetPostStepPoint()->GetGlobalTime()
132  << " ene=" << edep << G4endl;
133  }
134 
135  if ( m_nhits[ind] > 0 ) {
136  m_hit[ind]->add(edep,0.0,0.0);
137  } else {
138  // First hit in a cell
139  m_hit[ind] = new TileSimHit(m_id[ind],edep,0.0,0.0);
140  }
141 
142  ++m_nhits[ind];
143 
144  return true;
145 }

◆ StartOfAthenaEvent()

void MuonWallSD::StartOfAthenaEvent ( )

Definition at line 70 of file MuonWallSD.cxx.

70  {
71  if (verboseLevel >= 5) {
72  G4cout << "Initializing SD" << G4endl;
73  }
74 
75  memset(m_nhits, 0, sizeof(m_nhits));
76 }

Member Data Documentation

◆ m_hit

TileSimHit* MuonWallSD::m_hit[s_nCell]
private

Definition at line 52 of file MuonWallSD.h.

◆ m_HitColl

SG::WriteHandle<TileHitVector> MuonWallSD::m_HitColl
private

Definition at line 55 of file MuonWallSD.h.

◆ m_id

Identifier MuonWallSD::m_id[s_nCell]
private

Definition at line 53 of file MuonWallSD.h.

◆ m_nhits

int MuonWallSD::m_nhits[s_nCell]
private

Definition at line 51 of file MuonWallSD.h.

◆ m_tileTBID

const TileTBID* MuonWallSD::m_tileTBID
private

Definition at line 45 of file MuonWallSD.h.

◆ s_nCell

const int MuonWallSD::s_nCell = s_nCellMu+s_nCellS
staticprivate

Definition at line 49 of file MuonWallSD.h.

◆ s_nCellMu

const int MuonWallSD::s_nCellMu = 14
staticprivate

Definition at line 47 of file MuonWallSD.h.

◆ s_nCellS

const int MuonWallSD::s_nCellS = 4
staticprivate

Definition at line 48 of file MuonWallSD.h.


The documentation for this class was generated from the following files:
TileTBID::channel_id
Identifier channel_id(int type, int module, int channel) const
identifer for one channel of a Tile testbeam detector
Definition: TileTBID.cxx:196
TileTBID::to_string
std::string to_string(const Identifier &id, int level=0) const
extract all fields from TileTB identifier Identifier get_all_fields ( const Identifier & id,...
Definition: TileTBID.cxx:48
MuonWallSD::m_nhits
int m_nhits[s_nCell]
Definition: MuonWallSD.h:51
max
#define max(a, b)
Definition: cfImp.cxx:41
plotting.yearwise_efficiency.channel
channel
Definition: yearwise_efficiency.py:28
MuonWallSD::m_hit
TileSimHit * m_hit[s_nCell]
Definition: MuonWallSD.h:52
TileSimHit
Definition: TileSimHit.h:31
TileTBID::ADC_TYPE
@ ADC_TYPE
Definition: Calorimeter/CaloIdentifier/CaloIdentifier/TileTBID.h:84
MuonWallSD::m_tileTBID
const TileTBID * m_tileTBID
Definition: MuonWallSD.h:45
TileSimHit::add
int add(double energy, double time, double deltaT)
Add sub-hit to a given hit with time rounding to the center of nearest deltaT bin.
Definition: TileSimHit.cxx:41
python.PyAthena.module
module
Definition: PyAthena.py:134
MuonWallSD::m_id
Identifier m_id[s_nCell]
Definition: MuonWallSD.h:53
StoreGateSvc
The Athena Transient Store API.
Definition: StoreGateSvc.h:128
MuonWallSD::m_HitColl
SG::WriteHandle< TileHitVector > m_HitColl
Definition: MuonWallSD.h:55
MuonWallSD::s_nCell
static const int s_nCell
Definition: MuonWallSD.h:49
python.PyKernel.detStore
detStore
Definition: PyKernel.py:41
TileSimHit::energy
double energy(int ind=0) const
Return energy of ind-th sub-hit
Definition: TileSimHit.h:51
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
TileHit
Definition: TileSimEvent/TileSimEvent/TileHit.h:30
MuonWallSD::s_nCellMu
static const int s_nCellMu
Definition: MuonWallSD.h:47
TileTBID::CRACK_WALL
@ CRACK_WALL
Definition: Calorimeter/CaloIdentifier/CaloIdentifier/TileTBID.h:93
TileTBID::BACK_WALL
@ BACK_WALL
Definition: Calorimeter/CaloIdentifier/CaloIdentifier/TileTBID.h:91
MuonWallSD::s_nCellS
static const int s_nCellS
Definition: MuonWallSD.h:48
python.CaloScaleNoiseConfig.type
type
Definition: CaloScaleNoiseConfig.py:78
python.TriggerHandler.verbose
verbose
Definition: TriggerHandler.py:297
checkFileSG.ind
list ind
Definition: checkFileSG.py:118
description
std::string description
glabal timer - how long have I taken so far?
Definition: hcg.cxx:88