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 43 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, mergeVDriftHistosByStation::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.

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

References detNull_.

1111  {
1112 #ifdef printDebug
1113  int sum = detNull_[0] + detNull_[1] + detNull_[2];
1114  if (sum > 0)
1115  edm::LogVerbatim("HcalSim") << "NullDets " << detNull_[0] << " " << detNull_[1] << " " << detNull_[2] << " "
1116  << detNull_[3] << " " << (static_cast<float>(sum) / (sum + detNull_[3]));
1117 #endif
1118 }
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 305 of file HCalSD.cc.

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

Referenced by HCalSD().

307  {
308  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
309  const G4LogicalVolume* lv;
310  std::stringstream ss3;
311  ss3 << "HCalSD: " << lvnames.size() << " names to be tested for Volume <" << value << ">:";
312  for (unsigned int i = 0; i < lvnames.size(); ++i) {
314  lv = nullptr;
315  for (auto lvol : *lvs) {
316  if (DD4hep2DDDName::nameMatterLV(lvol->GetName(), dd4hep_) == namv) {
317  lv = lvol;
318  break;
319  }
320  }
321  lvvec.emplace_back(lv);
322  if (i / 10 * 10 == i) {
323  ss3 << "\n";
324  }
325  ss3 << " " << namv;
326  }
327  edm::LogVerbatim("HcalSim") << ss3.str();
328 }
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 575 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.

575  {
576  double threshold = 0;
577  DetId theId(aHit->getUnitID());
578  switch (theId.subdetId()) {
579  case HcalBarrel:
581  break;
582  case HcalEndcap:
584  break;
585  case HcalOuter:
587  break;
588  case HcalForward:
590  break;
591  default:
592  break;
593  }
594  return ((time <= tmaxHit) && (aHit->getEnergyDeposit() > threshold));
595 }
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 397 of file HCalSD.cc.

References agingFlagHB, agingFlagHE, birk1, birk2, birk3, HBHEDarkening::degradation(), deliveredLumi, hcalRecHitTable_cff::depth, depth_, ALCARECOPPSCalTrackBasedSel_cff::detid, 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().

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

◆ getFromHFLibrary()

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

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

737  {
738  std::vector<HFShowerLibrary::Hit> hits = showerLibrary->getHits(aStep, isKilled, weight_, false);
739  if (!isKilled || hits.empty()) {
740  return;
741  }
742 
743  int primaryID = setTrackID(aStep);
744 
745  // Reset entry point for new primary
746  resetForNewPrimary(aStep);
747 
748  auto const theTrack = aStep->GetTrack();
749  int det = 5;
750 
752  edepositEM = 1. * CLHEP::GeV;
753  edepositHAD = 0.;
754  } else {
755  edepositEM = 0.;
756  edepositHAD = 1. * CLHEP::GeV;
757  }
758 #ifdef EDM_ML_DEBUG
759  edm::LogVerbatim("HcalSim") << "HCalSD::getFromLibrary " << hits.size() << " hits for " << GetName() << " of "
760  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
761  << aStep->GetPreStepPoint()->GetKineticEnergy() / CLHEP::GeV << " GeV";
762 #endif
763  for (unsigned int i = 0; i < hits.size(); ++i) {
764  G4ThreeVector hitPoint = hits[i].position;
765  if (isItinFidVolume(hitPoint)) {
766  int depth = hits[i].depth;
767  double time = hits[i].time;
768  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
769  currentID[0].setID(unitID, time, primaryID, 0);
770 #ifdef plotDebug
771  plotProfile(aStep, hitPoint, 1.0 * CLHEP::GeV, time, depth);
772  bool emType = G4TrackToParticleID::isGammaElectronPositron(theTrack->GetDefinition()->GetPDGEncoding());
773  plotHF(hitPoint, emType);
774 #endif
775  processHit(aStep);
776  }
777  }
778 }
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:1018
void processHit(const G4Step *step)
Definition: CaloSD.h:117
float edepositHAD
Definition: CaloSD.h:144
void resetForNewPrimary(const G4Step *)
Definition: CaloSD.cc:715
void plotHF(const G4ThreeVector &pos, bool emType)
Definition: HCalSD.cc:1080
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:531
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:721
std::unique_ptr< HFShowerLibrary > showerLibrary
Definition: HCalSD.h:90
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:892
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 330 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_2024v14_cff::track, HFDarkening::upperZLimit, useParam, useShowerLibrary, weight_, and z.

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

