Belle II Software development
MaterialScanSpherical Class Reference

Specific implementation of MaterialScan to do Spherical scanning. More...

#include <MaterialScan.h>

Inheritance diagram for MaterialScanSpherical:
MaterialScan2D MaterialScanBase

Public Member Functions

 MaterialScanSpherical (TFile *rootFile, const G4ThreeVector &origin, const ScanParams &params, bool doCosTheta)
 Create a Spherical Scan object with the given parameters.
 
bool createNext (G4ThreeVector &origin, G4ThreeVector &direction) override
 Get the origin and direction for the next scan particle.
 
int getNRays () const override
 Return the number of rays in this scan.
 
void UserSteppingAction (const G4Step *step) override
 Record the material budget for each step of the particles.
 
std::string getName () const
 Return the name of the scan.
 

Protected Member Functions

void getRay (G4ThreeVector &origin, G4ThreeVector &direction) override
 Create a ray with the current parameter values according to a spherical distribution.
 
TH2D * getHistogram (const std::string &name)
 get histogram for a given name, create if needed.
 
void fillValue (const std::string &name, double value)
 Fill the recorded material budget into the corresponding histogram.
 
bool checkStep (const G4Step *step)
 check for stuck tracks by looking at the step length
 

Protected Attributes

G4ThreeVector m_origin
 Origin for the spherical scan.
 
bool m_doCosTheta
 Flag to indicate if polar-angular sampling is uniform in cos(theta) rather than theta.
 
ScanParams m_params
 Parameters for the scan.
 
double m_curU
 Current value of the parametetr u.
 
double m_stepU
 Stepsize for the parameter u.
 
double m_curV
 Current value of the parametetr v.
 
double m_stepV
 Stepsize for the parameter v.
 
double m_curDepth
 Tracklength of the current Ray.
 
std::map< std::string, std::unique_ptr< TH2D > > m_regions
 Map holding pointers to all created histograms.
 
TFile * m_rootFile
 Pointer to the root file for the histograms.
 
std::string m_name
 Name of the scan, will be prefixed to all histogram names.
 
std::string m_axisLabel
 Labels for the coordinate axes.
 

Private Attributes

int m_zeroSteps {0}
 Count the number of steps with (almost) zero length.
 

Static Private Attributes

static constexpr double c_zeroTolerance = 1e-6
 maximum Step length to be considered zero
 
static constexpr int c_maxZeroStepsNudge = 10
 maximum number of consecutive zero steps before nudging the track along
 
static constexpr int c_maxZeroStepsKill = 20
 maximum number of consecutive zero steps before killing the track
 

Detailed Description

Specific implementation of MaterialScan to do Spherical scanning.

That is shooting rays from the origin with varying aximuth and polar angle.

Definition at line 164 of file MaterialScan.h.

Constructor & Destructor Documentation

◆ MaterialScanSpherical()

MaterialScanSpherical ( TFile *  rootFile,
const G4ThreeVector &  origin,
const ScanParams params,
bool  doCosTheta 
)
inline

Create a Spherical Scan object with the given parameters.

Parameters
rootFilepointer to the ROOT File containing the histograms
originOrigin for the spherical scan
paramsParameters of the scan
doCosTheta

Definition at line 172 of file MaterialScan.h.

172 :
173 MaterialScan2D(rootFile, "Spherical", doCosTheta ? "cos(#theta);#phi [deg]" : "#theta [deg];#phi [deg]", params), m_origin(origin),
174 m_doCosTheta(doCosTheta)
175 {
176 if (doCosTheta) {
179 if (m_params.minU > m_params.maxU) std::swap(m_params.minU, m_params.maxU);
181 m_curU = m_params.minU - m_stepU / 2.;
182 }
183 }
double m_stepU
Stepsize for the parameter u.
Definition: MaterialScan.h:150
MaterialScan2D(TFile *rootFile, const std::string &name, const std::string &axisLabel, const ScanParams &params)
Constructor.
Definition: MaterialScan.cc:78
double m_curU
Current value of the parametetr u.
Definition: MaterialScan.h:148
ScanParams m_params
Parameters for the scan.
Definition: MaterialScan.h:146
G4ThreeVector m_origin
Origin for the spherical scan.
Definition: MaterialScan.h:189
bool m_doCosTheta
Flag to indicate if polar-angular sampling is uniform in cos(theta) rather than theta.
Definition: MaterialScan.h:192
static const double deg
degree to radians
Definition: Unit.h:109
double maxU
Maximum u value to scan.
Definition: MaterialScan.h:93
int nU
Number of rays along u coordinate.
Definition: MaterialScan.h:87
double minU
Minimum u value to scan.
Definition: MaterialScan.h:91

Member Function Documentation

◆ checkStep()

bool checkStep ( const G4Step *  step)
protectedinherited

check for stuck tracks by looking at the step length

Definition at line 44 of file MaterialScan.cc.

