CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes | Static Private Attributes
HCalSD Class Reference

#include <HCalSD.h>

Inheritance diagram for HCalSD:
CaloSD Observer< const BeginOfJob *> SensitiveCaloDetector Observer< const BeginOfRun *> Observer< const BeginOfEvent *> Observer< const BeginOfTrack *> Observer< const EndOfTrack *> Observer< const EndOfEvent *> SensitiveDetector

Public Member Functions

 HCalSD (const std::string &, const HcalDDDSimConstants *, const HcalDDDRecConstants *, const HcalSimulationConstants *, const HBHEDarkening *, const HBHEDarkening *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
uint32_t setDetUnitId (const G4Step *step) override
 
void setNumberingScheme (HcalNumberingScheme *)
 
 ~HCalSD () override=default
 
- Public Member Functions inherited from CaloSD
 CaloSD (const std::string &aSDname, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, float timeSlice=1., bool ignoreTkID=false, const std::string &newcolname="")
 
void clear () override
 
void clearHits () override
 
void DrawAll () override
 
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
void fillHits (edm::PCaloHitContainer &, const std::string &) override
 
void Initialize (G4HCofThisEvent *HCE) override
 
bool isItFineCalo (const G4VTouchable *touch)
 
void newCollection (const std::string &name, edm::ParameterSet const &p)
 
void PrintAll () override
 
G4bool ProcessHits (G4Step *step, G4TouchableHistory *) override
 
bool ProcessHits (G4GFlashSpot *aSpot, G4TouchableHistory *) override
 
void reset () override
 
 ~CaloSD () override
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (const std::string &iname, const SensitiveDetectorCatalog &clg, const std::string &newcollname="")
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const SensitiveDetectorCatalog &, bool calo, const std::string &newcollname="")
 
 ~SensitiveDetector () override
 
- Public Member Functions inherited from Observer< const BeginOfRun *>
 Observer ()
 
