CMS 3D CMS Logo

TrackingMaterialProducer.cc
Go to the documentation of this file.
1 #include <iostream> // FIXME: switch to MessagLogger & friends
2 #include <vector>
3 #include <string>
4 #include <cassert>
5 #include <exception>
6 #include <tuple>
7 
13 
18 
20 
21 // GEANT4
22 #include "G4Step.hh"
23 #include "G4Track.hh"
24 #include "G4VSolid.hh"
25 #include "G4LogicalVolumeStore.hh"
26 #include "G4TouchableHistory.hh"
27 #include "G4VPhysicalVolume.hh"
28 #include "G4AffineTransform.hh"
29 
31 
32 // Uncomment the following #define directive to have the full list of
33 // volumes known to G4 printed to LogInfo("TrackingMaterialProducer")
34 
35 #define DEBUG_G4_VOLUMES
36 
37 using namespace CLHEP;
38 using edm::LogInfo;
39 
40 // missing from GEANT4 < 9.0 : G4LogicalVolumeStore::GetVolume( name )
41 static
42 const G4LogicalVolume* GetVolume(const std::string& name) {
43  const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
44 
45 #ifdef DEBUG_G4_VOLUMES
46  for (G4LogicalVolumeStore::const_iterator volume = lvs->begin(); volume != lvs->end(); ++volume)
47  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer: G4 registered volumes "
48  << (*volume)->GetName() << std::endl;
49 #endif
50 
51  for (G4LogicalVolumeStore::const_iterator volume = lvs->begin(); volume != lvs->end(); ++volume) {
52  if ((const std::string&) (*volume)->GetName() == name)
53  return (*volume);
54  }
55  return 0;
56 }
57 
58 // missing from GEANT4 : G4TouchableHistory::GetTransform( depth )
59 static inline
60 const G4AffineTransform& GetTransform(const G4TouchableHistory* touchable, int depth)
61 {
62  return touchable->GetHistory()->GetTransform( touchable->GetHistory()->GetDepth() - depth );
63 }
64 
65 // navigate up the hierarchy of volumes until one with an attached sensitive detector is found
66 // return a tuple holding
67 // - pointer to the first (deepest) sensitive G4VPhysicalVolume
68 // - how may steps up in the hierarchy it is (0 is the starting volume)
69 // if no sensitive detector is found, return a NULL pointer and 0
70 
71 std::tuple<const G4VPhysicalVolume*, int> GetSensitiveVolume( const G4VTouchable* touchable )
72 {
73  int depth = touchable->GetHistoryDepth();
74  for (int level = 0; level < depth; ++level) { // 0 is self
75  const G4VPhysicalVolume* volume = touchable->GetVolume(level);
76  if (volume->GetLogicalVolume()->GetSensitiveDetector() != 0) {
77  return std::make_tuple(volume, level);
78  }
79  }
80  return std::tuple<const G4VPhysicalVolume*, int>(0, 0);
81 }
82 
83 //-------------------------------------------------------------------------
85 {
86  edm::ParameterSet config = iPSet.getParameter<edm::ParameterSet>("TrackingMaterialProducer");
87  m_selectedNames = config.getParameter< std::vector<std::string> >("SelectedVolumes");
88  m_primaryTracks = config.getParameter<bool>("PrimaryTracksOnly");
89  m_tracks = 0;
90 
91  produces< std::vector<MaterialAccountingTrack> >();
92  output_file_ = new TFile("radLen_vs_eta_fromProducer.root", "RECREATE");
93  output_file_->cd();
94  radLen_vs_eta_ = new TProfile("radLen", "radLen", 250., -5., 5., 0, 10.);
95 }
96 
97 //-------------------------------------------------------------------------
99 {
100 }
101 
102 //-------------------------------------------------------------------------
104 {
105  radLen_vs_eta_->Write();
106  output_file_->Close();
107 }
108 //-------------------------------------------------------------------------
110 {
111  // INFO
112  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer: List of the selected volumes: " << std::endl;
113  for (std::vector<std::string>::const_iterator volume_name = m_selectedNames.begin();
114  volume_name != m_selectedNames.end(); ++volume_name) {
115  const G4LogicalVolume* volume = GetVolume(*volume_name);
116  if (volume) {
117  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer: " << *volume_name << std::endl;
118  m_selectedVolumes.push_back( volume );
119  } else {
120  // FIXME: throw an exception ?
121  std::cerr << "TrackingMaterialProducer::update(const BeginOfJob*): WARNING: selected volume \"" << *volume_name << "\" not found in geometry " << std::endl;
122  }
123  }
124 }
125 
126 
127 //-------------------------------------------------------------------------
129 {
130  m_tracks = new std::vector<MaterialAccountingTrack>();
131 }
132 
133 
134 //-------------------------------------------------------------------------
136 {
137  m_track.reset();
138 
139  // prevent secondary tracks from propagating
140  G4Track* track = const_cast<G4Track*>((*event)());
141  if (m_primaryTracks and track->GetParentID() != 0) {
142  track->SetTrackStatus(fStopAndKill);
143  }
144 }
145 
146 bool TrackingMaterialProducer::isSelectedFast(const G4TouchableHistory* touchable) {
147  for (int d = touchable->GetHistoryDepth() -1; d >=0; --d) {
148  if (
149  std::find(
150  m_selectedNames.begin(),
151  m_selectedNames.end(),
152  touchable->GetVolume(d)->GetName())
153  != m_selectedNames.end())
154  return true;
155  }
156  return false;
157 }
158 
159 //-------------------------------------------------------------------------
161 {
162  const G4TouchableHistory* touchable = (G4TouchableHistory*)(step->GetTrack()->GetTouchable());
163  if (not isSelectedFast( touchable )) {
164  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer:\t[...] skipping "
165  << touchable->GetVolume()->GetName() << std::endl;
166  return;
167  }
168 
169  // material and step proterties
170  const G4Material* material = touchable->GetVolume()->GetLogicalVolume()->GetMaterial();
171  double length = step->GetStepLength() / cm; // mm -> cm
172  double X0 = material->GetRadlen() / cm; // mm -> cm
173  double Ne = material->GetElectronDensity() * cm3; // 1/mm3 -> 1/cm3
174  double Xi = Ne / 6.0221415e23 * 0.307075 / 2; // MeV / cm
175  double radiationLengths = length / X0; //
176  double energyLoss = length * Xi / 1000.; // GeV
177  //double energyLoss = step->GetDeltaEnergy()/MeV; should we use this??
178 
179  G4ThreeVector globalPosPre = step->GetPreStepPoint()->GetPosition();
180  G4ThreeVector globalPosPost = step->GetPostStepPoint()->GetPosition();
181  GlobalPoint globalPositionIn( globalPosPre.x() / cm, globalPosPre.y() / cm, globalPosPre.z() / cm ); // mm -> cm
182  GlobalPoint globalPositionOut( globalPosPost.x() / cm, globalPosPost.y() / cm, globalPosPost.z() / cm ); // mm -> cm
183 
184  // check for a sensitive detector
185  bool enter_sensitive = false;
186  bool leave_sensitive = false;
187  double cosThetaPre = 0.0;
188  double cosThetaPost = 0.0;
189  int level = 0;
190  const G4VPhysicalVolume* sensitive = 0;
192  std::tie(sensitive, level) = GetSensitiveVolume(touchable);
193  if (sensitive) {
194  const G4VSolid & solid = *touchable->GetSolid( level );
195  const G4AffineTransform & transform = GetTransform( touchable, level );
196  G4ThreeVector pos = transform.Inverse().TransformPoint( G4ThreeVector( 0., 0., 0. ) );
197  position = GlobalPoint( pos.x() / cm, pos.y() / cm, pos.z() / cm ); // mm -> cm
198 
199  G4ThreeVector localPosPre = transform.TransformPoint( globalPosPre );
200  EInside statusPre = solid.Inside( localPosPre );
201  if (statusPre == kSurface) {
202  enter_sensitive = true;
203  G4ThreeVector globalDirPre = step->GetPreStepPoint()->GetMomentumDirection();
204  G4ThreeVector localDirPre = transform.TransformAxis( globalDirPre );
205  G4ThreeVector normalPre = solid.SurfaceNormal( localPosPre );
206  cosThetaPre = normalPre.cosTheta( -localDirPre );
207  }
208 
209  G4ThreeVector localPosPost = transform.TransformPoint( globalPosPost );
210  EInside statusPost = solid.Inside( localPosPost );
211  if (statusPost == kSurface) {
212  leave_sensitive = true;
213  G4ThreeVector globalDirPost = step->GetPostStepPoint()->GetMomentumDirection();
214  G4ThreeVector localDirPost = transform.TransformAxis( globalDirPost );
215  G4ThreeVector normalPost = solid.SurfaceNormal( localPosPost );
216  cosThetaPost = normalPost.cosTheta( localDirPost );
217  }
218  }
219 
220  // update track accounting
221  if (enter_sensitive)
222  m_track.enterDetector( sensitive, position, cosThetaPre );
223  m_track.step(MaterialAccountingStep( length, radiationLengths, energyLoss, globalPositionIn, globalPositionOut ));
224  if (leave_sensitive)
225  m_track.leaveDetector( sensitive, cosThetaPost );
226 
227  if (sensitive)
228  LogInfo("TrackingMaterialProducer") << "Track was near sensitive volume "
229  << sensitive->GetName() << std::endl;
230  else
231  LogInfo("TrackingMaterialProducer") << "Track was near non-sensitive volume "
232  << touchable->GetVolume()->GetName() << std::endl;
233  LogInfo("TrackingMaterialProducer") << "Step length: "
234  << length << " cm\n"
235  << "globalPreStep(r,z): (" << globalPositionIn.perp()
236  << ", " << globalPositionIn.z() << ") cm\n"
237  << "globalPostStep(r,z): (" << globalPositionOut.perp()
238  << ", " << globalPositionOut.z() << ") cm\n"
239  << "position(r,z): ("
240  << position.perp()
241  << ", " << position.z() << ") cm\n"
242  << "Radiation lengths: "
243  << radiationLengths << " \t\t(X0: "
244  << X0 << " cm)\n"
245  << "Energy loss: "
246  << energyLoss << " MeV \t(Xi: "
247  << Xi << " MeV/cm)\n"
248  << "Track was " << (enter_sensitive ? "entering " : "in none ")
249  << "sensitive volume\n"
250  << "Track was " << (leave_sensitive ? "leaving " : "in none ")
251  << "sensitive volume\n";
252 
253 }
254 
255 
256 //-------------------------------------------------------------------------
258 {
259  const G4Track * track = (*event)();
260  if (m_primaryTracks and track->GetParentID() != 0)
261  return;
262 
263  radLen_vs_eta_->Fill(track->GetMomentum().eta(), m_track.summary().radiationLengths());
264  m_tracks->push_back(m_track);
265 
266  // LogInfo
267  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer: this track took "
268  << m_track.steps().size()
269  << " steps, and passed through "
270  << m_track.detectors().size()
271  << " sensitive detectors" << std::endl;
272  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer: track length: "
273  << m_track.summary().length()
274  << " cm" << std::endl;
275  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer: radiation lengths: "
276  << m_track.summary().radiationLengths() << std::endl;
277  LogInfo("TrackingMaterialProducer") << "TrackingMaterialProducer: energy loss: "
278  << m_track.summary().energyLoss()
279  << " MeV" << std::endl;
280 }
281 
282 //-------------------------------------------------------------------------
284 {
285  // transfer ownership to the Event
286  std::unique_ptr<std::vector<MaterialAccountingTrack> > tracks( m_tracks );
287  iEvent.put(std::move(tracks));
288  m_tracks = 0;
289 }
290 
291 //-------------------------------------------------------------------------
292 bool TrackingMaterialProducer::isSelected( const G4VTouchable* touchable )
293 {
294  for (size_t i = 0; i < m_selectedVolumes.size(); ++i)
295  if (m_selectedVolumes[i]->IsAncestor( touchable->GetVolume() )
296  or m_selectedVolumes[i] == touchable->GetVolume()->GetLogicalVolume())
297  return true;
298 
299  return false;
300 }
301 
302 //-------------------------------------------------------------------------
303 // define as a plugin
T getParameter(std::string const &) const
#define DEFINE_SIMWATCHER(type)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:122
T perp() const
Definition: PV3DBase.h:72
void produce(edm::Event &, const edm::EventSetup &)
Global3DPoint GlobalPoint
Definition: GlobalPoint.h:10
Definition: config.py:1
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:20
TrackingMaterialProducer(const edm::ParameterSet &)
static const G4AffineTransform & GetTransform(const G4TouchableHistory *touchable, int depth)
int iEvent
Definition: GenABIO.cc:230
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< void, edm::EventID const &, edm::Timestamp const & > We also list in braces which AR_WATCH_USING_METHOD_ is used for those or
Definition: Activities.doc:12
T z() const
Definition: PV3DBase.h:64
static const G4LogicalVolume * GetVolume(const std::string &name)
void update(const BeginOfJob *)
This routine will be called when the appropriate signal arrives.
static int position[264][3]
Definition: ReadPGInfo.cc:509
step
static const double X0
bool isSelected(const G4VTouchable *touch)
def move(src, dest)
Definition: eostools.py:510
bool isSelectedFast(const G4TouchableHistory *touch)
Definition: event.py:1
std::tuple< const G4VPhysicalVolume *, int > GetSensitiveVolume(const G4VTouchable *touchable)