CMS 3D CMS Logo

List of all members | Public Member Functions | Private Member Functions | Private Attributes
HGCalTBMBProducer Class Reference
Inheritance diagram for HGCalTBMBProducer:
SimProducer Observer< const BeginOfTrack *> Observer< const G4Step *> Observer< const EndOfTrack *> SimWatcher

Public Member Functions

 HGCalTBMBProducer (const edm::ParameterSet &)
 
 HGCalTBMBProducer (const HGCalTBMBProducer &)=delete
 
const HGCalTBMBProduceroperator= (const HGCalTBMBProducer &)=delete
 
void produce (edm::Event &, const edm::EventSetup &) override
 
 ~HGCalTBMBProducer () override=default
 
- Public Member Functions inherited from SimProducer
const SimProduceroperator= (const SimProducer &)=delete
 
void registerProducts (edm::ProducesCollector producesCollector)
 
 SimProducer ()
 
 SimProducer (const SimProducer &)=delete
 
- Public Member Functions inherited from SimWatcher
virtual void beginRun (edm::EventSetup const &)
 
bool isMT () const
 
const SimWatcheroperator= (const SimWatcher &)=delete
 
virtual void registerConsumes (edm::ConsumesCollector)
 
 SimWatcher ()
 
 SimWatcher (const SimWatcher &)=delete
 
virtual ~SimWatcher ()
 
- Public Member Functions inherited from Observer< const BeginOfTrack *>
 Observer ()
 
