CMS 3D CMS Logo

SimG4FluxProducer.cc
Go to the documentation of this file.
5 
7 
13 
14 #include "G4Step.hh"
15 #include "G4LogicalVolumeStore.hh"
16 #include "G4PhysicalVolumeStore.hh"
17 #include "G4Track.hh"
18 #include "G4TouchableHistory.hh"
19 #include "G4TransportationManager.hh"
20 
21 #include <CLHEP/Units/GlobalSystemOfUnits.h>
22 #include <CLHEP/Units/GlobalPhysicalConstants.h>
23 
24 #include <cmath>
25 #include <iostream>
26 #include <iomanip>
27 #include <map>
28 #include <memory>
29 #include <string>
30 #include <utility>
31 #include <vector>
32 
33 //#define EDM_ML_DEBUG
34 
36  public Observer<const BeginOfRun *>,
37  public Observer<const BeginOfEvent *>,
38  public Observer<const G4Step *> {
39 public:
41  SimG4FluxProducer(const SimG4FluxProducer &) = delete; // stop default
42  const SimG4FluxProducer &operator=(const SimG4FluxProducer &) = delete;
43  ~SimG4FluxProducer() override;
44 
45  void produce(edm::Event &, const edm::EventSetup &) override;
46 
47 private:
48  // observer classes
49  void update(const BeginOfRun *run) override;
50  void update(const BeginOfEvent *evt) override;
51  void update(const G4Step *step) override;
52 
53  void endOfEvent(ParticleFlux &pflx, unsigned int k);
54  G4VPhysicalVolume *getTopPV();
55  std::map<G4LogicalVolume *, std::pair<unsigned int, std::string>>::iterator findLV(G4LogicalVolume *plv);
56 
57  std::vector<std::string> LVNames_;
58  std::vector<int> LVTypes_;
59  G4VPhysicalVolume *topPV_;
60  std::map<G4LogicalVolume *, std::pair<unsigned int, std::string>> mapLV_;
61 
62  // some private members for ananlysis
63  unsigned int count_;
64  bool init_;
65  std::map<std::pair<G4LogicalVolume *, unsigned int>, ParticleFlux> store_;
66 };
67 
69  edm::ParameterSet m_FP = p.getParameter<edm::ParameterSet>("SimG4FluxProducer");
70  LVNames_ = m_FP.getUntrackedParameter<std::vector<std::string>>("LVNames");
71  LVTypes_ = m_FP.getUntrackedParameter<std::vector<int>>("LVTypes");
72 
73  for (unsigned int k = 0; k < LVNames_.size(); ++k) {
74  produces<ParticleFlux>(Form("%sParticleFlux", LVNames_[k].c_str()));
75 #ifdef EDM_ML_DEBUG
76  edm::LogVerbatim("SimG4FluxProducer")
77  << "SimG4FluxProducer::Collection name[" << k << "] ParticleFlux" << LVNames_[k] << " and type " << LVTypes_[k];
78 #endif
79  }
80 }
81 
83 
85 #ifdef EDM_ML_DEBUG
86  edm::LogVerbatim("SimG4FluxProducer") << "SimG4FluxProducer: enters produce with " << LVNames_.size() << " LV's";
87 #endif
88  for (unsigned int k = 0; k < LVNames_.size(); ++k) {
89  std::unique_ptr<ParticleFlux> pflux(new ParticleFlux);
90  endOfEvent(*pflux, k);
91  std::string name = LVNames_[k] + "ParticleFlux";
92  e.put(std::move(pflux), name);
93  }
94 }
95 
97 #ifdef EDM_ML_DEBUG
98  edm::LogVerbatim("SimG4FluxProducer") << "SimG4FluxProducer: enters BeginOfRun";
99 #endif
100  topPV_ = getTopPV();
101  if (topPV_ == nullptr) {
102  edm::LogWarning("SimG4FluxProducer") << "Cannot find top level volume\n";
103  } else {
104  init_ = true;
105  const G4LogicalVolumeStore *lvs = G4LogicalVolumeStore::GetInstance();
106  for (auto lvcite : *lvs) {
107  findLV(lvcite);
108  }
109 
110 #ifdef EDM_ML_DEBUG
111  edm::LogVerbatim("SimG4FluxProducer") << "SimG4FluxProducer::Finds " << mapLV_.size() << " logical volumes";
112  unsigned int k(0);
113  for (const auto &lvs : mapLV_) {
114  edm::LogVerbatim("SimG4FluxProducer")
115  << "Entry[" << k << "] " << lvs.first << ": (" << (lvs.second).first << ", " << (lvs.second).second << ")";
116  ++k;
117  }
118 #endif
119  }
120 }
121 
123  int iev = (*evt)()->GetEventID();
124  edm::LogVerbatim("SimG4FluxProducer") << "SimG4FluxProducer: =====> Begin event = " << iev;
125  ++count_;
126  store_.clear();
127 }
128 
129 void SimG4FluxProducer::update(const G4Step *aStep) {
130  if (aStep != nullptr) {
131  G4TouchableHistory *touchable = (G4TouchableHistory *)aStep->GetPreStepPoint()->GetTouchable();
132  G4LogicalVolume *plv = (G4LogicalVolume *)touchable->GetVolume()->GetLogicalVolume();
133  auto it = (init_) ? mapLV_.find(plv) : findLV(plv);
134  // edm::LogVerbatim("SimG4FluxProducer") << plv->GetName() << " Flag " << (it != mapLV_.end()) << " step " << aStep->IsFirstStepInVolume() << ":" << aStep->IsLastStepInVolume();
135  if (it != mapLV_.end()) {
136  int type = LVTypes_[(it->second).first];
137  if ((type == 0 && aStep->IsFirstStepInVolume()) || (type == 1 && aStep->IsLastStepInVolume())) {
138  unsigned int copy = (unsigned int)(touchable->GetReplicaNumber(0));
139  std::pair<G4LogicalVolume *, unsigned int> key(plv, copy);
140  auto itr = store_.find(key);
141  if (itr == store_.end()) {
142  store_[key] = ParticleFlux((it->second).second, copy);
143  itr = store_.find(key);
144  }
145  G4Track *track = aStep->GetTrack();
146  int pdgid = track->GetDefinition()->GetPDGEncoding();
147  int vxtyp = (track->GetCreatorProcess() == nullptr) ? 0 : 1;
148  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
149  const double mmTocm(0.1), MeVToGeV(0.001);
150  ParticleFlux::flux flx(pdgid, vxtyp, time);
151  flx.vertex = math::GlobalPoint(mmTocm * track->GetVertexPosition().x(),
152  mmTocm * track->GetVertexPosition().y(),
153  mmTocm * track->GetVertexPosition().z());
155  mmTocm * track->GetPosition().x(), mmTocm * track->GetPosition().y(), mmTocm * track->GetPosition().z());
156  flx.momentum = math::GlobalVector(MeVToGeV * track->GetMomentum().x(),
157  MeVToGeV * track->GetMomentum().y(),
158  MeVToGeV * track->GetMomentum().z());
159  (itr->second).addFlux(flx);
160 #ifdef EDM_ML_DEBUG
161  edm::LogVerbatim("SimG4FluxProducer")
162  << "SimG4FluxProducer: Element " << (it->second).first << ":" << (it->second).second << ":" << copy
163  << " ID " << pdgid << " VxType " << vxtyp << " TOF " << time << " Hit Point " << flx.hitPoint << " p "
164  << flx.momentum << " Vertex " << flx.vertex;
165 #endif
166  } //if(Step ok)
167  } //if( it != map.end() )
168  } //if (aStep != NULL)
169 }
170 
171 void SimG4FluxProducer::endOfEvent(ParticleFlux &flux, unsigned int k) {
172  bool done(false);
173  for (const auto &element : store_) {
174  G4LogicalVolume *lv = (element.first).first;
175  auto it = mapLV_.find(lv);
176  if (it != mapLV_.end()) {
177  if ((it->second).first == k) {
178  flux = element.second;
179  done = true;
180 #ifdef EDM_ML_DEBUG
181  edm::LogVerbatim("SimG4FluxProducer") << "SimG4FluxProducer[" << k << "] Flux " << flux.getName() << ":"
182  << flux.getId() << " with " << flux.getComponents() << " elements";
183  std::vector<ParticleFlux::flux> fluxes = flux.getFlux();
184  unsigned int k(0);
185  for (auto element : fluxes) {
186  edm::LogVerbatim("SimG4FluxProducer")
187  << "Flux[" << k << "] PDGId " << element.pdgId << " VT " << element.vxType << " ToF " << element.tof
188  << " Vertex " << element.vertex << " Hit " << element.hitPoint << " p " << element.momentum;
189  ++k;
190  }
191 #endif
192  }
193  }
194  if (done)
195  break;
196  }
197 }
198 
199 G4VPhysicalVolume *SimG4FluxProducer::getTopPV() {
200  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
201 }
202 
203 std::map<G4LogicalVolume *, std::pair<unsigned int, std::string>>::iterator SimG4FluxProducer::findLV(
204  G4LogicalVolume *plv) {
205  auto itr = mapLV_.find(plv);
206  if (itr == mapLV_.end()) {
207  std::string name = plv->GetName();
208  for (unsigned int k = 0; k < LVNames_.size(); ++k) {
209  if (name.find(LVNames_[k]) != std::string::npos) {
210  mapLV_[plv] = std::pair<unsigned int, std::string>(k, name);
211  itr = mapLV_.find(plv);
212  break;
213  }
214  }
215  }
216  return itr;
217 }
220 
Log< level::Info, true > LogVerbatim
unsigned int getComponents() const
Definition: ParticleFlux.h:26
#define DEFINE_SIMWATCHER(type)
std::map< G4LogicalVolume *, std::pair< unsigned int, std::string > > mapLV_
const SimG4FluxProducer & operator=(const SimG4FluxProducer &)=delete
SimG4FluxProducer(const edm::ParameterSet &p)
~SimG4FluxProducer() override
void produce(edm::Event &, const edm::EventSetup &) override
G4VPhysicalVolume * topPV_
std::vector< int > LVTypes_
TkSoAView< TrackerTraits > HitToTuple< TrackerTraits > const *__restrict__ int32_t int32_t int iev
void update(const BeginOfRun *run) override
This routine will be called when the appropriate signal arrives.
std::map< G4LogicalVolume *, std::pair< unsigned int, std::string > >::iterator findLV(G4LogicalVolume *plv)
std::string const & getName() const
Definition: ParticleFlux.h:24
T getUntrackedParameter(std::string const &, T const &) const
math::GlobalVector momentum
Definition: ParticleFlux.h:20
math::GlobalPoint vertex
Definition: ParticleFlux.h:19
G4VPhysicalVolume * getTopPV()
std::vector< std::string > LVNames_
void endOfEvent(ParticleFlux &pflx, unsigned int k)
ROOT::Math::DisplacementVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalVector
vector in glovbal coordinate system
Definition: Vector3D.h:28
math::GlobalPoint hitPoint
Definition: ParticleFlux.h:19
std::vector< ParticleFlux::flux > getFlux() const &
Definition: ParticleFlux.h:27
ROOT::Math::PositionVector3D< ROOT::Math::Cartesian3D< float >, ROOT::Math::GlobalCoordinateSystemTag > GlobalPoint
point in global coordinate system
Definition: Point3D.h:18
std::map< std::pair< G4LogicalVolume *, unsigned int >, ParticleFlux > store_
step
Definition: StallMonitor.cc:83
Log< level::Warning, false > LogWarning
def move(src, dest)
Definition: eostools.py:511
int getId() const
Definition: ParticleFlux.h:25