BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSScorerInfo.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSDebug.hh"
20#include "BDSException.hh"
21#include "BDSMaterials.hh"
22#include "BDSScorerInfo.hh"
23#include "BDSScorerType.hh"
24#include "BDSUtilities.hh"
25
26#include "parser/scorer.h"
27
28#include "G4ParticleTable.hh"
29#include "G4String.hh"
30#include "G4Types.hh"
31
32
33#include "CLHEP/Units/SystemOfUnits.h"
34
35#include <map>
36#include <string>
37#include <vector>
38
40 G4bool upgradeTo3D):
41 particle(nullptr)
42{
43 const std::map<std::string, std::string> replacements = {
44 {"cellcharge", "cellcharge3d"},
45 {"depositeddose", "depositeddose3d"},
46 {"depositedenergy", "depositedenergy3d"},
47 {"population", "population3d"},
48 {"cellflux", "cellflux3d"},
49 {"cellfluxscaled", "cellfluxscaled3d"},
50 {"cellfluxscaledperparticle", "cellfluxscaledperparticle3d"}
51 };
52
53 std::string scorerTypeNameOriginal = scorer.type;
54 std::string scorerTypeName = scorerTypeNameOriginal; // default copy
55 if (upgradeTo3D)
56 {
57 auto search = replacements.find(scorerTypeNameOriginal);
58 if (search != replacements.end())
59 {scorerTypeName = search->second;}
60 else if (!BDS::StrContains(scorerTypeNameOriginal, "3d") && !BDS::StrContains(scorerTypeNameOriginal, "4d") )
61 {
62 G4String msg = "scorer type \"" + scorerTypeNameOriginal + "\" specified which does not";
63 msg += "\ncontain the suffix \"3d\" or \"4d\", and there is no known mapping to the required";
64 msg += "\n3d one. Check the scorer name or specify it with the suffix explicitly.";
65 msg += "\nAvailable 3d scorers are:\n";
66 for (const auto& fs : replacements)
67 {msg += "\"" + fs.first + "\" -> \"" + fs.second + "\"\n";}
68 // now check if it's even a valid scorer type
69 try
70 {BDS::DetermineScorerType(G4String(scorerTypeNameOriginal));}
71 catch (const BDSException& e)
72 {
73 msg += "\n";
74 msg += e.what();
75 }
76 throw BDSException(__METHOD_NAME__, msg);
77 }
78 }
79
80 scorerType = BDS::DetermineScorerType(G4String(scorerTypeName));
81 name = scorer.name;
82 minimumKineticEnergy = scorer.minimumKineticEnergy * CLHEP::GeV;
83 maximumKineticEnergy = scorer.maximumKineticEnergy * CLHEP::GeV;
84 filename = scorer.conversionFactorFile;
85 pathname = scorer.conversionFactorPath;
86 minimumTime = scorer.minimumTime*CLHEP::second;
87 maximumTime = scorer.maximumTime*CLHEP::second;
88 worldVolumeOnly = scorer.scoreWorldVolumeOnly;
89 primariesOnly = scorer.scorePrimariesOnly;
90
91 if (scorer.particlePDGID != 0)
92 {
93 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
94 particle = particleTable->FindParticle(scorer.particlePDGID);
96
97 }
98 if (!(scorer.particleName.empty()))
99 {
100 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
101 particle = particleTable->FindParticle(scorer.particleName);
102 CheckParticle(particle, scorer.name);
103 }
104
105 for (const auto& mi : BDS::SplitOnWhiteSpace(scorer.materialToInclude))
106 {materialsToInclude.push_back(BDSMaterials::Instance()->GetMaterial(mi));}
107
108 for (const auto& me : BDS::SplitOnWhiteSpace(scorer.materialToExclude))
109 {materialsToExclude.push_back(BDSMaterials::Instance()->GetMaterial(me));}
110}
111
112void BDSScorerInfo::CheckParticle(G4ParticleDefinition* particleIn,
113 const G4String& nameIn)
114{
115 if (!particleIn)
116 {
117 G4cout << __METHOD_NAME__ << "Note, only 1 particle can be specified." << G4endl;
118 throw BDSException(__METHOD_NAME__, "Particle not found for scorer " + nameIn);
119 }
120}
General exception with possible name of object and message.
Definition: BDSException.hh:35
const char * what() const noexcept override
Override message in std::exception.
Definition: BDSException.hh:55
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:39
BDSScorerType scorerType
Scorer type.
G4bool worldVolumeOnly
Which materials to exclude for scoring.
std::vector< G4Material * > materialsToExclude
Which materials to include for scoring.
G4String name
Scorer name.
G4String filename
Name of the conversion factor file.
G4String pathname
Path of the conversion factor file (for ambient dose)
void CheckParticle(G4ParticleDefinition *particleIn, const G4String &nameIn)
Utility function to check valid pointer and throw exception if not.
BDSScorerInfo()=delete
No default constructor as unused.
G4ParticleDefinition * particle
Particle filter.
Scorer class for parser.
Definition: scorer.h:37
std::string name
Name.
Definition: scorer.h:39
int particlePDGID
PDG ID code for particle.
Definition: scorer.h:42
std::string type
Type of the scorer, ie forumula.
Definition: scorer.h:40
std::string particleName
Particle name as a string.
Definition: scorer.h:41
BDSScorerType DetermineScorerType(G4String scorerType)
function to determine the enum type of the cavity geometry (case-insensitive)
G4bool StrContains(const G4String &str, const G4String &test)
Utility function to simplify lots of syntax changes for pedantic g4 changes.
Definition: BDSUtilities.cc:66
std::vector< G4String > SplitOnWhiteSpace(const G4String &input)
Split a string on whitespace and return a vector of these 'words'.