9#include <simulation/modules/MaterialScan.h> 
   10#include <framework/core/ModuleParam.templateDetails.h> 
   11#include <framework/gearbox/Unit.h> 
   12#include <framework/utilities/Utils.h> 
   21#include <G4EventManager.hh> 
   22#include <G4UserEventAction.hh> 
   23#include <G4UserStackingAction.hh> 
   24#include <G4UserTrackingAction.hh> 
   25#include <G4UserSteppingAction.hh> 
   26#include <G4RayShooter.hh> 
   29#include <boost/format.hpp> 
   30#include <boost/algorithm/string.hpp> 
   46  double stlen = step->GetStepLength();
 
   47  G4StepPoint* preStepPoint = step->GetPreStepPoint();
 
   48  G4Region* region = preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
 
   56      B2ERROR(
"Track is stuck at " << preStepPoint->GetPosition() << 
" in volume '" 
   57              << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
 
   58              << 
" (" << region->GetName() << 
"): " 
   59              << 
m_zeroSteps << 
" consecutive steps with length less then " 
   61      step->GetTrack()->SetTrackStatus(fStopAndKill);
 
   63      B2WARNING(
"Track is stuck at " << preStepPoint->GetPosition() << 
" in volume '" 
   64                << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
 
   65                << 
" (" << region->GetName() << 
"): " 
   66                << 
m_zeroSteps << 
" consecutive steps with length less then " 
   68      G4ThreeVector pos = step->GetTrack()->GetPosition();
 
   69      G4ThreeVector dir = step->GetTrack()->GetMomentumDirection();
 
  119  getRay(origin, direction);
 
  127  std::unique_ptr<TH2D>& hist = 
m_regions[name];
 
  131    hist.reset(
new TH2D(name.c_str(), (name + 
";" + 
m_axisLabel).c_str(),
 
  147  G4StepPoint* preStepPoint = step->GetPreStepPoint();
 
  148  G4Region* region = preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
 
  149  G4Material* material = preStepPoint->GetMaterial();
 
  150  double stlen = step->GetStepLength();
 
  151  double x0 = stlen / (material->GetRadlen());
 
  152  double lambda = stlen / (material->GetNuclearInterLength());
 
  153  B2DEBUG(20, 
"Step in at " << preStepPoint->GetPosition() << 
" in volume '" 
  154          << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
 
  155          << 
" (" << region->GetName() << 
")" 
  156          << 
" with length=" << stlen << 
" mm");
 
  164      G4Track* track = step->GetTrack();
 
  165      track->SetTrackStatus(fStopAndKill);
 
  167      x0 *= remaining / stlen;
 
  168      lambda *= remaining / stlen;
 
  177  string x0_total = 
"All_Regions_x0";
 
  178  string lambda_total = 
"All_Regions_lambda";
 
  179  string x0_name = region->GetName() + 
"_x0";
 
  180  string lambda_name = region->GetName() + 
"_lambda";
 
  183    x0_total = 
"All_Materials_x0";
 
  184    lambda_total = 
"All_Materials_lambda";
 
  185    x0_name = material->GetName() + 
"_x0";
 
  186    lambda_name = material->GetName() + 
"_lambda";
 
  203  direction.set(sin(theta) * cos(phi),
 
  204                sin(theta) * sin(phi),
 
  228    B2INFO(
"Set maximum depth to " << 
m_maxDepth);
 
  243    double theta = acos(gRandom->Uniform(cos(
m_opening / 2.), 1));
 
  244    double phi = gRandom->Uniform(0, 2 * M_PI);
 
  245    auto orthogonal = direction.orthogonal();
 
  246    orthogonal.rotate(phi, direction);
 
  247    direction.rotate(theta, orthogonal);
 
  249  B2DEBUG(10, 
"Create Ray " << 
m_curRay << 
" from " << origin << 
" direction " <<
 
  250          direction << 
" [theta=" << (direction.theta() / M_PI * 180) << 
", phi=" << (direction.phi() / M_PI * 180) << 
"]");
 
  256  std::unique_ptr<TH1D>& hist = 
m_regions[name];
 
  274  const double dx = steplength / samples;
 
  275  const double dv = value / samples / 
m_count / h->GetBinWidth(1);
 
  276  for (
int i = 0; i < samples; ++i) {
 
  285  G4StepPoint* preStepPoint = step->GetPreStepPoint();
 
  286  G4Region* region = preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
 
  287  G4Material* material = preStepPoint->GetMaterial();
 
  289  double stlen = step->GetStepLength();
 
  290  double x0 = stlen / (material->GetRadlen());
 
  291  double lambda = stlen / (material->GetNuclearInterLength());
 
  292  B2DEBUG(20, 
"Step in at " << preStepPoint->GetPosition() << 
" in volume '" 
  293          << preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetName()
 
  294          << 
" (" << region->GetName() << 
")" 
  295          << 
" with length=" << stlen << 
" mm");
 
  302      G4Track* track = step->GetTrack();
 
  303      track->SetTrackStatus(fStopAndKill);
 
  305      x0 *= remaining / stlen;
 
  306      lambda *= remaining / stlen;
 
  313    string x0_total = 
"All_Regions_x0";
 
  314    string lambda_total = 
"All_Regions_lambda";
 
  315    string x0_name = region->GetName() + 
"_x0";
 
  316    string lambda_name = region->GetName() + 
"_lambda";
 
  319      x0_total = 
"All_Materials_x0";
 
  320      lambda_total = 
"All_Materials_lambda";
 
  321      x0_name = material->GetName() + 
"_x0";
 
  322      lambda_name = material->GetName() + 
"_lambda";
 
  336MaterialScanModule::MaterialScanModule(): m_rootFile(0), m_sphericalOrigin(3, 0), m_doCosTheta(false)
 
  339  setDescription(
"This Module is intended to scan the material budget of the " 
  340                 "geometry.  Currently, there are two different kinds of scans " 
  341                 "available: Spherical and Planar scan. \n" 
  342                 "Spherical scan will shoot rays from the origin of the " 
  343                 "detector and scan along the polar and azimuth angle.\n" 
  344                 "Planar scan will shoot rays perpendicular to a given " 
  348  m_spherical.ignoredMaterials.push_back(
"Vacuum");
 
  349  m_spherical.ignoredMaterials.push_back(
"Air");
 
  350  m_spherical.ignoredMaterials.push_back(
"G4_AIR");
 
  351  m_planar.ignoredMaterials = m_spherical.ignoredMaterials;
 
  352  m_rayIgnoredMaterials = m_spherical.ignoredMaterials;
 
  354  addParam(
"Filename",            m_filename,
 
  355           "The filename where the material scan will be stored",
 
  356           string(
"MaterialScan.root"));
 
  357  addParam(
"spherical",           m_doSpherical,
 
  358           "Do a spherical scan, that is shooting rays from the origin with " 
  359           "varying angles", 
true);
 
  360  addParam(
"spherical.origin",    m_sphericalOrigin,
 
  361           "Origin for the spherical scan", m_sphericalOrigin);
 
  362  addParam(
"spherical.nTheta",    m_spherical.nU,
 
  363           "Number of rays in theta", 200);
 
  364  addParam(
"spherical.minTheta",  m_spherical.minU,
 
  365           "Theta start angle", 17.);
 
  366  addParam(
"spherical.maxTheta",  m_spherical.maxU,
 
  367           "Theta stop angle", 150.);
 
  368  addParam(
"spherical.nPhi",      m_spherical.nV,
 
  369           "Number of rays in phi", 200);
 
  370  addParam(
"spherical.minPhi",    m_spherical.minV,
 
  371           "Phi start angle", 0.);
 
  372  addParam(
"spherical.maxPhi",    m_spherical.maxV,
 
  373           "Phi stop angle", 360.);
 
  374  addParam(
"spherical.maxDepth",  m_spherical.maxDepth,
 
  375           "Maximum scan depth in cm. The ray will be killed after having " 
  376           "reached the maximum Depth. <=0 means no Limit.", -1.0);
 
  377  addParam(
"spherical.ignored",      m_spherical.ignoredMaterials,
 
  378           "Names of Materials which should be ignored when doing the scan", m_spherical.ignoredMaterials);
 
  379  addParam(
"spherical.splitByMaterials", m_spherical.splitByMaterials,
 
  380           "If True, split output by material names instead of by regions", 
false);
 
  381  addParam(
"spherical.cosTheta", m_doCosTheta,
 
  382           "If True, perform the spherical scan uniform in cos(theta) instead of theta", 
false);
 
  384  addParam(
"planar",              m_doPlanar,
 
  385           "Do a plane scan, that is shooting parallel rays from a defined " 
  387  addParam(
"planar.plane",        m_planeName,
 
  388           "Plane to use for scanning, available are all two letter " 
  389           "combinations of X,Y and Z like XY and XZ and custom", 
string(
"ZX"));
 
  390  addParam(
"planar.nU",           m_planar.nU,
 
  391           "Number of rays in U", 200);
 
  392  addParam(
"planar.minU",         m_planar.minU,
 
  394  addParam(
"planar.maxU",         m_planar.maxU,
 
  396  addParam(
"planar.nV",           m_planar.nV,
 
  397           "Number of rays in V", 200);
 
  398  addParam(
"planar.minV",         m_planar.minV,
 
  400  addParam(
"planar.maxV",         m_planar.maxV,
 
  402  addParam(
"planar.maxDepth",     m_planar.maxDepth,
 
  403           "Maximum scan depth in cm. The ray will be killed after having " 
  404           "reached the maximum Depth. <=0 means no Limit.", -1.0);
 
  405  addParam(
"planar.custom",       m_customPlane,
 
  406           "Parameters of the plane when choosing custom. This is supposed to " 
  407           "be a list of 9 values: The first three are the coordinates of the " 
  408           "plane origin. The second three are the direction of U in the " 
  409           "r-phi coordinates. The last three are the direction of V in the " 
  410           "coordinates parallel to the detector axis (beamline). ", m_customPlane);
 
  411  addParam(
"planar.ignored",      m_planar.ignoredMaterials,
 
  412           "Names of Materials which should be ignored when doing the scan", m_planar.ignoredMaterials);
 
  413  addParam(
"planar.splitByMaterials", m_planar.splitByMaterials,
 
  414           "If True, split output by material names instead of by regions", 
false);
 
  416  addParam(
"ray", m_doRay, 
"Do a ray scan: Shoot from the given origin in a " 
  417           "given direction and record the material encountered along the way", 
false);
 
  418  addParam(
"ray.origin", m_rayOrigin, 
"Origin for the ray", m_rayOrigin);
 
  419  addParam(
"ray.theta", m_rayTheta, 
"Theta angle for the ray in degrees", m_rayTheta);
 
  420  addParam(
"ray.phi", m_rayPhi, 
"Phi angle of the ray in degrees", m_rayPhi);
 
  421  addParam(
"ray.direction", m_rayDirection, 
"Set the ray direction as vector " 
  422           "[x,y,z]. You cannot set this and theta/phi at the same time",
 
  424  addParam(
"ray.opening", m_rayOpening, 
"Opening angle for the ray in degrees. " 
  425           "If bigger than 0 then multiple rays (``ray.count``) will be shot " 
  426           "randomly inside the opening angle uniformly distributed in area " 
  427           "(cosine of the angle)", m_rayOpening);
 
  428  addParam(
"ray.count", m_rayCount, 
"Count of rays if the opening angle is >0",
 
  430  addParam(
"ray.sampleDepth", m_raySampleDepth, 
"Distance in which to sample the " 
  431           "material", m_raySampleDepth);
 
  432  addParam(
"ray.maxDepth", m_rayMaxDepth, 
"Max distance for the ray before it's " 
  433           "stopped, 0 for no limit. The actual depth can be smaller if the " 
  434           "simulation volume ends earlier", m_rayMaxDepth);
 
  435  addParam(
"ray.ignored", m_rayIgnoredMaterials, 
"Names of Materials which " 
  436           "should be ignored when doing the scan", m_rayIgnoredMaterials);
 
  443    B2ERROR(
"Cannot open rootfile for writing" << 
LogVar(
"rootfile", 
m_filename));
 
  452      B2ERROR(
"planar.custom: At least 9 values needed to define custom plane, only " << 
m_customPlane.size() << 
" provided");
 
  455    B2ERROR(
"planar.plane: Only custom or any two letter combinations of X,Y and Z are recognized");
 
  457    B2ERROR(
"planar.plane: " << 
m_planeName << 
" not valid, cannot use the same axis twice");
 
  462    B2ERROR(
"spherical.origin: Three values are needed to define a point, only " << 
m_sphericalOrigin.size() << 
" given.");
 
  470    if (getParam<double>(
"ray.theta").isSetInSteering() || getParam<double>(
"ray.phi").isSetInSteering()) {
 
  471      B2ERROR(
"Cannot set ray.theta/ray.phi and ray.direction at the same time, please only set one");
 
  480  G4EventManager* eventManager = G4EventManager::GetEventManager();
 
  483  G4UserEventAction*    vanillaEventAction    = eventManager->GetUserEventAction();
 
  484  G4UserStackingAction* vanillaStackingAction = eventManager->GetUserStackingAction();
 
  485  G4UserTrackingAction* vanillaTrackingAction = eventManager->GetUserTrackingAction();
 
  486  G4UserSteppingAction* vanillaSteppingAction = eventManager->GetUserSteppingAction();
 
  489  eventManager->SetUserAction((G4UserEventAction*)0);
 
  490  eventManager->SetUserAction((G4UserStackingAction*)0);
 
  491  eventManager->SetUserAction((G4UserTrackingAction*)0);
 
  496  vector<MaterialScanBase*> scans;
 
  502    G4ThreeVector origin(0, 0, 0);
 
  533    eventManager->SetUserAction(scan);
 
  535    G4RayShooter rayShooter;
 
  536    G4ThreeVector origin;
 
  537    G4ThreeVector direction;
 
  538    int maxRays = scan->getNRays();
 
  543    while (scan->createNext(origin, direction)) {
 
  544      G4Event* 
event = 
new G4Event(++curRay);
 
  545      rayShooter.Shoot(
event, origin, direction);
 
  546      eventManager->ProcessOneEvent(
event);
 
  550      int donePercent = 100 * curRay / maxRays;
 
  551      if (donePercent > lastPercent) {
 
  553        double eta = perRay * (maxRays - curRay);
 
  554        B2INFO(boost::format(
"%s Scan: %3d%%, %.3f ms per ray, ETA: %.2f seconds")
 
  555               % scan->getName() % donePercent
 
  557        lastPercent = donePercent;
 
  571  eventManager->SetUserAction(vanillaEventAction);
 
  572  eventManager->SetUserAction(vanillaStackingAction);
 
  573  eventManager->SetUserAction(vanillaTrackingAction);
 
  574  eventManager->SetUserAction(vanillaSteppingAction);
 
  580    case 'x': 
return G4ThreeVector(1, 0, 0);
 
  581    case 'y': 
return G4ThreeVector(0, 1, 0);
 
  582    case 'z': 
return G4ThreeVector(0, 0, 1);
 
  584      B2FATAL(
"Unknown axis: " << name);
 
  586  return G4ThreeVector(0, 0, 0);
 
void fillValue(const std::string &name, double value)
Fill the recorded material budget into the corresponding histogram.
virtual void getRay(G4ThreeVector &origin, G4ThreeVector &direction)=0
Get the origin and direction for the next scan particle.
double m_stepU
Stepsize for the parameter u.
bool createNext(G4ThreeVector &origin, G4ThreeVector &direction) override
Get the origin and direction for the next scan particle.
TH2D * getHistogram(const std::string &name)
get histogram for a given name, create if needed.
void UserSteppingAction(const G4Step *step) override
Record the material budget for each step of the particles.
MaterialScan2D(TFile *rootFile, const std::string &name, const std::string &axisLabel, const ScanParams ¶ms)
Constructor.
std::map< std::string, std::unique_ptr< TH2D > > m_regions
Map holding pointers to all created histograms.
double m_curDepth
Tracklength of the current Ray.
double m_curU
Current value of the parametetr u.
double m_curV
Current value of the parametetr v.
double m_stepV
Stepsize for the parameter v.
ScanParams m_params
Parameters for the scan.
Base class for Material Scans.
bool checkStep(const G4Step *step)
check for stuck tracks by looking at the step length
int m_zeroSteps
Count the number of steps with (almost) zero length.
static constexpr int c_maxZeroStepsNudge
maximum number of consecutive zero steps before nudging the track along
TFile * m_rootFile
Pointer to the root file for the histograms.
std::string m_axisLabel
Labels for the coordinate axes.
static constexpr int c_maxZeroStepsKill
maximum number of consecutive zero steps before killing the track
std::string m_name
Name of the scan, will be prefixed to all histogram names.
static constexpr double c_zeroTolerance
maximum Step length to be considered zero
double m_rayTheta
Theta direction of the ray if custom direction is not set.
G4ThreeVector getAxis(char name)
Return a vector along the axis with the given name.
double m_raySampleDepth
Sample depth of the ray: create one bin every x cm.
bool m_raySplitByMaterials
If true Split by materials instead of regions when doing ray scan.
double m_rayPhi
Phi direction of the ray if custom direction is not set.
void initialize() override
Initialize the output file and check the parameters.
TFile * m_rootFile
Pointer to the ROOT file for the histograms.
std::string m_planeName
Name of the plane to use for scanning.
std::vector< std::string > m_rayIgnoredMaterials
Materials ignored when ray scanning.
unsigned int m_rayCount
Number of rays to shoot (if opening angle is >0)
std::optional< std::vector< double > > m_rayDirection
Direction of the ray, cannot be set at the same time as Theta and Phi.
void beginRun() override
Do the actual scanning.
MaterialScan2D::ScanParams m_spherical
Scan parameters for the spherical scan.
bool m_doSpherical
Wether or not to do a spherical scan.
bool m_doPlanar
Wether or not to do a planar scan.
bool m_doRay
Perform a material ray scan: scan along a certain direction with a given opening angle and plot mater...
std::vector< double > m_sphericalOrigin
Origin for spherical scan.
std::string m_filename
Filename for the ROOT output.
MaterialScan2D::ScanParams m_planar
Scan parameters for the planar scan.
double m_rayOpening
Opening angle of the ray.
double m_rayMaxDepth
Max depth of the ray before stopping.
std::vector< double > m_customPlane
Custom plane definition if m_planName is "custom".
bool m_doCosTheta
Perform the spherical scan uniform in cos(theta) instead of theta.
std::vector< double > m_rayOrigin
Origin of the ray if a ray scan is performed
Specific implementaion of MaterialScan to scan parallel to a given plane.
G4ThreeVector m_origin
Origin of the scan plane.
G4ThreeVector m_dirV
v direction of the scan plane
G4ThreeVector m_dirU
u direction of the scan plane
G4ThreeVector m_dirW
direction perpendicluar to u and v
void getRay(G4ThreeVector &origin, G4ThreeVector &direction) override
Create a ray with the current parameter values according to a planar distribution.
MaterialScan implementation to shoot one ray along a defined direction and record the Material as a f...
void fillValue(const std::string &name, double value, double steplength)
Fill the recorded material budget into the corresponding histogram.
G4ThreeVector m_origin
Origin of the scan.
std::map< std::string, std::unique_ptr< TH1D > > m_regions
Map holding pointers to all created histograms.
bool m_splitByMaterials
If true Split by materials instead of regions.
double m_scanDepth
The first ray does not record any material but just checks for the maximum useful depth to not get a ...
bool createNext(G4ThreeVector &origin, G4ThreeVector &direction) override
Implement shooting along the ray.
int m_curRay
Current Ray number: 0 = scan for maximum depth, 1..N = record materials.
void UserSteppingAction(const G4Step *step) override
Record the material budget for each step of the particles.
G4ThreeVector m_dir
Direction of the ray.
double m_curDepth
Current depth of the current ray.
double m_maxDepth
Maximum depth for each ray after which it will be stopped.
double m_opening
Opening angle in radian.
std::set< std::string > m_ignoredMaterials
Materials ignored when scanning.
double m_sampleDepth
The ray length after which to sample.
TH1D * getHistogram(const std::string &name)
get histogram for a given name, create if needed.
int m_count
Amount of rays to shoot.
Specific implementation of MaterialScan to do Spherical scanning.
G4ThreeVector m_origin
Origin for the spherical scan.
void getRay(G4ThreeVector &origin, G4ThreeVector &direction) override
Create a ray with the current parameter values according to a spherical distribution.
bool m_doCosTheta
Flag to indicate if polar-angular sampling is uniform in cos(theta) rather than theta.
virtual void event()
This method is the core of the module.
static const double mm
[millimeters]
static const double deg
degree to radians
static const double ms
[millisecond]
static const double s
[second]
Class to store variables with their name which were sent to the logging service.
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
double getClock()
Return current value of the real-time clock.
Abstract base class for different kinds of events.
Helper struct to Store Parameters of a Scan.
double maxU
Maximum u value to scan.
bool splitByMaterials
If true, split output by Materials (otherwise by region)
double minV
Minimum v value to scan.
double maxDepth
Maximum depth of the scan.
int nU
Number of rays along u coordinate.
double minU
Minimum u value to scan.
double maxV
Maximum v value to scan.
std::vector< std::string > ignoredMaterials
Names of ignored Materials.
int nV
Number of rays along v coordinate.