|
DD4hep
1.30.0
Detector Description Toolkit for High Energy Physics
|
Go to the documentation of this file.
20 class G4LogicalVolume;
54 typedef std::vector<StepInfo*>
Steps;
67 virtual void operator()(
const G4Step* step, G4SteppingManager* mgr);
69 virtual void begin(
const G4Track* track);
71 virtual void end(
const G4Track* track);
93 #include <CLHEP/Units/SystemOfUnits.h>
94 #include <G4LogicalVolume.hh>
95 #include <G4Material.hh>
104 : pre(prePos), post(postPos), volume(vol)
110 : pre(c.pre), post(c.post), volume(c.volume)
143 std::string prePath = pre_handler.
path();
145 std::string postPath = post_handler.
path();
147 G4LogicalVolume* logVol = h.
logvol(h.
pre);
153 for_each(
m_steps.begin(),
m_steps.end(),detail::DestroyObject<StepInfo*>());
162 printP2(
"Starting tracking action for track ID=%d",track->GetTrackID());
163 for_each(
m_steps.begin(),
m_steps.end(),detail::DestroyObject<StepInfo*>());
173 constexpr
const char* line =
" +--------------------------------------------------------------------------------------------------------------------------------------------------\n";
174 constexpr
const char* fmt1 =
" | %5d %-20s %3.0f %8.3f %8.4f %11.4f %11.4f %10.3f %8.2f %11.6f %11.6f (%7.2f,%7.2f,%7.2f)\n";
175 constexpr
const char* fmt2 =
" | %5d %-20s %3.0f %8.3f %8.4f %11.6g %11.6g %10.3f %8.2f %11.6f %11.6f (%7.2f,%7.2f,%7.2f)\n";
179 ::printf(
"%s + Material scan between: x_0 = (%7.2f,%7.2f,%7.2f) [cm] and x_1 = (%7.2f,%7.2f,%7.2f) [cm] TrackID:%d: \n%s",
180 line,pre.X()/CLHEP::cm,pre.Y()/CLHEP::cm,pre.Z()/CLHEP::cm,post.X()/CLHEP::cm,post.Y()/CLHEP::cm,post.Z()/CLHEP::cm,track->GetTrackID(),line);
181 ::printf(
" | \\ %-11s Atomic Radiation Interaction Path Integrated Integrated Material\n",
"Material");
182 ::printf(
" | Num. \\ %-11s Number/Z Mass/A Density Length Length Thickness Length X0 Lambda Endpoint \n",
"Name");
183 ::printf(
" | Layer \\ %-11s [g/mole] [g/cm3] [cm] [cm] [cm] [cm] [cm] [cm] ( cm, cm, cm)\n",
"");
186 for(Steps::const_iterator i=
m_steps.begin(); i!=
m_steps.end(); ++i, ++count) {
187 const G4LogicalVolume* logVol = (*i)->volume;
188 G4Material* material = logVol->GetMaterial();
190 const Position& postPos = (*i)->post;
191 Position direction = postPos - prePos;
192 double length = direction.R()/CLHEP::cm;
193 double intLen = material->GetNuclearInterLength()/CLHEP::cm;
194 double radLen = material->GetRadlen()/CLHEP::cm;
195 double density = material->GetDensity()/(CLHEP::gram/CLHEP::cm3);
196 double nLambda = length / intLen;
197 double nx0 = length / radLen;
200 const char* fmt = radLen >= 1e5 ? fmt2 : fmt1;
201 const double* fractions = material->GetFractionVector();
202 for(
size_t j=0; j<material->GetNumberOfElements(); ++j) {
203 Zeff += fractions[j]*(material->GetElement(j)->GetZ());
204 Aeff += fractions[j]*(material->GetElement(j)->GetA())/CLHEP::gram;
209 ::printf(fmt,count,material->GetName().c_str(),
210 Zeff, Aeff, density, radLen, intLen, length,
212 postPos.X()/CLHEP::cm,postPos.Y()/CLHEP::cm,postPos.Z()/CLHEP::cm);
215 for_each(
m_steps.begin(),
m_steps.end(),detail::DestroyObject<StepInfo*>());
void printP2(const char *fmt,...) const
Support for messages with variable output level using output level+2.
~StepInfo()
Default destructor.
Helper class to ease the extraction of information from a G4Touchable object.
bool m_needsControl
Default property: Flag to create control instance.
Position prePos() const
Returns the pre-step position.
Concrete implementation of the Geant4 stepping action sequence.
Geant4MaterialScanner(Geant4Context *context, const std::string &name)
Standard constructor.
const G4LogicalVolume * volume
Reference to the logical volue.
static void increment(T *)
Increment count according to type information.
Structure to hold the information of one simulation step.
void callAtBegin(Q *p, void(T::*f)(const G4Track *), CallbackSequence::Location where=CallbackSequence::END)
Register Pre-track action callback.
#define DECLARE_GEANT4ACTION(name)
Plugin defintion to create Geant4Action objects.
Class to perform directional material scans using Geantinos.
Helper class to ease the extraction of information from a G4Step object.
virtual void end(const G4Track *track)
End-of-tracking callback.
Geant4EventActionSequence & eventAction() const
Access to the main event action sequence from the kernel object.
void callAtBegin(Q *p, void(T::*f)(const G4Event *))
Register begin-of-event callback.
void beginEvent(const G4Event *event)
Registered callback on Begin-event.
static void decrement(T *)
Decrement count according to type information.
const std::string & name() const
Access name of the action.
virtual void operator()(const G4Step *step, G4SteppingManager *mgr)
User stepping callback.
virtual ~Geant4MaterialScanner()
Default destructor.
Position postPos() const
Returns the post-step position.
virtual void begin(const G4Track *track)
Begin-of-tracking callback.
std::vector< StepInfo * > Steps
Position pre
Pre-step and Post-step position.
ROOT::Math::XYZVector Position
StepInfo & operator=(const StepInfo &c)
Assignment operator.
Namespace for the Geant4 based simulation part of the AIDA detector description toolkit.
std::string path() const
Helper: Access the placement path of a Geant4 touchable object as a string.
Namespace for the AIDA detector description toolkit.
StepInfo(const Position &pre, const Position &post, const G4LogicalVolume *volume)
Initializing constructor.
Generic context to extend user, run and event information.
G4LogicalVolume * logvol(const G4StepPoint *p) const
Geant4Context * context() const
Access the context.
Geant4TrackingActionSequence & trackingAction() const
Access to the main tracking action sequence from the kernel object.
void callAtEnd(Q *p, void(T::*f)(const G4Track *), CallbackSequence::Location where=CallbackSequence::END)
Register Post-track action callback.