Belle II Software development
MaterialScan.cc
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#include <simulation/modules/MaterialScan.h>
10#include <framework/core/ModuleParam.templateDetails.h>
11#include <framework/gearbox/Unit.h>
12#include <framework/utilities/Utils.h>
13#include <algorithm>
14
15#include <TFile.h>
16#include <TH1D.h>
17#include <TH2D.h>
18#include <TRandom.h>
19
20#include <G4Event.hh>
21#include <G4EventManager.hh>
22#include <G4UserEventAction.hh>
23#include <G4UserStackingAction.hh>
24#include <G4UserTrackingAction.hh>
25#include <G4UserSteppingAction.hh>
26#include <G4RayShooter.hh>
27#include <G4Track.hh>
28
29#include <boost/format.hpp>
30#include <boost/algorithm/string.hpp>
31
32using namespace std;
33using namespace Belle2;
34
35//-----------------------------------------------------------------
36// Register the Module
37//-----------------------------------------------------------------
38REG_MODULE(MaterialScan);
39
40//-----------------------------------------------------------------
41// Implementation
42//-----------------------------------------------------------------
43
44bool MaterialScanBase::checkStep(const G4Step* step)
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}
76
77
78MaterialScan2D::MaterialScan2D(TFile* rootFile, const std::string& name, const std::string& axisLabel, const ScanParams& params):
79 MaterialScanBase(rootFile, name, axisLabel), m_params(params), m_curDepth(0)
80{
81 //Sort the parameters accordingly
84 //Calculate step size of the parameters
87 //Set Start values
88 m_curU = m_params.minU - m_stepU / 2.;
89 m_curV = m_params.minV + m_stepV / 2.;
90 //Convert max depth to G4 units
92 //Sort the list of ignored materials so that we can use binary search
93 std::sort(m_params.ignoredMaterials.begin(), m_params.ignoredMaterials.end());
94}
95
96
97bool MaterialScan2D::createNext(G4ThreeVector& origin, G4ThreeVector& direction)
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}
124
125TH2D* MaterialScan2D::getHistogram(const std::string& name)
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}
137
138void MaterialScan2D::fillValue(const std::string& name, double value)
139{
140 TH2D* hist = getHistogram(name);
141 hist->Fill(m_curU, m_curV, value);
142}
143
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}
193
194void MaterialScanSpherical::getRay(G4ThreeVector& origin, G4ThreeVector& direction)
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}
207
208void MaterialScanPlanar::getRay(G4ThreeVector& origin, G4ThreeVector& direction)
209{
210 //We shoot perpendicular to the plane, so direction is always the same but the position varies.
211 origin = m_origin + m_curU / Unit::mm * m_dirU + m_curV / Unit::mm * m_dirV;
212 direction = m_dirW;
213}
214
215bool MaterialScanRay::createNext(G4ThreeVector& origin, G4ThreeVector& direction)
216{
217 if (m_curRay > m_count) return false;
218
219 origin = m_origin;
220 direction = m_dir;
221 ++m_curRay;
222 if (m_curRay == 0) {
223 // First ray is just to check what the maximum flight length is
224 return true;
225 } else if (m_curRay == 1) {
226 // this is the second ray, save the maximum depth for histograms
228 B2INFO("Set maximum depth to " << m_maxDepth);
229 // create summary histogram right now, not on demand, to make sure it is in the file
230 if (m_splitByMaterials) {
231 getHistogram("All_Materials_x0");
232 getHistogram("All_Materials_lambda");
233 } else {
234 getHistogram("All_Regions_x0");
235 getHistogram("All_Regions_lambda");
236 }
237 }
238 m_curDepth = 0;
239 if (m_opening > 0) {
240 // add small angle to the ray where we want the opening angle to be covered
241 // uniform by area so choose an angle uniform in cos(angle). We take half
242 // the opening angle because it can go in both sides
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);
248 }
249 B2DEBUG(10, "Create Ray " << m_curRay << " from " << origin << " direction " <<
250 direction << " [theta=" << (direction.theta() / M_PI * 180) << ", phi=" << (direction.phi() / M_PI * 180) << "]");
251 return true;
252}
253
254TH1D* MaterialScanRay::getHistogram(const std::string& name)
255{
256 std::unique_ptr<TH1D>& hist = m_regions[name];
257 if (!hist) {
258 //Create new histogram
259 m_rootFile->cd(m_name.c_str());
260 int bins = std::ceil(m_maxDepth / m_sampleDepth);
261 hist.reset(new TH1D(name.c_str(), (name + ";" + m_axisLabel).c_str(), bins, 0, m_maxDepth * Unit::mm));
262 }
263 return hist.get();
264}
265
266void MaterialScanRay::fillValue(const std::string& name, double value, double steplength)
267{
268 TH1D* h = getHistogram(name);
269 double x = m_curDepth;
270 //The steps are not aligned to the bin sizes so make sure we sub-sample the
271 //value and put it in the histogram in smaller steps. We choose 10 samples
272 //per bin here which should be good enough.
273 int samples = std::ceil(steplength / m_sampleDepth * 10);
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) {
277 h->Fill(x * Unit::mm, dv);
278 x += dx;
279 }
280}
281
283{
284 //Get information about radiation and interaction length
285 G4StepPoint* preStepPoint = step->GetPreStepPoint();
286 G4Region* region = preStepPoint->GetPhysicalVolume()->GetLogicalVolume()->GetRegion();
287 G4Material* material = preStepPoint->GetMaterial();
288 const bool ignored = m_ignoredMaterials.count(material->GetName()) > 0;
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");
296 checkStep(step);
297
298 //check if the depth is limited
299 if (m_maxDepth > 0) {
300 if (m_curDepth + stlen > m_maxDepth) {
301 //Depth reached, kill track and add remaining part of the material budget
302 G4Track* track = step->GetTrack();
303 track->SetTrackStatus(fStopAndKill);
304 double remaining = m_maxDepth - m_curDepth;
305 x0 *= remaining / stlen;
306 lambda *= remaining / stlen;
307 stlen = remaining;
308 }
309 }
310
311 if (m_curRay > 0 && !ignored) {
312 //Fill x0 and lambda in a histogram for each region
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";
317 //or for each Material
318 if (m_splitByMaterials) {
319 x0_total = "All_Materials_x0";
320 lambda_total = "All_Materials_lambda";
321 x0_name = material->GetName() + "_x0";
322 lambda_name = material->GetName() + "_lambda";
323 }
324 fillValue(x0_total, x0, stlen);
325 fillValue(lambda_total, lambda, stlen);
326 fillValue(x0_name, x0, stlen);
327 fillValue(lambda_name, lambda, stlen);
328 }
329 m_curDepth += stlen;
330 if (m_curRay == 0 && !ignored) {
331 // update max depth only for materials we do not ignore ...
333 }
334}
335
336MaterialScanModule::MaterialScanModule(): m_rootFile(0), m_sphericalOrigin(3, 0), m_doCosTheta(false)
337{
338 //Set module properties
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 "
345 "plane.");
346
347 //Set default ignored Materials
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;
353
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);
383
384 addParam("planar", m_doPlanar,
385 "Do a plane scan, that is shooting parallel rays from a defined "
386 "plane", true);
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,
393 "U start", -400.0);
394 addParam("planar.maxU", m_planar.maxU,
395 "U stop", 400.0);
396 addParam("planar.nV", m_planar.nV,
397 "Number of rays in V", 200);
398 addParam("planar.minV", m_planar.minV,
399 "V start", -400.0);
400 addParam("planar.maxV", m_planar.maxV,
401 "V stop", 400.0);
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);
415
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",
423 m_rayDirection);
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",
429 m_rayCount);
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);
437}
438
440{
441 m_rootFile = TFile::Open(m_filename.c_str(), "RECREATE");
442 if (!m_rootFile || m_rootFile->IsZombie()) {
443 B2ERROR("Cannot open rootfile for writing" << LogVar("rootfile", m_filename));
444 return;
445 }
446
447 boost::to_lower(m_planeName);
448 boost::trim(m_planeName);
449 //Check if plane definition makes sense
450 if (m_planeName == "custom") {
451 if (m_customPlane.size() < 9) {
452 B2ERROR("planar.custom: At least 9 values needed to define custom plane, only " << m_customPlane.size() << " provided");
453 }
454 } else if (m_planeName.size() != 2 || !boost::all(m_planeName, boost::is_any_of("xyz"))) {
455 B2ERROR("planar.plane: Only custom or any two letter combinations of X,Y and Z are recognized");
456 } else if (m_planeName[0] == m_planeName[1]) {
457 B2ERROR("planar.plane: " << m_planeName << " not valid, cannot use the same axis twice");
458 }
459
460 //Check if we have enough values to define the origin
461 if (m_sphericalOrigin.size() < 3) {
462 B2ERROR("spherical.origin: Three values are needed to define a point, only " << m_sphericalOrigin.size() << " given.");
463 }
464
465 //Convert plane definition to mm since Geant4 is of course using other units
466 for (double& value : m_customPlane) value /= Unit::mm;
467 for (double& value : m_sphericalOrigin) value /= Unit::mm;
468 for (double& value : m_rayOrigin) value /= Unit::mm;
469 if (m_rayDirection) {
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");
472 return;
473 }
474 for (double& value : *m_rayDirection) value /= Unit::mm;
475 }
476}
477
479{
480 G4EventManager* eventManager = G4EventManager::GetEventManager();
481
482 //First we save all user actions
483 G4UserEventAction* vanillaEventAction = eventManager->GetUserEventAction();
484 G4UserStackingAction* vanillaStackingAction = eventManager->GetUserStackingAction();
485 G4UserTrackingAction* vanillaTrackingAction = eventManager->GetUserTrackingAction();
486 G4UserSteppingAction* vanillaSteppingAction = eventManager->GetUserSteppingAction();
487
488 //Then we clear the user actions
489 eventManager->SetUserAction((G4UserEventAction*)0);
490 eventManager->SetUserAction((G4UserStackingAction*)0);
491 eventManager->SetUserAction((G4UserTrackingAction*)0);
492
493 //We create our own stepping actions which will
494 //- create the vectors for shooting rays
495 //- record the material budget for each ray
496 vector<MaterialScanBase*> scans;
497 if (m_doSpherical) {
498 G4ThreeVector origin(m_sphericalOrigin[0], m_sphericalOrigin[1], m_sphericalOrigin[2]);
499 scans.push_back(new MaterialScanSpherical(m_rootFile, origin, m_spherical, m_doCosTheta));
500 }
501 if (m_doPlanar) {
502 G4ThreeVector origin(0, 0, 0);
503 G4ThreeVector uDir;
504 G4ThreeVector vDir;
505 if (m_planeName == "custom") {
506 origin.set(m_customPlane[0], m_customPlane[1], m_customPlane[2]);
507 uDir.set(m_customPlane[3], m_customPlane[4], m_customPlane[5]);
508 vDir.set(m_customPlane[6], m_customPlane[7], m_customPlane[8]);
509 } else {
510 uDir = getAxis(m_planeName[0]);
511 vDir = getAxis(m_planeName[1]);
512 }
513 scans.push_back(new MaterialScanPlanar(m_rootFile, origin, uDir, vDir, m_planar));
514 }
515
516 if (m_doRay) {
517 G4ThreeVector origin(m_sphericalOrigin[0], m_sphericalOrigin[1], m_sphericalOrigin[2]);
518 G4ThreeVector dir;
519 if (m_rayDirection) {
520 for (int i = 0; i < 3; ++i) dir[0] = (*m_rayDirection)[0];
521 dir = dir.unit();
522 } else {
523 dir.setRThetaPhi(1., m_rayTheta * Unit::deg, m_rayPhi * Unit::deg);
524 }
525 scans.push_back(new MaterialScanRay(m_rootFile, origin, dir, m_rayOpening * Unit::deg, m_rayCount,
528 }
529
530 //Do each configured scan
531 for (MaterialScanBase* scan : scans) {
532 //Set the Scan as steppingaction to see material
533 eventManager->SetUserAction(scan);
534 //Now we can scan
535 G4RayShooter rayShooter;
536 G4ThreeVector origin;
537 G4ThreeVector direction;
538 int maxRays = scan->getNRays();
539 int curRay(0);
540 int lastPercent(-1);
541 double start = Utils::getClock();
542 //Create one event per ray and process it
543 while (scan->createNext(origin, direction)) {
544 G4Event* event = new G4Event(++curRay);
545 rayShooter.Shoot(event, origin, direction);
546 eventManager->ProcessOneEvent(event);
547 delete event;
548
549 //Show progress
550 int donePercent = 100 * curRay / maxRays;
551 if (donePercent > lastPercent) {
552 double perRay = (Utils::getClock() - start) / curRay;
553 double eta = perRay * (maxRays - curRay);
554 B2INFO(boost::format("%s Scan: %3d%%, %.3f ms per ray, ETA: %.2f seconds")
555 % scan->getName() % donePercent
556 % (perRay / Unit::ms) % (eta / Unit::s));
557 lastPercent = donePercent;
558 }
559 }
560 }
561 m_rootFile->Write();
562 // free the scan objects which also deletes the histograms ...
563 for (MaterialScanBase* scan : scans) {
564 delete scan;
565 }
566 // and delete the file, which we could not do before because that would also delete the histograms ...
567 m_rootFile->Close();
568 delete m_rootFile;
569 m_rootFile = nullptr;
570 //And now we reset the user actions
571 eventManager->SetUserAction(vanillaEventAction);
572 eventManager->SetUserAction(vanillaStackingAction);
573 eventManager->SetUserAction(vanillaTrackingAction);
574 eventManager->SetUserAction(vanillaSteppingAction);
575}
576
577G4ThreeVector MaterialScanModule::getAxis(char name)
578{
579 switch (name) {
580 case 'x': return G4ThreeVector(1, 0, 0);
581 case 'y': return G4ThreeVector(0, 1, 0);
582 case 'z': return G4ThreeVector(0, 0, 1);
583 default:
584 B2FATAL("Unknown axis: " << name);
585 }
586 return G4ThreeVector(0, 0, 0);
587}
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:150
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:158
double m_curDepth
Tracklength of the current Ray.
Definition: MaterialScan.h:156
double m_curU
Current value of the parametetr u.
Definition: MaterialScan.h:148
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
ScanParams m_params
Parameters for the scan.
Definition: MaterialScan.h:146
Base class for Material Scans.
Definition: MaterialScan.h:33
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:73
static constexpr int c_maxZeroStepsNudge
maximum number of consecutive zero steps before nudging the track along
Definition: MaterialScan.h:69
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
static constexpr int c_maxZeroStepsKill
maximum number of consecutive zero steps before killing the track
Definition: MaterialScan.h:71
std::string m_name
Name of the scan, will be prefixed to all histogram names.
Definition: MaterialScan.h:62
static constexpr double c_zeroTolerance
maximum Step length to be considered zero
Definition: MaterialScan.h:67
double m_rayTheta
Theta direction of the ray if custom direction is not set.
Definition: MaterialScan.h:351
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:357
bool m_raySplitByMaterials
If true Split by materials instead of regions when doing ray scan.
Definition: MaterialScan.h:372
double m_rayPhi
Phi direction of the ray if custom direction is not set.
Definition: MaterialScan.h:353
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:331
std::string m_planeName
Name of the plane to use for scanning.
Definition: MaterialScan.h:335
std::vector< std::string > m_rayIgnoredMaterials
Materials ignored when ray scanning.
Definition: MaterialScan.h:347
unsigned int m_rayCount
Number of rays to shoot (if opening angle is >0)
Definition: MaterialScan.h:361
std::optional< std::vector< double > > m_rayDirection
Direction of the ray, cannot be set at the same time as Theta and Phi.
Definition: MaterialScan.h:349
void beginRun() override
Do the actual scanning.
MaterialScan2D::ScanParams m_spherical
Scan parameters for the spherical scan.
Definition: MaterialScan.h:337
bool m_doSpherical
Wether or not to do a spherical scan.
Definition: MaterialScan.h:363
bool m_doPlanar
Wether or not to do a planar scan.
Definition: MaterialScan.h:365
bool m_doRay
Perform a material ray scan: scan along a certain direction with a given opening angle and plot mater...
Definition: MaterialScan.h:370
std::vector< double > m_sphericalOrigin
Origin for spherical scan.
Definition: MaterialScan.h:343
std::string m_filename
Filename for the ROOT output.
Definition: MaterialScan.h:333
MaterialScan2D::ScanParams m_planar
Scan parameters for the planar scan.
Definition: MaterialScan.h:339
double m_rayOpening
Opening angle of the ray.
Definition: MaterialScan.h:355
double m_rayMaxDepth
Max depth of the ray before stopping.
Definition: MaterialScan.h:359
std::vector< double > m_customPlane
Custom plane definition if m_planName is "custom".
Definition: MaterialScan.h:341
bool m_doCosTheta
Perform the spherical scan uniform in cos(theta) instead of theta.
Definition: MaterialScan.h:367
std::vector< double > m_rayOrigin
Origin of the ray if a ray scan is performed
Definition: MaterialScan.h:345
Specific implementaion of MaterialScan to scan parallel to a given plane.
Definition: MaterialScan.h:202
G4ThreeVector m_origin
Origin of the scan plane.
Definition: MaterialScan.h:222
G4ThreeVector m_dirV
v direction of the scan plane
Definition: MaterialScan.h:226
G4ThreeVector m_dirU
u direction of the scan plane
Definition: MaterialScan.h:224
G4ThreeVector m_dirW
direction perpendicluar to u and v
Definition: MaterialScan.h:228
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:239
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:276
std::map< std::string, std::unique_ptr< TH1D > > m_regions
Map holding pointers to all created histograms.
Definition: MaterialScan.h:274
bool m_splitByMaterials
If true Split by materials instead of regions.
Definition: MaterialScan.h:298
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:292
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:296
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:278
double m_curDepth
Current depth of the current ray.
Definition: MaterialScan.h:286
double m_maxDepth
Maximum depth for each ray after which it will be stopped.
Definition: MaterialScan.h:284
double m_opening
Opening angle in radian.
Definition: MaterialScan.h:280
std::set< std::string > m_ignoredMaterials
Materials ignored when scanning.
Definition: MaterialScan.h:272
double m_sampleDepth
The ray length after which to sample.
Definition: MaterialScan.h:282
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:294
Specific implementation of MaterialScan to do Spherical scanning.
Definition: MaterialScan.h:164
G4ThreeVector m_origin
Origin for the spherical scan.
Definition: MaterialScan.h:189
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:192
virtual void event()
This method is the core of the module.
Definition: Module.h:157
static const double mm
[millimeters]
Definition: Unit.h:70
static const double deg
degree to radians
Definition: Unit.h:109
static const double ms
[millisecond]
Definition: Unit.h:96
static const double s
[second]
Definition: Unit.h:95
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.
Definition: Module.h:650
double getClock()
Return current value of the real-time clock.
Definition: Utils.cc:66
Abstract base class for different kinds of events.
STL namespace.
Helper struct to Store Parameters of a Scan.
Definition: MaterialScan.h:83
double maxU
Maximum u value to scan.
Definition: MaterialScan.h:93
bool splitByMaterials
If true, split output by Materials (otherwise by region)
Definition: MaterialScan.h:103
double minV
Minimum v value to scan.
Definition: MaterialScan.h:95
double maxDepth
Maximum depth of the scan.
Definition: MaterialScan.h:99
int nU
Number of rays along u coordinate.
Definition: MaterialScan.h:87
double minU
Minimum u value to scan.
Definition: MaterialScan.h:91
double maxV
Maximum v value to scan.
Definition: MaterialScan.h:97
std::vector< std::string > ignoredMaterials
Names of ignored Materials.
Definition: MaterialScan.h:101
int nV
Number of rays along v coordinate.
Definition: MaterialScan.h:89