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)
 
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 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)
 
- 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)
 
 ~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 ()
 
CaloG4HitcreateNewHit (const G4Step *, const G4Track *)
 
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 ()
 
double getResponseWt (const G4Track *)
 
virtual int getTrackID (const G4Track *)
 
bool hitExists (const G4Step *)
 
void ignoreRejection ()
 
void printDetectorLevels (const G4VTouchable *) const
 
void processHit (const G4Step *step)
 
void resetForNewPrimary (const G4Step *)
 
void setNumberCheckedHits (int val)
 
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 *)
 
- 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 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
 
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
CaloG4HitcurrentHit
 
CaloHitID currentID
 
float edepositEM
 
float edepositHAD
 
double eminHit
 
double energyCut
 
G4ThreeVector entranceLocal
 
G4ThreeVector entrancePoint
 
bool forceSave
 
float incidentEnergy
 
double kmaxIon
 
double kmaxNeutron
 
double kmaxProton
 
G4ThreeVector posGlobal
 
CaloHitID previousID
 
bool suppressHeavy
 
double tmaxHit
 

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 43 of file HCalSD.cc.

References agingFlagHB, agingFlagHE, applyFidCut, betaThr, birk1, birk2, birk3, 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, 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.

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

◆ ~HCalSD()

HCalSD::~HCalSD ( )
overridedefault

Member Function Documentation

◆ endEvent()

void HCalSD::endEvent ( )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 1091 of file HCalSD.cc.

References detNull_.

1091  {
1092 #ifdef printDebug
1093  int sum = detNull_[0] + detNull_[1] + detNull_[2];
1094  if (sum > 0)
1095  edm::LogVerbatim("HcalSim") << "NullDets " << detNull_[0] << " " << detNull_[1] << " " << detNull_[2] << " "
1096  << detNull_[3] << " " << (static_cast<float>(sum) / (sum + detNull_[3]));
1097 #endif
1098 }
Log< level::Info, true > LogVerbatim
std::vector< int > detNull_
Definition: HCalSD.h:125

◆ fillLogVolumeVector()

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

Definition at line 310 of file HCalSD.cc.

References mps_fire::i.

Referenced by HCalSD().

312  {
313  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
314  const G4LogicalVolume* lv;
315  std::stringstream ss3;
316  ss3 << "HCalSD: " << lvnames.size() << " names to be tested for Volume <" << value << ">:";
317  for (unsigned int i = 0; i < lvnames.size(); ++i) {
318  G4String namv(static_cast<std::string>(dd4hep::dd::noNamespace(lvnames[i])));
319  lv = nullptr;
320  for (auto lvol : *lvs) {
321  if (dd4hep::dd::noNamespace(lvol->GetName()) == namv) {
322  lv = lvol;
323  break;
324  }
325  }
326  lvvec.emplace_back(lv);
327  if (i / 10 * 10 == i) {
328  ss3 << "\n";
329  }
330  ss3 << " " << namv;
331  }
332  edm::LogVerbatim("HcalSim") << ss3.str();
333 }
Log< level::Info, true > LogVerbatim
Definition: value.py:1

◆ filterHit()

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

Reimplemented from CaloSD.

Definition at line 556 of file HCalSD.cc.

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

556  {
557  double threshold = 0;
558  DetId theId(aHit->getUnitID());
559  switch (theId.subdetId()) {
560  case HcalBarrel:
562  break;
563  case HcalEndcap:
565  break;
566  case HcalOuter:
568  break;
569  case HcalForward:
571  break;
572  default:
573  break;
574  }
575  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > threshold));
576 }
double eminHitHE
Definition: HCalSD.h:112
double eminHitHB
Definition: HCalSD.h:112
double eminHitHF
Definition: HCalSD.h:112
double eminHitHO
Definition: HCalSD.h:112
double tmaxHit
Definition: CaloSD.h:144
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 399 of file HCalSD.cc.

References agingFlagHB, agingFlagHE, birk1, birk2, birk3, HBHEDarkening::degradation(), deliveredLumi, LEDCalibrationChannels::depth, depth_, CaloSD::getAttenuation(), getHitFibreBundle(), getHitPMT(), HcalDDDSimConstants::getLayer0Wt(), CaloSD::getResponseWt(), hcalConstants_, hitForFibre(), LEDCalibrationChannels::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().

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

◆ getFromHFLibrary()

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

Definition at line 718 of file HCalSD.cc.

References CaloSD::currentID, LEDCalibrationChannels::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, protons_cff::time, and weight_.

Referenced by getFromLibrary().

718  {
719  std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, isKilled, weight_, false);
720  if (!isKilled || hits.empty()) {
721  return;
722  }
723 
724  int primaryID = setTrackID(aStep);
725 
726  // Reset entry point for new primary
727  resetForNewPrimary(aStep);
728 
729  auto const theTrack = aStep->GetTrack();
730  int det = 5;
731 
733  edepositEM = 1. * GeV;
734  edepositHAD = 0.;
735  } else {
736  edepositEM = 0.;
737  edepositHAD = 1. * GeV;
738  }
739 #ifdef EDM_ML_DEBUG
740  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary " << hits.size() << " hits for " << GetName() << " of "
741  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
742  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV";
743 #endif
744  for (unsigned int i = 0; i < hits.size(); ++i) {
745  G4ThreeVector hitPoint = hits[i].position;
746  if (isItinFidVolume(hitPoint)) {
747  int depth = hits[i].depth;
748  double time = hits[i].time;
749  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
750  currentID.setID(unitID, time, primaryID, 0);
751 #ifdef plotDebug
752  plotProfile(aStep, hitPoint, 1.0 * GeV, time, depth);
753  bool emType = G4TrackToParticleID::isGammaElectronPositron(theTrack->GetDefinition()->GetPDGEncoding());
754  plotHF(hitPoint, emType);
755 #endif
756  processHit(aStep);
757  }
758  }
759 }
float edepositEM
Definition: CaloSD.h:140
Log< level::Info, true > LogVerbatim
double weight_
Definition: HCalSD.h:114
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:998
void processHit(const G4Step *step)
Definition: CaloSD.h:113
float edepositHAD
Definition: CaloSD.h:140
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:653
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1060
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:533
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
CaloHitID currentID
Definition: CaloSD.h:142
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:702
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:89
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:826
static bool isGammaElectronPositron(int pdgCode)