◆ getFromParam()

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

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

823  {
824  std::vector<HFShowerParam::Hit> hits = showerParam->getHits(aStep, weight_, isKilled);
825  if (!isKilled || hits.empty()) {
826  return;
827  }
828 
829  int primaryID = setTrackID(aStep);
830  int det = 5;
831 
832 #ifdef EDM_ML_DEBUG
833  edm::LogVerbatim("HcalSim") << "HCalSD::getFromParam " << hits.size() << " hits for " << GetName() << " of "
834  << primaryID << " with " << aStep->GetTrack()->GetDefinition()->GetParticleName()
835  << " of " << aStep->GetPreStepPoint()->GetKineticEnergy() / CLHEP::GeV
836  << " GeV in detector type " << det;
837 #endif
838  for (unsigned int i = 0; i < hits.size(); ++i) {
839  G4ThreeVector hitPoint = hits[i].position;
840  int depth = hits[i].depth;
841  double time = hits[i].time;
842  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
843  currentID[0].setID(unitID, time, primaryID, 0);
844  edepositEM = hits[i].edep * CLHEP::GeV;
845  edepositHAD = 0.;
846 #ifdef plotDebug
847  plotProfile(aStep, hitPoint, edepositEM, time, depth);
848 #endif
849  processHit(aStep);
850  }
851 }
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:1018
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:531
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:892
CaloHitID currentID[2]
Definition: CaloSD.h:146

◆ getHitFibreBundle()

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

Definition at line 912 of file HCalSD.cc.

References Matriplex::atan2(), HLT_2024v14_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().

912  {
913  auto const preStepPoint = aStep->GetPreStepPoint();
914  auto const theTrack = aStep->GetTrack();
915  double edep = showerBundle->getHits(aStep, type);
916 
917  if (edep >= 0.0) {
918  edep *= CLHEP::GeV;
919  double etrack = preStepPoint->GetKineticEnergy();
920  int primaryID = 0;
921  if (etrack >= energyCut) {
922  primaryID = theTrack->GetTrackID();
923  } else {
924  primaryID = theTrack->GetParentID();
925  if (primaryID == 0)
926  primaryID = theTrack->GetTrackID();
927  }
928  // Reset entry point for new primary
929  resetForNewPrimary(aStep);
930  //
931  int det = static_cast<int>(HcalForward);
932  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
933  double rr = hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y();
934  double phi = rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x());
935  double etaR = showerBundle->getRadius();
936  int depth = 1;
937  if (etaR < 0.) {
938  depth = 2;
939  etaR = -etaR;
940  }
941  if (hitPoint.z() < 0.)
942  etaR = -etaR;
943 #ifdef EDM_ML_DEBUG
944  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / CLHEP::deg
945  << " depth " << depth;
946 #endif
947  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
948  uint32_t unitID = 0;
949  if (numberingFromDDD) {
950  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
951  unitID = setDetUnitId(tmp);
952  }
953  if (type)
954  currentID[0].setID(unitID, time, primaryID, 3);
955  else
956  currentID[0].setID(unitID, time, primaryID, 2);
957 
958  edepositHAD = aStep->GetTotalEnergyDeposit();
959  edepositEM = -edepositHAD + edep;
960 #ifdef plotDebug
961  plotProfile(aStep, hitPoint, edep, time, depth);
962 #endif
963 #ifdef EDM_ML_DEBUG
964  double beta = preStepPoint->GetBeta();
965  edm::LogVerbatim("HcalSim") << "HCalSD::getHitFibreBundle 1 hit for " << GetName() << " of " << primaryID
966  << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
967  << preStepPoint->GetKineticEnergy() / CLHEP::GeV << " GeV with velocity " << beta
968  << " UnitID " << std::hex << unitID << std::dec;
969 #endif
970  processHit(aStep);
971  } // non-zero energy deposit
972 }
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:1018
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:715
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:531
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
MPlex< T, D1, D2, N > atan2(const MPlex< T, D1, D2, N > &y, const MPlex< T, D1, D2, N > &x)
Definition: Matriplex.h:648

◆ getHitPMT()

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

Definition at line 853 of file HCalSD.cc.

References Matriplex::atan2(), HLT_2024v14_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().

