BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSPSCellFluxScaledPerParticle3D.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 "BDSHistBinMapper.hh"
22#include "BDSScorerConversionLoader.hh"
23#include "BDSPSCellFluxScaledPerParticle3D.hh"
24#include "BDSUtilities.hh"
25
26#include "globals.hh"
27#include "G4PhysicsVector.hh"
28#include "G4SystemOfUnits.hh"
29#include "G4VSolid.hh"
30#include "G4VPhysicalVolume.hh"
31#include "G4VPVParameterisation.hh"
32
33#include <fstream>
34#include <map>
35#include <string>
36
37#ifdef USE_GZSTREAM
38#include "src-external/gzstream/gzstream.h"
39#endif
40
41BDSPSCellFluxScaledPerParticle3D::BDSPSCellFluxScaledPerParticle3D(const G4String& scorerName,
42 const BDSHistBinMapper* mapperIn,
43 const G4String& pathname,
44 const G4String& unitIn,
45 G4int ni,
46 G4int nj,
47 G4int nk,
48 G4int depi,
49 G4int depj,
50 G4int depk):
51 BDSPSCellFluxScaled3D(scorerName, mapperIn, unitIn, ni, nj, nk, depi, depj, depk)
52{
53 G4String filePath = BDS::GetFullPath(pathname);
54 if (filePath.back() != '/')
55 {filePath += '/';}
56
57 const std::map<std::string, G4int> files = {{"protons.dat", 2212},
58 {"neutrons.dat", 2112},
59 {"photons.dat", 22},
60 {"electrons.dat", 11},
61 {"positrons.dat", -11}};
62
63 G4cout << "Scorer \"" << GetName() << "\" - adding conversionFiles:" << G4endl;
65 for (const auto& filePDG : files)
66 {
67 G4String uncompressedFile = filePath + filePDG.first;
68 G4String compressedFile = filePath + filePDG.first+".gz";
69 if (BDS::FileExists(uncompressedFile))
70 {
71 conversionFactors[filePDG.second] = loader.Load(uncompressedFile);
72 G4cout << "Adding: " << uncompressedFile << G4endl;
73 }
74 else if (BDS::FileExists(compressedFile))
75 {
76#ifdef USE_GZSTREAM
78 conversionFactors[filePDG.second] = loaderC.Load(compressedFile);
79 G4cout << "Adding: " << compressedFile << G4endl;
80#else
81 throw BDSException(__METHOD_NAME__, "Compressed file loading - but BDSIM not compiled with ZLIB.");
82#endif
83 }
84 }
85 if (conversionFactors.empty())
86 {
87 G4String message = "no conversion files for scorer \"" + scorerName + "\" found. Please specify one of:\n";
88 for (const auto &kv : files)
89 {message += kv.first + "\n";}
90 message += "In the specified directory \"" + pathname + "\"\n";
91 message += "Would result in a factor of 0 for all particles and all energies.";
92 throw BDSException(__METHOD_NAME__, message);
93 }
94}
95
96BDSPSCellFluxScaledPerParticle3D::~BDSPSCellFluxScaledPerParticle3D()
97{
98 for (auto cf : conversionFactors)
99 {delete cf.second;}
100}
101
102G4double BDSPSCellFluxScaledPerParticle3D::GetConversionFactor(G4int particleID, G4double kineticEnergy) const
103{
104 auto search = conversionFactors.find(particleID);
105 if (search != conversionFactors.end())
106 {return search->second->Value(kineticEnergy);}
107 else
108 {return 0;}
109}
General exception with possible name of object and message.
Definition: BDSException.hh:35
Mapping from axis indices to 1D index.
Primitive scorer for cell flux in a 3D mesh with a conversion factor.
Loader for scoring conversion tables as function of energy.
G4PhysicsVector * Load(const G4String &fileName, G4bool silent=false)
Load the file.
G4String GetFullPath(G4String filename, bool excludeNameFromPath=false, bool useCWDForPrefix=false)
G4bool FileExists(const G4String &filename)
Checks if filename exists.
const char * GetName(int)
Name of element.
Definition: python.cc:38