void slotForUpdate (const BeginOfRun * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfEvent *>
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfTrack *>
 Observer ()
 
void slotForUpdate (const BeginOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfTrack *>
 Observer ()
 
void slotForUpdate (const EndOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent *>
 Observer ()
 
void slotForUpdate (const EndOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfJob *>
 Observer ()
 
void slotForUpdate (const BeginOfJob * iT)
 
virtual ~Observer ()
 

Protected Member Functions

void endEvent () override
 
bool filterHit (CaloG4Hit *, double) override
 
double getEnergyDeposit (const G4Step *) override
 
bool getFromLibrary (const G4Step *) override
 
void initEvent (const BeginOfEvent *) override
 
void initRun () override
 
void update (const BeginOfRun *) override
 
void update (const ::EndOfEvent *) override
 
void update (const EndOfTrack *trk) override
 
void update (const BeginOfTrack *trk) override
 
void update (const BeginOfEvent *) override
 
void update (const BeginOfJob *) override
 This routine will be called when the appropriate signal arrives. More...
 
- Protected Member Functions inherited from CaloSD
bool checkHit (int k=0)
 
CaloG4HitcreateNewHit (const G4Step *, const G4Track *, int k)
 
virtual double EnergyCorrected (const G4Step &step, const G4Track *)
 
unsigned int findBoundaryCrossingParent (const G4Track *track, bool markParentAsSaveable=true)
 
double getAttenuation (const G4Step *aStep, double birk1, double birk2, double birk3) const
 
virtual uint16_t getDepth (const G4Step *)
 
int getNumberOfHits (int k=0)
 
double getResponseWt (const G4Track *, int k=0)
 
virtual int getTrackID (const G4Track *)
 
bool hitExists (const G4Step *, int k)
 
void ignoreRejection ()
 
void printDetectorLevels (const G4VTouchable *) const
 
void processHit (const G4Step *step)
 
virtual void processSecondHit (const G4Step *, const G4Track *)
 
void resetForNewPrimary (const G4Step *)
 
void setNumberCheckedHits (int val, int k=0)
 
void setParameterized (bool val)
 
G4ThreeVector setToGlobal (const G4ThreeVector &, const G4VTouchable *) const
 
G4ThreeVector setToLocal (const G4ThreeVector &, const G4VTouchable *) const
 
virtual int setTrackID (const G4Step *)
 
void setUseMap (bool val)
 
std::string shortreprID (const CaloHitID &ID)
 
std::string shortreprID (const CaloG4Hit *hit)
 
void update (const BeginOfRun *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfEvent *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const EndOfTrack *trk) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const ::EndOfEvent *) override
 
void updateHit (CaloG4Hit *, int k)
 
- Protected Member Functions inherited from SensitiveDetector
TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 
- Protected Member Functions inherited from Observer< const EndOfEvent *>
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 

Private Member Functions

void fillLogVolumeVector (const std::string &, const std::vector< std::string > &, std::vector< const G4LogicalVolume *> &)
 
void getFromHFLibrary (const G4Step *step, bool &isKilled)
 
void getFromParam (const G4Step *step, bool &isKilled)
 
void getHitFibreBundle (const G4Step *step, bool type)
 
void getHitPMT (const G4Step *step)
 
void hitForFibre (const G4Step *step)
 
bool isItConicalBundle (const G4LogicalVolume *)
 
bool isItFibre (const G4LogicalVolume *)
 
bool isItFibre (const G4String &)
 
bool isItHF (const G4Step *)
 
bool isItHF (const G4String &)
 
bool isItinFidVolume (const G4ThreeVector &)
 
bool isItPMT (const G4LogicalVolume *)
 
bool isItScintillator (const G4Material *)
 
bool isItStraightBundle (const G4LogicalVolume *)
 
double layerWeight (int, const G4ThreeVector &, int, int)
 
void modifyDepth (HcalNumberingFromDDD::HcalID &id)
 
void plotHF (const G4ThreeVector &pos, bool emType)
 
void plotProfile (const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
 
void printVolume (const G4VTouchable *touch) const
 
void readWeightFromFile (const std::string &)
 
uint32_t setDetUnitId (int, const G4ThreeVector &, int, int)
 
uint32_t setDetUnitId (HcalNumberingFromDDD::HcalID &tmp)
 

Private Attributes

bool agingFlagHB
 
bool agingFlagHE
 
bool applyFidCut
 
double betaThr
 
double birk1
 
double birk2
 
double birk3
 
bool dd4hep_
 
double deliveredLumi
 
int depth_
 
std::vector< int > detNull_
 
TH1F * dist_ [9]
 
double eminHitHB
 
double eminHitHE
 
double eminHitHF
 
double eminHitHO
 
std::vector< const G4LogicalVolume * > fibre1LV
 
std::vector< const G4LogicalVolume * > fibre2LV
 
std::vector< const G4LogicalVolume * > fibreLV
 
std::vector< std::string > fibreNames
 
std::vector< double > gpar
 
const HcalDDDSimConstantshcalConstants_
 
const HcalSimulationConstantshcalSimConstants_
 
std::vector< int > hfLevels
 
std::vector< const G4LogicalVolume * > hfLV
 
std::vector< std::string > hfNames
 
std::unique_ptr< HFShowerhfshower
 
TH1F * hit_ [9]
 
TH1F * hzvem
 
TH1F * hzvhad
 
bool isHF
 
std::map< uint32_t, double > layerWeights
 
const HBHEDarkeningm_HBDarkening
 
std::unique_ptr< HcalTestNSm_HcalTestNS
 
const HBHEDarkeningm_HEDarkening
 
std::unique_ptr< HFDarkeningm_HFDarkening
 
std::vector< const G4Material * > materials
 
std::vector< std::string > matNames
 
bool neutralDensity
 
std::unique_ptr< HcalNumberingFromDDDnumberingFromDDD
 
std::unique_ptr< HcalNumberingSchemenumberingScheme
 
std::vector< const G4LogicalVolume * > pmtLV
 
std::unique_ptr< HFShowerFibreBundleshowerBundle
 
std::unique_ptr< HFShowerLibraryshowerLibrary
 
std::unique_ptr< HFShowerParamshowerParam
 
std::unique_ptr< HFShowerPMTshowerPMT
 
bool testNS_
 
bool testNumber
 
TH1F * time_ [9]
 
bool useBirk
 
bool useFibreBundle
 
bool useHF
 
bool useLayerWt
 
bool useParam
 
bool usePMTHit
 
bool useShowerLibrary
 
double weight_
 

Static Private Attributes

static constexpr double maxRoff_ = 450.0
 
static constexpr double maxZ_ = 10000.0
 
static constexpr double minRoff_ = -1500.0
 
static constexpr double slopeHE_ = 0.4
 

Additional Inherited Members

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 
- Static Protected Member Functions inherited from CaloSD
static std::string printableDecayChain (const std::vector< unsigned int > &decayChain)
 
- Protected Attributes inherited from CaloSD
std::string collName_ [2]
 
CaloG4HitcurrentHit [2]
 
CaloHitID currentID [2]
 
std::string detName_
 
float edepositEM
 
float edepositHAD
 
double eminHit
 
double energyCut
 
G4ThreeVector entranceLocal
 
G4ThreeVector entrancePoint
 
bool forceSave
 
std::vector< std::string > hcn_
 
float incidentEnergy
 
double kmaxIon
 
double kmaxNeutron
 
double kmaxProton
 
int nHC_
 
G4ThreeVector posGlobal
 
CaloHitID previousID [2]
 
bool suppressHeavy
 
double tmaxHit
 
std::vector< int > useResMap_
 

Detailed Description

Definition at line 38 of file HCalSD.h.

Constructor & Destructor Documentation

◆ HCalSD()

HCalSD::HCalSD ( const std::string &  name,
const HcalDDDSimConstants hcns,
const HcalDDDRecConstants hcnr,
const HcalSimulationConstants hscs,
const HBHEDarkening hbd,
const HBHEDarkening hed,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 42 of file HCalSD.cc.

References agingFlagHB, agingFlagHE, applyFidCut, betaThr, birk1, birk2, birk3, dd4hep_, deliveredLumi, dist_, eminHitHB, eminHitHE, eminHitHF, eminHitHO, fibre1LV, fibre2LV, fibreLV, fibreNames, geometryDiff::file, fillLogVolumeVector(), g, relativeConstraints::geom, HcalDDDSimConstants::getGparHF(), HcalDDDSimConstants::getMaxDepth(), gpar, hcalConstants_, HcalSimulationParameters::hcalMaterialNames_, hcalSimConstants_, HcalSimulationConstants::hcalsimpar(), HcalSimulationParameters::hfFibreConicalNames_, HcalSimulationParameters::hfFibreNames_, HcalSimulationParameters::hfFibreStraightNames_, hfLevels, HcalSimulationParameters::hfLevels_, hfLV, hfNames, HcalSimulationParameters::hfNames_, HcalSimulationParameters::hfPMTNames_, hfshower, hit_, hzvem, hzvhad, mps_fire::i, CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, SummaryClient_cfi::labels, SensitiveDetectorCatalog::logicalNames(), materialBudgetVolume_cfi::lvNames, m_HcalTestNS, m_HFDarkening, TFileDirectory::make(), materials, matNames, Skims_PA_cff::name, DD4hep2DDDName::nameMatterLV(), neutralDensity, numberingFromDDD, numberingScheme, AlCaHLTBitMon_ParallelJobs::p, pmtLV, readWeightFromFile(), generator_cfi::scheme, setNumberingScheme(), CaloSD::setParameterized(), showerBundle, showerLibrary, showerParam, showerPMT, AlCaHLTBitMon_QueryRunRegistry::string, CaloSD::suppressHeavy, groupFilesInBlocks::temp, testNS_, testNumber, compare::tfile, time_, runGCPTkAlMap::title, useBirk, useFibreBundle, useHF, useLayerWt, useParam, usePMTHit, useShowerLibrary, and relativeConstraints::value.

51  : CaloSD(name,
52  clg,
53  p,
54  manager,
55  (float)(p.getParameter<edm::ParameterSet>("HCalSD").getParameter<double>("TimeSliceUnit")),
56  p.getParameter<edm::ParameterSet>("HCalSD").getParameter<bool>("IgnoreTrackID")),
57  hcalConstants_(hcns),
58  hcalSimConstants_(hscs),
59  m_HBDarkening(hbd),
60  m_HEDarkening(hed),
61  isHF(false),
62  weight_(1.0),
63  depth_(1) {
64  numberingFromDDD.reset(nullptr);
65  numberingScheme.reset(nullptr);
66  showerLibrary.reset(nullptr);
67  hfshower.reset(nullptr);
68  showerParam.reset(nullptr);
69  showerPMT.reset(nullptr);
70  showerBundle.reset(nullptr);
71  m_HFDarkening.reset(nullptr);
72  m_HcalTestNS.reset(nullptr);
73 
74  dd4hep_ = p.getParameter<bool>("g4GeometryDD4hepSource");
75  edm::ParameterSet m_HC = p.getParameter<edm::ParameterSet>("HCalSD");
76  useBirk = m_HC.getParameter<bool>("UseBirkLaw");
77  double bunit = (CLHEP::g / (CLHEP::MeV * CLHEP::cm2));
78  birk1 = m_HC.getParameter<double>("BirkC1") * bunit;
79  birk2 = m_HC.getParameter<double>("BirkC2");
80  birk3 = m_HC.getParameter<double>("BirkC3");
81  useShowerLibrary = m_HC.getParameter<bool>("UseShowerLibrary");
82  useParam = m_HC.getParameter<bool>("UseParametrize");
83  testNumber = m_HC.getParameter<bool>("TestNumberingScheme");
84  neutralDensity = m_HC.getParameter<bool>("doNeutralDensityFilter");
85  usePMTHit = m_HC.getParameter<bool>("UsePMTHits");
86  betaThr = m_HC.getParameter<double>("BetaThreshold");
87  eminHitHB = m_HC.getParameter<double>("EminHitHB") * MeV;
88  eminHitHE = m_HC.getParameter<double>("EminHitHE") * MeV;
89  eminHitHO = m_HC.getParameter<double>("EminHitHO") * MeV;
90  eminHitHF = m_HC.getParameter<double>("EminHitHF") * MeV;
91  useFibreBundle = m_HC.getParameter<bool>("UseFibreBundleHits");
92  deliveredLumi = m_HC.getParameter<double>("DelivLuminosity");
93  agingFlagHB = m_HC.getParameter<bool>("HBDarkening");
94  agingFlagHE = m_HC.getParameter<bool>("HEDarkening");
95  bool agingFlagHF = m_HC.getParameter<bool>("HFDarkening");
96  useHF = m_HC.getUntrackedParameter<bool>("UseHF", true);
97  bool forTBHC = m_HC.getUntrackedParameter<bool>("ForTBHCAL", false);
98  bool forTBH2 = m_HC.getUntrackedParameter<bool>("ForTBH2", false);
99  useLayerWt = m_HC.getUntrackedParameter<bool>("UseLayerWt", false);
100  std::string file = m_HC.getUntrackedParameter<std::string>("WtFile", "None");
101  testNS_ = m_HC.getUntrackedParameter<bool>("TestNS", false);
102  edm::ParameterSet m_HF = p.getParameter<edm::ParameterSet>("HFShower");
103  applyFidCut = m_HF.getParameter<bool>("ApplyFiducialCut");
104  bool dumpGeom = m_HC.getUntrackedParameter<bool>("DumpGeometry", false);
105 
106 #ifdef EDM_ML_DEBUG
107  edm::LogVerbatim("HcalSim") << "***************************************************"
108  << "\n"
109  << "* Constructing a HCalSD with name " << name << "\n"
110  << "\n"
111  << "***************************************************";
112 #endif
113  edm::LogVerbatim("HcalSim") << "HCalSD:: Use of HF code is set to " << useHF
114  << "\n Use of shower parametrization set to " << useParam
115  << "\n Use of shower library is set to " << useShowerLibrary
116  << "\n Use PMT Hit is set to " << usePMTHit << " with beta Threshold " << betaThr
117  << "\n USe of FibreBundle Hit set to " << useFibreBundle
118  << "\n Use of Birks law is set to " << useBirk
119  << " with three constants kB = " << birk1 / bunit << ", C1 = " << birk2
120  << ", C2 = " << birk3;
121  edm::LogVerbatim("HcalSim") << "HCalSD:: Suppression Flag " << suppressHeavy << " protons below " << kmaxProton
122  << " MeV,"
123  << " neutrons below " << kmaxNeutron << " MeV and"
124  << " ions below " << kmaxIon << " MeV\n"
125  << " Threshold for storing hits in HB: " << eminHitHB << " HE: " << eminHitHE
126  << " HO: " << eminHitHO << " HF: " << eminHitHF << "\n"
127  << "Delivered luminosity for Darkening " << deliveredLumi << " Flag (HE) " << agingFlagHE
128  << " Flag (HB) " << agingFlagHB << " Flag (HF) " << agingFlagHF << "\n"
129  << "Application of Fiducial Cut " << applyFidCut
130  << "Flag for test number|neutral density filter " << testNumber << " " << neutralDensity;
131 
132  if (forTBHC) {
133  useHF = false;
134  matNames.emplace_back("Scintillator");
135  } else {
137  }
138 
140  if (testNumber || forTBH2) {
141  scheme = dynamic_cast<HcalNumberingScheme*>(new HcalTestNumberingScheme(forTBH2));
142  } else {
143  scheme = new HcalNumberingScheme();
144  }
146 
147  // always call getFromLibrary() method to identify HF region
148  setParameterized(true);
149 
150  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
151  // std::vector<G4LogicalVolume *>::const_iterator lvcite;
152  const G4LogicalVolume* lv;
153  std::string attribute, value;
154 
155  if (useHF) {
156  if (useParam) {
157  showerParam = std::make_unique<HFShowerParam>(name, hcalConstants_, hcalSimConstants_->hcalsimpar(), p);
158  } else {
159  if (useShowerLibrary) {
160  showerLibrary = std::make_unique<HFShowerLibrary>(name, hcalConstants_, hcalSimConstants_->hcalsimpar(), p);
161  }
162  hfshower = std::make_unique<HFShower>(name, hcalConstants_, hcalSimConstants_->hcalsimpar(), p, 0);
163  }
164 
165  // HF volume names
167  const std::vector<int>& temp = hcalSimConstants_->hcalsimpar()->hfLevels_;
168 #ifdef EDM_ML_DEBUG
169  std::stringstream ss0;
170  ss0 << "HCalSD: Names to be tested for Volume = HF has " << hfNames.size() << " elements";
171 #endif
172  int addlevel = dd4hep_ ? 1 : 0;
173  for (unsigned int i = 0; i < hfNames.size(); ++i) {
175  lv = nullptr;
176  for (auto lvol : *lvs) {
177  if (DD4hep2DDDName::nameMatterLV(lvol->GetName(), dd4hep_) == namv) {
178  lv = lvol;
179  break;
180  }
181  }
182  hfLV.emplace_back(lv);
183  hfLevels.emplace_back(temp[i] + addlevel);
184 #ifdef EDM_ML_DEBUG
185  ss0 << "\n HF[" << i << "] = " << namv << " LV " << lv << " at level " << (temp[i] + addlevel);
186 #endif
187  }
188 #ifdef EDM_ML_DEBUG
189  edm::LogVerbatim("HcalSim") << ss0.str();
190 #endif
191  // HF Fibre volume names
194  const std::vector<std::string>& pmtNames = hcalSimConstants_->hcalsimpar()->hfPMTNames_;
195  fillLogVolumeVector("HFPMT", pmtNames, pmtLV);
196  const std::vector<std::string>& straightNames = hcalSimConstants_->hcalsimpar()->hfFibreStraightNames_;
197  fillLogVolumeVector("HFFibreBundleStraight", straightNames, fibre1LV);
198  const std::vector<std::string>& conicalNames = hcalSimConstants_->hcalsimpar()->hfFibreConicalNames_;
199  fillLogVolumeVector("HFFibreBundleConical", conicalNames, fibre2LV);
200  }
201 
202  //Material list for HB/HE/HO sensitive detectors
203  const G4MaterialTable* matTab = G4Material::GetMaterialTable();
204  std::vector<G4Material*>::const_iterator matite;
205  for (auto const& namx : matNames) {
206  const G4Material* mat = nullptr;
207  for (matite = matTab->begin(); matite != matTab->end(); ++matite) {
208  if (DD4hep2DDDName::nameMatterLV((*matite)->GetName(), dd4hep_) == namx) {
209  mat = (*matite);
210  break;
211  }
212  }
213  materials.emplace_back(mat);
214  }
215 #ifdef EDM_ML_DEBUG
216  std::stringstream ss1;
217  for (unsigned int i = 0; i < matNames.size(); ++i) {
218  if (i / 10 * 10 == i) {
219  ss1 << "\n";
220  }
221  ss1 << " " << matNames[i];
222  }
223  edm::LogVerbatim("HcalSim") << "HCalSD: Material names for HCAL: " << ss1.str();
224 #endif
225  if (useLayerWt) {
227  }
228  numberingFromDDD = std::make_unique<HcalNumberingFromDDD>(hcalConstants_);
229 
230  //Special Geometry parameters
232 #ifdef EDM_ML_DEBUG
233  std::stringstream ss2;
234  for (unsigned int ig = 0; ig < gpar.size(); ig++) {
235  ss2 << "\n gpar[" << ig << "] = " << gpar[ig] / cm << " cm";
236  }
237  edm::LogVerbatim("HcalSim") << "Maximum depth for HF " << hcalConstants_->getMaxDepth(2) << gpar.size()
238  << " gpar (cm)" << ss2.str();
239 #endif
240 
241  //Test Hcal Numbering Scheme
242  if (testNS_)
243  m_HcalTestNS = std::make_unique<HcalTestNS>(hcnr);
244 
245  for (int i = 0; i < 9; ++i) {
246  hit_[i] = time_[i] = dist_[i] = nullptr;
247  }
248  hzvem = hzvhad = nullptr;
249 
250  if (agingFlagHF) {
251  m_HFDarkening = std::make_unique<HFDarkening>(m_HC.getParameter<edm::ParameterSet>("HFDarkeningParameterBlock"));
252  }
253 #ifdef plotDebug
255 
256  if (tfile.isAvailable()) {
257  static const char* const labels[] = {"HB",
258  "HE",
259  "HO",
260  "HF Absorber",
261  "HF PMT",
262  "HF Absorber Long",
263  "HF Absorber Short",
264  "HF PMT Long",
265  "HF PMT Short"};
266  TFileDirectory hcDir = tfile->mkdir("ProfileFromHCalSD");
267  char name[20], title[60];
268  for (int i = 0; i < 9; ++i) {
269  sprintf(title, "Hit energy in %s", labels[i]);
270  sprintf(name, "HCalSDHit%d", i);
271  hit_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
272  sprintf(title, "Energy (MeV)");
273  hit_[i]->GetXaxis()->SetTitle(title);
274  hit_[i]->GetYaxis()->SetTitle("Hits");
275  sprintf(title, "Time of the hit in %s", labels[i]);
276  sprintf(name, "HCalSDTime%d", i);
277  time_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
278  sprintf(title, "Time (ns)");
279  time_[i]->GetXaxis()->SetTitle(title);
280  time_[i]->GetYaxis()->SetTitle("Hits");
281  sprintf(title, "Longitudinal profile in %s", labels[i]);
282  sprintf(name, "HCalSDDist%d", i);
283  dist_[i] = hcDir.make<TH1F>(name, title, 2000, 0., 2000.);
284  sprintf(title, "Distance (mm)");
285  dist_[i]->GetXaxis()->SetTitle(title);
286  dist_[i]->GetYaxis()->SetTitle("Hits");
287  }
288  if (useHF && (!useParam)) {
289  hzvem = hcDir.make<TH1F>("hzvem", "Longitudinal Profile (EM Part)", 330, 0.0, 1650.0);
290  hzvem->GetXaxis()->SetTitle("Longitudinal Profile (EM Part)");
291  hzvhad = hcDir.make<TH1F>("hzvhad", "Longitudinal Profile (Had Part)", 330, 0.0, 1650.0);
292  hzvhad->GetXaxis()->SetTitle("Longitudinal Profile (Hadronic Part)");
293  }
294  }
295 #endif
296  if (dumpGeom) {
297  std::unique_ptr<HcalNumberingFromDDD> hcn = std::make_unique<HcalNumberingFromDDD>(hcalConstants_);
298  const auto& lvNames = clg.logicalNames(name);
299  HcalDumpGeometry geom(lvNames, hcn.get(), testNumber, false);
300  geom.update();
301  }
302 }
Log< level::Info, true > LogVerbatim
std::vector< std::string > matNames
Definition: HCalSD.h:122
void readWeightFromFile(const std::string &)
Definition: HCalSD.cc:972
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:124
bool useParam
Definition: HCalSD.h:112
double eminHitHE
Definition: HCalSD.h:113
TH1F * time_[9]
Definition: HCalSD.h:126
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
double kmaxNeutron
Definition: CaloSD.h:155
bool useLayerWt
Definition: HCalSD.h:109
std::vector< std::string > fibreNames
Definition: HCalSD.h:121
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:92
double weight_
Definition: HCalSD.h:115
std::vector< double > gpar
Definition: HCalSD.h:118
bool useFibreBundle
Definition: HCalSD.h:109
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:96
double betaThr
Definition: HCalSD.h:111
std::vector< std::string > hfNames_
double deliveredLumi
Definition: HCalSD.h:114
std::vector< std::string > hfFibreNames_
double eminHitHB
Definition: HCalSD.h:113
bool useShowerLibrary
Definition: HCalSD.h:112
std::vector< std::string > hfNames
Definition: HCalSD.h:120
double birk2
Definition: HCalSD.h:111
bool usePMTHit
Definition: HCalSD.h:109
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:124
void setNumberingScheme(HcalNumberingScheme *)
Definition: HCalSD.cc:563
bool agingFlagHE
Definition: HCalSD.h:108
double birk1
Definition: HCalSD.h:111
double kmaxProton
Definition: CaloSD.h:155
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:93
const std::vector< double > & getGparHF() const
The Signals That Services Can Subscribe To This is based on ActivityRegistry and is current per Services can connect to the signals distributed by the ActivityRegistry in order to monitor the activity of the application Each possible callback has some defined which we here list in angle e g
Definition: Activities.doc:4
TH1F * hzvhad
Definition: HCalSD.h:126
double eminHitHF
Definition: HCalSD.h:113
const HcalSimulationParameters * hcalsimpar() const
T * make(const Args &...args) const
make new ROOT object
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
TH1F * hit_[9]
Definition: HCalSD.h:126
bool agingFlagHB
Definition: HCalSD.h:108
bool testNumber
Definition: HCalSD.h:110
double kmaxIon
Definition: CaloSD.h:155
bool suppressHeavy
Definition: CaloSD.h:154
bool dd4hep_
Definition: HCalSD.h:117
std::vector< std::string > hfFibreConicalNames_
bool useBirk
Definition: HCalSD.h:109
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:100
double birk3
Definition: HCalSD.h:111
Definition: tfile.py:1
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:99
bool applyFidCut
Definition: HCalSD.h:112
TH1F * dist_[9]
Definition: HCalSD.h:126
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:101
std::vector< int > hfLevels
Definition: HCalSD.h:119
bool testNS_
Definition: HCalSD.h:110
double eminHitHO
Definition: HCalSD.h:113
TH1F * hzvem
Definition: HCalSD.h:126
CaloSD(const std::string &aSDname, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *, float timeSlice=1., bool ignoreTkID=false, const std::string &newcolname="")
Definition: CaloSD.cc:33
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:124
bool neutralDensity
Definition: HCalSD.h:110
bool isHF
Definition: HCalSD.h:107
void fillLogVolumeVector(const std::string &, const std::vector< std::string > &, std::vector< const G4LogicalVolume *> &)
Definition: HCalSD.cc:304
std::vector< std::string > hcalMaterialNames_
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:94
std::vector< std::string > hfPMTNames_
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:90
const HcalSimulationConstants * hcalSimConstants_
Definition: HCalSD.h:97
int depth_
Definition: HCalSD.h:116
std::vector< std::string > hfFibreStraightNames_
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:89
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:91
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:124
std::string nameMatterLV(const std::string &name, bool dd4hep)
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:124
std::vector< const G4Material * > materials
Definition: HCalSD.h:123
int getMaxDepth(const int &type) const
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:98
bool useHF
Definition: HCalSD.h:112
void setParameterized(bool val)
Definition: CaloSD.h:114

◆ ~HCalSD()

HCalSD::~HCalSD ( )
overridedefault

Member Function Documentation

◆ endEvent()

void HCalSD::endEvent ( )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 1109 of file HCalSD.cc.

References detNull_.

1109  {
1110 #ifdef printDebug
1111  int sum = detNull_[0] + detNull_[1] + detNull_[2];
1112  if (sum > 0)
1113  edm::LogVerbatim("HcalSim") << "NullDets " << detNull_[0] << " " << detNull_[1] << " " << detNull_[2] << " "
1114  << detNull_[3] << " " << (static_cast<float>(sum) / (sum + detNull_[3]));
1115 #endif
1116 }
Log< level::Info, true > LogVerbatim
std::vector< int > detNull_
Definition: HCalSD.h:127

◆ fillLogVolumeVector()

void HCalSD::fillLogVolumeVector ( const std::string &  value,
const std::vector< std::string > &  lvnames,
std::vector< const G4LogicalVolume *> &  lvvec 
)
private

Definition at line 304 of file HCalSD.cc.

References dd4hep_, mps_fire::i, DD4hep2DDDName::nameMatterLV(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by HCalSD().

306  {
307  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
308  const G4LogicalVolume* lv;
309  std::stringstream ss3;
310  ss3 << "HCalSD: " << lvnames.size() << " names to be tested for Volume <" << value << ">:";
311  for (unsigned int i = 0; i < lvnames.size(); ++i) {
313  lv = nullptr;
314  for (auto lvol : *lvs) {
315  if (DD4hep2DDDName::nameMatterLV(lvol->GetName(), dd4hep_) == namv) {
316  lv = lvol;
317  break;
318  }
319  }
320  lvvec.emplace_back(lv);
321  if (i / 10 * 10 == i) {
322  ss3 << "\n";
323  }
324  ss3 << " " << namv;
325  }
326  edm::LogVerbatim("HcalSim") << ss3.str();
327 }
Log< level::Info, true > LogVerbatim
bool dd4hep_
Definition: HCalSD.h:117
Definition: value.py:1
std::string nameMatterLV(const std::string &name, bool dd4hep)

◆ filterHit()

bool HCalSD::filterHit ( CaloG4Hit aHit,
double  time 
)
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 574 of file HCalSD.cc.

References eminHitHB, eminHitHE, eminHitHF, eminHitHO, CaloG4Hit::getEnergyDeposit(), CaloG4Hit::getUnitID(), HcalBarrel, HcalEndcap, HcalForward, HcalOuter, DiMuonV_cfg::threshold, hcalRecHitTable_cff::time, and CaloSD::tmaxHit.

574  {
575  double threshold = 0;
576  DetId theId(aHit->getUnitID());
577  switch (theId.subdetId()) {
578  case HcalBarrel:
580  break;
581  case HcalEndcap:
583  break;
584  case HcalOuter:
586  break;
587  case HcalForward:
589  break;
590  default:
591  break;
592  }
593  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > threshold));
594 }
double eminHitHE
Definition: HCalSD.h:113
double eminHitHB
Definition: HCalSD.h:113
double eminHitHF
Definition: HCalSD.h:113
double eminHitHO
Definition: HCalSD.h:113
double tmaxHit
Definition: CaloSD.h:148
Definition: DetId.h:17
double getEnergyDeposit() const
Definition: CaloG4Hit.h:79
uint32_t getUnitID() const
Definition: CaloG4Hit.h:66

◆ getEnergyDeposit()

double HCalSD::getEnergyDeposit ( const G4Step *  aStep)
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 396 of file HCalSD.cc.

References agingFlagHB, agingFlagHE, birk1, birk2, birk3, HBHEDarkening::degradation(), deliveredLumi, hcalRecHitTable_cff::depth, depth_, CaloSD::getAttenuation(), getHitFibreBundle(), getHitPMT(), HcalDDDSimConstants::getLayer0Wt(), CaloSD::getResponseWt(), hcalConstants_, hitForFibre(), hcalRecHitTable_cff::ieta, HcalDetId::ietaAbs(), HcalDetId::iphi(), isHF, isItConicalBundle(), isItFibre(), isItPMT(), isItScintillator(), isItStraightBundle(), G4TrackToParticleID::isMuon(), TrackInformation::isPrimary(), CaloSD::kmaxIon, CaloSD::kmaxNeutron, CaloSD::kmaxProton, layerWeight(), m_HBDarkening, m_HEDarkening, neutralDensity, phi, setDetUnitId(), showerBundle, showerPMT, DetId::subdetId(), CaloSD::suppressHeavy, testNumber, HcalTestNumbering::unpackHcalIndex(), useBirk, useFibreBundle, useLayerWt, usePMTHit, useShowerLibrary, weight_, z, and HcalDetId::zside().

396  {
397  double destep(0.0);
398  auto const lv = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume();
399  auto const theTrack = aStep->GetTrack();
400 
401  if (isHF) {
402  if (useShowerLibrary && G4TrackToParticleID::isMuon(theTrack) && isItFibre(lv)) {
403 #ifdef EDM_ML_DEBUG
404  edm::LogVerbatim("HcalSim") << "HCalSD: Hit at Fibre in LV " << lv->GetName() << " for track "
405  << aStep->GetTrack()->GetTrackID() << " ("
406  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
407 #endif
408  hitForFibre(aStep);
409  }
410  return destep;
411  }
412 
413  if (isItPMT(lv)) {
414  if (usePMTHit && showerPMT) {
415  getHitPMT(aStep);
416  }
417 #ifdef EDM_ML_DEBUG
418  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from PMT parametrization in LV " << lv->GetName() << " for Track "
419  << aStep->GetTrack()->GetTrackID() << " ("
420  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
421 #endif
422  return destep;
423 
424  } else if (isItStraightBundle(lv)) {
425  if (useFibreBundle && showerBundle) {
426  getHitFibreBundle(aStep, false);
427  }
428 #ifdef EDM_ML_DEBUG
429  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from straight FibreBundle in LV: " << lv->GetName() << " for track "
430  << aStep->GetTrack()->GetTrackID() << " ("
431  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
432 #endif
433  return destep;
434 
435  } else if (isItConicalBundle(lv)) {
436  if (useFibreBundle && showerBundle) {
437  getHitFibreBundle(aStep, true);
438  }
439 #ifdef EDM_ML_DEBUG
440  edm::LogVerbatim("HcalSim") << "HCalSD: Hit from conical FibreBundle PV: " << lv->GetName() << " for track "
441  << aStep->GetTrack()->GetTrackID() << " ("
442  << aStep->GetTrack()->GetDefinition()->GetParticleName() << ")";
443 #endif
444  return destep;
445  }
446 
447  // normal hit
448  destep = aStep->GetTotalEnergyDeposit();
449 
450  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
451  uint32_t detid = setDetUnitId(aStep);
452  int det(0), ieta(0), phi(0), z(0), lay, depth(-1);
453  if (testNumber) {
454  HcalTestNumbering::unpackHcalIndex(detid, det, z, depth, ieta, phi, lay);
455  if (z == 0) {
456  z = -1;
457  }
458  } else {
459  HcalDetId hcid(detid);
460  det = hcid.subdetId();
461  ieta = hcid.ietaAbs();
462  phi = hcid.iphi();
463  z = hcid.zside();
464  }
465  lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
466 #ifdef EDM_ML_DEBUG
467  edm::LogVerbatim("HcalSim") << "HCalSD: det: " << det << " ieta: " << ieta << " iphi: " << phi << " zside " << z
468  << " lay: " << lay - 2;
469 #endif
470  if (depth_ == 0 && (det == 1 || det == 2) && ((!testNumber) || neutralDensity))
472  if (useLayerWt) {
473  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
474  weight_ = layerWeight(det + 2, hitPoint, depth_, lay);
475  }
476 
477  if (agingFlagHB && m_HBDarkening && det == 1) {
478  double dweight = m_HBDarkening->degradation(deliveredLumi, ieta, lay);
479  weight_ *= dweight;
480 #ifdef EDM_ML_DEBUG
481  edm::LogVerbatim("HcalSim") << "HCalSD: HB Lumi: " << deliveredLumi << " coefficient = " << dweight
482  << " Weight= " << weight_;
483 #endif
484  }
485 
486  if (agingFlagHE && m_HEDarkening && det == 2) {
487  double dweight = m_HEDarkening->degradation(deliveredLumi, ieta, lay);
488  weight_ *= dweight;
489 #ifdef EDM_ML_DEBUG
490  edm::LogVerbatim("HcalSim") << "HCalSD: HB Lumi: " << deliveredLumi << " coefficient = " << dweight
491  << " Weight= " << weight_;
492 #endif
493  }
494 
495  if (suppressHeavy) {
496  TrackInformation* trkInfo = (TrackInformation*)(theTrack->GetUserInformation());
497  if (trkInfo) {
498  int pdg = theTrack->GetDefinition()->GetPDGEncoding();
499  if (!(trkInfo->isPrimary())) { // Only secondary particles
500  double ke = theTrack->GetKineticEnergy();
501  if (pdg / 1000000000 == 1 && (pdg / 10000) % 100 > 0 && (pdg / 10) % 100 > 0 && ke < kmaxIon)
502  weight_ = 0;
503  if ((pdg == 2212) && (ke < kmaxProton))
504  weight_ = 0;
505  if ((pdg == 2112) && (ke < kmaxNeutron))
506  weight_ = 0;
507  }
508  }
509  }
510  double wt0(1.0);
511  if (useBirk) {
512  const G4Material* mat = aStep->GetPreStepPoint()->GetMaterial();
513  if (isItScintillator(mat))
514  wt0 = getAttenuation(aStep, birk1, birk2, birk3);
515  }
516  weight_ *= wt0;
517  double wt1 = getResponseWt(theTrack);
518  double wt2 = theTrack->GetWeight();
519  double edep = weight_ * wt1 * destep;
520  if (wt2 > 0.0) {
521  edep *= wt2;
522  }
523 #ifdef EDM_ML_DEBUG
524  edm::LogVerbatim("HcalSim") << "HCalSD: edep= " << edep << " Det: " << det + 2 << " depth= " << depth_
525  << " weight= " << weight_ << " wt0= " << wt0 << " wt1= " << wt1 << " wt2= " << wt2;
526 #endif
527  return edep;
528 }
Log< level::Info, true > LogVerbatim
bool isPrimary() const
double kmaxNeutron
Definition: CaloSD.h:155
void hitForFibre(const G4Step *step)
Definition: HCalSD.cc:779
static bool isMuon(int pdgCode)
bool useLayerWt
Definition: HCalSD.h:109
double weight_
Definition: HCalSD.h:115
bool useFibreBundle
Definition: HCalSD.h:109
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:96
double deliveredLumi
Definition: HCalSD.h:114
float degradation(float intlumi, int ieta, int lay) const
bool useShowerLibrary
Definition: HCalSD.h:112
double birk2
Definition: HCalSD.h:111
bool usePMTHit
Definition: HCalSD.h:109
double getLayer0Wt(const int &det, const int &phi, const int &zside) const
bool agingFlagHE
Definition: HCalSD.h:108
double birk1
Definition: HCalSD.h:111
double kmaxProton
Definition: CaloSD.h:155
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:93
bool isItConicalBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:704
double getResponseWt(const G4Track *, int k=0)
Definition: CaloSD.cc:928
bool agingFlagHB
Definition: HCalSD.h:108
bool testNumber
Definition: HCalSD.h:110
double kmaxIon
Definition: CaloSD.h:155
bool suppressHeavy
Definition: CaloSD.h:154
bool useBirk
Definition: HCalSD.h:109
double birk3
Definition: HCalSD.h:111
const HBHEDarkening * m_HEDarkening
Definition: HCalSD.h:99
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
bool isItStraightBundle(const G4LogicalVolume *)
Definition: HCalSD.cc:696
bool neutralDensity
Definition: HCalSD.h:110
bool isHF
Definition: HCalSD.h:107
void getHitFibreBundle(const G4Step *step, bool type)
Definition: HCalSD.cc:910
double getAttenuation(const G4Step *aStep, double birk1, double birk2, double birk3) const
Definition: CaloSD.cc:727
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:94
bool isItPMT(const G4LogicalVolume *)
Definition: HCalSD.cc:688
void getHitPMT(const G4Step *step)
Definition: HCalSD.cc:851
double layerWeight(int, const G4ThreeVector &, int, int)
Definition: HCalSD.cc:996
int depth_
Definition: HCalSD.h:116
bool isItScintillator(const G4Material *)
Definition: HCalSD.cc:712
const HBHEDarkening * m_HBDarkening
Definition: HCalSD.h:98
bool isItFibre(const G4LogicalVolume *)
Definition: HCalSD.cc:672

◆ getFromHFLibrary()

void HCalSD::getFromHFLibrary ( const G4Step *  step,
bool &  isKilled 
)
private

Definition at line 736 of file HCalSD.cc.

References CaloSD::currentID, hcalRecHitTable_cff::depth, CaloSD::edepositEM, CaloSD::edepositHAD, hfClusterShapes_cfi::hits, mps_fire::i, G4TrackToParticleID::isGammaElectronPositron(), isItinFidVolume(), plotHF(), plotProfile(), CaloSD::processHit(), CaloSD::resetForNewPrimary(), setDetUnitId(), CaloHitID::setID(), CaloSD::setTrackID(), showerLibrary, hcalRecHitTable_cff::time, and weight_.

Referenced by getFromLibrary().

736  {
737  std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, isKilled, weight_, false);
738  if (!isKilled || hits.empty()) {
739  return;
740  }
741 
742  int primaryID = setTrackID(aStep);
743 
744  // Reset entry point for new primary
745  resetForNewPrimary(aStep);
746 
747  auto const theTrack = aStep->GetTrack();
748  int det = 5;
749 
751  edepositEM = 1. * GeV;
752  edepositHAD = 0.;
753  } else {
754  edepositEM = 0.;
755  edepositHAD = 1. * GeV;
756  }
757 #ifdef EDM_ML_DEBUG
758  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary " << hits.size() << " hits for " << GetName() << " of "
759  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
760  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV";
761 #endif
762  for (unsigned int i = 0; i < hits.size(); ++i) {
763  G4ThreeVector hitPoint = hits[i].position;
764  if (isItinFidVolume(hitPoint)) {
765  int depth = hits[i].depth;
766  double time = hits[i].time;
767  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
768  currentID[0].setID(unitID, time, primaryID, 0);
769 #ifdef plotDebug
770  plotProfile(aStep, hitPoint, 1.0 * GeV, time, depth);
771  bool emType = G4TrackToParticleID::isGammaElectronPositron(theTrack->GetDefinition()->GetPDGEncoding());
772  plotHF(hitPoint, emType);
773 #endif
774  processHit(aStep);
775  }
776  }
777 }
float edepositEM
Definition: CaloSD.h:144
Log< level::Info, true > LogVerbatim
double weight_
Definition: HCalSD.h:115
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1016
void processHit(const G4Step *step)
Definition: CaloSD.h:117
float edepositHAD
Definition: CaloSD.h:144
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:714
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1078
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:720
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:90
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:891
static bool isGammaElectronPositron(int pdgCode)
CaloHitID currentID[2]
Definition: CaloSD.h:146

◆ getFromLibrary()

bool HCalSD::getFromLibrary ( const G4Step *  aStep)
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 329 of file HCalSD.cc.

References funct::abs(), deliveredLumi, depth_, getFromHFLibrary(), getFromParam(), CaloSD::getNumberOfHits(), mps_fire::i, createfilelist::int, G4TrackToParticleID::isGammaElectronPositron(), isHF, isItHF(), G4TrackToParticleID::isMuon(), G4TrackToParticleID::isStableHadronIon(), HFDarkening::lowZLimit, m_HFDarkening, HFDarkening::numberOfZLayers, alignCSCRings::r, HLT_2024v11_cff::track, HFDarkening::upperZLimit, useParam, useShowerLibrary, weight_, and z.

329  {
330  auto const track = aStep->GetTrack();
331  depth_ = (aStep->GetPreStepPoint()->GetTouchable()->GetReplicaNumber(0)) % 10;
332  weight_ = 1.0;
333  bool kill(false);
334  isHF = isItHF(aStep);
335 #ifdef EDM_ML_DEBUG
336  edm::LogVerbatim("HcalSim") << "GetFromLibrary: "
337  << (aStep->GetPreStepPoint()->GetTouchable()->GetVolume(0)->GetName())
338 
339  << " isHF " << isHF << " darken " << (m_HFDarkening != nullptr) << " useParam "
340  << useParam << " useShowerLibrary " << useShowerLibrary << " Muon? "
341  << G4TrackToParticleID::isMuon(track) << " electron? "
342  << G4TrackToParticleID::isGammaElectronPositron(track) << " Stable Hadron? "
344 #endif
345  if (isHF) {
346  if (m_HFDarkening) {
347  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
348  const double invcm = 1. / CLHEP::cm;
349  double r = hitPoint.perp() * invcm;
350  double z = std::abs(hitPoint.z()) * invcm;
351  double dose_acquired = 0.;
353  unsigned int hfZLayer = (unsigned int)((z - HFDarkening::lowZLimit) / 5);
354  if (hfZLayer >= HFDarkening::upperZLimit)
355  hfZLayer = (HFDarkening::upperZLimit - 1);
356  float normalized_lumi = m_HFDarkening->int_lumi(deliveredLumi);
357  for (int i = hfZLayer; i != HFDarkening::numberOfZLayers; ++i) {
358  dose_acquired = m_HFDarkening->dose(i, r);
359  weight_ *= m_HFDarkening->degradation(normalized_lumi * dose_acquired);
360  }
361  }
362 #ifdef EDM_ML_DEBUG
363  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary: HFLumiDarkening at "
364  << "r= " << r << ", z= " << z << " Dose= " << dose_acquired << " weight= " << weight_;
365 #endif
366  }
367 
368  if (useParam) {
369  getFromParam(aStep, kill);
370 #ifdef EDM_ML_DEBUG
371  G4String nameVolume = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
372  edm::LogVerbatim("HcalSim") << "HCalSD: " << getNumberOfHits() << " hits from parametrization in " << nameVolume
373  << " for Track " << track->GetTrackID() << " ("
374  << track->GetDefinition()->GetParticleName() << ")";
375 #endif
378 #ifdef EDM_ML_DEBUG
379  auto nameVolume = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
380  edm::LogVerbatim("HcalSim") << "HCalSD: Starts shower library from " << nameVolume << " for Track "
381  << track->GetTrackID() << " (" << track->GetDefinition()->GetParticleName() << ")";
382 
383 #endif
384  getFromHFLibrary(aStep, kill);
385  }
386  }
387  }
388 #ifdef EDM_ML_DEBUG
389  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary ID= " << track->GetTrackID() << " ("
390  << track->GetDefinition()->GetParticleName() << ") kill= " << kill
391  << " weight= " << weight_ << " depth= " << depth_ << " isHF: " << isHF;
392 #endif
393  return kill;
394 }
int getNumberOfHits(int k=0)
Definition: CaloSD.cc:521
Log< level::Info, true > LogVerbatim
bool useParam
Definition: HCalSD.h:112
static bool isMuon(int pdgCode)
double weight_
Definition: HCalSD.h:115
double deliveredLumi
Definition: HCalSD.h:114
bool useShowerLibrary
Definition: HCalSD.h:112
static const unsigned int numberOfZLayers
Definition: HFDarkening.h:24
std::unique_ptr< HFDarkening > m_HFDarkening
Definition: HCalSD.h:100
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void getFromParam(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:821
void getFromHFLibrary(const G4Step *step, bool &isKilled)
Definition: HCalSD.cc:736
static bool isStableHadronIon(const G4Track *)
static const unsigned int lowZLimit
Definition: HFDarkening.h:27
static const unsigned int upperZLimit
Definition: HFDarkening.h:28
bool isHF
Definition: HCalSD.h:107
int depth_
Definition: HCalSD.h:116
bool isItHF(const G4Step *)
Definition: HCalSD.cc:651
static bool isGammaElectronPositron(int pdgCode)

◆ getFromParam()

void HCalSD::getFromParam ( const G4Step *  step,
bool &  isKilled 
)
private

Definition at line 821 of file HCalSD.cc.

References CaloSD::currentID, hcalRecHitTable_cff::depth, CaloSD::edepositEM, CaloSD::edepositHAD, hfClusterShapes_cfi::hits, mps_fire::i, plotProfile(), CaloSD::processHit(), setDetUnitId(), CaloHitID::setID(), CaloSD::setTrackID(), showerParam, hcalRecHitTable_cff::time, and weight_.

Referenced by getFromLibrary().

821  {
822  std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep, weight_, isKilled);
823  if (!isKilled || hits.empty()) {
824  return;
825  }
826 
827  int primaryID = setTrackID(aStep);
828  int det = 5;
829 
830 #ifdef EDM_ML_DEBUG
831  edm::LogVerbatim("HcalSim") << "HCalSD::getFromParam " << hits.size() << " hits for " << GetName() << " of "
832  << primaryID << " with " << aStep->GetTrack()->GetDefinition()->GetParticleName()
833  << " of " << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV
834  << " GeV in detector type " << det;
835 #endif
836  for (unsigned int i = 0; i < hits.size(); ++i) {
837  G4ThreeVector hitPoint = hits[i].position;
838  int depth = hits[i].depth;
839  double time = hits[i].time;
840  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
841  currentID[0].setID(unitID, time, primaryID, 0);
842  edepositEM = hits[i].edep * GeV;
843  edepositHAD = 0.;
844 #ifdef plotDebug
845  plotProfile(aStep, hitPoint, edepositEM, time, depth);
846 #endif
847  processHit(aStep);
848  }
849 }
float edepositEM
Definition: CaloSD.h:144
Log< level::Info, true > LogVerbatim
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:92
double weight_
Definition: HCalSD.h:115
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1016
void processHit(const G4Step *step)
Definition: CaloSD.h:117
float edepositHAD
Definition: CaloSD.h:144
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:891
CaloHitID currentID[2]
Definition: CaloSD.h:146

◆ getHitFibreBundle()

void HCalSD::getHitFibreBundle ( const G4Step *  step,
bool  type 
)
private

Definition at line 910 of file HCalSD.cc.

References HLT_2024v11_cff::beta, CaloSD::currentID, TauDecayModes::dec, hcalRecHitTable_cff::depth, CaloSD::edepositEM, CaloSD::edepositHAD, CaloSD::energyCut, HcalForward, numberingFromDDD, phi, plotProfile(), CaloSD::processHit(), CaloSD::resetForNewPrimary(), findQualityFiles::rr, setDetUnitId(), CaloHitID::setID(), showerBundle, hcalRecHitTable_cff::time, and createJobs::tmp.

Referenced by getEnergyDeposit().

910  {
911  auto const preStepPoint = aStep->GetPreStepPoint();
912  auto const theTrack = aStep->GetTrack();
913  double edep = showerBundle->getHits(aStep, type);
914 
915  if (edep >= 0.0) {
916  edep *= GeV;
917  double etrack = preStepPoint->GetKineticEnergy();
918  int primaryID = 0;
919  if (etrack >= energyCut) {
920  primaryID = theTrack->GetTrackID();
921  } else {
922  primaryID = theTrack->GetParentID();
923  if (primaryID == 0)
924  primaryID = theTrack->GetTrackID();
925  }
926  // Reset entry point for new primary
927  resetForNewPrimary(aStep);
928  //
929  int det = static_cast<int>(HcalForward);
930  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
931  double rr = hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y();
932  double phi = rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x());
933  double etaR = showerBundle->getRadius();
934  int depth = 1;
935  if (etaR < 0.) {
936  depth = 2;
937  etaR = -etaR;
938  }
939  if (hitPoint.z() < 0.)
940  etaR = -etaR;
941 #ifdef EDM_ML_DEBUG
942  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
943  << " depth " << depth;
944 #endif
945  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
946  uint32_t unitID = 0;
947  if (numberingFromDDD) {
948  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
949  unitID = setDetUnitId(tmp);
950  }
951  if (type)
952  currentID[0].setID(unitID, time, primaryID, 3);
953  else
954  currentID[0].setID(unitID, time, primaryID, 2);
955 
956  edepositHAD = aStep->GetTotalEnergyDeposit();
957  edepositEM = -edepositHAD + edep;
958 #ifdef plotDebug
959  plotProfile(aStep, hitPoint, edep, time, depth);
960 #endif
961 #ifdef EDM_ML_DEBUG
962  double beta = preStepPoint->GetBeta();
963  edm::LogVerbatim("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName() << " of " << primaryID
964  << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
965  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
966  << std::hex << unitID << std::dec;
967 #endif
968  processHit(aStep);
969  } // non-zero energy deposit
970 }
float edepositEM
Definition: CaloSD.h:144
double energyCut
Definition: CaloSD.h:148
Log< level::Info, true > LogVerbatim
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1016
void processHit(const G4Step *step)
Definition: CaloSD.h:117
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
float edepositHAD
Definition: CaloSD.h:144
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:714
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:94
CaloHitID currentID[2]
Definition: CaloSD.h:146
tmp
align.sh
Definition: createJobs.py:716

◆ getHitPMT()

void HCalSD::getHitPMT ( const G4Step *  step)
private

Definition at line 851 of file HCalSD.cc.

References HLT_2024v11_cff::beta, CaloSD::currentID, TauDecayModes::dec, hcalRecHitTable_cff::depth, CaloSD::edepositEM, CaloSD::edepositHAD, CaloSD::energyCut, HcalForward, numberingFromDDD, phi, plotProfile(), CaloSD::processHit(), CaloSD::resetForNewPrimary(), findQualityFiles::rr, setDetUnitId(), CaloHitID::setID(), showerPMT, hcalRecHitTable_cff::time, and createJobs::tmp.

Referenced by getEnergyDeposit().

851  {
852  auto const preStepPoint = aStep->GetPreStepPoint();
853  auto const theTrack = aStep->GetTrack();
854  double edep = showerPMT->getHits(aStep);
855 
856  if (edep >= 0.) {
857  edep *= GeV;
858  double etrack = preStepPoint->GetKineticEnergy();
859  int primaryID = 0;
860  if (etrack >= energyCut) {
861  primaryID = theTrack->GetTrackID();
862  } else {
863  primaryID = theTrack->GetParentID();
864  if (primaryID == 0)
865  primaryID = theTrack->GetTrackID();
866  }
867  // Reset entry point for new primary
868  resetForNewPrimary(aStep);
869  //
870  int det = static_cast<int>(HcalForward);
871  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
872  double rr = (hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y());
873  double phi = (rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x()));
874  double etaR = showerPMT->getRadius();
875  int depth = 1;
876  if (etaR < 0) {
877  depth = 2;
878  etaR = -etaR;
879  }
880  if (hitPoint.z() < 0)
881  etaR = -etaR;
882 #ifdef EDM_ML_DEBUG
883  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
884  << " depth " << depth;
885 #endif
886  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
887  uint32_t unitID = 0;
888  if (numberingFromDDD) {
889  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
890  unitID = setDetUnitId(tmp);
891  }
892  currentID[0].setID(unitID, time, primaryID, 1);
893 
894  edepositHAD = aStep->GetTotalEnergyDeposit();
895  edepositEM = -edepositHAD + edep;
896 #ifdef plotDebug
897  plotProfile(aStep, hitPoint, edep, time, depth);
898 #endif
899 #ifdef EDM_ML_DEBUG
900  double beta = preStepPoint->GetBeta();
901  edm::LogVerbatim("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() << " of " << primaryID << " with "
902  << theTrack->GetDefinition()->GetParticleName() << " of "
903  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
904  << std::hex << unitID << std::dec;
905 #endif
906  processHit(aStep);
907  }
908 }
float edepositEM
Definition: CaloSD.h:144
double energyCut
Definition: CaloSD.h:148
Log< level::Info, true > LogVerbatim
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1016
void processHit(const G4Step *step)
Definition: CaloSD.h:117
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:93
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
float edepositHAD
Definition: CaloSD.h:144
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:714
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
CaloHitID currentID[2]
Definition: CaloSD.h:146
tmp
align.sh
Definition: createJobs.py:716

◆ hitForFibre()

void HCalSD::hitForFibre ( const G4Step *  step)
private

Definition at line 779 of file HCalSD.cc.

References CaloSD::currentID, hcalRecHitTable_cff::depth, CaloSD::edepositEM, CaloSD::edepositHAD, hfshower, hfClusterShapes_cfi::hits, mps_fire::i, G4TrackToParticleID::isGammaElectronPositron(), isItinFidVolume(), plotHF(), plotProfile(), CaloSD::processHit(), setDetUnitId(), CaloHitID::setID(), CaloSD::setTrackID(), hcalRecHitTable_cff::time, funct::true, and weight_.

Referenced by getEnergyDeposit().

779  { // if not ParamShower
780 
781  std::vector<HFShower::Hit> hits = hfshower->getHits(aStep, weight_);
782  if (hits.empty()) {
783  return;
784  }
785 
786  auto const theTrack = aStep->GetTrack();
787  int primaryID = setTrackID(aStep);
788  int det = 5;
789 
791  edepositEM = 1. * GeV;
792  edepositHAD = 0.;
793  } else {
794  edepositEM = 0.;
795  edepositHAD = 1. * GeV;
796  }
797 
798 #ifdef EDM_ML_DEBUG
799  edm::LogVerbatim("HcalSim") << "HCalSD::hitForFibre " << hits.size() << " hits for " << GetName() << " of "
800  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
801  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV in detector type " << det;
802 #endif
803 
804  for (unsigned int i = 0; i < hits.size(); ++i) {
805  G4ThreeVector hitPoint = hits[i].position;
806  if (isItinFidVolume(hitPoint)) {
807  int depth = hits[i].depth;
808  double time = hits[i].time;
809  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
810  currentID[0].setID(unitID, time, primaryID, 0);
811 #ifdef plotDebug
812  plotProfile(aStep, hitPoint, edepositEM, time, depth);
813  bool emType = (edepositEM > 0.) ? true : false;
814  plotHF(hitPoint, emType);
815 #endif
816  processHit(aStep);
817  }
818  }
819 }
float edepositEM
Definition: CaloSD.h:144
Log< level::Info, true > LogVerbatim
double weight_
Definition: HCalSD.h:115
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:1016
void processHit(const G4Step *step)
Definition: CaloSD.h:117
float edepositHAD
Definition: CaloSD.h:144
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1078
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:720
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:891
static bool isGammaElectronPositron(int pdgCode)
CaloHitID currentID[2]
Definition: CaloSD.h:146
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:91