853  {
854  auto const preStepPoint = aStep->GetPreStepPoint();
855  auto const theTrack = aStep->GetTrack();
856  double edep = showerPMT->getHits(aStep);
857 
858  if (edep >= 0.) {
859  edep *= CLHEP::GeV;
860  double etrack = preStepPoint->GetKineticEnergy();
861  int primaryID = 0;
862  if (etrack >= energyCut) {
863  primaryID = theTrack->GetTrackID();
864  } else {
865  primaryID = theTrack->GetParentID();
866  if (primaryID == 0)
867  primaryID = theTrack->GetTrackID();
868  }
869  // Reset entry point for new primary
870  resetForNewPrimary(aStep);
871  //
872  int det = static_cast<int>(HcalForward);
873  const G4ThreeVector& hitPoint = preStepPoint->GetPosition();
874  double rr = (hitPoint.x() * hitPoint.x() + hitPoint.y() * hitPoint.y());
875  double phi = (rr == 0. ? 0. : atan2(hitPoint.y(), hitPoint.x()));
876  double etaR = showerPMT->getRadius();
877  int depth = 1;
878  if (etaR < 0) {
879  depth = 2;
880  etaR = -etaR;
881  }
882  if (hitPoint.z() < 0)
883  etaR = -etaR;
884 #ifdef EDM_ML_DEBUG
885  edm::LogVerbatim("HcalSim") << "HCalSD::Hit for Detector " << det << " etaR " << etaR << " phi " << phi / CLHEP::deg
886  << " depth " << depth;
887 #endif
888  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
889  uint32_t unitID = 0;
890  if (numberingFromDDD) {
891  HcalNumberingFromDDD::HcalID tmp = numberingFromDDD->unitID(det, etaR, phi, depth, 1);
892  unitID = setDetUnitId(tmp);
893  }
894  currentID[0].setID(unitID, time, primaryID, 1);
895 
896  edepositHAD = aStep->GetTotalEnergyDeposit();
897  edepositEM = -edepositHAD + edep;
898 #ifdef plotDebug
899  plotProfile(aStep, hitPoint, edep, time, depth);
900 #endif
901 #ifdef EDM_ML_DEBUG
902  double beta = preStepPoint->GetBeta();
903  edm::LogVerbatim("HcalSim") << "HCalSD::getHitPMT 1 hit for " << GetName() << " of " << primaryID << " with "
904  << theTrack->GetDefinition()->GetParticleName() << " of "
905  << preStepPoint->GetKineticEnergy() / CLHEP::GeV << " GeV with velocity " << beta
906  << " UnitID " << std::hex << unitID << std::dec;
907 #endif
908  processHit(aStep);
909  }
910 }
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:1018
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:715
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:531
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
MPlex< T, D1, D2, N > atan2(const MPlex< T, D1, D2, N > &y, const MPlex< T, D1, D2, N > &x)
Definition: Matriplex.h:648

◆ hitForFibre()

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

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

780  { // if not ParamShower
781 
782  std::vector<HFShower::Hit> hits = hfshower->getHits(aStep, weight_);
783  if (hits.empty()) {
784  return;
785  }
786 
787  auto const theTrack = aStep->GetTrack();
788  int primaryID = setTrackID(aStep);
789  int det = 5;
790 
792  edepositEM = 1. * CLHEP::GeV;
793  edepositHAD = 0.;
794  } else {
795  edepositEM = 0.;
796  edepositHAD = 1. * CLHEP::GeV;
797  }
798 
799 #ifdef EDM_ML_DEBUG
800  edm::LogVerbatim("HcalSim") << "HCalSD::hitForFibre " << hits.size() << " hits for " << GetName() << " of "
801  << primaryID << " with " << theTrack->GetDefinition()->GetParticleName() << " of "
802  << aStep->GetPreStepPoint()->GetKineticEnergy() / CLHEP::GeV << " GeV in detector type "
803  << det;
804 #endif
805 
806  for (unsigned int i = 0; i < hits.size(); ++i) {
807  G4ThreeVector hitPoint = hits[i].position;
808  if (isItinFidVolume(hitPoint)) {
809  int depth = hits[i].depth;
810  double time = hits[i].time;
811  unsigned int unitID = setDetUnitId(det, hitPoint, depth);
812  currentID[0].setID(unitID, time, primaryID, 0);
813 #ifdef plotDebug
814  plotProfile(aStep, hitPoint, edepositEM, time, depth);
815  bool emType = (edepositEM > 0.) ? true : false;
816  plotHF(hitPoint, emType);
817 #endif
818  processHit(aStep);
819  }
820  }
821 }
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:1018
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:1080
uint32_t setDetUnitId(const G4Step *step) override
Definition: HCalSD.cc:531
void setID(uint32_t unitID, double timeSlice, int trackID, uint16_t depth=0)
Definition: CaloHitID.cc:41
bool isItinFidVolume(const G4ThreeVector &)
Definition: HCalSD.cc:721
virtual int setTrackID(const G4Step *)
Definition: CaloSD.cc:892
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 1105 of file HCalSD.cc.

