24 #include "G4VSolid.hh"
25 #include "G4LogicalVolumeStore.hh"
26 #include "G4TouchableHistory.hh"
27 #include "G4VPhysicalVolume.hh"
28 #include "G4AffineTransform.hh"
35 #define DEBUG_G4_VOLUMES
37 using namespace CLHEP;
42 const G4LogicalVolumeStore* lvs = G4LogicalVolumeStore::GetInstance();
44 #ifdef DEBUG_G4_VOLUMES
45 for (G4LogicalVolumeStore::const_iterator volume = lvs->begin(); volume != lvs->end(); ++volume)
46 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer: G4 registered volumes " << (*volume)->GetName()
50 for (G4LogicalVolumeStore::const_iterator volume = lvs->begin(); volume != lvs->end(); ++volume) {
58 static inline const G4AffineTransform&
GetTransform(
const G4TouchableHistory* touchable,
int depth) {
59 return touchable->GetHistory()->GetTransform(touchable->GetHistory()->GetDepth() -
depth);
69 int depth = touchable->GetHistoryDepth();
71 const G4VPhysicalVolume* volume = touchable->GetVolume(
level);
72 if (volume->GetLogicalVolume()->GetSensitiveDetector() !=
nullptr) {
73 return std::make_tuple(volume,
level);
76 return std::tuple<const G4VPhysicalVolume*, int>(
nullptr, 0);
82 m_selectedNames =
config.getParameter<std::vector<std::string> >(
"SelectedVolumes");
83 m_primaryTracks =
config.getParameter<
bool>(
"PrimaryTracksOnly");
85 m_hgcalzfront =
config.getParameter<
double>(
"hgcalzfront");
88 produces<std::vector<MaterialAccountingTrack> >();
89 output_file_ =
new TFile(
"radLen_vs_eta_fromProducer.root",
"RECREATE");
91 radLen_vs_eta_ =
new TProfile(
"radLen",
"radLen", 250., -5., 5., 0, 10.);
95 if (
std::find(m_selectedNames.begin(), m_selectedNames.end(),
"HGCal") != m_selectedNames.end()) {
100 if (
std::find(m_selectedNames.begin(), m_selectedNames.end(),
"HFNose") != m_selectedNames.end()) {
104 outVolumeZpositionTxt.open(m_txtOutFile.c_str(),
std::ios::out);
113 radLen_vs_eta_->Write();
114 output_file_->Close();
119 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer: List of the selected volumes: " << std::endl;
120 for (std::vector<std::string>::const_iterator volume_name = m_selectedNames.begin();
121 volume_name != m_selectedNames.end();
123 const G4LogicalVolume* volume =
GetVolume(*volume_name);
125 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer: " << *volume_name << std::endl;
126 m_selectedVolumes.push_back(volume);
129 std::cerr <<
"TrackingMaterialProducer::update(const BeginOfJob*): WARNING: selected volume \"" << *volume_name
130 <<
"\" not found in geometry " << std::endl;
137 m_tracks =
new std::vector<MaterialAccountingTrack>();
145 G4Track*
track = const_cast<G4Track*>((*
event)());
146 if (m_primaryTracks and
track->GetParentID() != 0) {
147 track->SetTrackStatus(fStopAndKill);
157 if (
isHGCal &&
track->GetTrackStatus() != fStopAndKill && fabs(
track->GetMomentum().eta()) > outerHGCalEta &&
158 fabs(
track->GetMomentum().eta()) < innerHGCalEta) {
159 if (
track->GetMomentum().eta() > 0.) {
160 outVolumeZpositionTxt <<
"StainlessSteel " << m_hgcalzfront <<
" " << 0 <<
" " << 0 <<
" " << 0 <<
" " << 0
162 }
else if (
track->GetMomentum().eta() <= 0.) {
163 outVolumeZpositionTxt <<
"StainlessSteel " << -m_hgcalzfront <<
" " << 0 <<
" " << 0 <<
" " << 0 <<
" " << 0
172 if (
isHFNose &&
track->GetTrackStatus() != fStopAndKill && fabs(
track->GetMomentum().eta()) > outerHFnoseEta &&
173 fabs(
track->GetMomentum().eta()) < innerHFnoseEta) {
174 if (
track->GetMomentum().eta() > 0.) {
175 outVolumeZpositionTxt <<
"Polyethylene " << m_hgcalzfront <<
" " << 0 <<
" " << 0 <<
" " << 0 <<
" " << 0
177 }
else if (
track->GetMomentum().eta() <= 0.) {
178 outVolumeZpositionTxt <<
"Polyethylene " << -m_hgcalzfront <<
" " << 0 <<
" " << 0 <<
" " << 0 <<
" " << 0
185 for (
int d = touchable->GetHistoryDepth() - 1;
d >= 0; --
d) {
186 if (
std::find(m_selectedNames.begin(), m_selectedNames.end(), touchable->GetVolume(
d)->GetName()) !=
187 m_selectedNames.end())
195 const G4TouchableHistory* touchable = static_cast<const G4TouchableHistory*>(
step->GetTrack()->GetTouchable());
196 if (not isSelectedFast(touchable)) {
197 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer:\t[...] skipping "
198 << touchable->GetVolume()->GetName() << std::endl;
203 const G4Material* material = touchable->GetVolume()->GetLogicalVolume()->GetMaterial();
204 double length =
step->GetStepLength() / cm;
205 double X0 = material->GetRadlen() / cm;
206 double Ne = material->GetElectronDensity() * cm3;
207 double Xi = Ne / 6.0221415e23 * 0.307075 / 2;
208 double radiationLengths = length /
X0;
212 G4ThreeVector globalPosPre =
step->GetPreStepPoint()->GetPosition();
213 G4ThreeVector globalPosPost =
step->GetPostStepPoint()->GetPosition();
214 GlobalPoint globalPositionIn(globalPosPre.x() / cm, globalPosPre.y() / cm, globalPosPre.z() / cm);
215 GlobalPoint globalPositionOut(globalPosPost.x() / cm, globalPosPost.y() / cm, globalPosPost.z() / cm);
217 G4StepPoint* prePoint =
step->GetPreStepPoint();
218 G4StepPoint* postPoint =
step->GetPostStepPoint();
219 const CLHEP::Hep3Vector& postPos = postPoint->GetPosition();
225 if (postPoint->GetStepStatus() == fGeomBoundary && fabs(postPoint->GetMomentum().eta()) > outerHGCalEta &&
226 fabs(postPoint->GetMomentum().eta()) < innerHGCalEta) {
231 outVolumeZpositionTxt << prePoint->GetMaterial()->GetName() <<
" " << postPos.
z() <<
" "
232 << postPoint->GetMomentum().eta() <<
" "
233 <<
sqrt(postPos.x() * postPos.x() + postPos.y() * postPos.y()) <<
" " << 0 <<
" " << 0
237 if (postPoint->GetStepStatus() == fGeomBoundary && fabs(postPoint->GetMomentum().eta()) > outerHFnoseEta &&
238 fabs(postPoint->GetMomentum().eta()) < innerHFnoseEta) {
239 outVolumeZpositionTxt << prePoint->GetMaterial()->GetName() <<
" " << postPos.z() <<
" "
240 << postPoint->GetMomentum().eta() <<
" "
241 <<
sqrt(postPos.x() * postPos.x() + postPos.y() * postPos.y()) <<
" " << 0 <<
" " << 0
248 bool enter_sensitive =
false;
249 bool leave_sensitive =
false;
250 double cosThetaPre = 0.0;
251 double cosThetaPost = 0.0;
253 const G4VPhysicalVolume* sensitive =
nullptr;
257 const G4VSolid& solid = *touchable->GetSolid(
level);
259 G4ThreeVector
pos =
transform.Inverse().TransformPoint(G4ThreeVector(0., 0., 0.));
262 G4ThreeVector localPosPre =
transform.TransformPoint(globalPosPre);
263 EInside statusPre = solid.Inside(localPosPre);
264 if (statusPre == kSurface) {
265 enter_sensitive =
true;
266 G4ThreeVector globalDirPre =
step->GetPreStepPoint()->GetMomentumDirection();
267 G4ThreeVector localDirPre =
transform.TransformAxis(globalDirPre);
268 G4ThreeVector normalPre = solid.SurfaceNormal(localPosPre);
269 cosThetaPre = normalPre.cosTheta(-localDirPre);
272 G4ThreeVector localPosPost =
transform.TransformPoint(globalPosPost);
273 EInside statusPost = solid.Inside(localPosPost);
274 if (statusPost == kSurface) {
275 leave_sensitive =
true;
276 G4ThreeVector globalDirPost =
step->GetPostStepPoint()->GetMomentumDirection();
277 G4ThreeVector localDirPost =
transform.TransformAxis(globalDirPost);
278 G4ThreeVector normalPost = solid.SurfaceNormal(localPosPost);
279 cosThetaPost = normalPost.cosTheta(localDirPost);
285 m_track.enterDetector(sensitive,
position, cosThetaPre);
288 m_track.leaveDetector(sensitive, cosThetaPost);
291 LogInfo(
"TrackingMaterialProducer") <<
"Track was near sensitive volume " << sensitive->GetName() << std::endl;
293 LogInfo(
"TrackingMaterialProducer") <<
"Track was near non-sensitive volume " << touchable->GetVolume()->GetName()
295 LogInfo(
"TrackingMaterialProducer") <<
"Step length: " << length <<
" cm\n"
296 <<
"globalPreStep(r,z): (" << globalPositionIn.perp() <<
", "
297 << globalPositionIn.z() <<
") cm\n"
298 <<
"globalPostStep(r,z): (" << globalPositionOut.perp() <<
", "
299 << globalPositionOut.z() <<
") cm\n"
300 <<
"position(r,z): (" <<
position.perp() <<
", " <<
position.z() <<
") cm\n"
301 <<
"Radiation lengths: " << radiationLengths <<
" \t\t(X0: " <<
X0
303 <<
"Energy loss: " <<
energyLoss <<
" MeV \t(Xi: " << Xi
305 <<
"Track was " << (enter_sensitive ?
"entering " :
"in none ")
306 <<
"sensitive volume\n"
307 <<
"Track was " << (leave_sensitive ?
"leaving " :
"in none ")
308 <<
"sensitive volume\n";
313 const G4Track*
track = (*event)();
314 if (m_primaryTracks and
track->GetParentID() != 0)
317 radLen_vs_eta_->Fill(
track->GetMomentum().eta(), m_track.summary().radiationLengths());
318 m_tracks->push_back(m_track);
321 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer: this track took " << m_track.steps().size()
322 <<
" steps, and passed through " << m_track.detectors().size()
323 <<
" sensitive detectors" << std::endl;
324 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer: track length: " << m_track.summary().length()
325 <<
" cm" << std::endl;
326 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer: radiation lengths: "
327 << m_track.summary().radiationLengths() << std::endl;
328 LogInfo(
"TrackingMaterialProducer") <<
"TrackingMaterialProducer: energy loss: "
329 << m_track.summary().energyLoss() <<
" MeV" << std::endl;
335 std::unique_ptr<std::vector<MaterialAccountingTrack> >
tracks(m_tracks);
342 for (
size_t i = 0;
i < m_selectedVolumes.size(); ++
i)
343 if (m_selectedVolumes[
i]->IsAncestor(touchable->GetVolume())
or
344 m_selectedVolumes[
i] == touchable->GetVolume()->GetLogicalVolume())