◆ initEvent()

void HCalSD::initEvent ( const BeginOfEvent )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 1103 of file HCalSD.cc.

References detNull_.

1103  {
1104 #ifdef printDebug
1105  detNull_ = {0, 0, 0, 0};
1106 #endif
1107 }
std::vector< int > detNull_
Definition: HCalSD.h:127

◆ initRun()

void HCalSD::initRun ( )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 572 of file HCalSD.cc.

572 {}

◆ isItConicalBundle()

bool HCalSD::isItConicalBundle ( const G4LogicalVolume *  lv)
private

Definition at line 704 of file HCalSD.cc.

References fibre2LV.

Referenced by getEnergyDeposit().

704  {
705  for (auto lvol : fibre2LV)
706  if (lv == lvol) {
707  return true;
708  }
709  return false;
710 }
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:124

◆ isItFibre() [1/2]

bool HCalSD::isItFibre ( const G4LogicalVolume *  lv)
private

Definition at line 672 of file HCalSD.cc.

References fibreLV.

Referenced by getEnergyDeposit().

672  {
673  for (auto lvol : fibreLV)
674  if (lv == lvol) {
675  return true;
676  }
677  return false;
678 }
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:124

◆ isItFibre() [2/2]

bool HCalSD::isItFibre ( const G4String &  name)
private