References detNull_.

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

◆ initRun()

void HCalSD::initRun ( )
overrideprotectedvirtual

Reimplemented from CaloSD.

Definition at line 573 of file HCalSD.cc.

573 {}

◆ isItConicalBundle()

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

Definition at line 705 of file HCalSD.cc.

References fibre2LV.

Referenced by getEnergyDeposit().

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

◆ isItFibre() [1/2]

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

Definition at line 673 of file HCalSD.cc.

References fibreLV.

Referenced by getEnergyDeposit().

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

◆ isItFibre() [2/2]

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

Definition at line 681 of file HCalSD.cc.

References fibreNames, and mergeVDriftHistosByStation::name.

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

◆ isItHF() [1/2]

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

Definition at line 652 of file HCalSD.cc.

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

Referenced by getFromLibrary().

652  {
653  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
654  int levels = (touch->GetHistoryDepth()) + 1;
655  for (unsigned int it = 0; it < hfNames.size(); ++it) {
656  if (levels >= hfLevels[it]) {
657  const G4LogicalVolume* lv = touch->GetVolume(levels - hfLevels[it])->GetLogicalVolume();
658  if (lv == hfLV[it])
659  return true;
660  }
661  }
662  return false;
663 }
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 665 of file HCalSD.cc.

References hfNames, and mergeVDriftHistosByStation::name.

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

◆ isItinFidVolume()

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

Definition at line 721 of file HCalSD.cc.

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

Referenced by getFromHFLibrary(), and hitForFibre().

721  {
722  bool flag = true;
723  if (applyFidCut) {
724  int npmt = HFFibreFiducial::PMTNumber(hitPoint);
725 #ifdef EDM_ML_DEBUG
726  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume:#PMT= " << npmt << " for hit point " << hitPoint;
727 #endif
728  if (npmt <= 0)
729  flag = false;
730  }
731 #ifdef EDM_ML_DEBUG
732  edm::LogVerbatim("HcalSim") << "HCalSD::isItinFidVolume: point " << hitPoint << " return flag " << flag;
733 #endif
734  return flag;
735 }
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 689 of file HCalSD.cc.

References pmtLV.

Referenced by getEnergyDeposit().

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

◆ isItScintillator()

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

Definition at line 713 of file HCalSD.cc.

References materials.

Referenced by getEnergyDeposit().

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

◆ isItStraightBundle()

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

Definition at line 697 of file HCalSD.cc.

References fibre1LV.

Referenced by getEnergyDeposit().

697  {
698  for (auto lvol : fibre1LV)
699  if (lv == lvol) {
700  return true;
701  }
702  return false;
703 }
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 998 of file HCalSD.cc.

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

Referenced by getEnergyDeposit().

998  {
999  double wt = 1.;
1000  if (numberingFromDDD) {
1001  //get the ID's as eta, phi, depth, ... indices
1003  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
1004  modifyDepth(tmp);
1005  uint32_t id = HcalTestNumbering::packHcalIndex(tmp.subdet, tmp.zside, 1, tmp.etaR, tmp.phis, tmp.lay);
1006  std::map<uint32_t, double>::const_iterator ite = layerWeights.find(id);
1007  if (ite != layerWeights.end())
1008  wt = ite->second;
1009 #ifdef EDM_ML_DEBUG
1010  edm::LogVerbatim("HcalSim") << "HCalSD::layerWeight: ID " << std::hex << id << std::dec << " (" << tmp.subdet << "/"
1011  << tmp.zside << "/1/" << tmp.etaR << "/" << tmp.phis << "/" << tmp.lay << ") Weight "
1012  << wt;
1013 #endif
1014  }
1015  return wt;
1016 }
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1091
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 1091 of file HCalSD.cc.

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

Referenced by layerWeight(), and setDetUnitId().