void slotForUpdate (const BeginOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const G4Step *>
 Observer ()
 
void slotForUpdate (const G4Step * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfTrack *>
 Observer ()
 
void slotForUpdate (const EndOfTrack * iT)
 
virtual ~Observer ()
 

Private Member Functions

int findVolume (const G4VTouchable *touch, bool stop) const
 
bool stopAfter (const G4Step *)
 
void update (const BeginOfTrack *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const G4Step *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const EndOfTrack *) override
 This routine will be called when the appropriate signal arrives. More...
 

Private Attributes

std::vector< double > intLen_
 
const std::vector< std::string > listNames_
 
const edm::ParameterSet m_p
 
MaterialAccountingCaloCollection matcoll_
 
const unsigned int nList_
 
std::vector< double > radLen_
 
std::vector< double > stepLen_
 
const std::string stopName_
 
const double stopZ_
 

Additional Inherited Members

- Protected Member Functions inherited from SimProducer
template<class T >
void produces ()
 
template<class T >
void produces (const std::string &instanceName)
 
- Protected Member Functions inherited from SimWatcher
void setMT (bool val)
 

Detailed Description

Definition at line 23 of file HGCalTBMBProducer.cc.

Constructor & Destructor Documentation

◆ HGCalTBMBProducer() [1/2]

HGCalTBMBProducer::HGCalTBMBProducer ( const edm::ParameterSet p)

Definition at line 52 of file HGCalTBMBProducer.cc.

References dqmdumpme::k, listNames_, nList_, stopName_, and stopZ_.

53  : m_p(p.getParameter<edm::ParameterSet>("HGCalTBMB")),
54  listNames_(m_p.getParameter<std::vector<std::string> >("DetectorNames")),
55  stopName_(m_p.getParameter<std::string>("StopName")),
56  stopZ_(m_p.getParameter<double>("MaximumZ")),
57  nList_(listNames_.size()) {
58  edm::LogVerbatim("HGCSim") << "HGCalTBMBProducer initialized for " << nList_ << " volumes";
59  for (unsigned int k = 0; k < nList_; ++k)
60  edm::LogVerbatim("HGCSim") << " [" << k << "] " << listNames_[k];
61  edm::LogVerbatim("HGCSim") << "Stop after " << stopZ_ << " or reaching volume " << stopName_;
62 
63  produces<MaterialAccountingCaloCollection>("HGCalTBMB");
64 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
const unsigned int nList_
const std::string stopName_
const edm::ParameterSet m_p
const std::vector< std::string > listNames_

◆ HGCalTBMBProducer() [2/2]

HGCalTBMBProducer::HGCalTBMBProducer ( const HGCalTBMBProducer )
delete

◆ ~HGCalTBMBProducer()

HGCalTBMBProducer::~HGCalTBMBProducer ( )
overridedefault

Member Function Documentation

◆ findVolume()

int HGCalTBMBProducer::findVolume ( const G4VTouchable *  touch,
bool  stop 
) const
private

Definition at line 154 of file HGCalTBMBProducer.cc.

References cuy::ii, dqmdumpme::k, personalPlayback::level, listNames_, Skims_PA_cff::name, nList_, stopName_, and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by stopAfter(), and update().

154  {
155  int ivol = -1;
156  int level = (touch->GetHistoryDepth()) + 1;
157  for (int ii = 0; ii < level; ii++) {
158  std::string name = touch->GetVolume(ii)->GetName();
159  if (stop) {
160  if (strcmp(name.c_str(), stopName_.c_str()) == 0)
161  ivol = 0;
162  } else {
163  for (unsigned int k = 0; k < nList_; ++k) {
164  if (strcmp(name.c_str(), listNames_[k].c_str()) == 0) {
165  ivol = k;
166  break;
167  }
168  }
169  }
170  if (ivol >= 0)
171  break;
172  }
173  return ivol;
174 }
const unsigned int nList_
const std::string stopName_
ii
Definition: cuy.py:589
const std::vector< std::string > listNames_

◆ operator=()

const HGCalTBMBProducer& HGCalTBMBProducer::operator= ( const HGCalTBMBProducer )
delete

◆ produce()

void HGCalTBMBProducer::produce ( edm::Event e,
const edm::EventSetup  
)
overridevirtual

Implements SimProducer.

Definition at line 66 of file HGCalTBMBProducer.cc.

References MillePedeFileConverter_cfg::e, caloTruthProducer_cfi::hgc, matcoll_, and eostools::move().

66  {
67  std::unique_ptr<MaterialAccountingCaloCollection> hgc(new MaterialAccountingCaloCollection);
68  for (auto const& mbc : matcoll_) {
69  hgc->emplace_back(mbc);
70  }
71  e.put(std::move(hgc), "HGCalTBMB");
72 }
MaterialAccountingCaloCollection matcoll_
std::vector< MaterialAccountingCalo > MaterialAccountingCaloCollection
def move(src, dest)
Definition: eostools.py:511

◆ stopAfter()

bool HGCalTBMBProducer::stopAfter ( const G4Step *  aStep)
private

Definition at line 137 of file HGCalTBMBProducer.cc.

References findVolume(), RemoveAddSevLevel::flag, stopZ_, and geometryCSVtoXML::zz.

Referenced by update().

137  {
138  bool flag(false);
139  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
140  G4ThreeVector hitPoint = aStep->GetPreStepPoint()->GetPosition();
141  if (aStep->GetPostStepPoint() != nullptr)
142  hitPoint = aStep->GetPostStepPoint()->GetPosition();
143  double zz = hitPoint.z();
144 
145  if ((findVolume(touch, true) == 0) || (zz > stopZ_))
146  flag = true;
147 #ifdef EDM_ML_DEBUG
148  edm::LogVerbatim("HGCSim") << " HGCalTBMBProducer::Name " << touch->GetVolume(0)->GetName() << " z " << zz << " Flag"
149  << flag;
150 #endif
151  return flag;
152 }
Log< level::Info, true > LogVerbatim
int findVolume(const G4VTouchable *touch, bool stop) const

◆ update() [1/3]

void HGCalTBMBProducer::update ( const BeginOfTrack )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack *>.

Definition at line 74 of file HGCalTBMBProducer.cc.

References createfilelist::int, intLen_, nList_, radLen_, and stepLen_.

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

74  {
75  radLen_ = std::vector<double>(nList_ + 1, 0);
76  intLen_ = std::vector<double>(nList_ + 1, 0);
77  stepLen_ = std::vector<double>(nList_ + 1, 0);
78 
79 #ifdef EDM_ML_DEBUG
80  const G4Track* aTrack = (*trk)(); // recover G4 pointer if wanted
81  const G4ThreeVector& mom = aTrack->GetMomentum();
82  double theEnergy = aTrack->GetTotalEnergy();
83  int theID = (int)(aTrack->GetDefinition()->GetPDGEncoding());
84  edm::LogVerbatim("HGCSim") << "HGCalTBMBProducer: Track " << aTrack->GetTrackID() << " Code " << theID << " Energy "
85  << theEnergy / CLHEP::GeV << " GeV; Momentum " << mom;
86 #endif
87 }
Log< level::Info, true > LogVerbatim
const unsigned int nList_
std::vector< double > radLen_
std::vector< double > stepLen_
std::vector< double > intLen_

◆ update() [2/3]

void HGCalTBMBProducer::update ( const G4Step *  )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const G4Step *>.

Definition at line 89 of file HGCalTBMBProducer.cc.

References findVolume(), intLen_, nList_, radLen_, stepLen_, stopAfter(), and HLT_2022v12_cff::track.

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

89  {
90  G4Material* material = aStep->GetPreStepPoint()->GetMaterial();
91  double step = aStep->GetStepLength();
92  double radl = material->GetRadlen();
93  double intl = material->GetNuclearInterLength();
94 
95  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
96  int indx = findVolume(touch, false);
97 
98  if (indx >= 0) {
99  stepLen_[indx] += step;
100  radLen_[indx] += (step / radl);
101  intLen_[indx] += (step / intl);
102  }
103  stepLen_[nList_] += step;
104  radLen_[nList_] += (step / radl);
105  intLen_[nList_] += (step / intl);
106 #ifdef EDM_ML_DEBUG
107  edm::LogVerbatim("HGCSim") << "HGCalTBMBProducer::Step in " << touch->GetVolume(0)->GetLogicalVolume()->GetName()
108  << " Index " << indx << " Step " << step << " RadL " << step / radl << " IntL "
109  << step / intl;
110 #endif
111 
112  if (stopAfter(aStep)) {
113  G4Track* track = aStep->GetTrack();
114  track->SetTrackStatus(fStopAndKill);
115  }
116 }
Log< level::Info, true > LogVerbatim
const unsigned int nList_
std::vector< double > radLen_
std::vector< double > stepLen_
std::vector< double > intLen_
step
Definition: StallMonitor.cc:98
bool stopAfter(const G4Step *)
int findVolume(const G4VTouchable *touch, bool stop) const

◆ update() [3/3]

void HGCalTBMBProducer::update ( const EndOfTrack )
overrideprivatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const EndOfTrack *>.

Definition at line 118 of file HGCalTBMBProducer.cc.

References cuy::ii, intLen_, listNames_, MaterialAccountingCalo::m_eta, MaterialAccountingCalo::m_intLen, MaterialAccountingCalo::m_phi, MaterialAccountingCalo::m_radLen, MaterialAccountingCalo::m_stepLen, matcoll_, Skims_PA_cff::name, nList_, radLen_, stepLen_, and AlCaHLTBitMon_QueryRunRegistry::string.

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

118  {
119  MaterialAccountingCalo matCalo;
120  matCalo.m_eta = 0;
121  matCalo.m_phi = 0;
122  for (unsigned int ii = 0; ii <= nList_; ++ii) {
123  matCalo.m_stepLen.emplace_back(stepLen_[ii]);
124  matCalo.m_radLen.emplace_back(radLen_[ii]);
125  matCalo.m_intLen.emplace_back(intLen_[ii]);
126 #ifdef EDM_ML_DEBUG
127  std::string name("Total");
128  if (ii < nList_)
129  name = listNames_[ii];
130  edm::LogVerbatim("HGCSim") << "HGCalTBMBProducer::Volume[" << ii << "]: " << name << " == Step " << stepLen_[ii]
131  << " RadL " << radLen_[ii] << " IntL " << intLen_[ii];
132 #endif
133  }
134  matcoll_.emplace_back(matCalo);
135 }
Log< level::Info, true > LogVerbatim
MaterialAccountingCaloCollection matcoll_
const unsigned int nList_
std::vector< double > m_radLen
std::vector< double > radLen_
std::vector< double > stepLen_
ii
Definition: cuy.py:589
const std::vector< std::string > listNames_
std::vector< double > m_intLen
std::vector< double > intLen_
std::vector< double > m_stepLen

Member Data Documentation

◆ intLen_

std::vector<double> HGCalTBMBProducer::intLen_
private

Definition at line 49 of file HGCalTBMBProducer.cc.

Referenced by update().

◆ listNames_

const std::vector<std::string> HGCalTBMBProducer::listNames_
private

Definition at line 44 of file HGCalTBMBProducer.cc.

Referenced by findVolume(), HGCalTBMBProducer(), and update().

◆ m_p

const edm::ParameterSet HGCalTBMBProducer::m_p
private

Definition at line 43 of file HGCalTBMBProducer.cc.

◆ matcoll_

MaterialAccountingCaloCollection HGCalTBMBProducer::matcoll_
private

Definition at line 48 of file HGCalTBMBProducer.cc.

Referenced by produce(), and update().

◆ nList_

const unsigned int HGCalTBMBProducer::nList_
private

Definition at line 47 of file HGCalTBMBProducer.cc.

Referenced by findVolume(), HGCalTBMBProducer(), and update().

◆ radLen_

std::vector<double> HGCalTBMBProducer::radLen_
private

Definition at line 49 of file HGCalTBMBProducer.cc.

Referenced by update().

◆ stepLen_

std::vector<double> HGCalTBMBProducer::stepLen_
private

Definition at line 49 of file HGCalTBMBProducer.cc.

Referenced by update().

◆ stopName_

const std::string HGCalTBMBProducer::stopName_
private

Definition at line 45 of file HGCalTBMBProducer.cc.

Referenced by findVolume(), and HGCalTBMBProducer().

◆ stopZ_

const double HGCalTBMBProducer::stopZ_
private

Definition at line 46 of file HGCalTBMBProducer.cc.

Referenced by HGCalTBMBProducer(), and stopAfter().