Definition at line 680 of file HCalSD.cc.

References fibreNames, and Skims_PA_cff::name.

680  {
681  for (const auto& nam : fibreNames)
682  if (name == static_cast<G4String>(nam)) {
683  return true;
684  }
685  return false;
686 }
std::vector< std::string > fibreNames
Definition: HCalSD.h:121

◆ isItHF() [1/2]

bool HCalSD::isItHF ( const G4Step *  aStep)
private

Definition at line 651 of file HCalSD.cc.

References hfLevels, hfLV, hfNames, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, and GCP_tree_cfg::levels.

Referenced by getFromLibrary().

651  {
652  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
653  int levels = (touch->GetHistoryDepth()) + 1;
654  for (unsigned int it = 0; it < hfNames.size(); ++it) {
655  if (levels >= hfLevels[it]) {
656  const G4LogicalVolume* lv = touch->GetVolume(levels - hfLevels[it])->GetLogicalVolume();
657  if (lv == hfLV[it])
658  return true;
659  }
660  }
661  return false;
662 }
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:124
std::vector< std::string > hfNames
Definition: HCalSD.h:120
std::vector< int > hfLevels
Definition: HCalSD.h:119

◆ isItHF() [2/2]

bool HCalSD::isItHF ( const G4String &  name)
private