45{
46 double stlen = step->GetStepLength();
47 G4StepPoint* preStepPoint = step->GetPreStepPoint();
48 G4Region* region = preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
49 if (stlen < c_zeroTolerance) {
51 } else {
52 m_zeroSteps = 0;
53 }
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 "
60 << c_zeroTolerance << " mm, killing it");
61 step->GetTrack()->SetTrackStatus(fStopAndKill);
62 } else {
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 "
67 << c_zeroTolerance << " mm, nudging it along");
68 G4ThreeVector pos = step->GetTrack()->GetPosition();
69 G4ThreeVector dir = step->GetTrack()->GetMomentumDirection();
70 step->GetTrack()->SetPosition(pos + c_zeroTolerance * dir);
71 }
72 return false;
73 }
74 return true;
75}
int m_zeroSteps
Count the number of steps with (almost) zero length.
Definition: MaterialScan.h:73
static constexpr int c_maxZeroStepsNudge
maximum number of consecutive zero steps before nudging the track along
Definition: MaterialScan.h:69
static constexpr int c_maxZeroStepsKill
maximum number of consecutive zero steps before killing the track
Definition: MaterialScan.h:71
static constexpr double c_zeroTolerance
maximum Step length to be considered zero
Definition: MaterialScan.h:67

◆ createNext()

bool createNext ( G4ThreeVector &  origin,
G4ThreeVector &  direction 
)
overridevirtualinherited

Get the origin and direction for the next scan particle.

Parameters
[out]originOrigin of the next scan particle
[out]directionDirection of the next scan particle
Returns
false if the scan is finished

Implements MaterialScanBase.

Definition at line 97 of file MaterialScan.cc.

98{
99 if (m_regions.empty()) {
100 // create summary histogram right now, not on demand, to make sure it is in the file
102 getHistogram("All_Materials_x0");
103 getHistogram("All_Materials_lambda");
104 } else {
105 getHistogram("All_Regions_x0");
106 getHistogram("All_Regions_lambda");
107 }
108 }
109 //Increase the internal coordinates
110 m_curU += m_stepU;
111 if (m_curU >= m_params.maxU) {
112 m_curU = m_params.minU + m_stepU / 2.;
113 m_curV += m_stepV;
114 }
115 //Reset depth counter
116 m_curDepth = 0;
117
118 //Get the origin and direction of the ray
119 getRay(origin, direction);
120
121 //Check wether we are finished
122 return (m_curV <= m_params.maxV);
123}
virtual void getRay(G4ThreeVector &origin, G4ThreeVector &direction)=0
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.
std::map< std::string, std::unique_ptr< TH2D > > m_regions
Map holding pointers to all created histograms.
Definition: MaterialScan.h:158
double m_curDepth
Tracklength of the current Ray.
Definition: MaterialScan.h:156
double m_curV
Current value of the parametetr v.
Definition: MaterialScan.h:152
double m_stepV
Stepsize for the parameter v.
Definition: MaterialScan.h:154
bool splitByMaterials
If true, split output by Materials (otherwise by region)
Definition: MaterialScan.h:103
double maxV
Maximum v value to scan.
Definition: MaterialScan.h:97

◆ fillValue()

void fillValue ( const std::string &  name,
double  value 
)
protectedinherited

Fill the recorded material budget into the corresponding histogram.

Parameters
nameName of the histogram
valueValue to store

Definition at line 138 of file MaterialScan.cc.

139{
140 TH2D* hist = getHistogram(name);
141 hist->Fill(m_curU, m_curV, value);
142}

◆ getHistogram()

TH2D * getHistogram ( const std::string &  name)
protectedinherited

get histogram for a given name, create if needed.

Parameters
nameName of the histogram

Definition at line 125 of file MaterialScan.cc.