◆ getFromLibrary()

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

Reimplemented from CaloSD.

Definition at line 335 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_2023v12_cff::track, HFDarkening::upperZLimit, useParam, useShowerLibrary, weight_, and z.

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

◆ getFromParam()

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

Definition at line 803 of file HCalSD.cc.

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

Referenced by getFromLibrary().

803  {
804  std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep, weight_, isKilled);
805  if (!isKilled || hits.empty()) {
806  return;
807  }
808 
809  int primaryID = setTrackID(aStep);
810  int det = 5;
811 
812 #ifdef EDM_ML_DEBUG
813  edm::LogVerbatim("HcalSim") << "HCalSD::getFromParam " << hits.size() << " hits for " << GetName() << " of "
814  << primaryID << " with " << aStep->GetTrack()->GetDefinition()->GetParticleName()
815  << " of " << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV
816  << " GeV in detector type " << det;
817 #endif
818  for (unsigned int i = 0; i < hits.size(); ++i) {
819  G4ThreeVector hitPoint = hits[i].position;
820  int depth = hits[i].depth;
821  double time = hits[i].time;
822  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
823  currentID.setID(unitID, time, primaryID, 0);
824  edepositEM = hits[i].edep * GeV;
825  edepositHAD = 0.;
826 #ifdef plotDebug
827  plotProfile(aStep, hitPoint, edepositEM, time, depth);
828 #endif
829  processHit(aStep);
830  }
831 }
float edepositEM
Definition: CaloSD.h:140
Log< level::Info, true > LogVerbatim
std::unique_ptr< HFShowerParam > showerParam
Definition: HCalSD.h:91
double weight_
Definition: HCalSD.h:114
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:998
void processHit(const G4Step *step)
Definition: CaloSD.h:113
float edepositHAD
Definition: CaloSD.h:140
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:533
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
CaloHitID currentID
Definition: CaloSD.h:142
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:826

◆ getHitFibreBundle()

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

Definition at line 892 of file HCalSD.cc.

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

Referenced by getEnergyDeposit().

892  {
893  auto const preStepPoint = aStep->GetPreStepPoint();
894  auto const theTrack = aStep->GetTrack();
895  double edep = showerBundle->getHits(aStep, type);
896 
897  if (edep >= 0.0) {
898  edep *= GeV;
899  double etrack = preStepPoint->GetKineticEnergy();
900  int primaryID = 0;
901  if (etrack >= energyCut) {
902  primaryID = theTrack->GetTrackID();
903  } else {
904  primaryID = theTrack->GetParentID();
905  if (primaryID == 0)
906  primaryID = theTrack->GetTrackID();
907  }
908  // Reset entry point for new primary
909  resetForNewPrimary(aStep);
910  //
911  int det = static_cast<int>(HcalForward);
912  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
913  double rr = hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y();
914  double phi = rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x());
915  double etaR = showerBundle->getRadius();
916  int depth = 1;
917  if (etaR < 0.) {
918  depth = 2;
919  etaR = -etaR;
920  }
921  if (hitPoint.z() < 0.)
922  etaR = -etaR;
923 #ifdef EDM_ML_DEBUG
924  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
925  << " depth " << depth;
926 #endif
927  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
928  uint32_t unitID = 0;
929  if (numberingFromDDD) {
930  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
931  unitID = setDetUnitId(tmp);
932  }
933  if (type)
934  currentID.setID(unitID, time, primaryID, 3);
935  else
936  currentID.setID(unitID, time, primaryID, 2);
937 
938  edepositHAD = aStep->GetTotalEnergyDeposit();
939  edepositEM = -edepositHAD + edep;
940 #ifdef plotDebug
941  plotProfile(aStep, hitPoint, edep, time, depth);
942 #endif
943 #ifdef EDM_ML_DEBUG
944  double beta = preStepPoint->GetBeta();
945  edm::LogVerbatim("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName() << " of " << primaryID
946  << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
947  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
948  << std::hex << unitID << std::dec;
949 #endif
950  processHit(aStep);
951  } // non-zero energy deposit
952 }
float edepositEM
Definition: CaloSD.h:140
double energyCut
Definition: CaloSD.h:144
Log< level::Info, true > LogVerbatim
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:998
void processHit(const G4Step *step)
Definition: CaloSD.h:113
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:87
float edepositHAD
Definition: CaloSD.h:140
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:653
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:533
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
CaloHitID currentID
Definition: CaloSD.h:142
std::unique_ptr< HFShowerFibreBundle > showerBundle
Definition: HCalSD.h:93
tmp
align.sh
Definition: createJobs.py:716

◆ getHitPMT()

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

Definition at line 833 of file HCalSD.cc.

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

Referenced by getEnergyDeposit().