Definition at line 664 of file HCalSD.cc.

References hfNames, and Skims_PA_cff::name.

664  {
665  for (const auto& nam : hfNames)
666  if (name == static_cast<G4String>(nam)) {
667  return true;
668  }
669  return false;
670 }
std::vector< std::string > hfNames
Definition: HCalSD.h:120

◆ isItinFidVolume()

bool HCalSD::isItinFidVolume ( const G4ThreeVector &  hitPoint)
private

Definition at line 720 of file HCalSD.cc.

References applyFidCut, RemoveAddSevLevel::flag, and HFFibreFiducial::PMTNumber().

Referenced by getFromHFLibrary(), and hitForFibre().

720  {
721  bool flag = true;
722  if (applyFidCut) {
723  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
724 #ifdef EDM_ML_DEBUG
725  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume:#PMT= " << npmt << " for hit point " << hitPoint;
726 #endif
727  if (npmt <= 0)
728  flag = false;
729  }
730 #ifdef EDM_ML_DEBUG
731  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume: point " << hitPoint << " return flag " << flag;
732 #endif
733  return flag;
734 }
Log< level::Info, true > LogVerbatim
bool applyFidCut
Definition: HCalSD.h:112
int PMTNumber(const G4ThreeVector &pe_effect)

◆ isItPMT()

