ATLAS Offline Software
MaterialStepRecorder.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 // MaterialStepRecorder.cxx, (c) ATLAS Detector software
8 
11 #include "GaudiKernel/IDataProviderSvc.h"
12 #include "GaudiKernel/INTupleSvc.h"
13 #include "GaudiKernel/NTuple.h"
14 #include "GaudiKernel/SmartDataPtr.h"
15 #include "GaudiKernel/MsgStream.h"
16 #include "GaudiKernel/Bootstrap.h"
17 #include "GaudiKernel/ISvcLocator.h"
18 #include "GaudiKernel/IMessageSvc.h"
20 #include "StoreGate/StoreGateSvc.h"
21 #include "G4Step.hh"
22 #include "G4Material.hh"
23 #include "G4Element.hh"
24 #include "G4StepPoint.hh"
25 #include "G4TouchableHistory.hh"
26 #include "G4LogicalVolume.hh"
27 #include <climits>
28 #include <cmath>
29 
30 namespace G4UA
31 {
32 
34  AthMessaging(Gaudi::svcLocator()->service< IMessageSvc >( "MessageSvc" ),"MaterialStepRecorder"),
35  m_evtStore("StoreGateSvc/StoreGateSvc","MaterialStepRecorder"), //FIXME name should be passed in via a Config struct rather than hardcoded.
36  m_matStepCollection(nullptr),
37  m_matStepCollectionName("MaterialStepRecords"), //FIXME should be passed in via a Config struct rather than hardcoded.
38  m_recordComposition(true), //FIXME should be passed in via a Config struct rather than hardcoded.
39  m_totalSteps(0),
40  m_eventID(0),
41  m_elementTable(nullptr),
42  m_elementTableName("ElementTable"),
43  m_runElementTable(nullptr)
44  {
45  }
46 
48  {
49  ATH_MSG_DEBUG(" BeginOfEventAction");
50 
51  // create a new Collection
53  // m_eventStepLength = 0;
54  }
55 
57  {
58  ATH_MSG_DEBUG(" EndOfEventAction");
59 
60  ++m_eventID;
61  // write the Collection to StoreGate
62  if(m_evtStore->record(m_matStepCollection, m_matStepCollectionName, false).isFailure())
63  {
64  ATH_MSG_ERROR("cannot record step collection to StoreGate");
65  }
66  // write out the ElementTable of this event
67  if ( m_evtStore->record(m_elementTable, m_elementTableName, false).isFailure() )
68  {
69  ATH_MSG_ERROR("EndOfEventAction : recording ElementTable in StoreGate FAILED");
70  delete m_elementTable;
71  }
72  m_elementTable = nullptr;
73  }
74 
76  {
77  ATH_MSG_DEBUG(" BeginOfRunAction");
78 
79  // initialize
80  m_totalSteps = 0;
81  m_eventID = 0;
82  }
83 
84  void MaterialStepRecorder::UserSteppingAction(const G4Step* aStep)
85  {
86  // kill secondaries
87  if (aStep->GetTrack()->GetParentID()) {
88  aStep->GetTrack()->SetTrackStatus(fStopAndKill);
89  return;
90  }
91 
92  // ElementTable preparation
96  }
97 
98  // the material information
99  const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
100  // G4LogicalVolume
101  const G4LogicalVolume *lv= touchHist ? touchHist->GetVolume()->GetLogicalVolume() : nullptr;
102  const G4Material *mat = lv ? lv->GetMaterial() : nullptr;
103 
104  std::vector<unsigned char> elements;
105  std::vector<unsigned char> fractions;
106 
107  // log the information // cut off air
108  if (mat && mat->GetRadlen() < 200000.) {
109 
110  ++m_totalSteps;
111 
112  // the step information
113  double steplength = aStep->GetStepLength();
114 
115  // the position information
116  G4ThreeVector pos = aStep->GetPreStepPoint()->GetPosition();
117 
118  double X0 = mat->GetRadlen();
119  double L0 = mat->GetNuclearInterLength();
120  double A = 0.;
121  double Z = 0.;
122  // double totAtoms = 0.;
123 
124  double rho = mat->GetDensity()*CLHEP::mm3/CLHEP::gram;
125 
126  // get the number of Elements
127  size_t elNumber = mat->GetNumberOfElements();
128  const G4ElementVector* elVector = mat->GetElementVector();
129  double totAtoms = mat->GetTotNbOfAtomsPerVolume();
130 
131  // reserve the right number of elements
132  elements.reserve(elNumber);
133  fractions.reserve(elNumber);
134 
135  if (1 == elNumber) {
136 
137  A = mat->GetA()/CLHEP::gram;
138  Z = mat->GetZ();
139 
140  unsigned int Zint = (unsigned int)Z;
141  // the element and fraction vector
142  elements.push_back(static_cast<unsigned char>(Z));
143  fractions.push_back(UCHAR_MAX);
144  // record the Element
146  // the element Material
147  Trk::Material elMat(X0,L0,A,Z,rho);
148  G4String elName = (*elVector)[0]->GetName();
149  // add it to the table
150  m_elementTable->addElement(elMat, elName);
151  m_runElementTable->addElement(elMat,elName);
152  }
153 
154 
155  } else {
156 
157  const G4double* atVector = mat->GetVecNbOfAtomsPerVolume();
158  double totalFrac = 0.;
159 
160  for (size_t iel = 0; iel < elNumber; ++iel) {
161 
162  const G4Element* currentEl = (*elVector)[iel];
163  double currentNum = atVector ? atVector[iel] : 1.;
164  double relNbAtoms = currentNum/totAtoms;
165 
166  double currentZ = currentEl->GetZ();
167 
168  A += relNbAtoms*currentEl->GetA()/CLHEP::gram;
169  Z += relNbAtoms*currentEl->GetZ();
170  unsigned int Zint = (unsigned int)(currentEl->GetZ());
171 
172  // the element and fraction vector
173  elements.push_back(int(currentZ));
174  // calculate the fraction with a accuracy of 1/256.
175  unsigned int relNbAtomsChar = (unsigned int)(relNbAtoms*(UCHAR_MAX));
176  relNbAtomsChar = relNbAtomsChar > UCHAR_MAX ? UCHAR_MAX : relNbAtomsChar;
177 
178  // smaller components than 0.5 % are automatically ignored
179  totalFrac += relNbAtoms;
180  if (relNbAtomsChar) {
181  fractions.push_back(relNbAtomsChar);
182  // record composition
184  double curA = currentEl->GetA()/CLHEP::gram;
185  double curZ = currentEl->GetZ();
186  // approximate formulas for X0 and L0
187  // X0 from : PH-EP-Tech-Note-2010-013 g/cm3 -> g/mm3
188  // L0 from :
189  double curX0 = 1432.8*curA/(curZ*(curZ+1)*(11.319-std::log(curZ)));
190  double curL0 = 0.;
191  double curRho = rho*relNbAtoms;
192  Trk::Material elMat(curX0,curL0,curA,curZ,curRho);
193  const G4String& elName = currentEl->GetName();
194  // add it to the table
195  m_elementTable->addElement(elMat, elName);
196  m_runElementTable->addElement(elMat,elName);
197  }
198  }
199  }
200  if ((totalFrac-1.)*(totalFrac-1.) > 10e-4 )
201  ATH_MSG_DEBUG("Total fractions do not add up to one at INPUT (" << totalFrac << ") !");
202  }
203 
204  // is it a Geantino?
205  if (aStep->GetTrack()->GetParticleDefinition()->GetPDGEncoding() == 0) {
207  m_matStepCollection->push_back(new Trk::MaterialStep(pos.x(), pos.y(), pos.z(), steplength, X0, L0, A, Z, rho, elements, fractions));
208  else
209  m_matStepCollection->push_back(new Trk::MaterialStep(pos.x(), pos.y(), pos.z(), steplength, X0, L0, A, Z, rho));
210  }
211  }
212  }
213 
214 } // namespace G4UA
G4UA::MaterialStepRecorder::m_eventID
size_t m_eventID
Definition: MaterialStepRecorder.h:59
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
G4UA
for nSW
Definition: CalibrationDefaultProcessing.h:19
G4UA::MaterialStepRecorder::m_matStepCollectionName
std::string m_matStepCollectionName
Definition: MaterialStepRecorder.h:54
mat
GeoMaterial * mat
Definition: LArDetectorConstructionTBEC.cxx:55
Trk::ElementTable::contains
bool contains(unsigned int Z) const
quick check
Definition: ElementTable.h:93
ServiceAccessor.h
G4UA::MaterialStepRecorder::MaterialStepRecorder
MaterialStepRecorder()
Definition: MaterialStepRecorder.cxx:33
Trk::MaterialStepCollection
DataVector< Trk::MaterialStep > MaterialStepCollection
Definition: MaterialStepCollection.h:17
python.SystemOfUnits.gram
int gram
Definition: SystemOfUnits.py:165
G4UA::MaterialStepRecorder::m_matStepCollection
Trk::MaterialStepCollection * m_matStepCollection
Definition: MaterialStepRecorder.h:53
G4UA::MaterialStepRecorder::UserSteppingAction
virtual void UserSteppingAction(const G4Step *) override
Definition: MaterialStepRecorder.cxx:84
python.SystemOfUnits.mm3
int mm3
Definition: SystemOfUnits.py:85
pdg_comparison.X0
X0
Definition: pdg_comparison.py:314
Trk::MaterialStep
Definition: MaterialStep.h:34
A
ATH_MSG_ERROR
#define ATH_MSG_ERROR(x)
Definition: AthMsgStreamMacros.h:33
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
G4UA::MaterialStepRecorder::m_recordComposition
bool m_recordComposition
Definition: MaterialStepRecorder.h:56
AthMessaging
Class to provide easy MsgStream access and capabilities.
Definition: AthMessaging.h:55
G4UA::MaterialStepRecorder::m_elementTable
Trk::ElementTable * m_elementTable
Definition: MaterialStepRecorder.h:61
G4UA::MaterialStepRecorder::m_runElementTable
Trk::ElementTable * m_runElementTable
Definition: MaterialStepRecorder.h:64
MaterialStepRecorder.h
DataVector::push_back
value_type push_back(value_type pElem)
Add an element to the end of the collection.
python.LumiBlobConversion.pos
pos
Definition: LumiBlobConversion.py:18
G4UA::MaterialStepRecorder::EndOfEventAction
virtual void EndOfEventAction(const G4Event *) override
Definition: MaterialStepRecorder.cxx:56
G4UA::MaterialStepRecorder::m_evtStore
ServiceHandle< StoreGateSvc > m_evtStore
Pointer to StoreGate (event store by default)
Definition: MaterialStepRecorder.h:51
Trk::ElementTable
Definition: ElementTable.h:24
G4UA::MaterialStepRecorder::m_elementTableName
std::string m_elementTableName
Definition: MaterialStepRecorder.h:62
G4UA::MaterialStepRecorder::BeginOfRunAction
virtual void BeginOfRunAction(const G4Run *) override
Definition: MaterialStepRecorder.cxx:75
python.CaloCondTools.log
log
Definition: CaloCondTools.py:20
Gaudi
=============================================================================
Definition: CaloGPUClusterAndCellDataMonitorOptions.h:273
Trk::Material
Definition: Material.h:116
G4UA::MaterialStepRecorder::BeginOfEventAction
virtual void BeginOfEventAction(const G4Event *) override
Definition: MaterialStepRecorder.cxx:47
StoreGateSvc.h
fitman.rho
rho
Definition: fitman.py:532
Trk::L0
@ L0
Definition: AlignModuleList.h:32
G4UA::MaterialStepRecorder::m_totalSteps
size_t m_totalSteps
Definition: MaterialStepRecorder.h:58
Trk::ElementTable::addElement
void addElement(const Material &mat, const std::string &mname="")
Add material to the Table - if the elment is already filled ignore.
Definition: ElementTable.h:75
MaterialStep.h