ATLAS Offline Software
Loading...
Searching...
No Matches
FCS_StepInfoSD Class Reference

Common sensitive detector class for LAr systems. More...

#include <FCS_StepInfoSD.h>

Inheritance diagram for FCS_StepInfoSD:
Collaboration diagram for FCS_StepInfoSD:

Public Member Functions

 FCS_StepInfoSD (G4String a_name, const FCS_Param::Config &config)
 Constructor.
virtual G4bool ProcessHits (G4Step *a_step, G4TouchableHistory *) override
 Main processing method.
void EndOfAthenaEvent (ISF_FCS_Parametrization::FCS_StepInfoCollection *hitContnainer)
 End of athena event processing.
void setupHelpers (const LArEM_ID *EM, const LArFCAL_ID *FCAL, const LArHEC_ID *HEC, const TileID *tile)
 Sets the ID helper pointers.

Protected Member Functions

void getCaloDDManager ()
 Keep a map instead of trying to keep the full vector.
void update_map (const CLHEP::Hep3Vector &l_vec, const Identifier &l_identifier, double l_energy, double l_time, bool l_valid, int l_detector)

Protected Attributes

FCS_Param::Config m_config
const LArEM_IDm_larEmID {nullptr}
 Pointers to the identifier helpers.
const LArFCAL_IDm_larFcalID {nullptr}
const LArHEC_IDm_larHecID {nullptr}
const TileIDm_tileID {nullptr}
CxxUtils::CachedPointer< const CaloDetDescrManagerm_calo_dd_man
std::map< Identifier, std::vector< ISF_FCS_Parametrization::FCS_StepInfo * > * > m_hit_map

Private Member Functions

double getMaxTime (const CaloCell_ID::CaloSample &layer) const

Detailed Description

Common sensitive detector class for LAr systems.

This SD implementation saves the standard LArHits. See LArG4CalibSD for an SD that handles calibration hits.

Definition at line 71 of file FCS_StepInfoSD.h.

Constructor & Destructor Documentation

◆ FCS_StepInfoSD()

FCS_StepInfoSD::FCS_StepInfoSD ( G4String a_name,
const FCS_Param::Config & config )

Constructor.

Definition at line 21 of file FCS_StepInfoSD.cxx.

22 : G4VSensitiveDetector(std::move(a_name)),
23 m_config(config),
24 m_calo_dd_man(nullptr) {}
FCS_Param::Config m_config
CxxUtils::CachedPointer< const CaloDetDescrManager > m_calo_dd_man

Member Function Documentation

◆ EndOfAthenaEvent()

void FCS_StepInfoSD::EndOfAthenaEvent ( ISF_FCS_Parametrization::FCS_StepInfoCollection * hitContnainer)

End of athena event processing.

Definition at line 146 of file FCS_StepInfoSD.cxx.

147 {
148 // Unpack map into vector
149 for (auto it : m_hit_map) {
150 for (auto* a_s : *it.second) {
151 // Giving away ownership of the objects!
152 hitContainer->push_back(a_s);
153 }
154 it.second->clear();
155 delete it.second;
156 } // Vector of IDs in the map
157 m_hit_map.clear();
158 if (m_config.verboseLevel > 4) {
159 G4cout << this->GetName()
160 << " DEBUG EndOfAthenaEvent: After initial cleanup, N="
161 << hitContainer->size() << G4endl;
162 }
163 return;
164}
std::map< Identifier, std::vector< ISF_FCS_Parametrization::FCS_StepInfo * > * > m_hit_map

◆ getCaloDDManager()

void FCS_StepInfoSD::getCaloDDManager ( )
protected

Keep a map instead of trying to keep the full vector.

At the end of the event we'll push the map back into the FCS_StepInfoCollection in StoreGate.

Definition at line 49 of file FCS_StepInfoSD.cxx.

