BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSScorerFactory.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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 "BDSScorerFactory.hh"
22#include "BDSScorerInfo.hh"
23#include "BDSPSCellFluxScaledPerParticle3D.hh"
24#include "BDSPSPopulationScaled.hh"
25#include "BDSPSBLMEnergyDeposit.hh"
26#include "BDSPSCellFluxScaled3D.hh"
27#include "BDSSDFilterAnd.hh"
28#include "BDSSDFilterLogicalVolume.hh"
29#include "BDSSDFilterMaterial.hh"
30#include "BDSSDFilterPrimary.hh"
31#include "BDSSDFilterTime.hh"
32#include "BDSScorerType.hh"
33#include "BDSUtilities.hh"
34
35#include "G4PSCellCharge.hh"
36#include "G4PSCellCharge3D.hh"
37#include "G4PSCellFlux.hh"
38#include "G4PSCellFlux3D.hh"
39#include "BDSPSCellFlux4D.hh"
40#include "G4PSDoseDeposit.hh"
41#include "G4PSDoseDeposit3D.hh"
42#include "G4PSEnergyDeposit.hh"
43#include "G4PSEnergyDeposit3D.hh"
44#include "G4PSPopulation.hh"
45#include "G4PSPopulation3D.hh"
46#include "G4SDParticleWithEnergyFilter.hh"
47#include "G4String.hh"
48#include "G4Types.hh"
49#include "G4VPrimitiveScorer.hh"
50
51BDSScorerFactory::BDSScorerFactory()
52{;}
53
54G4VPrimitiveScorer* BDSScorerFactory::CreateScorer(const BDSScorerInfo* info,
55 const BDSHistBinMapper* mapper,
56 G4double* unit,
57 G4LogicalVolume* worldLV)
58{
59 // here we create the scorer with the information from BDSScorerInfo.
60 G4VPrimitiveScorer* primitiveScorer = GetAppropriateScorer(*info, mapper, unit);
61
62 BDSSDFilterAnd* filter = CreateFilter(info->name + "_scorer_filter", info, worldLV);
63 if (filter)
64 {primitiveScorer->SetFilter(filter);}
65
66 return primitiveScorer;
67}
68
70 const BDSHistBinMapper* mapper,
71 G4double* unit)
72{
73 G4VPrimitiveScorer* result = nullptr;
74 switch (info.scorerType.underlying())
75 {// if adding new 3D ones, remember to update the 3D upgrade mapping in BDSScorerInfo
76 case BDSScorerType::cellcharge:
77 {result = new G4PSCellCharge(info.name); break;}
78 case BDSScorerType::cellcharge3d:
79 {result = new G4PSCellCharge3D(info.name); break;}
80 case BDSScorerType::depositeddose:
81 {result = new G4PSDoseDeposit(info.name); break;}
82 case BDSScorerType::depositeddose3d:
83 {result = new G4PSDoseDeposit3D(info.name); break;}
84 case BDSScorerType::depositedenergy:
85 {result = new G4PSEnergyDeposit(info.name, "GeV"); break;}
86 case BDSScorerType::depositedenergyblm:
87 {result = new BDSPSBLMEnergyDeposit(info.name, "GeV"); break;}
88 case BDSScorerType::depositedenergy3d:
89 {result = new G4PSEnergyDeposit3D(info.name, "GeV"); break;}
90 case BDSScorerType::population:
91 {
92 G4PSPopulation* scorer = new G4PSPopulation(info.name);
93 scorer->Weighted(true);
94 result = scorer;
95 break;
96 }
97 case BDSScorerType::population3d:
98 {
99 G4PSPopulation3D* scorer = new G4PSPopulation3D(info.name);
100 scorer->Weighted(true);
101 result = scorer;
102 break;
103 }
104 case BDSScorerType::populationscaled:
105 {result = new BDSPSPopulationScaled(info.name, info.pathname); break;}
106 case BDSScorerType::cellflux:
107 {
108 G4PSCellFlux* scorer = new G4PSCellFlux(info.name, "percm2");
109 scorer->Weighted(true);
110 result = scorer;
111 break;
112 }
113 case BDSScorerType::cellflux3d:
114 {
115 G4PSCellFlux3D* scorer = new G4PSCellFlux3D(info.name, "percm2");
116 scorer->Weighted(true);
117 result = scorer;
118 break;
119 }
120 case BDSScorerType::cellflux4d:
121 {
122 BDSPSCellFlux4D* scorer = new BDSPSCellFlux4D(info.name, mapper, "percm2");
123 scorer->Weighted(true);
124 result = scorer;
125 break;
126 }
127 case BDSScorerType::cellfluxscaledperparticle3d:
128 {result = new BDSPSCellFluxScaledPerParticle3D(info.name, mapper, info.pathname); break;}
129 case BDSScorerType::cellfluxscaled3d:
130 {result = new BDSPSCellFluxScaled3D(info.name, mapper, info.filename, "percm2"); break;}
131 case BDSScorerType::cellfluxscaled:
132 case BDSScorerType::cellfluxscaledperparticle:
133 {
134 throw BDSException(__METHOD_NAME__, "unimplemented scorer \"" + info.scorerType.ToString() + "\"");
135 break;
136 }
137 default:
138 {
139 throw BDSException(__METHOD_NAME__, "unknown scorer type \"" + info.scorerType.ToString() + "\"");
140 break;
141 }
142 }
143 if (unit)
144 {*unit = result->GetUnitValue();}
145 return result;
146}
147
149 const BDSScorerInfo* info,
150 G4LogicalVolume* worldLV) const
151{
152 BDSSDFilterAnd* result = new BDSSDFilterAnd(name, /*ownsFilters=*/true);
153
154 if (info->particle)
155 {
156 G4String particleName = info->particle->GetParticleName();
157 auto pwkef = new G4SDParticleWithEnergyFilter("particle_filter",
158 info->minimumKineticEnergy,
159 info->maximumKineticEnergy);
160 pwkef->add(particleName);
161 result->RegisterFilter(pwkef);
162 }
163 if (BDS::IsFinite(info->maximumTime) || BDS::IsFinite(info->minimumTime))
164 {
165 auto timeFilter = new BDSSDFilterTime("time_filter",
166 info->minimumTime,
167 info->maximumTime);
168 result->RegisterFilter(timeFilter);
169 }
170 if (!(info->materialsToInclude.empty()))
171 {
172 auto materialFilterInc = new BDSSDFilterMaterial("material_filter_include",
173 info->materialsToInclude,
174 /*inclusive*/true);
175 result->RegisterFilter(materialFilterInc);
176 }
177 if (!(info->materialsToExclude.empty()))
178 {
179 auto materialFilterExc = new BDSSDFilterMaterial("material_filter_exclude",
180 info->materialsToExclude,
181 /*inclusive*/false);
182 result->RegisterFilter(materialFilterExc);
183 }
184 if (info->worldVolumeOnly)
185 {
186 auto worldLVFilter = new BDSSDFilterLogicalVolume("world_lv_only", worldLV);
187 result->RegisterFilter(worldLVFilter);
188 }
189 if (info->primariesOnly)
190 {
191 auto primaryFilter = new BDSSDFilterPrimary("primary_filter");
192 result->RegisterFilter(primaryFilter);
193 }
194
195 // if we didn't register any filters, just delete it
196 if (result->size() == 0)
197 {
198 delete result;
199 return nullptr;
200 }
201 return result;
202}
General exception with possible name of object and message.
Definition: BDSException.hh:35
Mapping from axis indices to 1D index.
Primitive scorer for energy deposition in BLMs.
Primitive scorer for cell flux in a 4D mesh.
Primitive scorer for cell flux in a 3D mesh with a conversion factor.
Primitive scorer for a 3D mesh with a conversion factor.
Primitive scorer for population in a volume with a conversion factor based on angle and kinetic energ...
Filter that applies AND to a vector of filters.
void RegisterFilter(G4VSDFilter *filterIn)
Register the filter.
size_t size() const
Accessor.
SD filter for a particular volume.
SD filter for a particular volume.
Filter for only primary particles.
Filter for time value for a sensitive detector.
BDSSDFilterAnd * CreateFilter(const G4String &name, const BDSScorerInfo *info, G4LogicalVolume *worldLV=nullptr) const
Create a combined filter with AND logic for the scorer.
G4VPrimitiveScorer * CreateScorer(const BDSScorerInfo *info, const BDSHistBinMapper *mapper, G4double *unit=nullptr, G4LogicalVolume *worldLV=nullptr)
Main function to create a scorer.
G4VPrimitiveScorer * GetAppropriateScorer(const BDSScorerInfo &info, const BDSHistBinMapper *mapper, G4double *unit=nullptr)
Construct the primitive scorer required.
Recipe class for scorer. Checks values.
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)
G4ParticleDefinition * particle
Particle filter.
type underlying() const
return underlying value (can be used in switch statement)
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())