1091  {
1092  if (id.subdet == 4) {
1093  int ieta = (id.zside == 0) ? -id.etaR : id.etaR;
1094  if (hcalConstants_->maxHFDepth(ieta, id.phis) > 2) {
1095  if (id.depth <= 2) {
1096  if (G4UniformRand() > 0.5)
1097  id.depth += 2;
1098  }
1099  }
1100  } else if ((id.subdet == 1 || id.subdet == 2) && testNumber) {
1101  id.depth = (depth_ == 0) ? 1 : 2;
1102  }
1103 }
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 1080 of file HCalSD.cc.

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

Referenced by getFromHFLibrary(), and hitForFibre().

1080  {
1081  double zv = std::abs(hitPoint.z()) - gpar[4];
1082  if (emType) {
1083  if (hzvem != nullptr)
1084  hzvem->Fill(zv);
1085  } else {
1086  if (hzvhad != nullptr)
1087  hzvhad->Fill(zv);
1088  }
1089 }
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 1018 of file HCalSD.cc.

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

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

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

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

Referenced by setDetUnitId().

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

◆ readWeightFromFile()

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

Definition at line 974 of file HCalSD.cc.

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

Referenced by HCalSD().

974  {
975  std::ifstream infile;
976  int entry = 0;
977  infile.open(fName.c_str(), std::ios::in);
978  if (infile) {
979  int det, zside, etaR, phi, lay;
980  double wt;
981  while (infile >> det >> zside >> etaR >> phi >> lay >> wt) {
982  uint32_t id = HcalTestNumbering::packHcalIndex(det, zside, 1, etaR, phi, lay);
983  layerWeights.insert(std::pair<uint32_t, double>(id, wt));
984  ++entry;
985 #ifdef EDM_ML_DEBUG
986  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile:Entry " << entry << " ID " << std::hex << id
987  << std::dec << " (" << det << "/" << zside << "/1/" << etaR << "/" << phi << "/"
988  << lay << ") Weight " << wt;
989 #endif
990  }
991  infile.close();
992  }
993  edm::LogVerbatim("HcalSim") << "HCalSD::readWeightFromFile: reads " << entry << " weights from " << fName;
994  if (entry <= 0)
995  useLayerWt = false;
996 }
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 531 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(), ppsFastLocalSimulation_cfi::z0, and HcalDetId::zside().

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

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

◆ setDetUnitId() [2/3]

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

Definition at line 597 of file HCalSD.cc.

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

597  {
598  uint32_t id = 0;
599  if (det == 0) {
600 #ifdef printDebug
601  double eta = std::abs(pos.eta());
602 #endif
603  if (std::abs(pos.z()) > maxZ_) {
604  det = 5;
605 #ifdef printDebug
606  if (eta < 2.868)
607  ++detNull_[2];
608 #endif
609  } else if (!(hcalConstants_->isHE())) {
610  det = 3;
611 #ifdef printDebug
612  ++detNull_[0];
613 #endif
614  } else {
615  double minR = minRoff_ + slopeHE_ * std::abs(pos.z());
616  double maxR = maxRoff_ + slopeHE_ * std::abs(pos.z());
617  det = ((pos.perp() > minR) && (pos.perp() < maxR)) ? 4 : 3;
618 #ifdef printDebug
619  ++detNull_[det - 3];
620 #endif
621  }
622 #ifdef printDEBUG
623  edm::LogVerbatim("HcalSim") << "Position " << pos.perp() << ":" << std::abs(pos.z()) << " Limits "
624  << !(hcalConstants_->isHE()) << ":" << maxZ_ << " det " << det;
625  } else {
626  ++detNull_[3];
627 #endif
628  }
629 
630  if (numberingFromDDD.get()) {
631  //get the ID's as eta, phi, depth, ... indices
633  numberingFromDDD->unitID(det, math::XYZVectorD(pos.x(), pos.y(), pos.z()), depth, lay);
634  id = setDetUnitId(tmp);
635  }
636  return id;
637 }
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:531
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 639 of file HCalSD.cc.

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

639  {
640  modifyDepth(tmp);
641  uint32_t id = (numberingScheme.get()) ? numberingScheme->getUnitID(tmp) : 0;
642  if ((!testNumber) && m_HcalTestNS.get()) {
643  bool ok = m_HcalTestNS->compare(tmp, id);
644  if (!ok)
645  edm::LogWarning("HcalSim") << "Det ID from HCalSD " << HcalDetId(id) << " " << std::hex << id << std::dec
646  << " does not match one from relabller for " << tmp.subdet << ":" << tmp.etaR << ":"
647  << tmp.phi << ":" << tmp.phis << ":" << tmp.depth << ":" << tmp.lay << std::endl;
648  }
649  return id;
650 }
void modifyDepth(HcalNumberingFromDDD::HcalID &id)
Definition: HCalSD.cc:1091
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 564 of file HCalSD.cc.