833  {
834  auto const preStepPoint = aStep->GetPreStepPoint();
835  auto const theTrack = aStep->GetTrack();
836  double edep = showerPMT->getHits(aStep);
837 
838  if (edep >= 0.) {
839  edep *= GeV;
840  double etrack = preStepPoint->GetKineticEnergy();
841  int primaryID = 0;
842  if (etrack >= energyCut) {
843  primaryID = theTrack->GetTrackID();
844  } else {
845  primaryID = theTrack->GetParentID();
846  if (primaryID == 0)
847  primaryID = theTrack->GetTrackID();
848  }
849  // Reset entry point for new primary
850  resetForNewPrimary(aStep);
851  //
852  int det = static_cast<int>(HcalForward);
853  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
854  double rr = (hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y());
855  double phi = (rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x()));
856  double etaR = showerPMT->getRadius();
857  int depth = 1;
858  if (etaR < 0) {
859  depth = 2;
860  etaR = -etaR;
861  }
862  if (hitPoint.z() < 0)
863  etaR = -etaR;
864 #ifdef EDM_ML_DEBUG
865  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / deg
866  << " depth " << depth;
867 #endif
868  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
869  uint32_t unitID = 0;
870  if (numberingFromDDD) {
871  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
872  unitID = setDetUnitId(tmp);
873  }
874  currentID.setID(unitID, time, primaryID, 1);
875 
876  edepositHAD = aStep->GetTotalEnergyDeposit();
877  edepositEM = -edepositHAD + edep;
878 #ifdef plotDebug
879  plotProfile(aStep, hitPoint, edep, time, depth);
880 #endif
881 #ifdef EDM_ML_DEBUG
882  double beta = preStepPoint->GetBeta();
883  edm::LogVerbatim("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() << " of " << primaryID << " with "
884  << theTrack->GetDefinition()->GetParticleName() << " of "
885  << preStepPoint->GetKineticEnergy() / GeV << " GeV with velocity " << beta << " UnitID "
886  << std::hex << unitID << std::dec;
887 #endif
888  processHit(aStep);
889  }
890 }
float edepositEM
Definition: CaloSD.h:140
double energyCut
Definition: CaloSD.h:144
Log< level::Info, true > LogVerbatim
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:998
void processHit(const G4Step *step)
Definition: CaloSD.h:113
std::unique_ptr< HFShowerPMT > showerPMT
Definition: HCalSD.h:92
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:87
float edepositHAD
Definition: CaloSD.h:140
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:653
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:533
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
CaloHitID currentID
Definition: CaloSD.h:142
tmp
align.sh
Definition: createJobs.py:716

◆ hitForFibre()

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

Definition at line 761 of file HCalSD.cc.

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

Referenced by getEnergyDeposit().

761  { // if not ParamShower
762 
763  std::vector<HFShower::Hit> hits = hfshower->getHits(aStep, weight_);
764  if (hits.empty()) {
765  return;
766  }
767 
768  auto const theTrack = aStep->GetTrack();
769  int primaryID = setTrackID(aStep);
770  int det = 5;
771 
773  edepositEM = 1. * GeV;
774  edepositHAD = 0.;
775  } else {
776  edepositEM = 0.;
777  edepositHAD = 1. * GeV;
778  }
779 
780 #ifdef EDM_ML_DEBUG
781  edm::LogVerbatim("HcalSim") << "HCalSD::hitForFibre " << hits.size() << " hits for " << GetName() << " of "
782  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
783  << aStep->GetPreStepPoint()->GetKineticEnergy() / GeV << " GeV in detector type " << det;
784 #endif
785 
786  for (unsigned int i = 0; i < hits.size(); ++i) {
787  G4ThreeVector hitPoint = hits[i].position;
788  if (isItinFidVolume(hitPoint)) {
789  int depth = hits[i].depth;
790  double time = hits[i].time;
791  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
792  currentID.setID(unitID, time, primaryID, 0);
793 #ifdef plotDebug
794  plotProfile(aStep, hitPoint, edepositEM, time, depth);
795  bool emType = (edepositEM > 0.) ? true : false;
796  plotHF(hitPoint, emType);
797 #endif
798  processHit(aStep);
799  }
800  }
801 }
float edepositEM
Definition: CaloSD.h:140
Log< level::Info, true > LogVerbatim
double weight_
Definition: HCalSD.h:114
void plotProfile(const G4Step *step, const G4ThreeVector &pos, double edep, double time, int id)
Definition: HCalSD.cc:998
void processHit(const G4Step *step)
Definition: CaloSD.h:113
float edepositHAD
Definition: CaloSD.h:140
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1060
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:533
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
CaloHitID currentID
Definition: CaloSD.h:142
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:702
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:826
static bool isGammaElectronPositron(int pdgCode)
std::unique_ptr< HFShower > hfshower
Definition: HCalSD.h:90

◆ initEvent()

void HCalSD::initEvent ( const BeginOfEvent )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 1085 of file HCalSD.cc.

References detNull_.

1085  {
1086 #ifdef printDebug
1087  detNull_ = {0, 0, 0, 0};
1088 #endif
1089 }
std::vector< int > detNull_
Definition: HCalSD.h:125

◆ initRun()

void HCalSD::initRun ( )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 554 of file HCalSD.cc.

554 {}

◆ isItConicalBundle()

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

Definition at line 686 of file HCalSD.cc.

References fibre2LV.

Referenced by getEnergyDeposit().

686  {
687  for (auto lvol : fibre2LV)
688  if (lv == lvol) {
689  return true;
690  }
691  return false;
692 }
std::vector< const G4LogicalVolume * > fibre2LV
Definition: HCalSD.h:122

◆ isItFibre() [1/2]

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

Definition at line 654 of file HCalSD.cc.

References fibreLV.

Referenced by getEnergyDeposit().

654  {
655  for (auto lvol : fibreLV)
656  if (lv == lvol) {
657  return true;
658  }
659  return false;
660 }
std::vector< const G4LogicalVolume * > fibreLV
Definition: HCalSD.h:122

◆ isItFibre() [2/2]

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

Definition at line 662 of file HCalSD.cc.

References fibreNames, and Skims_PA_cff::name.

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

◆ isItHF() [1/2]

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

Definition at line 633 of file HCalSD.cc.