126{
127 std::unique_ptr<TH2D>& hist = m_regions[name];
128 if (!hist) {
129 //Create new histogram
130 m_rootFile->cd(m_name.c_str());
131 hist.reset(new TH2D(name.c_str(), (name + ";" + m_axisLabel).c_str(),
134 }
135 return hist.get();
136}
TFile * m_rootFile
Pointer to the root file for the histograms.
Definition: MaterialScan.h:60
std::string m_axisLabel
Labels for the coordinate axes.
Definition: MaterialScan.h:64
std::string m_name
Name of the scan, will be prefixed to all histogram names.
Definition: MaterialScan.h:62
double minV
Minimum v value to scan.
Definition: MaterialScan.h:95
int nV
Number of rays along v coordinate.
Definition: MaterialScan.h:89

◆ getName()

std::string getName ( ) const
inlineinherited

Return the name of the scan.

Definition at line 43 of file MaterialScan.h.

43{ return m_name; }

◆ getNRays()

int getNRays ( ) const
inlineoverridevirtualinherited

Return the number of rays in this scan.

Implements MaterialScanBase.

Definition at line 122 of file MaterialScan.h.

122{ return m_params.nU * m_params.nV; }

◆ getRay()

void getRay ( G4ThreeVector &  origin,
G4ThreeVector &  direction 
)
overrideprotectedvirtual

Create a ray with the current parameter values according to a spherical distribution.

Implements MaterialScan2D.

Definition at line 194 of file MaterialScan.cc.

195{
196 //We always shoot from the origin
197 origin = m_origin;
198 double theta = m_curU * Unit::deg;
199 double phi = m_curV * Unit::deg;
200 if (m_doCosTheta) {
201 theta = acos(m_curU);
202 }
203 direction.set(sin(theta) * cos(phi),
204 sin(theta) * sin(phi),
205 cos(theta));
206}

◆ UserSteppingAction()

void UserSteppingAction ( const G4Step *  step)
overrideinherited

Record the material budget for each step of the particles.

Definition at line 144 of file MaterialScan.cc.

145{
146 //Get information about radiation and interaction length
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");
157 checkStep(step);
158
159 //check if the depth is limited
160 if (m_params.maxDepth > 0) {
161 m_curDepth += stlen;
163 //Depth reached, kill track and add remaining part of the material budget
164 G4Track* track = step->GetTrack();
165 track->SetTrackStatus(fStopAndKill);
166 double remaining = stlen - m_curDepth + m_params.maxDepth;
167 x0 *= remaining / stlen;
168 lambda *= remaining / stlen;
169 }
170 }
171
172 if (std::binary_search(m_params.ignoredMaterials.begin(), m_params.ignoredMaterials.end(), material->GetName())) {
173 return;
174 }
175
176 //Fill x0 and lambda in a histogram for each region
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";
181 //or for each Material
183 x0_total = "All_Materials_x0";
184 lambda_total = "All_Materials_lambda";
185 x0_name = material->GetName() + "_x0";
186 lambda_name = material->GetName() + "_lambda";
187 }
188 fillValue(x0_total, x0);
189 fillValue(lambda_total, lambda);
190 fillValue(x0_name, x0);
191 fillValue(lambda_name, lambda);
192}
void fillValue(const std::string &name, double value)
Fill the recorded material budget into the corresponding histogram.
bool checkStep(const G4Step *step)
check for stuck tracks by looking at the step length
Definition: MaterialScan.cc:44
double maxDepth
Maximum depth of the scan.
Definition: MaterialScan.h:99
std::vector< std::string > ignoredMaterials
Names of ignored Materials.
Definition: MaterialScan.h:101

Member Data Documentation

◆ c_maxZeroStepsKill

constexpr int c_maxZeroStepsKill = 20
staticconstexprprivateinherited

maximum number of consecutive zero steps before killing the track

Definition at line 71 of file MaterialScan.h.

◆ c_maxZeroStepsNudge

constexpr int c_maxZeroStepsNudge = 10
staticconstexprprivateinherited

maximum number of consecutive zero steps before nudging the track along

Definition at line 69 of file MaterialScan.h.

◆ c_zeroTolerance

constexpr double c_zeroTolerance = 1e-6
staticconstexprprivateinherited

maximum Step length to be considered zero

Definition at line 67 of file MaterialScan.h.

◆ m_axisLabel

std::string m_axisLabel
protectedinherited

Labels for the coordinate axes.

Definition at line 64 of file MaterialScan.h.

◆ m_curDepth

double m_curDepth
protectedinherited

Tracklength of the current Ray.

Definition at line 156 of file MaterialScan.h.

◆ m_curU

double m_curU
protectedinherited

Current value of the parametetr u.

Definition at line 148 of file MaterialScan.h.

◆ m_curV

double m_curV
protectedinherited

Current value of the parametetr v.

Definition at line 152 of file MaterialScan.h.

◆ m_doCosTheta

bool m_doCosTheta
protected

Flag to indicate if polar-angular sampling is uniform in cos(theta) rather than theta.

Definition at line 192 of file MaterialScan.h.

◆ m_name

std::string m_name
protectedinherited

Name of the scan, will be prefixed to all histogram names.

Definition at line 62 of file MaterialScan.h.

◆ m_origin

G4ThreeVector m_origin
protected

Origin for the spherical scan.

Definition at line 189 of file MaterialScan.h.

◆ m_params

ScanParams m_params
protectedinherited

Parameters for the scan.

Definition at line 146 of file MaterialScan.h.

◆ m_regions

std::map<std::string, std::unique_ptr<TH2D> > m_regions
protectedinherited

Map holding pointers to all created histograms.

Definition at line 158 of file MaterialScan.h.

◆ m_rootFile

TFile* m_rootFile
protectedinherited

Pointer to the root file for the histograms.

Definition at line 60 of file MaterialScan.h.

◆ m_stepU

double m_stepU
protectedinherited

Stepsize for the parameter u.

Definition at line 150 of file MaterialScan.h.

◆ m_stepV

double m_stepV
protectedinherited

Stepsize for the parameter v.

Definition at line 154 of file MaterialScan.h.

◆ m_zeroSteps

int m_zeroSteps {0}
privateinherited

Count the number of steps with (almost) zero length.

Definition at line 73 of file MaterialScan.h.


The documentation for this class was generated from the following files: