CMS 3D CMS Logo

SimG4FluxProducer.cc
Go to the documentation of this file.
1 #include "SimG4FluxProducer.h"
3 
4 #include "G4Track.hh"
5 #include "G4TouchableHistory.hh"
6 #include "G4TransportationManager.hh"
7 #include "CLHEP/Units/GlobalSystemOfUnits.h"
8 #include "CLHEP/Units/GlobalPhysicalConstants.h"
9 
10 #include <cmath>
11 #include <iostream>
12 #include <iomanip>
13 #include <memory>
14 #include <utility>
15 
16 //#define EDM_ML_DEBUG
17 
19  edm::ParameterSet m_FP = p.getParameter<edm::ParameterSet>("SimG4FluxProducer");
20  LVNames_ = m_FP.getUntrackedParameter<std::vector<std::string>>("LVNames");
21  LVTypes_ = m_FP.getUntrackedParameter<std::vector<int>>("LVTypes");
22 
23  for (unsigned int k = 0; k < LVNames_.size(); ++k) {
24  produces<ParticleFlux>(Form("%sParticleFlux", LVNames_[k].c_str()));
25 #ifdef EDM_ML_DEBUG
26  std::cout << "Collection name[" << k << "] ParticleFlux" << LVNames_[k] << " and type " << LVTypes_[k] << std::endl;
27 #endif
28  }
29 }
30 
32 
34  for (unsigned int k = 0; k < LVNames_.size(); ++k) {
35  std::unique_ptr<ParticleFlux> pflux(new ParticleFlux);
36  endOfEvent(*pflux, k);
37  std::string name = LVNames_[k] + "ParticleFlux";
38  e.put(std::move(pflux), name);
39  }
40 }
41 
43  topPV_ = getTopPV();
44  if (topPV_ == nullptr) {
45  edm::LogWarning("SimG4FluxProducer") << "Cannot find top level volume\n";
46  } else {
47  init_ = true;
48  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
49  for (auto lvcite : *lvs) {
50  findLV(lvcite);
51  }
52 
53 #ifdef EDM_ML_DEBUG
54  std::cout << "SimG4FluxProducer::Finds " << mapLV_.size() << " logical volumes\n";
55  unsigned int k(0);
56  for (const auto& lvs : mapLV_) {
57  std::cout << "Entry[" << k << "] " << lvs.first << ": (" << (lvs.second).first << ", " << (lvs.second).second
58  << ")\n";
59  ++k;
60  }
61 #endif
62  }
63 }
64 
66  int iev = (*evt)()->GetEventID();
67  edm::LogInfo("ValidHGCal") << "SimG4FluxProducer: =====> Begin event = " << iev << std::endl;
68  ++count_;
69  store_.clear();
70 }
71 
72 void SimG4FluxProducer::update(const G4Step* aStep) {
73  if (aStep != nullptr) {
74  G4TouchableHistory* touchable = (G4TouchableHistory*)aStep->GetPreStepPoint()->GetTouchable();
75  G4LogicalVolume* plv = (G4LogicalVolume*)touchable->GetVolume()->GetLogicalVolume();
76  auto it = (init_) ? mapLV_.find(plv) : findLV(plv);
77  // std::cout << plv->GetName() << " Flag " << (it != mapLV_.end()) << " step " << aStep->IsFirstStepInVolume() << ":" << aStep->IsLastStepInVolume() << std::endl;
78  if (it != mapLV_.end()) {
79  int type = LVTypes_[(it->second).first];
80  if ((type == 0 && aStep->IsFirstStepInVolume()) || (type == 1 && aStep->IsLastStepInVolume())) {
81  unsigned int copy = (unsigned int)(touchable->GetReplicaNumber(0));
82  std::pair<G4LogicalVolume*, unsigned int> key(plv, copy);
83  auto itr = store_.find(key);
84  if (itr == store_.end()) {
85  store_[key] = ParticleFlux((it->second).second, copy);
86  itr = store_.find(key);
87  }
88  G4Track* track = aStep->GetTrack();
89  int pdgid = track->GetDefinition()->GetPDGEncoding();
90  int vxtyp = (track->GetCreatorProcess() == nullptr) ? 0 : 1;
91  double time = (aStep->GetPostStepPoint()->GetGlobalTime());
92  const double mmTocm(0.1), MeVToGeV(0.001);
93  ParticleFlux::flux flx(pdgid, vxtyp, time);
94  flx.vertex = math::GlobalPoint(mmTocm * track->GetVertexPosition().x(),
95  mmTocm * track->GetVertexPosition().y(),
96  mmTocm * track->GetVertexPosition().z());
98  mmTocm * track->GetPosition().x(), mmTocm * track->GetPosition().y(), mmTocm * track->GetPosition().z());
99  flx.momentum = math::GlobalVector(MeVToGeV * track->GetMomentum().x(),
100  MeVToGeV * track->GetMomentum().y(),
101  MeVToGeV * track->GetMomentum().z());
102  (itr->second).addFlux(flx);
103 #ifdef EDM_ML_DEBUG
104  std::cout << "SimG4FluxProducer: Element " << (it->second).first << ":" << (it->second).second << ":" << copy
105  << " ID " << pdgid << " VxType " << vxtyp << " TOF " << time << " Hit Point " << flx.hitPoint << " p "
106  << flx.momentum << " Vertex " << flx.vertex << std::endl;
107 #endif
108  } //if(Step ok)
109  } //if( it != map.end() )
110  } //if (aStep != NULL)
111 }
112 
113 void SimG4FluxProducer::endOfEvent(ParticleFlux& flux, unsigned int k) {
114  bool done(false);
115  for (const auto& element : store_) {
116  G4LogicalVolume* lv = (element.first).first;
117  auto it = mapLV_.find(lv);
118  if (it != mapLV_.end()) {
119  if ((it->second).first == k) {
120  flux = element.second;
121  done = true;
122 #ifdef EDM_ML_DEBUG
123  std::cout << "SimG4FluxProducer[" << k << "] Flux " << flux.getName() << ":" << flux.getId() << " with "
124  << flux.getComponents() << " elements" << std::endl;
125  std::vector<ParticleFlux::flux> fluxes = flux.getFlux();
126  unsigned int k(0);
127  for (auto element : fluxes) {
128  std::cout << "Flux[" << k << "] PDGId " << element.pdgId << " VT " << element.vxType << " ToF " << element.tof
129  << " Vertex " << element.vertex << " Hit " << element.hitPoint << " p " << element.momentum
130  << std::endl;
131  ++k;
132  }
133 #endif
134  }
135  }
136  if (done)
137  break;
138  }
139 }
140 
141 G4VPhysicalVolume* SimG4FluxProducer::getTopPV() {
142  return G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume();
143 }
144 
145 std::map<G4LogicalVolume*, std::pair<unsigned int, std::string>>::iterator SimG4FluxProducer::findLV(
146  G4LogicalVolume* plv) {
147  auto itr = mapLV_.find(plv);
148  if (itr == mapLV_.end()) {
149  std::string name = plv->GetName();
150  for (unsigned int k = 0; k < LVNames_.size(); ++k) {
151  if (name.find(LVNames_[k]) != std::string::npos) {
152  mapLV_[plv] = std::pair<unsigned int, std::string>(k, name);
153  itr = mapLV_.find(plv);
154  break;
155  }
156  }
157  }
158  return itr;
159 }
type
Definition: HCALResponse.h:21
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
std::map< G4LogicalVolume *, std::pair< unsigned int, std::string > > mapLV_
SimG4FluxProducer(const edm::ParameterSet &p)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
~SimG4FluxProducer() override
void produce(edm::Event &, const edm::EventSetup &) override
G4VPhysicalVolume * topPV_
std::vector< int > LVTypes_
unsigned int getComponents() const
Definition: ParticleFlux.h:26
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)
int getId() const
Definition: ParticleFlux.h:25
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
std::string const & getName() const
Definition: ParticleFlux.h:24
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_
def move(src, dest)
Definition: eostools.py:511