bool HCalSD::isItPMT ( const G4LogicalVolume *  lv)
private

Definition at line 688 of file HCalSD.cc.

References pmtLV.

Referenced by getEnergyDeposit().

688  {
689  for (auto lvol : pmtLV)
690  if (lv == lvol) {
691  return true;
692  }
693  return false;
694 }
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:124

◆ isItScintillator()

bool HCalSD::isItScintillator ( const G4Material *  mat)
private

Definition at line 712 of file HCalSD.cc.

References materials.

Referenced by getEnergyDeposit().

712  {
713  for (auto amat : materials)
714  if (amat == mat) {
715  return true;
716  }
717  return false;
718 }
std::vector< const G4Material * > materials
Definition: HCalSD.h:123

◆ isItStraightBundle()

bool HCalSD::isItStraightBundle ( const G4LogicalVolume *  lv)
private

Definition at line 696 of file HCalSD.cc.

References fibre1LV.

Referenced by getEnergyDeposit().

696  {
697  for (auto lvol : fibre1LV)
698  if (lv == lvol) {
699  return true;
700  }
701  return false;
702 }
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:124

◆ layerWeight()

double HCalSD::layerWeight ( int  det,
const G4ThreeVector &  pos,
int  depth,
int  lay 
)
private

Definition at line 996 of file HCalSD.cc.

References TauDecayModes::dec, hcalRecHitTable_cff::depth, EDM_ML_DEBUG, layerWeights, modifyDepth(), numberingFromDDD, HcalTestNumbering::packHcalIndex(), and createJobs::tmp.

Referenced by getEnergyDeposit().

996  {
997  double wt = 1.;
998  if (numberingFromDDD) {
999  //get the ID's as eta, phi, depth, ... indices
1001  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
1002  modifyDepth(tmp);
1003  uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1, tmp.etaR, tmp.phis, tmp.lay);
1004  std::map<uint32_t, double>::const_iterator ite = layerWeights.find(id);
1005  if (ite != layerWeights.end())
1006  wt = ite->second;
1007 #ifdef EDM_ML_DEBUG
1008  edm::LogVerbatim("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id << std::dec << " (" << tmp.subdet << "/"
1009  << tmp.zside << "/1/" << tmp.etaR << "/" << tmp.phis << "/" << tmp.lay << ") Weight "
1010  << wt;
1011 #endif
1012  }
1013  return wt;
1014 }
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1089
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
#define EDM_ML_DEBUG
Definition: MPUnpacker.cc:1
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:125
tmp
align.sh
Definition: createJobs.py:716

◆ modifyDepth()

void HCalSD::modifyDepth ( HcalNumberingFromDDD::HcalID id)
private

Definition at line 1089 of file HCalSD.cc.

References hcalRecHitTable_cff::depth, depth_, hcalConstants_, hcalRecHitTable_cff::ieta, HcalDDDSimConstants::maxHFDepth(), and testNumber.

Referenced by layerWeight(), and setDetUnitId().

1089  {
1090  if (id.subdet == 4) {
1091  int ieta = (id.zside == 0) ? -id.etaR : id.etaR;
1092  if (hcalConstants_->maxHFDepth(ieta, id.phis) > 2) {
1093  if (id.depth <= 2) {
1094  if (G4UniformRand() > 0.5)
1095  id.depth += 2;
1096  }
1097  }
1098  } else if ((id.subdet == 1 || id.subdet == 2) && testNumber) {
1099  id.depth = (depth_ == 0) ? 1 : 2;
1100  }
1101 }
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:96
int maxHFDepth(const int &ieta, const int &iphi) const
bool testNumber
Definition: HCalSD.h:110
int depth_
Definition: HCalSD.h:116

◆ plotHF()

void HCalSD::plotHF ( const G4ThreeVector &  pos,
bool  emType 
)
private

Definition at line 1078 of file HCalSD.cc.

References funct::abs(), gpar, hzvem, hzvhad, and ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::zv.

Referenced by getFromHFLibrary(), and hitForFibre().

1078  {
1079  double zv = std::abs(hitPoint.z()) - gpar[4];
1080  if (emType) {
1081  if (hzvem != nullptr)
1082  hzvem->Fill(zv);
1083  } else {
1084  if (hzvhad != nullptr)
1085  hzvhad->Fill(zv);
1086  }
1087 }
std::vector< double > gpar
Definition: HCalSD.h:118
TH1F * hzvhad
Definition: HCalSD.h:126
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TH1F * hzvem
Definition: HCalSD.h:126

◆ plotProfile()

void HCalSD::plotProfile ( const G4Step *  step,
const G4ThreeVector &  pos,
double  edep,
double  time,
int  id 
)
private

Definition at line 1016 of file HCalSD.cc.

References funct::abs(), dd4hep_, hcalRecHitTable_cff::depth, dist_, newFWLiteAna::found, hit_, heavyIonCSV_trainingSettings::idx, cuy::ii, DTRecHitClients_cfi::local, dqmiodumpmetadata::n, Skims_PA_cff::name, DD4hep2DDDName::nameMatterLV(), names, hcalRecHitTable_cff::time, and time_.

Referenced by getFromHFLibrary(), getFromParam(), getHitFibreBundle(), getHitPMT(), and hitForFibre().

1016  {
1017  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
1018  static const unsigned int names = 10;
1019  static const G4String modName[names] = {
1020  "HEModule", "HVQF", "HBModule", "MBAT", "MBBT", "MBBTC", "MBBT_R1P", "MBBT_R1M", "MBBT_R1PX", "MBBT_R1MX"};
1021  G4ThreeVector local;
1022  bool found = false;
1023  double depth = -2000;
1024  int idx = 4;
1025  for (int n = 0; n < touch->GetHistoryDepth(); ++n) {
1026  G4String name(static_cast<G4String>(DD4hep2DDDName::nameMatterLV(touch->GetVolume(n)->GetName(), dd4hep_)));
1027 #ifdef EDM_ML_DEBUG
1028  edm::LogVerbatim("HcalSim") << "plotProfile Depth " << n << " Name " << name;
1029 #endif
1030  for (unsigned int ii = 0; ii < names; ++ii) {
1031  if (name == modName[ii]) {
1032  found = true;
1033  int dn = touch->GetHistoryDepth() - n;
1034  local = touch->GetHistory()->GetTransform(dn).TransformPoint(global);
1035  if (ii == 0) {
1036  depth = local.z() - 4006.5;
1037  idx = 1;
1038  } else if (ii == 1) {
1039  depth = local.z() + 825.0;
1040  idx = 3;
1041  } else if (ii == 2) {
1042  depth = local.x() - 1775.;
1043  idx = 0;
1044  } else {
1045  depth = local.y() + 15.;
1046  idx = 2;
1047  }
1048  break;
1049  }
1050  }
1051  if (found)
1052  break;
1053  }
1054  if (!found)
1055  depth = std::abs(global.z()) - 11500;
1056 #ifdef EDM_ML_DEBUG
1057  edm::LogVerbatim("HcalSim") << "plotProfile Found " << found << " Global " << global << " Local " << local
1058  << " depth " << depth << " ID " << id << " EDEP " << edep << " Time " << time;
1059 #endif
1060  if (hit_[idx] != nullptr)
1061  hit_[idx]->Fill(edep);
1062  if (time_[idx] != nullptr)
1063  time_[idx]->Fill(time, edep);
1064  if (dist_[idx] != nullptr)
1065  dist_[idx]->Fill(depth, edep);
1066  int jd = 2 * idx + id - 7;
1067  if (jd >= 0 && jd < 4) {
1068  jd += 5;
1069  if (hit_[jd] != nullptr)
1070  hit_[jd]->Fill(edep);
1071  if (time_[jd] != nullptr)
1072  time_[jd]->Fill(time, edep);
1073  if (dist_[jd] != nullptr)
1074  dist_[jd]->Fill(depth, edep);
1075  }
1076 }
Log< level::Info, true > LogVerbatim
TH1F * time_[9]
Definition: HCalSD.h:126
const std::string names[nVars_]
TH1F * hit_[9]
Definition: HCalSD.h:126
bool dd4hep_
Definition: HCalSD.h:117
TH1F * dist_[9]
Definition: HCalSD.h:126
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ii
Definition: cuy.py:589
std::string nameMatterLV(const std::string &name, bool dd4hep)

◆ printVolume()

void HCalSD::printVolume ( const G4VTouchable *  touch) const
private

Definition at line 1118 of file HCalSD.cc.

References mps_fire::i, cuy::ii, personalPlayback::level, Skims_PA_cff::name, AlCaHLTBitMon_QueryRunRegistry::string, and susybsm::HSCParticleType::unknown.

Referenced by setDetUnitId().

1118  {
1119  if (touch) {
1120 #ifdef EDM_ML_DEBUG
1121  int level = ((touch->GetHistoryDepth()) + 1);
1122  edm::LogVerbatim("CaloSimX") << "HCalSD::printVolume with " << level << " levels";
1123  static const std::string unknown("Unknown");
1124  //Get name and copy numbers
1125  for (int ii = 0; ii < level; ii++) {
1126  int i = level - ii - 1;
1127  G4VPhysicalVolume* pv = touch->GetVolume(i);
1128  G4String name = (pv != nullptr) ? pv->GetName() : unknown;
1129  G4int copyno = touch->GetReplicaNumber(i);
1130  edm::LogVerbatim("HcalSim") << "[" << ii << "] " << name << ":" << copyno;
1131  }
1132 #endif
1133  }
1134 }
Log< level::Info, true > LogVerbatim
ii
Definition: cuy.py:589

◆ readWeightFromFile()

void HCalSD::readWeightFromFile ( const std::string &  fName)
private

Definition at line 972 of file HCalSD.cc.

References TauDecayModes::dec, mps_splice::entry, fileCollector2::fName, recoMuon::in, timingPdfMaker::infile, layerWeights, HcalTestNumbering::packHcalIndex(), phi, useLayerWt, and ecaldqm::zside().

Referenced by HCalSD().

972  {
973  std::ifstream infile;
974  int entry = 0;
975  infile.open(fName.c_str(), std::ios::in);
976  if (infile) {
977  int det, zside, etaR, phi, lay;
978  double wt;
979  while (infile >> det >> zside >> etaR >> phi >> lay >> wt) {
980  uint32_t id = HcalTestNumbering::packHcalIndex(det, zside, 1, etaR, phi, lay);
981  layerWeights.insert(std::pair<uint32_t, double>(id, wt));
982  ++entry;
983 #ifdef EDM_ML_DEBUG
984  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry << " ID " << std::hex << id
985  << std::dec << " (" << det << "/" << zside << "/1/" << etaR << "/" << phi << "/"
986  << lay << ") Weight " << wt;
987 #endif
988  }
989  infile.close();
990  }
991  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry << " weights from " << fName;
992  if (entry <= 0)
993  useLayerWt = false;
994 }
Log< level::Info, true > LogVerbatim
bool useLayerWt
Definition: HCalSD.h:109
int zside(DetId const &)
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
std::map< uint32_t, double > layerWeights
Definition: HCalSD.h:125

◆ setDetUnitId() [1/3]

uint32_t HCalSD::setDetUnitId ( const G4Step *  step)
overridevirtual

Implements CaloSD.

Definition at line 530 of file HCalSD.cc.

References TauDecayModes::dec, hcalRecHitTable_cff::depth, HcalDetId::depth(), heavyIonCSV_trainingSettings::idx, HcalDetId::ietaAbs(), HcalDetId::iphi(), printVolume(), DetId::subdetId(), testNumber, HcalTestNumbering::unpackHcalIndex(), HLTMuonOfflineAnalyzer_cfi::z0, and HcalDetId::zside().

Referenced by getEnergyDeposit(), getFromHFLibrary(), getFromParam(), getHitFibreBundle(), getHitPMT(), hitForFibre(), and setDetUnitId().

530  {
531  auto const prePoint = aStep->GetPreStepPoint();
532  auto const touch = prePoint->GetTouchable();
533  const G4ThreeVector& hitPoint = prePoint->GetPosition();
534 
535  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
536  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
537  int det = (touch->GetReplicaNumber(1)) / 1000;
538 
539  uint32_t idx = setDetUnitId(det, hitPoint, depth, lay);
540 #ifdef EDM_ML_DEBUG
541  if (depth == 1) {
542  edm::LogVerbatim("HcalSim") << "HCalSD: Check for " << det << ":" << depth << ":" << lay << " ID " << std::hex
543  << idx << std::dec;
544  int det0, z0, depth0, eta0, phi0, lay0(-1);
545  if (testNumber) {
546  HcalTestNumbering::unpackHcalIndex(idx, det0, z0, depth0, eta0, phi0, lay0);
547  } else {
548  HcalDetId hcid0(idx);
549  det0 = hcid0.subdetId();
550  eta0 = hcid0.ietaAbs();
551  phi0 = hcid0.iphi();
552  z0 = hcid0.zside();
553  depth0 = hcid0.depth();
554  }
555  edm::LogVerbatim("HcalSim") << "HCalSD: det|z|depth|eta|phi|lay " << det0 << ":" << z0 << ":" << depth0 << ":"
556  << eta0 << ":" << phi0 << ":" << lay0;
557  printVolume(touch);
558  }
559 #endif
560  return idx;
561 }
Log< level::Info, true > LogVerbatim
bool testNumber
Definition: HCalSD.h:110
static void unpackHcalIndex(const uint32_t &idx, int &det, int &z, int &depth, int &eta, int &phi, int &lay)
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
void printVolume(const G4VTouchable *touch) const
Definition: HCalSD.cc:1118

◆ setDetUnitId() [2/3]

uint32_t HCalSD::setDetUnitId ( int  det,
const G4ThreeVector &  pos,
int  depth,
int  lay = 1 
)
private

Definition at line 596 of file HCalSD.cc.

References funct::abs(), hcalRecHitTable_cff::depth, detNull_, PVValHelper::eta, hcalConstants_, l1ctLayer2EG_cff::id, HcalDDDSimConstants::isHE(), maxRoff_, maxZ_, minRoff_, numberingFromDDD, setDetUnitId(), slopeHE_, and createJobs::tmp.

596  {
597  uint32_t id = 0;
598  if (det == 0) {
599 #ifdef printDebug
600  double eta = std::abs(pos.eta());
601 #endif
602  if (std::abs(pos.z()) > maxZ_) {
603  det = 5;
604 #ifdef printDebug
605  if (eta < 2.868)
606  ++detNull_[2];
607 #endif
608  } else if (!(hcalConstants_->isHE())) {
609  det = 3;
610 #ifdef printDebug
611  ++detNull_[0];
612 #endif
613  } else {
614  double minR = minRoff_ + slopeHE_ * std::abs(pos.z());
615  double maxR = maxRoff_ + slopeHE_ * std::abs(pos.z());
616  det = ((pos.perp() > minR) && (pos.perp() < maxR)) ? 4 : 3;
617 #ifdef printDebug
618  ++detNull_[det - 3];
619 #endif
620  }
621 #ifdef printDEBUG
622  edm::LogVerbatim("HcalSim") << "Position " << pos.perp() << ":" << std::abs(pos.z()) << " Limits "
623  << !(hcalConstants_->isHE()) << ":" << maxZ_ << " det " << det;
624  } else {
625  ++detNull_[3];
626 #endif
627  }
628 
629  if (numberingFromDDD.get()) {
630  //get the ID's as eta, phi, depth, ... indices
632  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
633  id = setDetUnitId(tmp);
634  }
635  return id;
636 }
Log< level::Info, true > LogVerbatim
static constexpr double slopeHE_
Definition: HCalSD.h:106
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:96
static constexpr double minRoff_
Definition: HCalSD.h:104
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:88
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< double > > XYZVectorD
spatial vector with cartesian internal representation
Definition: Vector3D.h:8
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:530
static constexpr double maxZ_
Definition: HCalSD.h:103
std::vector< int > detNull_
Definition: HCalSD.h:127
tmp
align.sh
Definition: createJobs.py:716
static constexpr double maxRoff_
Definition: HCalSD.h:105

◆ setDetUnitId() [3/3]

uint32_t HCalSD::setDetUnitId ( HcalNumberingFromDDD::HcalID tmp)
private

Definition at line 638 of file HCalSD.cc.

References TauDecayModes::dec, l1ctLayer2EG_cff::id, m_HcalTestNS, modifyDepth(), numberingScheme, convertSQLiteXML::ok, testNumber, and createJobs::tmp.

638  {
639  modifyDepth(tmp);
640  uint32_t id = (numberingScheme.get()) ? numberingScheme->getUnitID(tmp) : 0;
641  if ((!testNumber) && m_HcalTestNS.get()) {
642  bool ok = m_HcalTestNS->compare(tmp, id);
643  if (!ok)
644  edm::LogWarning("HcalSim") << "Det ID from HCalSD " << HcalDetId(id) << " " << std::hex << id << std::dec
645  << " does not match one from relabller for " << tmp.subdet << ":" << tmp.etaR << ":"
646  << tmp.phi << ":" << tmp.phis << ":" << tmp.depth << ":" << tmp.lay << std::endl;
647  }
648  return id;
649 }
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1089
bool testNumber
Definition: HCalSD.h:110
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:101
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:89
Log< level::Warning, false > LogWarning
tmp
align.sh
Definition: createJobs.py:716

◆ setNumberingScheme()

void HCalSD::setNumberingScheme ( HcalNumberingScheme scheme)

Definition at line 563 of file HCalSD.cc.

References numberingScheme, and generator_cfi::scheme.

Referenced by HCalSD(), HcalTestAnalysis::update(), SimG4HcalValidation::update(), and HcalTB04Analysis::update().

563  {
564  if (scheme != nullptr) {
565  edm::LogVerbatim("HcalSim") << "HCalSD: updates numbering scheme for " << GetName();
566  numberingScheme.reset(scheme);
567  }
568 }
Log< level::Info, true > LogVerbatim
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:89

◆ update() [1/6]

void CaloSD::update
overrideprotected

Definition at line 751 of file CaloSD.cc.

Referenced by MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), and MatrixUtil.Steps::overwrite().