References hfLevels, hfLV, hfNames, and GCP_tree_cfg::levels.

Referenced by getFromLibrary().

633  {
634  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
635  int levels = (touch->GetHistoryDepth()) + 1;
636  for (unsigned int it = 0; it < hfNames.size(); ++it) {
637  if (levels >= hfLevels[it]) {
638  const G4LogicalVolume* lv = touch->GetVolume(levels - hfLevels[it])->GetLogicalVolume();
639  if (lv == hfLV[it])
640  return true;
641  }
642  }
643  return false;
644 }
std::vector< const G4LogicalVolume * > hfLV
Definition: HCalSD.h:122
std::vector< std::string > hfNames
Definition: HCalSD.h:118
std::vector< int > hfLevels
Definition: HCalSD.h:117

◆ isItHF() [2/2]

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

Definition at line 646 of file HCalSD.cc.

References hfNames, and Skims_PA_cff::name.

646  {
647  for (const auto& nam : hfNames)
648  if (name == static_cast<G4String>(nam)) {
649  return true;
650  }
651  return false;
652 }
std::vector< std::string > hfNames
Definition: HCalSD.h:118

◆ isItinFidVolume()

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

Definition at line 702 of file HCalSD.cc.

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

Referenced by getFromHFLibrary(), and hitForFibre().

702  {
703  bool flag = true;
704  if (applyFidCut) {
705  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
706 #ifdef EDM_ML_DEBUG
707  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume:#PMT= " << npmt << " for hit point " << hitPoint;
708 #endif
709  if (npmt <= 0)
710  flag = false;
711  }
712 #ifdef EDM_ML_DEBUG
713  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume: point " << hitPoint << " return flag " << flag;
714 #endif
715  return flag;
716 }
Log< level::Info, true > LogVerbatim
bool applyFidCut
Definition: HCalSD.h:111
int PMTNumber(const G4ThreeVector &pe_effect)

◆ isItPMT()

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

Definition at line 670 of file HCalSD.cc.

References pmtLV.

Referenced by getEnergyDeposit().

670  {
671  for (auto lvol : pmtLV)
672  if (lv == lvol) {
673  return true;
674  }
675  return false;
676 }
std::vector< const G4LogicalVolume * > pmtLV
Definition: HCalSD.h:122

◆ isItScintillator()

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

Definition at line 694 of file HCalSD.cc.

References materials.

Referenced by getEnergyDeposit().

694  {
695  for (auto amat : materials)
696  if (amat == mat) {
697  return true;
698  }
699  return false;
700 }
std::vector< const G4Material * > materials
Definition: HCalSD.h:121

◆ isItStraightBundle()

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

Definition at line 678 of file HCalSD.cc.

References fibre1LV.

Referenced by getEnergyDeposit().

678  {
679  for (auto lvol : fibre1LV)
680  if (lv == lvol) {
681  return true;
682  }
683  return false;
684 }
std::vector< const G4LogicalVolume * > fibre1LV
Definition: HCalSD.h:122

◆ layerWeight()

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

Definition at line 978 of file HCalSD.cc.

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

Referenced by getEnergyDeposit().

