Belle II Software  release-06-02-00
MaterialScan.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 #include <framework/core/Module.h>
12 #include <framework/gearbox/Unit.h>
13 
14 #include <string>
15 #include <vector>
16 #include <set>
17 
18 #include <G4UserSteppingAction.hh>
19 #include <G4ThreeVector.hh>
20 
21 #include <TFile.h>
22 class TH1D;
23 class TH2D;
24 
25 #include <boost/optional.hpp>
26 
27 namespace Belle2 {
34  class MaterialScanBase: public G4UserSteppingAction {
35  public:
38  MaterialScanBase(TFile* rootFile, const std::string& name, const std::string& axisLabel):
39  m_rootFile(rootFile), m_name(name), m_axisLabel(axisLabel)
40  {
41  m_rootFile->mkdir(name.c_str());
42  }
44  std::string getName() const { return m_name; }
46  virtual int getNRays() const = 0;
47 
55  virtual bool createNext(G4ThreeVector& origin, G4ThreeVector& direction) = 0;
56 
57  protected:
59  bool checkStep(const G4Step* step);
61  TFile* m_rootFile;
63  std::string m_name;
65  std::string m_axisLabel;
66  private:
68  static constexpr double c_zeroTolerance = 1e-6;
70  static constexpr int c_maxZeroStepsNudge = 10;
72  static constexpr int c_maxZeroStepsKill = 20;
74  int m_zeroSteps {0};
75  };
76 
82  public:
84  struct ScanParams {
86  ScanParams(): nU(0), nV(0), minU(0), maxU(0), minV(0), maxV(0), maxDepth(-1), splitByMaterials(false) {}
88  int nU;
90  int nV;
92  double minU;
94  double maxU;
96  double minV;
98  double maxV;
100  double maxDepth;
102  std::vector<std::string> ignoredMaterials;
105  };
106 
113  MaterialScan2D(TFile* rootFile, const std::string& name, const std::string& axisLabel, const ScanParams& params);
114 
120  bool createNext(G4ThreeVector& origin, G4ThreeVector& direction) override;
121 
123  int getNRays() const override { return m_params.nU * m_params.nV; }
124 
126  void UserSteppingAction(const G4Step* step) override;
127  protected:
133  virtual void getRay(G4ThreeVector& origin, G4ThreeVector& direction) = 0;
134 
138  TH2D* getHistogram(const std::string& name);
139 
144  void fillValue(const std::string& name, double value);
145 
149  double m_curU;
151  double m_stepU;
153  double m_curV;
155  double m_stepV;
157  double m_curDepth;
159  std::map<std::string, std::unique_ptr<TH2D>> m_regions;
160  };
161 
166  public:
173  MaterialScanSpherical(TFile* rootFile, const G4ThreeVector& origin, const ScanParams& params, bool doCosTheta):
174  MaterialScan2D(rootFile, "Spherical", doCosTheta ? "cos(#theta);#phi [deg]" : "#theta [deg];#phi [deg]", params), m_origin(origin),
175  m_doCosTheta(doCosTheta)
176  {
177  if (doCosTheta) {
180  if (m_params.minU > m_params.maxU) std::swap(m_params.minU, m_params.maxU);
182  m_curU = m_params.minU - m_stepU / 2.;
183  }
184  }
185  protected:
187  void getRay(G4ThreeVector& origin, G4ThreeVector& direction) override;
188 
190  G4ThreeVector m_origin;
191 
194  };
195 
204  public:
212  MaterialScanPlanar(TFile* rootFile, const G4ThreeVector& origin, const G4ThreeVector& dirU, const G4ThreeVector& dirV,
213  const ScanParams& params):
214  MaterialScan2D(rootFile, "Planar", "u [cm];v [cm]", params), m_origin(origin), m_dirU(dirU.unit()), m_dirV(dirV.unit()),
215  m_dirW(m_dirU.cross(m_dirV))
216  {
217  }
218  protected:
220  void getRay(G4ThreeVector& origin, G4ThreeVector& direction) override;
221  protected:
223  G4ThreeVector m_origin;
225  G4ThreeVector m_dirU;
227  G4ThreeVector m_dirV;
229  G4ThreeVector m_dirW;
230  };
231 
241  public:
243  MaterialScanRay(TFile* rootFile, const G4ThreeVector& origin, const G4ThreeVector& dir, double opening,
244  int count, double sampleDepth, double maxDepth, bool splitByMaterials,
245  const std::vector<std::string>& ignoredMaterials):
246  MaterialScanBase(rootFile, "Ray", "ray length [cm]; material budget [X_0/cm]"),
247  m_ignoredMaterials(ignoredMaterials.begin(), ignoredMaterials.end()),
248  m_origin(origin), m_dir(dir), m_opening(opening), m_sampleDepth(sampleDepth),
249  m_maxDepth(maxDepth), m_count(std::max(1, count)), m_splitByMaterials(splitByMaterials)
250  {
251  if (m_opening <= 0) m_count = 1;
252  }
253 
255  int getNRays() const override { return m_count + 1; }
257  void UserSteppingAction(const G4Step* step) override;
259  bool createNext(G4ThreeVector& origin, G4ThreeVector& direction) override;
260  private:
264  TH1D* getHistogram(const std::string& name);
270  void fillValue(const std::string& name, double value, double steplength);
271 
273  std::set<std::string> m_ignoredMaterials;
275  std::map<std::string, std::unique_ptr<TH1D>> m_regions;
277  G4ThreeVector m_origin;
279  G4ThreeVector m_dir;
281  double m_opening;
285  double m_maxDepth{0};
287  double m_curDepth{0};
293  double m_scanDepth{0};
295  int m_count;
297  int m_curRay{ -1};
300  };
301 
312  class MaterialScanModule : public Module {
313 
314  public:
315 
316  /* The constructor of the module.
317  * Sets the description and the parameters of the module.
318  */
320 
322  void initialize() override;
323 
325  void beginRun() override;
326  private:
329  G4ThreeVector getAxis(char name);
330 
332  TFile* m_rootFile;
334  std::string m_filename;
336  std::string m_planeName;
342  std::vector<double> m_customPlane;
344  std::vector<double> m_sphericalOrigin;
346  std::vector<double> m_rayOrigin{0, 0, 0};
348  std::vector<std::string> m_rayIgnoredMaterials;
350  boost::optional<std::vector<double>> m_rayDirection;
352  double m_rayTheta{0};
354  double m_rayPhi{0};
356  double m_rayOpening{0};
358  double m_raySampleDepth{0.1};
360  double m_rayMaxDepth{1000};
362  unsigned int m_rayCount{0};
371  bool m_doRay;
374  };
376 }
Base class to create a Material Scan of the detector geometry.
Definition: MaterialScan.h:81
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.
Definition: MaterialScan.h:151
bool createNext(G4ThreeVector &origin, G4ThreeVector &direction) override
Get the origin and direction for the next scan particle.
Definition: MaterialScan.cc:97
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 &params)
Constructor.
Definition: MaterialScan.cc:78
std::map< std::string, std::unique_ptr< TH2D > > m_regions
Map holding pointers to all created histograms.
Definition: MaterialScan.h:159
double m_curDepth
Tracklength of the current Ray.
Definition: MaterialScan.h:157
int getNRays() const override
Return the number of rays in this scan.
Definition: MaterialScan.h:123
double m_curU
Current value of the parametetr u.
Definition: MaterialScan.h:149
double m_curV
Current value of the parametetr v.
Definition: MaterialScan.h:153
double m_stepV
Stepsize for the parameter v.
Definition: MaterialScan.h:155
ScanParams m_params
Parameters for the scan.
Definition: MaterialScan.h:147
Base class for Material Scans.
Definition: MaterialScan.h:34
virtual bool createNext(G4ThreeVector &origin, G4ThreeVector &direction)=0
Belle2::MaterialScanBase::createNext() implemention Get the origin and direction for the next scan pa...
bool checkStep(const G4Step *step)
check for stuck tracks by looking at the step length
Definition: MaterialScan.cc:44
int m_zeroSteps
Count the number of steps with (almost) zero length.
Definition: MaterialScan.h:74
static constexpr int c_maxZeroStepsNudge
maximum number of consecutive zero steps before nudging the track along
Definition: MaterialScan.h:70
virtual int getNRays() const =0
Return the number of rays necessary to perform the scan.
TFile * m_rootFile
Pointer to the root file for the histograms.
Definition: MaterialScan.h:61
std::string m_axisLabel
Labels for the coordinate axes.
Definition: MaterialScan.h:65
MaterialScanBase(TFile *rootFile, const std::string &name, const std::string &axisLabel)
This is indeed the constructor and it creates a TDirectory in the output root file and sets all varia...
Definition: MaterialScan.h:38
std::string getName() const
Return the name of the scan.
Definition: MaterialScan.h:44
static constexpr int c_maxZeroStepsKill
maximum number of consecutive zero steps before killing the track
Definition: MaterialScan.h:72
std::string m_name
Name of the scan, will be prefixed to all histogram names.
Definition: MaterialScan.h:63
static constexpr double c_zeroTolerance
maximum Step length to be considered zero
Definition: MaterialScan.h:68
The MaterialScan module.
Definition: MaterialScan.h:312
double m_rayTheta
Theta direction of the ray if custom direction is not set.
Definition: MaterialScan.h:352
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.
Definition: MaterialScan.h:358
bool m_raySplitByMaterials
If true Split by materials instead of regions when doing ray scan.
Definition: MaterialScan.h:373
double m_rayPhi
Phi direction of the ray if custom direction is not set.
Definition: MaterialScan.h:354
void initialize() override
Initialize the output file and check the parameters.
TFile * m_rootFile
Pointer to the ROOT file for the histograms.
Definition: MaterialScan.h:332
std::string m_planeName
Name of the plane to use for scanning.
Definition: MaterialScan.h:336
std::vector< std::string > m_rayIgnoredMaterials
Materials ignored when ray scanning.
Definition: MaterialScan.h:348
unsigned int m_rayCount
Number of rays to shoot (if opening angle is >0)
Definition: MaterialScan.h:362
void beginRun() override
Do the actual scanning.
MaterialScan2D::ScanParams m_spherical
Scan parameters for the spherical scan.
Definition: MaterialScan.h:338
bool m_doSpherical
Wether or not to do a spherical scan.
Definition: MaterialScan.h:364
bool m_doPlanar
Wether or not to do a planar scan.
Definition: MaterialScan.h:366
bool m_doRay
Perform a material ray scan: scan along a certain direction with a given opening angle and plot mater...
Definition: MaterialScan.h:371
std::vector< double > m_sphericalOrigin
Origin for spherical scan.
Definition: MaterialScan.h:344
boost::optional< std::vector< double > > m_rayDirection
Direction of the ray, cannot be set at the same time as Theta and Phi.
Definition: MaterialScan.h:350
std::string m_filename
Filename for the ROOT output.
Definition: MaterialScan.h:334
MaterialScan2D::ScanParams m_planar
Scan parameters for the planar scan.
Definition: MaterialScan.h:340
double m_rayOpening
Opening angle of the ray.
Definition: MaterialScan.h:356
double m_rayMaxDepth
Max depth of the ray before stopping.
Definition: MaterialScan.h:360
std::vector< double > m_customPlane
Custom plane definition if m_planName is "custom".
Definition: MaterialScan.h:342
bool m_doCosTheta
Perform the spherical scan uniform in cos(theta) instead of theta.
Definition: MaterialScan.h:368
std::vector< double > m_rayOrigin
Origin of the ray if a ray scan is performed
Definition: MaterialScan.h:346
Specific implementaion of MaterialScan to scan parallel to a given plane.
Definition: MaterialScan.h:203
G4ThreeVector m_origin
Origin of the scan plane.
Definition: MaterialScan.h:223
G4ThreeVector m_dirV
v direction of the scan plane
Definition: MaterialScan.h:227
G4ThreeVector m_dirU
u direction of the scan plane
Definition: MaterialScan.h:225
MaterialScanPlanar(TFile *rootFile, const G4ThreeVector &origin, const G4ThreeVector &dirU, const G4ThreeVector &dirV, const ScanParams &params)
Create a Planar Scan object with the given parameters.
Definition: MaterialScan.h:212
G4ThreeVector m_dirW
direction perpendicluar to u and v
Definition: MaterialScan.h:229
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...
Definition: MaterialScan.h:240
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.
Definition: MaterialScan.h:277
std::map< std::string, std::unique_ptr< TH1D > > m_regions
Map holding pointers to all created histograms.
Definition: MaterialScan.h:275
MaterialScanRay(TFile *rootFile, const G4ThreeVector &origin, const G4ThreeVector &dir, double opening, int count, double sampleDepth, double maxDepth, bool splitByMaterials, const std::vector< std::string > &ignoredMaterials)
Construct a new instance and set all parameters.
Definition: MaterialScan.h:243
bool m_splitByMaterials
If true Split by materials instead of regions.
Definition: MaterialScan.h:299
double m_scanDepth
The first ray does not record any material but just checks for the maximum useful depth to not get a ...
Definition: MaterialScan.h:293
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.
Definition: MaterialScan.h:297
void UserSteppingAction(const G4Step *step) override
Record the material budget for each step of the particles.
G4ThreeVector m_dir
Direction of the ray.
Definition: MaterialScan.h:279
double m_curDepth
Current depth of the current ray.
Definition: MaterialScan.h:287
double m_maxDepth
Maximum depth for each ray after which it will be stopped.
Definition: MaterialScan.h:285
int getNRays() const override
Return the number of rays.
Definition: MaterialScan.h:255
double m_opening
Opening angle in radian.
Definition: MaterialScan.h:281
std::set< std::string > m_ignoredMaterials
Materials ignored when scanning.
Definition: MaterialScan.h:273
double m_sampleDepth
The ray length after which to sample.
Definition: MaterialScan.h:283
TH1D * getHistogram(const std::string &name)
get histogram for a given name, create if needed.
int m_count
Amount of rays to shoot.
Definition: MaterialScan.h:295
Specific implementation of MaterialScan to do Spherical scanning.
Definition: MaterialScan.h:165
G4ThreeVector m_origin
Origin for the spherical scan.
Definition: MaterialScan.h:190
MaterialScanSpherical(TFile *rootFile, const G4ThreeVector &origin, const ScanParams &params, bool doCosTheta)
Create a Spherical Scan object with the given parameters.
Definition: MaterialScan.h:173
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.
Definition: MaterialScan.h:193
Base class for Modules.
Definition: Module.h:72
static const double deg
degree to radians
Definition: Unit.h:109
Abstract base class for different kinds of events.
Helper struct to Store Parameters of a Scan.
Definition: MaterialScan.h:84
double maxU
Maximum u value to scan.
Definition: MaterialScan.h:94
bool splitByMaterials
If true, split output by Materials (otherwise by region)
Definition: MaterialScan.h:104
ScanParams()
Default Constructor.
Definition: MaterialScan.h:86
double minV
Minimum v value to scan.
Definition: MaterialScan.h:96
double maxDepth
Maximum depth of the scan.
Definition: MaterialScan.h:100
int nU
Number of rays along u coordinate.
Definition: MaterialScan.h:88
double minU
Minimum u value to scan.
Definition: MaterialScan.h:92
double maxV
Maximum v value to scan.
Definition: MaterialScan.h:98
std::vector< std::string > ignoredMaterials
Names of ignored Materials.
Definition: MaterialScan.h:102
int nV
Number of rays along v coordinate.
Definition: MaterialScan.h:90