49 {
50 SG::ReadCondHandleKey<CaloDetDescrManager> caloMgrKey{"CaloDetDescrManager"};
51 if (caloMgrKey.initialize().isFailure()) {
52 G4ExceptionDescription description;
53 description << "Failed to get CaloDetDescrManager!";
54 G4Exception("FCS_StepInfoSD", "FCSBadCall", FatalException, description);
55 abort();
56 }
57 SG::ReadCondHandle<CaloDetDescrManager> caloMgr(
58 caloMgrKey, Gaudi::Hive::currentContext());
59 m_calo_dd_man.set(*caloMgr);
60}
StatusCode initialize(bool used=true)
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:91

◆ getMaxTime()

double FCS_StepInfoSD::getMaxTime ( const CaloCell_ID::CaloSample & layer) const
inlineprivate

NB The result of this function should actually be constant for each SD

Definition at line 34 of file FCS_StepInfoSD.cxx.

35 {
37 if (layer >= CaloCell_ID::PreSamplerB && layer <= CaloCell_ID::EME3) {
38 return m_config.m_maxTimeLAr;
39 }
40 if (layer >= CaloCell_ID::HEC0 && layer <= CaloCell_ID::HEC3) {
41 return m_config.m_maxTimeHEC;
42 }
43 if (layer >= CaloCell_ID::FCAL0 && layer <= CaloCell_ID::FCAL2) {
44 return m_config.m_maxTimeFCAL;
45 }
46 return m_config.m_maxTime;
47}

◆ ProcessHits()

G4bool FCS_StepInfoSD::ProcessHits ( G4Step * a_step,
G4TouchableHistory *  )
overridevirtual

Main processing method.

Reimplemented in LArFCS_StepInfoSD, and TileFCS_StepInfoSD.

Definition at line 26 of file FCS_StepInfoSD.cxx.

26 {
27 G4ExceptionDescription description;
28 description << "ProcessHits: Base class method should not be called!!!";
29 G4Exception("FCS_StepInfoSD", "FCSBadCall", FatalException, description);
30 abort();
31 return false;
32}

◆ setupHelpers()

void FCS_StepInfoSD::setupHelpers ( const LArEM_ID * EM,
const LArFCAL_ID * FCAL,
const LArHEC_ID * HEC,
const TileID * tile )
inline

Sets the ID helper pointers.

Definition at line 84 of file FCS_StepInfoSD.h.

85 {
86 m_larEmID = EM;
89 m_tileID = tile;
90 }
const LArHEC_ID * m_larHecID
const LArFCAL_ID * m_larFcalID
const LArEM_ID * m_larEmID
Pointers to the identifier helpers.
const TileID * m_tileID

◆ update_map()

void FCS_StepInfoSD::update_map ( const CLHEP::Hep3Vector & l_vec,
const Identifier & l_identifier,
double l_energy,
double l_time,
bool l_valid,
int l_detector )
protected

Definition at line 62 of file FCS_StepInfoSD.cxx.