978  {
979  double wt = 1.;
980  if (numberingFromDDD) {
981  //get the ID's as eta, phi, depth, ... indices
983  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
984  modifyDepth(tmp);
985  uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1, tmp.etaR, tmp.phis, tmp.lay);
986  std::map<uint32_t, double>::const_iterator ite = layerWeights.find(id);
987  if (ite != layerWeights.end())
988  wt = ite->second;
989 #ifdef EDM_ML_DEBUG
990  edm::LogVerbatim("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id << std::dec << " (" << tmp.subdet << "/"
991  << tmp.zside << "/1/" << tmp.etaR << "/" << tmp.phis << "/" << tmp.lay << ") Weight "
992  << wt;
993 #endif
994  }
995  return wt;
996 }
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1071
static uint32_t packHcalIndex(int det, int z, int depth, int eta, int phi, int lay)
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:87
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:123
tmp
align.sh
Definition: createJobs.py:716

◆ modifyDepth()

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

Definition at line 1071 of file HCalSD.cc.

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

Referenced by layerWeight(), and setDetUnitId().

1071  {
1072  if (id.subdet == 4) {
1073  int ieta = (id.zside == 0) ? -id.etaR : id.etaR;
1074  if (hcalConstants_->maxHFDepth(ieta, id.phis) > 2) {
1075  if (id.depth <= 2) {
1076  if (G4UniformRand() > 0.5)
1077  id.depth += 2;
1078  }
1079  }
1080  } else if ((id.subdet == 1 || id.subdet == 2) && testNumber) {
1081  id.depth = (depth_ == 0) ? 1 : 2;
1082  }
1083 }
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:95
int maxHFDepth(const int &ieta, const int &iphi) const
bool testNumber
Definition: HCalSD.h:109
int depth_
Definition: HCalSD.h:115

◆ plotHF()

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

Definition at line 1060 of file HCalSD.cc.

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

Referenced by getFromHFLibrary(), and hitForFibre().

1060  {
1061  double zv = std::abs(hitPoint.z()) - gpar[4];
1062  if (emType) {
1063  if (hzvem != nullptr)
1064  hzvem->Fill(zv);
1065  } else {
1066  if (hzvhad != nullptr)
1067  hzvhad->Fill(zv);
1068  }
1069 }
std::vector< double > gpar
Definition: HCalSD.h:116
float *__restrict__ zv
TH1F * hzvhad
Definition: HCalSD.h:124
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TH1F * hzvem
Definition: HCalSD.h:124

◆ plotProfile()

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

Definition at line 998 of file HCalSD.cc.

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

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

998  {
999  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
1000  static const unsigned int names = 10;
1001  static const G4String modName[names] = {
1002  "HEModule", "HVQF", "HBModule", "MBAT", "MBBT", "MBBTC", "MBBT_R1P", "MBBT_R1M", "MBBT_R1PX", "MBBT_R1MX"};
1003  G4ThreeVector local;
1004  bool found = false;
1005  double depth = -2000;
1006  int idx = 4;
1007  for (int n = 0; n < touch->GetHistoryDepth(); ++n) {
1008  G4String name(static_cast<std::string>(dd4hep::dd::noNamespace(touch->GetVolume(n)->GetName())));
1009 #ifdef EDM_ML_DEBUG
1010  edm::LogVerbatim("HcalSim") << "plotProfile Depth " << n << " Name " << name;
1011 #endif
1012  for (unsigned int ii = 0; ii < names; ++ii) {
1013  if (name == modName[ii]) {
1014  found = true;
1015  int dn = touch->GetHistoryDepth() - n;
1016  local = touch->GetHistory()->GetTransform(dn).TransformPoint(global);
1017  if (ii == 0) {
1018  depth = local.z() - 4006.5;
1019  idx = 1;
1020  } else if (ii == 1) {
1021  depth = local.z() + 825.0;
1022  idx = 3;
1023  } else if (ii == 2) {
1024  depth = local.x() - 1775.;
1025  idx = 0;
1026  } else {
1027  depth = local.y() + 15.;
1028  idx = 2;
1029  }
1030  break;
1031  }
1032  }
1033  if (found)
1034  break;
1035  }
1036  if (!found)
1037  depth = std::abs(global.z()) - 11500;
1038 #ifdef EDM_ML_DEBUG
1039  edm::LogVerbatim("HcalSim") << "plotProfile Found " << found << " Global " << global << " Local " << local
1040  << " depth " << depth << " ID " << id << " EDEP " << edep << " Time " << time;
1041 #endif
1042  if (hit_[idx] != nullptr)
1043  hit_[idx]->Fill(edep);
1044  if (time_[idx] != nullptr)
1045  time_[idx]->Fill(time, edep);
1046  if (dist_[idx] != nullptr)
1047  dist_[idx]->Fill(depth, edep);
1048  int jd = 2 * idx + id - 7;
1049  if (jd >= 0 && jd < 4) {
1050  jd += 5;
1051  if (hit_[jd] != nullptr)
1052  hit_[jd]->Fill(edep);
1053  if (time_[jd] != nullptr)
1054  time_[jd]->Fill(time, edep);
1055  if (dist_[jd] != nullptr)
1056  dist_[jd]->Fill(depth, edep);
1057  }
1058 }
Log< level::Info, true > LogVerbatim
TH1F * time_[9]
Definition: HCalSD.h:124
const std::string names[nVars_]
TH1F * hit_[9]
Definition: HCalSD.h:124
TH1F * dist_[9]
Definition: HCalSD.h:124
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
ii
Definition: cuy.py:589

◆ readWeightFromFile()

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

Definition at line 954 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().

954  {
955  std::ifstream infile;
956  int entry = 0;
957  infile.open(fName.c_str(), std::ios::in);
958  if (infile) {
959  int det, zside, etaR, phi, lay;
960  double wt;
961  while (infile >> det >> zside >> etaR >> phi >> lay >> wt) {
962  uint32_t id = HcalTestNumbering::packHcalIndex(det, zside, 1, etaR, phi, lay);
963  layerWeights.insert(std::pair<uint32_t, double>(id, wt));
964  ++entry;
965 #ifdef EDM_ML_DEBUG
966  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry << " ID " << std::hex << id
967  << std::dec << " (" << det << "/" << zside << "/1/" << etaR << "/" << phi << "/"
968  << lay << ") Weight " << wt;
969 #endif
970  }
971  infile.close();
972  }
973  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry << " weights from " << fName;
974  if (entry <= 0)
975  useLayerWt = false;
976 }
Log< level::Info, true > LogVerbatim
bool useLayerWt
Definition: HCalSD.h:108
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:123

◆ setDetUnitId() [1/3]

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

Implements CaloSD.

Definition at line 533 of file HCalSD.cc.

References LEDCalibrationChannels::depth.

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

533  {
534  auto const prePoint = aStep->GetPreStepPoint();
535  auto const touch = prePoint->GetTouchable();
536  const G4ThreeVector& hitPoint = prePoint->GetPosition();
537 
538  int depth = (touch->GetReplicaNumber(0)) % 10 + 1;
539  int lay = (touch->GetReplicaNumber(0) / 10) % 100 + 1;
540  int det = (touch->GetReplicaNumber(1)) / 1000;
541 
542  return setDetUnitId(det, hitPoint, depth, lay);
543 }
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:533

◆ setDetUnitId() [2/3]

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

Definition at line 578 of file HCalSD.cc.

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

578  {
579  uint32_t id = 0;
580  if (det == 0) {
581 #ifdef printDebug
582  double eta = std::abs(pos.eta());
583 #endif
584  if (std::abs(pos.z()) > maxZ_) {
585  det = 5;
586 #ifdef printDebug
587  if (eta < 2.868)
588  ++detNull_[2];
589 #endif
590  } else if (!(hcalConstants_->isHE())) {
591  det = 3;
592 #ifdef printDebug
593  ++detNull_[0];
594 #endif
595  } else {
596  double minR = minRoff_ + slopeHE_ * std::abs(pos.z());
597  double maxR = maxRoff_ + slopeHE_ * std::abs(pos.z());
598  det = ((pos.perp() > minR) && (pos.perp() < maxR)) ? 4 : 3;
599 #ifdef printDebug
600  ++detNull_[det - 3];
601 #endif
602  }
603 #ifdef printDEBUG
604  edm::LogVerbatim("HcalSim") << "Position " << pos.perp() << ":" << std::abs(pos.z()) << " Limits "
605  << !(hcalConstants_->isHE()) << ":" << maxZ_ << " det " << det;
606  } else {
607  ++detNull_[3];
608 #endif
609  }
610 
611  if (numberingFromDDD.get()) {
612  //get the ID's as eta, phi, depth, ... indices
614  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
615  id = setDetUnitId(tmp);
616  }
617  return id;
618 }
Log< level::Info, true > LogVerbatim
static constexpr double slopeHE_
Definition: HCalSD.h:105
const HcalDDDSimConstants * hcalConstants_
Definition: HCalSD.h:95
static constexpr double minRoff_
Definition: HCalSD.h:103
std::unique_ptr< HcalNumberingFromDDD > numberingFromDDD
Definition: HCalSD.h:87
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:533
static constexpr double maxZ_
Definition: HCalSD.h:102
std::vector< int > detNull_
Definition: HCalSD.h:125
tmp
align.sh
Definition: createJobs.py:716
static constexpr double maxRoff_
Definition: HCalSD.h:104

◆ setDetUnitId() [3/3]

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

Definition at line 620 of file HCalSD.cc.

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

620  {
621  modifyDepth(tmp);
622  uint32_t id = (numberingScheme.get()) ? numberingScheme->getUnitID(tmp) : 0;
623  if ((!testNumber) && m_HcalTestNS.get()) {
624  bool ok = m_HcalTestNS->compare(tmp, id);
625  if (!ok)
626  edm::LogWarning("HcalSim") << "Det ID from HCalSD " << HcalDetId(id) << " " << std::hex << id << std::dec
627  << " does not match one from relabller for " << tmp.subdet << ":" << tmp.etaR << ":"
628  << tmp.phi << ":" << tmp.phis << ":" << tmp.depth << ":" << tmp.lay << std::endl;
629  }
630  return id;
631 }
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1071
bool testNumber
Definition: HCalSD.h:109
std::unique_ptr< HcalTestNS > m_HcalTestNS
Definition: HCalSD.h:100
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:88
Log< level::Warning, false > LogWarning
tmp
align.sh
Definition: createJobs.py:716

◆ setNumberingScheme()

void HCalSD::setNumberingScheme ( HcalNumberingScheme scheme)

Definition at line 545 of file HCalSD.cc.

References numberingScheme, and generator_cfi::scheme.

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

545  {
546  if (scheme != nullptr) {
547  edm::LogVerbatim("HcalSim") << "HCalSD: updates numbering scheme for " << GetName();
548  numberingScheme.reset(scheme);
549  }
550 }
Log< level::Info, true > LogVerbatim
std::unique_ptr< HcalNumberingScheme > numberingScheme
Definition: HCalSD.h:88

◆ update() [1/6]

void CaloSD::update
overrideprotected

Definition at line 688 of file CaloSD.cc.

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

688 { initRun(); }
void initRun() override
Definition: HCalSD.cc:554

◆ update() [2/6]

void CaloSD::update
overrideprotected

Definition at line 945 of file CaloSD.cc.

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

945  {
946  int primary = -1;
947  TrackInformation* trkInfo = cmsTrackInformation((*trk)());
948  if (trkInfo->isPrimary())
949  primary = (*trk)()->GetTrackID();
950 
951 #ifdef EDM_ML_DEBUG
952  edm::LogVerbatim("CaloSim") << "New track: isPrimary " << trkInfo->isPrimary() << " primary ID = " << primary
953  << " primary ancestor ID " << primAncestor;
954 #endif
955 
956  // update the information if a different primary track ID
957 
958  if (primary > 0 && primary != primAncestor) {
959  primAncestor = primary;
960 
961  // clean the hits information
962 
963  if (theHC->entries() > 0)
965  }
966 }
Log< level::Info, true > LogVerbatim
bool isPrimary() const
int primAncestor
Definition: CaloSD.h:175
void cleanHitCollection()
Definition: CaloSD.cc:968
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
CaloG4HitCollection * theHC
Definition: CaloSD.h:166

◆ update() [3/6]

void CaloSD::update
overrideprotected

Definition at line 698 of file CaloSD.cc.

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

698  {
699  int id = (*trk)()->GetTrackID();
700  TrackInformation* trkI = cmsTrackInformation((*trk)());
701  int lastTrackID = -1;
702  if (trkI)
703  lastTrackID = trkI->getIDonCaloSurface();
704  if (id == lastTrackID) {
705  const TrackContainer* trksForThisEvent = m_trackManager->trackContainer();
706  if (trksForThisEvent != nullptr) {
707  int it = (int)(trksForThisEvent->size()) - 1;
708  if (it >= 0) {
709  TrackWithHistory* trkH = (*trksForThisEvent)[it];
710  if (trkH->trackID() == (unsigned int)(id))
711  tkMap[id] = trkH;
712 #ifdef EDM_ML_DEBUG
713  edm::LogVerbatim("CaloSim") << "CaloSD: get track " << it << " from Container of size "
714  << trksForThisEvent->size() << " with ID " << trkH->trackID();
715  } else {
716  edm::LogVerbatim("CaloSim") << "CaloSD: get track " << it << " from Container of size "
717  << trksForThisEvent->size() << " with no ID";
718 #endif
719  }
720  }
721  }
722 }
Log< level::Info, true > LogVerbatim
std::vector< TrackWithHistory * > TrackContainer
Definition: TrackContainer.h:8
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
int getIDonCaloSurface() const
const TrackContainer * trackContainer() const
std::map< int, TrackWithHistory * > tkMap
Definition: CaloSD.h:188
const SimTrackManager * m_trackManager
Definition: CaloSD.h:161
unsigned int trackID() const

◆ update() [4/6]

void CaloSD::update
overrideprotected

Definition at line 724 of file CaloSD.cc.

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

724  {
725  endEvent();
726  slave.get()->ReserveMemory(theHC->entries());
727 
728  int count(0);
729  int wrong(0);
730  double eEM(0.0);
731  double eHAD(0.0);
732  double eEM2(0.0);
733  double eHAD2(0.0);
734 #ifdef EDM_ML_DEBUG
735  double tt(0.0);
736  double zloc(0.0);
737  double zglob(0.0);
738  double ee(0.0);
739 #endif
740  int hc_entries = theHC->entries();
741  for (int i = 0; i < hc_entries; ++i) {
742  if (!saveHit((*theHC)[i])) {
743  ++wrong;
744  }
745  ++count;
746  double x = (*theHC)[i]->getEM();
747  eEM += x;
748  eEM2 += x * x;
749  x = (*theHC)[i]->getHadr();
750  eHAD += x;
751  eHAD2 += x * x;
752 #ifdef EDM_ML_DEBUG
753  tt += (*theHC)[i]->getTimeSlice();
754  ee += (*theHC)[i]->getIncidentEnergy();
755  zglob += std::abs((*theHC)[i]->getEntry().z());
756  zloc += std::abs((*theHC)[i]->getEntryLocal().z());
757 #endif
758  }
759 
760  double norm = (count > 0) ? 1.0 / count : 0.0;
761  eEM *= norm;
762  eEM2 *= norm;
763  eHAD *= norm;
764  eHAD2 *= norm;
765  eEM2 = std::sqrt(eEM2 - eEM * eEM);
766  eHAD2 = std::sqrt(eHAD2 - eHAD * eHAD);
767 #ifdef EDM_ML_DEBUG
768  tt *= norm;
769  ee *= norm;
770  zglob *= norm;
771  zloc *= norm;
772  edm::LogVerbatim("CaloSim") << "CaloSD: " << GetName() << " store " << count << " hits; " << wrong
773  << " track IDs not given properly and " << totalHits - count
774  << " hits not passing cuts\n EmeanEM= " << eEM << " ErmsEM= " << eEM2
775  << "\n EmeanHAD= " << eHAD << " ErmsHAD= " << eHAD2 << " TimeMean= " << tt
776  << " E0mean= " << ee << " Zglob= " << zglob << " Zloc= " << zloc << " ";
777 #endif
778  tkMap.erase(tkMap.begin(), tkMap.end());
779  std::vector<std::unique_ptr<CaloG4Hit>>().swap(reusehit);
780  if (useMap)
781  hitMap.erase(hitMap.begin(), hitMap.end());
783 }
Int_t getEntry(TBranch *branch, EntryNumber entryNumber)
Definition: RootTree.cc:527
Log< level::Info, true > LogVerbatim
int totalHits
Definition: CaloSD.h:177
bool useMap
Definition: CaloSD.h:171
void endEvent() override
Definition: HCalSD.cc:1091
void swap(Association< C > &lhs, Association< C > &rhs)
Definition: Association.h:117
std::unordered_map< unsigned int, unsigned int > boundaryCrossingParentMap_
Definition: CaloSD.h:189
Definition: TTTypes.h:54
std::unique_ptr< CaloSlaveSD > slave
Definition: CaloSD.h:163
T sqrt(T t)
Definition: SSEVec.h:19
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
std::map< CaloHitID, CaloG4Hit * > hitMap
Definition: CaloSD.h:187
CaloG4HitCollection * theHC
Definition: CaloSD.h:166
std::map< int, TrackWithHistory * > tkMap
Definition: CaloSD.h:188
std::vector< std::unique_ptr< CaloG4Hit > > reusehit
Definition: CaloSD.h:190
bool saveHit(CaloG4Hit *)
Definition: CaloSD.cc:882

◆ update() [5/6]

void CaloSD::update
overrideprotected

Definition at line 690 of file CaloSD.cc.

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

690  {
691 #ifdef EDM_ML_DEBUG
692  edm::LogVerbatim("CaloSim") << "CaloSD: Dispatched BeginOfEvent for " << GetName() << " !";
693 #endif
694  clearHits();
695  initEvent(g4Event);
696 }
Log< level::Info, true > LogVerbatim
void initEvent(const BeginOfEvent *) override
Definition: HCalSD.cc:1085
void clearHits() override
Definition: CaloSD.cc:785

◆ 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 552 of file HCalSD.cc.

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

552 {}

Member Data Documentation

◆ agingFlagHB

bool HCalSD::agingFlagHB
private

Definition at line 107 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ agingFlagHE

bool HCalSD::agingFlagHE
private

Definition at line 107 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ applyFidCut

bool HCalSD::applyFidCut
private

Definition at line 111 of file HCalSD.h.

Referenced by HCalSD(), and isItinFidVolume().

◆ betaThr

double HCalSD::betaThr
private

Definition at line 110 of file HCalSD.h.

Referenced by HCalSD().

◆ birk1

double HCalSD::birk1
private

Definition at line 110 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ birk2

double HCalSD::birk2
private

Definition at line 110 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ birk3

double HCalSD::birk3
private

Definition at line 110 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ deliveredLumi

double HCalSD::deliveredLumi
private

Definition at line 113 of file HCalSD.h.

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

◆ depth_

int HCalSD::depth_
private

Definition at line 115 of file HCalSD.h.

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

◆ detNull_

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

Definition at line 125 of file HCalSD.h.

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

◆ dist_

TH1F * HCalSD::dist_[9]
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and plotProfile().

◆ eminHitHB

double HCalSD::eminHitHB
private

Definition at line 112 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ eminHitHE

double HCalSD::eminHitHE
private

Definition at line 112 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ eminHitHF

double HCalSD::eminHitHF
private

Definition at line 112 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ eminHitHO

double HCalSD::eminHitHO
private

Definition at line 112 of file HCalSD.h.

Referenced by filterHit(), and HCalSD().

◆ fibre1LV

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

Definition at line 122 of file HCalSD.h.

Referenced by HCalSD(), and isItStraightBundle().

◆ fibre2LV

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

Definition at line 122 of file HCalSD.h.

Referenced by HCalSD(), and isItConicalBundle().

◆ fibreLV

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

Definition at line 122 of file HCalSD.h.

Referenced by HCalSD(), and isItFibre().

◆ fibreNames

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

Definition at line 119 of file HCalSD.h.

Referenced by HCalSD(), and isItFibre().

◆ gpar

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

Definition at line 116 of file HCalSD.h.

Referenced by HCalSD(), and plotHF().

◆ hcalConstants_

const HcalDDDSimConstants* HCalSD::hcalConstants_
private

Definition at line 95 of file HCalSD.h.

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

◆ hcalSimConstants_

const HcalSimulationConstants* HCalSD::hcalSimConstants_
private

Definition at line 96 of file HCalSD.h.

Referenced by HCalSD().

◆ hfLevels

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

Definition at line 117 of file HCalSD.h.

Referenced by HCalSD(), and isItHF().

◆ hfLV

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

Definition at line 122 of file HCalSD.h.

Referenced by HCalSD(), and isItHF().

◆ hfNames

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

Definition at line 118 of file HCalSD.h.

Referenced by HCalSD(), and isItHF().

◆ hfshower

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

Definition at line 90 of file HCalSD.h.

Referenced by HCalSD(), and hitForFibre().

◆ hit_

TH1F* HCalSD::hit_[9]
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and plotProfile().

◆ hzvem

TH1F * HCalSD::hzvem
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and plotHF().

◆ hzvhad

TH1F * HCalSD::hzvhad
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and plotHF().

◆ isHF

bool HCalSD::isHF
private

Definition at line 106 of file HCalSD.h.

Referenced by getEnergyDeposit(), and getFromLibrary().

◆ layerWeights

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

Definition at line 123 of file HCalSD.h.

Referenced by layerWeight(), and readWeightFromFile().

◆ m_HBDarkening

const HBHEDarkening* HCalSD::m_HBDarkening
private

Definition at line 97 of file HCalSD.h.

Referenced by getEnergyDeposit().

◆ m_HcalTestNS

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

Definition at line 100 of file HCalSD.h.

Referenced by HCalSD(), and setDetUnitId().

◆ m_HEDarkening

const HBHEDarkening* HCalSD::m_HEDarkening
private

Definition at line 98 of file HCalSD.h.

Referenced by getEnergyDeposit().

◆ m_HFDarkening

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

Definition at line 99 of file HCalSD.h.

Referenced by getFromLibrary(), and HCalSD().

◆ materials

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

Definition at line 121 of file HCalSD.h.

Referenced by HCalSD(), and isItScintillator().

◆ matNames

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

Definition at line 120 of file HCalSD.h.

Referenced by HCalSD().

◆ maxRoff_

constexpr double HCalSD::maxRoff_ = 450.0
staticprivate

Definition at line 104 of file HCalSD.h.

Referenced by setDetUnitId().

◆ maxZ_

constexpr double HCalSD::maxZ_ = 10000.0
staticprivate

Definition at line 102 of file HCalSD.h.

Referenced by setDetUnitId().

◆ minRoff_

constexpr double HCalSD::minRoff_ = -1500.0
staticprivate

Definition at line 103 of file HCalSD.h.

Referenced by setDetUnitId().

◆ neutralDensity

bool HCalSD::neutralDensity
private

Definition at line 109 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ numberingFromDDD

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

Definition at line 87 of file HCalSD.h.

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

◆ numberingScheme

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

Definition at line 88 of file HCalSD.h.

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

◆ pmtLV

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

Definition at line 122 of file HCalSD.h.

Referenced by HCalSD(), and isItPMT().

◆ showerBundle

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

Definition at line 93 of file HCalSD.h.

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

◆ showerLibrary

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

Definition at line 89 of file HCalSD.h.

Referenced by getFromHFLibrary(), and HCalSD().

◆ showerParam

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

Definition at line 91 of file HCalSD.h.

Referenced by getFromParam(), and HCalSD().

◆ showerPMT

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

Definition at line 92 of file HCalSD.h.

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

◆ slopeHE_

constexpr double HCalSD::slopeHE_ = 0.4
staticprivate

Definition at line 105 of file HCalSD.h.

Referenced by setDetUnitId().

◆ testNS_

bool HCalSD::testNS_
private

Definition at line 109 of file HCalSD.h.

Referenced by HCalSD().

◆ testNumber

bool HCalSD::testNumber
private

Definition at line 109 of file HCalSD.h.

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

◆ time_

TH1F * HCalSD::time_[9]
private

Definition at line 124 of file HCalSD.h.

Referenced by HCalSD(), and plotProfile().

◆ useBirk

bool HCalSD::useBirk
private

Definition at line 108 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ useFibreBundle

bool HCalSD::useFibreBundle
private

Definition at line 108 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ useHF

bool HCalSD::useHF
private

Definition at line 111 of file HCalSD.h.

Referenced by HCalSD().

◆ useLayerWt

bool HCalSD::useLayerWt
private

Definition at line 108 of file HCalSD.h.

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

◆ useParam

bool HCalSD::useParam
private

Definition at line 111 of file HCalSD.h.

Referenced by getFromLibrary(), and HCalSD().

◆ usePMTHit

bool HCalSD::usePMTHit
private

Definition at line 108 of file HCalSD.h.

Referenced by getEnergyDeposit(), and HCalSD().

◆ useShowerLibrary

bool HCalSD::useShowerLibrary
private

Definition at line 111 of file HCalSD.h.

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

◆ weight_

double HCalSD::weight_
private