CMS 3D CMS Logo

HGCalTBMBProducer.cc
Go to the documentation of this file.
6 
12 
13 #include "G4LogicalVolumeStore.hh"
14 #include "G4Step.hh"
15 #include "G4Track.hh"
16 
17 #include <iostream>
18 #include <string>
19 #include <vector>
20 
21 //#define EDM_ML_DEBUG
22 
24  public Observer<const BeginOfTrack*>,
25  public Observer<const G4Step*>,
26  public Observer<const EndOfTrack*> {
27 public:
29  HGCalTBMBProducer(const HGCalTBMBProducer&) = delete; // stop default
30  const HGCalTBMBProducer& operator=(const HGCalTBMBProducer&) = delete; // ...
31  ~HGCalTBMBProducer() override = default;
32 
33  void produce(edm::Event&, const edm::EventSetup&) override;
34 
35 private:
36  void update(const BeginOfTrack*) override;
37  void update(const G4Step*) override;
38  void update(const EndOfTrack*) override;
39 
40  bool stopAfter(const G4Step*);
41  int findVolume(const G4VTouchable* touch, bool stop) const;
42 
44  const std::vector<std::string> listNames_;
46  const double stopZ_;
47  const unsigned int nList_;
49  std::vector<double> radLen_, intLen_, stepLen_;
50 };
51 
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 }
65 
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 }
73 
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 }
88 
89 void HGCalTBMBProducer::update(const G4Step* aStep) {
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 }
117 
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 }
136 
137 bool HGCalTBMBProducer::stopAfter(const G4Step* aStep) {
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 }
153 
154 int HGCalTBMBProducer::findVolume(const G4VTouchable* touch, bool stop) const {
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 }
175 
179 
size
Write out results.
Log< level::Info, true > LogVerbatim
#define DEFINE_SIMWATCHER(type)
MaterialAccountingCaloCollection matcoll_
const HGCalTBMBProducer & operator=(const HGCalTBMBProducer &)=delete
std::vector< MaterialAccountingCalo > MaterialAccountingCaloCollection
HGCalTBMBProducer(const edm::ParameterSet &)
const unsigned int nList_
std::vector< double > m_radLen
const std::string stopName_
void produce(edm::Event &, const edm::EventSetup &) override
std::vector< double > radLen_
void update(const BeginOfTrack *) override
This routine will be called when the appropriate signal arrives.
std::vector< double > stepLen_
const edm::ParameterSet m_p
ii
Definition: cuy.py:589
const std::vector< std::string > listNames_
std::vector< double > m_intLen
~HGCalTBMBProducer() override=default
HLT enums.
std::vector< double > intLen_
step
Definition: StallMonitor.cc:98
std::vector< double > m_stepLen
def move(src, dest)
Definition: eostools.py:511
bool stopAfter(const G4Step *)
int findVolume(const G4VTouchable *touch, bool stop) const