751  {
752 #ifdef EDM_ML_DEBUG
753  edm::LogVerbatim("CaloSim") << "CaloSD: Dispatched BeginOfEvent for " << GetName() << " !";
754 #endif
755  clearHits();
756  initEvent(g4Event);
757 }
Log< level::Info, true > LogVerbatim
void initEvent(const BeginOfEvent *) override
Definition: HCalSD.cc:1103
void clearHits() override
Definition: CaloSD.cc:848

◆ update() [2/6]

void CaloSD::update
overrideprotected

Definition at line 780 of file CaloSD.cc.

Referenced by MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), and MatrixUtil.Steps::overwrite().

780  {
781  endEvent();
782  for (int k = 0; k < nHC_; ++k) {
783  slave[k].get()->ReserveMemory(theHC[k]->entries());
784 
785  int count(0);
786  double eEM(0.0);
787  double eHAD(0.0);
788  double eEM2(0.0);
789  double eHAD2(0.0);
790 #ifdef EDM_ML_DEBUG
791  int wrong(0);
792  double tt(0.0);
793  double zloc(0.0);
794  double zglob(0.0);
795  double ee(0.0);
796 #endif
797  int hc_entries = theHC[k]->entries();
798  for (int i = 0; i < hc_entries; ++i) {
799 #ifdef EDM_ML_DEBUG
800  if (!saveHit((*theHC[k])[i], k)) {
801  ++wrong;
802  }
803 #else
804  saveHit((*theHC[k])[i], k);
805 #endif
806 
807  ++count;
808  double x = (*theHC[k])[i]->getEM();
809  eEM += x;
810  eEM2 += x * x;
811  x = (*theHC[k])[i]->getHadr();
812  eHAD += x;
813  eHAD2 += x * x;
814 #ifdef EDM_ML_DEBUG
815  tt += (*theHC[k])[i]->getTimeSlice();
816  ee += (*theHC[k])[i]->getIncidentEnergy();
817  zglob += std::abs((*theHC[k])[i]->getEntry().z());
818  zloc += std::abs((*theHC[k])[i]->getEntryLocal().z());
819 #endif
820  }
821 
822  double norm = (count > 0) ? 1.0 / count : 0.0;
823  eEM *= norm;
824  eEM2 *= norm;
825  eHAD *= norm;
826  eHAD2 *= norm;
827  eEM2 = std::sqrt(eEM2 - eEM * eEM);
828  eHAD2 = std::sqrt(eHAD2 - eHAD * eHAD);
829 #ifdef EDM_ML_DEBUG
830  tt *= norm;
831  ee *= norm;
832  zglob *= norm;
833  zloc *= norm;
834  edm::LogVerbatim("CaloSim") << "CaloSD: " << GetName() << " store " << count << " hits; " << wrong
835  << " track IDs not given properly and " << totalHits - count
836  << " hits not passing cuts\n EmeanEM= " << eEM << " ErmsEM= " << eEM2
837  << "\n EmeanHAD= " << eHAD << " ErmsHAD= " << eHAD2 << " TimeMean= " << tt
838  << " E0mean= " << ee << " Zglob= " << zglob << " Zloc= " << zloc << " ";
839 #endif
840  std::vector<std::unique_ptr<CaloG4Hit>>().swap(reusehit[k]);
841  if (useMap)
842  hitMap[k].erase(hitMap[k].begin(), hitMap[k].end());
843  }
844  tkMap.erase(tkMap.begin(), tkMap.end());
846 }
Int_t getEntry(TBranch *branch, EntryNumber entryNumber)
Definition: RootTree.cc:527
Log< level::Info, true > LogVerbatim
std::map< CaloHitID, CaloG4Hit * > hitMap[2]
Definition: CaloSD.h:195
bool useMap
Definition: CaloSD.h:179
void endEvent() override
Definition: HCalSD.cc:1109
CaloG4HitCollection * theHC[2]
Definition: CaloSD.h:174
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:112
bool saveHit(CaloG4Hit *, int k=0)
Definition: CaloSD.cc:949
std::unordered_map< unsigned int, unsigned int > boundaryCrossingParentMap_
Definition: CaloSD.h:197
int nHC_
Definition: CaloSD.h:158
Definition: TTTypes.h:54
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::vector< std::unique_ptr< CaloG4Hit > > reusehit[2]
Definition: CaloSD.h:198
int totalHits[2]
Definition: CaloSD.h:185
std::map< int, TrackWithHistory * > tkMap
Definition: CaloSD.h:196
std::unique_ptr< CaloSlaveSD > slave[2]
Definition: CaloSD.h:171

◆ update() [3/6]

void CaloSD::update
overrideprotected

Definition at line 1013 of file CaloSD.cc.

Referenced by MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), and MatrixUtil.Steps::overwrite().