66{
67 // NB l_identifier refers to:
68 // - the cell identifier for LAr
69 // - the PMT identifier for Tile
70
71 // Drop any hits that don't have a good identifier attached
72 if (!m_calo_dd_man.get()->get_element(l_identifier)) {
73 if (m_config.verboseLevel > 4) {
74 G4cout << this->GetName() << " DEBUG update_map: bad identifier: "
75 << l_identifier.getString() << " skipping this hit." << G4endl;
76 }
77 return;
78 }
79
80 auto map_item = m_hit_map.find(l_identifier);
81 if (map_item == m_hit_map.end()) {
82 m_hit_map[l_identifier] =
83 new std::vector<ISF_FCS_Parametrization::FCS_StepInfo*>;
84 m_hit_map[l_identifier]->reserve(200);
85 m_hit_map[l_identifier]->push_back(
86 new ISF_FCS_Parametrization::FCS_StepInfo(l_vec, l_identifier, l_energy,
87 l_time, l_valid, l_detector));
88 } else {
89
90 // Get the appropriate merging limits
92 m_calo_dd_man.get()->get_element(l_identifier)->getSampling();
93
94 double timeWindow = m_config.m_maxTime;
95 const double distWinLong = m_config.m_maxRadiusLongitudinal.at(layer);
96 const double distWinLat = m_config.m_maxRadiusLateral.at(layer);
97
98 const double tsame(this->getMaxTime(layer));
99 bool match = false;
100 for (auto* map_it : *map_item->second) {
101 // Time check ... both a global flag and a check on the layer
102 const double delta_t = std::fabs(map_it->time() - l_time);
103 if (delta_t >= tsame) {
104 continue;
105 }
106 if (delta_t >= timeWindow) {
107 continue;
108 }
109
110 // Distance check
111 const CLHEP::Hep3Vector& currentPosition = map_it->position();
112 const double currentPosition_mag = currentPosition.mag();
113 const double proj_longitudinal =
114 currentPosition.dot(l_vec) / currentPosition_mag;
115 const double delta_longitudinal = currentPosition_mag - proj_longitudinal;
116 if (std::fabs(delta_longitudinal) >= distWinLong) {
117 continue;
118 }
119
120 // Lateral distance check
121 double delta_lateral_2 = l_vec.mag2() - proj_longitudinal * proj_longitudinal;
122 if (delta_lateral_2 < 0) {
123 delta_lateral_2 = 0; // Avoid negative square root
124 }
125 const double delta_lateral =
126 std::sqrt(delta_lateral_2);
127 if (delta_lateral >= distWinLat) {
128 continue;
129 }
130
131 // Found a match. Make a temporary that will be deleted!
132 const ISF_FCS_Parametrization::FCS_StepInfo my_info(
133 l_vec, l_identifier, l_energy, l_time, l_valid, l_detector);
134 *map_it += my_info;
135 match = true;
136 break;
137 } // End of search for match in time and space
138 if (!match) {
139 map_item->second->push_back(new ISF_FCS_Parametrization::FCS_StepInfo(
140 l_vec, l_identifier, l_energy, l_time, l_valid, l_detector));
141 } // Didn't match
142 } // ID already in the map
143 return;
144} // That's it for updating the map!
CaloSampling::CaloSample CaloSample
Definition CaloCell_ID.h:53
double getMaxTime(const CaloCell_ID::CaloSample &layer) const
std::string getString() const
Provide a string form of the identifier - hexadecimal.
bool match(std::string s1, std::string s2)
match the individual directories of two strings
Definition hcg.cxx:357
@ layer
Definition HitInfo.h:79

Member Data Documentation

◆ m_calo_dd_man

CxxUtils::CachedPointer<const CaloDetDescrManager> FCS_StepInfoSD::m_calo_dd_man
protected

Definition at line 106 of file FCS_StepInfoSD.h.

◆ m_config

FCS_Param::Config FCS_StepInfoSD::m_config
protected

Definition at line 100 of file FCS_StepInfoSD.h.

◆ m_hit_map

std::map<Identifier, std::vector<ISF_FCS_Parametrization::FCS_StepInfo*>*> FCS_StepInfoSD::m_hit_map
protected

Definition at line 108 of file FCS_StepInfoSD.h.

◆ m_larEmID

const LArEM_ID* FCS_StepInfoSD::m_larEmID {nullptr}
protected

Pointers to the identifier helpers.

Definition at line 102 of file FCS_StepInfoSD.h.

102{nullptr};

◆ m_larFcalID

const LArFCAL_ID* FCS_StepInfoSD::m_larFcalID {nullptr}
protected

Definition at line 103 of file FCS_StepInfoSD.h.

103{nullptr};

◆ m_larHecID

const LArHEC_ID* FCS_StepInfoSD::m_larHecID {nullptr}
protected

Definition at line 104 of file FCS_StepInfoSD.h.

104{nullptr};

◆ m_tileID

const TileID* FCS_StepInfoSD::m_tileID {nullptr}
protected

Definition at line 105 of file FCS_StepInfoSD.h.

105{nullptr};

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