References numberingScheme, and generator_cfi::scheme.

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

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

◆ update() [1/6]

void CaloSD::update
overrideprotected

Definition at line 752 of file CaloSD.cc.

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

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

◆ update() [2/6]

void CaloSD::update
overrideprotected

Definition at line 781 of file CaloSD.cc.

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

781  {
782  endEvent();
783  for (int k = 0; k < nHC_; ++k) {
784  slave[k].get()->ReserveMemory(theHC[k]->entries());
785 
786  int count(0);
787  double eEM(0.0);
788  double eHAD(0.0);
789  double eEM2(0.0);
790  double eHAD2(0.0);
791 #ifdef EDM_ML_DEBUG
792  int wrong(0);
793  double tt(0.0);
794  double zloc(0.0);
795  double zglob(0.0);
796  double ee(0.0);
797 #endif
798  int hc_entries = theHC[k]->entries();
799  for (int i = 0; i < hc_entries; ++i) {
800 #ifdef EDM_ML_DEBUG
801  if (!saveHit((*theHC[k])[i], k)) {
802  ++wrong;
803  }
804 #else
805  saveHit((*theHC[k])[i], k);
806 #endif
807 
808  ++count;
809  double x = (*theHC[k])[i]->getEM();
810  eEM += x;
811  eEM2 += x * x;
812  x = (*theHC[k])[i]->getHadr();
813  eHAD += x;
814  eHAD2 += x * x;
815 #ifdef EDM_ML_DEBUG
816  tt += (*theHC[k])[i]->getTimeSlice();
817  ee += (*theHC[k])[i]->getIncidentEnergy();
818  zglob += std::abs((*theHC[k])[i]->getEntry().z());
819  zloc += std::abs((*theHC[k])[i]->getEntryLocal().z());
820 #endif
821  }
822 
823  double norm = (count > 0) ? 1.0 / count : 0.0;
824  eEM *= norm;
825  eEM2 *= norm;
826  eHAD *= norm;
827  eHAD2 *= norm;
828  eEM2 = std::sqrt(eEM2 - eEM * eEM);
829  eHAD2 = std::sqrt(eHAD2 - eHAD * eHAD);
830 #ifdef EDM_ML_DEBUG
831  tt *= norm;
832  ee *= norm;
833  zglob *= norm;
834  zloc *= norm;
835  edm::LogVerbatim("CaloSim") << "CaloSD: " << GetName() << " store " << count << " hits; " << wrong
836  << " track IDs not given properly and " << totalHits - count
837  << " hits not passing cuts\n EmeanEM= " << eEM << " ErmsEM= " << eEM2
838  << "\n EmeanHAD= " << eHAD << " ErmsHAD= " << eHAD2 << " TimeMean= " << tt
839  << " E0mean= " << ee << " Zglob= " << zglob << " Zloc= " << zloc << " ";
840 #endif
841  std::vector<std::unique_ptr<CaloG4Hit>>().swap(reusehit[k]);
842  if (useMap)
843  hitMap[k].erase(hitMap[k].begin(), hitMap[k].end());
844  }
845  tkMap.erase(tkMap.begin(), tkMap.end());
847 }
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:1111
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:950
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:23
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 1014 of file CaloSD.cc.

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

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

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

760  {
761  int id = (*trk)()->GetTrackID();
763  int lastTrackID = -1;
764  if (trkI)
765  lastTrackID = trkI->getIDonCaloSurface();
766  if (id == lastTrackID) {
767  auto trksForThisEvent = m_trackManager->trackContainer();
768  if (!trksForThisEvent->empty()) {
769  TrackWithHistory* trkH = trksForThisEvent->back();
770  if (trkH->trackID() == id) {
771  tkMap[id] = trkH;
772 #ifdef EDM_ML_DEBUG
773  edm::LogVerbatim("CaloSim") << "CaloSD: get track " << id << " from Container of size "
774  << trksForThisEvent->size() << " with ID " << trkH->trackID();
775 #endif
776  }
777  }
778  }
779 }
Log< level::Info, true > LogVerbatim
Trktree trk
Definition: Trktree.cc:2
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 750 of file CaloSD.cc.

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

750 { initRun(); }
void initRun() override
Definition: HCalSD.cc:573

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

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

571 {}

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