1013  {
1014  int primary = -1;
1015  TrackInformation* trkInfo = cmsTrackInformation((*trk)());
1016  if (trkInfo->isPrimary())
1017  primary = (*trk)()->GetTrackID();
1018 
1019 #ifdef EDM_ML_DEBUG
1020  edm::LogVerbatim("CaloSim") << "New track: isPrimary " << trkInfo->isPrimary() << " primary ID = " << primary
1021  << " primary ancestor ID " << primAncestor;
1022 #endif
1023 
1024  // update the information if a different primary track ID
1025 
1026  if (primary > 0 && primary != primAncestor) {
1027  primAncestor = primary;
1028 
1029  // clean the hits information
1030 
1031  bool clean(false);
1032  for (int k = 0; k < nHC_; ++k)
1033  if (theHC[k]->entries() > 0)
1034  clean = true;
1035  if (clean)
1037  }
1038 }
Log< level::Info, true > LogVerbatim
bool isPrimary() const
CaloG4HitCollection * theHC[2]
Definition: CaloSD.h:174
int primAncestor
Definition: CaloSD.h:183
void cleanHitCollection()
Definition: CaloSD.cc:1040
static void clean(char *s)
int nHC_
Definition: CaloSD.h:158
TrackInformation * cmsTrackInformation(const G4Track *aTrack)

◆ update() [4/6]

void CaloSD::update
overrideprotected

Definition at line 759 of file CaloSD.cc.

Referenced by MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), and MatrixUtil.Steps::overwrite().

759  {
760  int id = (*trk)()->GetTrackID();
761  TrackInformation* trkI = cmsTrackInformation((*trk)());
762  int lastTrackID = -1;
763  if (trkI)
764  lastTrackID = trkI->getIDonCaloSurface();
765  if (id == lastTrackID) {
766  auto trksForThisEvent = m_trackManager->trackContainer();
767  if (!trksForThisEvent->empty()) {
768  TrackWithHistory* trkH = trksForThisEvent->back();
769  if (trkH->trackID() == id) {
770  tkMap[id] = trkH;
771 #ifdef EDM_ML_DEBUG
772  edm::LogVerbatim("CaloSim") << "CaloSD: get track " << id << " from Container of size "
773  << trksForThisEvent->size() << " with ID " << trkH->trackID();
774 #endif
775  }
776  }
777  }
778 }
Log< level::Info, true > LogVerbatim
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
int getIDonCaloSurface() const
const std::vector< TrackWithHistory * > * trackContainer() const
int trackID() const
std::map< int, TrackWithHistory * > tkMap
Definition: CaloSD.h:196
const SimTrackManager * m_trackManager
Definition: CaloSD.h:169

◆ update() [5/6]

void CaloSD::update
overrideprotected

Definition at line 749 of file CaloSD.cc.

Referenced by MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), and MatrixUtil.Steps::overwrite().

749 { initRun(); }
void initRun() override
Definition: HCalSD.cc:572

◆ update() [6/6]

void HCalSD::update ( const BeginOfJob )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfJob *>.

Definition at line 570 of file HCalSD.cc.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

570 {}

Member Data Documentation

◆ agingFlagHB

bool HCalSD::agingFlagHB
private

Definition at line 108 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ agingFlagHE

bool HCalSD::agingFlagHE
private

Definition at line 108 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ applyFidCut

bool HCalSD::applyFidCut
private

Definition at line 112 of file HCalSD.h.

Referenced by HCalSD(), and isItinFidVolume().

◆ betaThr

double HCalSD::betaThr
private

Definition at line 111 of file HCalSD.h.

Referenced by HCalSD().

◆ birk1

double HCalSD::birk1
private

Definition at line 111 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ birk2

double HCalSD::birk2
private

Definition at line 111 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ birk3

double HCalSD::birk3
private

Definition at line 111 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ dd4hep_

bool HCalSD::dd4hep_
private

Definition at line 117 of file HCalSD.h.

Referenced by fillLogVolumeVector(), HCalSD(), and plotProfile().

◆ deliveredLumi

double HCalSD::deliveredLumi
private

Definition at line 114 of file HCalSD.h.

Referenced by getEnergyDeposit(), getFromLibrary(), and HCalSD().

◆ depth_

int HCalSD::depth_
private

Definition at line 116 of file HCalSD.h.

Referenced by getEnergyDeposit(), getFromLibrary(), and modifyDepth().

◆ detNull_

std::vector<int> HCalSD::detNull_
private

Definition at line 127 of file HCalSD.h.

Referenced by endEvent(), initEvent(), and setDetUnitId().

◆ dist_

TH1F * HCalSD::dist_[9]
private

Definition at line 126 of file HCalSD.h.

Referenced by HCalSD(), and plotProfile().

◆ eminHitHB

double HCalSD::eminHitHB
private

Definition at line 113 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ eminHitHE

double HCalSD::eminHitHE
private

Definition at line 113 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ eminHitHF

double HCalSD::eminHitHF
private

Definition at line 113 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ eminHitHO

double HCalSD::eminHitHO
private

Definition at line 113 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ fibre1LV

std::vector<const G4LogicalVolume*> HCalSD::fibre1LV
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and isItStraightBundle().

◆ fibre2LV

std::vector<const G4LogicalVolume*> HCalSD::fibre2LV
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and isItConicalBundle().

◆ fibreLV

std::vector<const G4LogicalVolume*> HCalSD::fibreLV
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and isItFibre().

◆ fibreNames

std::vector<std::string> HCalSD::fibreNames
private

Definition at line 121 of file HCalSD.h.

Referenced by HCalSD(), and isItFibre().

◆ gpar

std::vector<double> HCalSD::gpar
private

Definition at line 118 of file HCalSD.h.

Referenced by HCalSD(), and plotHF().

◆ hcalConstants_

const HcalDDDSimConstants* HCalSD::hcalConstants_
private

Definition at line 96 of file HCalSD.h.

Referenced by getEnergyDeposit(), HCalSD(), modifyDepth(), and setDetUnitId().

◆ hcalSimConstants_

const HcalSimulationConstants* HCalSD::hcalSimConstants_
private

Definition at line 97 of file HCalSD.h.

Referenced by HCalSD().

◆ hfLevels

std::vector<int> HCalSD::hfLevels
private

Definition at line 119 of file HCalSD.h.

Referenced by HCalSD(), and isItHF().

◆ hfLV

std::vector<const G4LogicalVolume*> HCalSD::hfLV
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and isItHF().

◆ hfNames

std::vector<std::string> HCalSD::hfNames
private

Definition at line 120 of file HCalSD.h.

Referenced by HCalSD(), and isItHF().

◆ hfshower

std::unique_ptr<HFShower> HCalSD::hfshower
private

Definition at line 91 of file HCalSD.h.

Referenced by HCalSD(), and hitForFibre().

◆ hit_

TH1F* HCalSD::hit_[9]
private

Definition at line 126 of file HCalSD.h.

Referenced by HCalSD(), and plotProfile().

◆ hzvem

TH1F * HCalSD::hzvem
private

Definition at line 126 of file HCalSD.h.

Referenced by HCalSD(), and plotHF().

◆ hzvhad

TH1F * HCalSD::hzvhad
private

Definition at line 126 of file HCalSD.h.

Referenced by HCalSD(), and plotHF().

◆ isHF

bool HCalSD::isHF
private

Definition at line 107 of file HCalSD.h.

Referenced by getEnergyDeposit(), and getFromLibrary().

◆ layerWeights

std::map<uint32_t, double> HCalSD::layerWeights
private

Definition at line 125 of file HCalSD.h.

Referenced by layerWeight(), and readWeightFromFile().

◆ m_HBDarkening

const HBHEDarkening* HCalSD::m_HBDarkening
private

Definition at line 98 of file HCalSD.h.

Referenced by getEnergyDeposit().

◆ m_HcalTestNS

std::unique_ptr<HcalTestNS> HCalSD::m_HcalTestNS
private

Definition at line 101 of file HCalSD.h.

Referenced by HCalSD(), and setDetUnitId().

◆ m_HEDarkening

const HBHEDarkening* HCalSD::m_HEDarkening
private

Definition at line 99 of file HCalSD.h.

Referenced by getEnergyDeposit().

◆ m_HFDarkening

std::unique_ptr<HFDarkening> HCalSD::m_HFDarkening
private

Definition at line 100 of file HCalSD.h.

Referenced by getFromLibrary(), and HCalSD().

◆ materials

std::vector<const G4Material*> HCalSD::materials
private

Definition at line 123 of file HCalSD.h.

Referenced by HCalSD(), and isItScintillator().

◆ matNames

std::vector<std::string> HCalSD::matNames
private

Definition at line 122 of file HCalSD.h.

Referenced by HCalSD().

◆ maxRoff_

constexpr double HCalSD::maxRoff_ = 450.0
staticprivate

Definition at line 105 of file HCalSD.h.

Referenced by setDetUnitId().

◆ maxZ_

constexpr double HCalSD::maxZ_ = 10000.0
staticprivate

Definition at line 103 of file HCalSD.h.

Referenced by setDetUnitId().

◆ minRoff_

constexpr double HCalSD::minRoff_ = -1500.0
staticprivate

Definition at line 104 of file HCalSD.h.

Referenced by setDetUnitId().

◆ neutralDensity

bool HCalSD::neutralDensity
private

Definition at line 110 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ numberingFromDDD

std::unique_ptr<HcalNumberingFromDDD> HCalSD::numberingFromDDD
private

Definition at line 88 of file HCalSD.h.

Referenced by getHitFibreBundle(), getHitPMT(), HCalSD(), layerWeight(), and setDetUnitId().

◆ numberingScheme

std::unique_ptr<HcalNumberingScheme> HCalSD::numberingScheme
private

Definition at line 89 of file HCalSD.h.

Referenced by HCalSD(), setDetUnitId(), and setNumberingScheme().

◆ pmtLV

std::vector<const G4LogicalVolume*> HCalSD::pmtLV
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and isItPMT().

◆ showerBundle

std::unique_ptr<HFShowerFibreBundle> HCalSD::showerBundle
private

Definition at line 94 of file HCalSD.h.

Referenced by getEnergyDeposit(), getHitFibreBundle(), and HCalSD().

◆ showerLibrary

std::unique_ptr<HFShowerLibrary> HCalSD::showerLibrary
private

Definition at line 90 of file HCalSD.h.

Referenced by getFromHFLibrary(), and HCalSD().

◆ showerParam

std::unique_ptr<HFShowerParam> HCalSD::showerParam
private

Definition at line 92 of file HCalSD.h.

Referenced by getFromParam(), and HCalSD().

◆ showerPMT

std::unique_ptr<HFShowerPMT> HCalSD::showerPMT
private

Definition at line 93 of file HCalSD.h.

Referenced by getEnergyDeposit(), getHitPMT(), and HCalSD().

◆ slopeHE_

constexpr double HCalSD::slopeHE_ = 0.4
staticprivate

Definition at line 106 of file HCalSD.h.

Referenced by setDetUnitId().

◆ testNS_

bool HCalSD::testNS_
private

Definition at line 110 of file HCalSD.h.

Referenced by HCalSD().

◆ testNumber

bool HCalSD::testNumber
private

Definition at line 110 of file HCalSD.h.

Referenced by getEnergyDeposit(), HCalSD(), modifyDepth(), and setDetUnitId().

◆ time_

TH1F * HCalSD::time_[9]
private

Definition at line 126 of file HCalSD.h.

Referenced by HCalSD(), and plotProfile().

◆ useBirk

bool HCalSD::useBirk
private

Definition at line 109 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ useFibreBundle

bool HCalSD::useFibreBundle
private

Definition at line 109 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ useHF

bool HCalSD::useHF
private

Definition at line 112 of file HCalSD.h.

Referenced by HCalSD().

◆ useLayerWt

bool HCalSD::useLayerWt
private

Definition at line 109 of file HCalSD.h.

Referenced by getEnergyDeposit(), HCalSD(), and readWeightFromFile().

◆ useParam

bool HCalSD::useParam
private

Definition at line 112 of file HCalSD.h.

Referenced by getFromLibrary(), and HCalSD().

◆ usePMTHit

bool HCalSD::usePMTHit
private

Definition at line 109 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ useShowerLibrary

bool HCalSD::useShowerLibrary
private

Definition at line 112 of file HCalSD.h.

Referenced by getEnergyDeposit(), getFromLibrary(), and HCalSD().

◆ weight_

double HCalSD::weight_
private