CMS 3D CMS Logo

NuclearInteractionFTFSimulator.cc
Go to the documentation of this file.
1 // system headers
2 #include <mutex>
3 
4 // Framework Headers
6 
7 // Fast Sim headers
12 
13 // Geant4 headers
14 #include "G4ParticleDefinition.hh"
15 #include "G4DynamicParticle.hh"
16 #include "G4Track.hh"
17 #include "G4Step.hh"
18 #include "G4StepPoint.hh"
19 #include "G4TheoFSGenerator.hh"
20 #include "G4FTFModel.hh"
21 #include "G4ExcitedStringDecay.hh"
22 #include "G4LundStringFragmentation.hh"
23 #include "G4GeneratorPrecompoundInterface.hh"
24 #include "G4CascadeInterface.hh"
25 #include "G4DiffuseElastic.hh"
26 
27 #include "G4Proton.hh"
28 #include "G4Neutron.hh"
29 #include "G4PionPlus.hh"
30 #include "G4PionMinus.hh"
31 #include "G4AntiProton.hh"
32 #include "G4AntiNeutron.hh"
33 #include "G4KaonPlus.hh"
34 #include "G4KaonMinus.hh"
35 #include "G4KaonZeroLong.hh"
36 #include "G4KaonZeroShort.hh"
37 #include "G4KaonZero.hh"
38 #include "G4AntiKaonZero.hh"
39 #include "G4GenericIon.hh"
40 
41 #include "G4Lambda.hh"
42 #include "G4OmegaMinus.hh"
43 #include "G4SigmaMinus.hh"
44 #include "G4SigmaPlus.hh"
45 #include "G4SigmaZero.hh"
46 #include "G4XiMinus.hh"
47 #include "G4XiZero.hh"
48 #include "G4AntiLambda.hh"
49 #include "G4AntiOmegaMinus.hh"
50 #include "G4AntiSigmaMinus.hh"
51 #include "G4AntiSigmaPlus.hh"
52 #include "G4AntiSigmaZero.hh"
53 #include "G4AntiXiMinus.hh"
54 #include "G4AntiXiZero.hh"
55 #include "G4AntiAlpha.hh"
56 #include "G4AntiDeuteron.hh"
57 #include "G4AntiTriton.hh"
58 #include "G4AntiHe3.hh"
59 
60 #include "G4Material.hh"
61 #include "G4DecayPhysics.hh"
62 #include "G4ParticleTable.hh"
63 #include "G4IonTable.hh"
64 #include "G4ProcessManager.hh"
65 #include "G4PhysicsLogVector.hh"
66 #include "G4SystemOfUnits.hh"
67 
68 static std::once_flag initializeOnce;
69 [[cms::thread_guard("initializeOnce")]] const G4ParticleDefinition* NuclearInteractionFTFSimulator::theG4Hadron[] = {nullptr};
70 [[cms::thread_guard("initializeOnce")]] int NuclearInteractionFTFSimulator::theId[] = {0};
71 
72 const double fact = 1.0/CLHEP::GeV;
73 
74 // inelastic interaction length corrections per particle and energy
76  {1.0872, 1.1026, 1.111, 1.111, 1.0105, 0.97622, 0.9511, 0.9526, 0.97591, 0.99277, 1.0099, 1.015, 1.0217, 1.0305, 1.0391, 1.0438, 1.0397, 1.0328, 1.0232, 1.0123, 1.0},
77  {1.0416, 1.1044, 1.1467, 1.1273, 1.026, 0.99085, 0.96572, 0.96724, 0.99091, 1.008, 1.0247, 1.0306, 1.0378, 1.0427, 1.0448, 1.0438, 1.0397, 1.0328, 1.0232, 1.0123, 1.0},
78  {0.5308, 0.53589, 0.67059, 0.80253, 0.82341, 0.79083, 0.85967, 0.90248, 0.93792, 0.9673, 1.0034, 1.022, 1.0418, 1.0596, 1.0749, 1.079, 1.0704, 1.0576, 1.0408, 1.0214, 1.0},
79  {0.49107, 0.50571, 0.64149, 0.77209, 0.80472, 0.78166, 0.83509, 0.8971, 0.93234, 0.96154, 0.99744, 1.0159, 1.0355, 1.0533, 1.0685, 1.0732, 1.0675, 1.0485, 1.0355, 1.0191, 1.0},
80  {1.9746, 1.7887, 1.5645, 1.2817, 1.0187, 0.95216, 0.9998, 1.035, 1.0498, 1.0535, 1.0524, 1.0495, 1.0461, 1.0424, 1.0383, 1.0338, 1.0287, 1.0228, 1.0161, 1.0085, 1.0},
81  {0.46028, 0.59514, 0.70355, 0.70698, 0.62461, 0.65103, 0.71945, 0.77753, 0.83582, 0.88422, 0.92117, 0.94889, 0.96963, 0.98497, 0.99596, 1.0033, 1.0075, 1.0091, 1.0081, 1.005, 1.0},
82  {0.75016, 0.89607, 0.97185, 0.91083, 0.77425, 0.77412, 0.8374, 0.88848, 0.93104, 0.96174, 0.98262, 0.99684, 1.0065, 1.0129, 1.0168, 1.0184, 1.018, 1.0159, 1.0121, 1.0068, 1.0},
83  {0.75016, 0.89607, 0.97185, 0.91083, 0.77425, 0.77412, 0.8374, 0.88848, 0.93104, 0.96174, 0.98262, 0.99684, 1.0065, 1.0129, 1.0168, 1.0184, 1.018, 1.0159, 1.0121, 1.0068, 1.0},
84  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
85  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
86  {1.1006, 1.1332, 1.121, 1.1008, 1.086, 1.077, 1.0717, 1.0679, 1.0643, 1.0608, 1.057, 1.053, 1.0487, 1.0441, 1.0392, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
87  {1.1318, 1.1255, 1.1062, 1.0904, 1.0802, 1.0742, 1.0701, 1.0668, 1.0636, 1.0602, 1.0566, 1.0527, 1.0485, 1.044, 1.0391, 1.0337, 1.028, 1.0217, 1.015, 1.0078, 1.0},
88  {1.1094, 1.1332, 1.1184, 1.0988, 1.0848, 1.0765, 1.0714, 1.0677, 1.0642, 1.0607, 1.0569, 1.053, 1.0487, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
89  {1.1087, 1.1332, 1.1187, 1.099, 1.0849, 1.0765, 1.0715, 1.0677, 1.0642, 1.0607, 1.057, 1.053, 1.0487, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
90  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
91  {1.1192, 1.132, 1.1147, 1.0961, 1.0834, 1.0758, 1.0711, 1.0674, 1.064, 1.0606, 1.0569, 1.0529, 1.0486, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
92  {1.1188, 1.1321, 1.1149, 1.0963, 1.0834, 1.0758, 1.0711, 1.0675, 1.0641, 1.0606, 1.0569, 1.0529, 1.0486, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
93  {0.50776, 0.5463, 0.5833, 0.61873, 0.65355, 0.68954, 0.72837, 0.7701, 0.81267, 0.85332, 0.89037, 0.92329, 0.95177, 0.97539, 0.99373, 1.0066, 1.014, 1.0164, 1.0144, 1.0087, 1.0},
94  {0.50787, 0.5464, 0.58338, 0.6188, 0.65361, 0.6896, 0.72841, 0.77013, 0.8127, 0.85333, 0.89038, 0.92329, 0.95178, 0.9754, 0.99373, 1.0066, 1.014, 1.0164, 1.0144, 1.0087, 1.0},
95  {1.1006, 1.1332, 1.121, 1.1008, 1.086, 1.077, 1.0717, 1.0679, 1.0643, 1.0608, 1.057, 1.053, 1.0487, 1.0441, 1.0392, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
96  {1.1318, 1.1255, 1.1062, 1.0904, 1.0802, 1.0742, 1.0701, 1.0668, 1.0636, 1.0602, 1.0566, 1.0527, 1.0485, 1.044, 1.0391, 1.0337, 1.028, 1.0217, 1.015, 1.0078, 1.0},
97  {1.1094, 1.1332, 1.1184, 1.0988, 1.0848, 1.0765, 1.0714, 1.0677, 1.0642, 1.0607, 1.0569, 1.053, 1.0487, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
98  {1.1087, 1.1332, 1.1187, 1.099, 1.0849, 1.0765, 1.0715, 1.0677, 1.0642, 1.0607, 1.057, 1.053, 1.0487, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
99  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
100  {1.1192, 1.132, 1.1147, 1.0961, 1.0834, 1.0758, 1.0711, 1.0674, 1.064, 1.0606, 1.0569, 1.0529, 1.0486, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
101  {1.1188, 1.1321, 1.1149, 1.0963, 1.0834, 1.0758, 1.0711, 1.0675, 1.0641, 1.0606, 1.0569, 1.0529, 1.0486, 1.0441, 1.0391, 1.0338, 1.028, 1.0218, 1.015, 1.0078, 1.0},
102  {0.47677, 0.51941, 0.56129, 0.60176, 0.64014, 0.67589, 0.70891, 0.73991, 0.77025, 0.80104, 0.83222, 0.86236, 0.8901, 0.91518, 0.9377, 0.95733, 0.97351, 0.98584, 0.9942, 0.99879, 1.0},
103  {0.49361, 0.53221, 0.56976, 0.60563, 0.63954, 0.67193, 0.70411, 0.73777, 0.77378, 0.81114, 0.84754, 0.88109, 0.91113, 0.93745, 0.95974, 0.97762, 0.99081, 0.99929, 1.0033, 1.0034, 1.0},
104  {0.4873, 0.52744, 0.56669, 0.60443, 0.64007, 0.67337, 0.70482, 0.73572, 0.76755, 0.80086, 0.83456, 0.86665, 0.8959, 0.92208, 0.94503, 0.96437, 0.97967, 0.99072, 0.99756, 1.0005, 1.0},
105  {0.48729, 0.52742, 0.56668, 0.60442, 0.64006, 0.67336, 0.70482, 0.73571, 0.76754, 0.80086, 0.83455, 0.86665, 0.8959, 0.92208, 0.94503, 0.96437, 0.97967, 0.99072, 0.99756, 1.0005, 1.0},
106 };
107 
108 // elastic interaction length corrections per particle and energy
110  {0.58834, 1.1238, 1.7896, 1.4409, 0.93175, 0.80587, 0.80937, 0.83954, 0.87453, 0.91082, 0.94713, 0.98195, 1.0134, 1.0397, 1.0593, 1.071, 1.0739, 1.0678, 1.053, 1.03, 1.0},
111  {0.40938, 0.92337, 1.3365, 1.1607, 1.008, 0.82206, 0.81163, 0.79489, 0.82919, 0.91812, 0.96688, 1.0225, 1.0734, 1.0833, 1.0874, 1.0854, 1.0773, 1.0637, 1.0448, 1.0235, 1.0},
112  {0.43699, 0.42165, 0.46594, 0.64917, 0.85314, 0.80782, 0.83204, 0.91162, 1.0155, 1.0665, 1.0967, 1.1125, 1.1275, 1.1376, 1.1464, 1.1477, 1.1312, 1.1067, 1.0751, 1.039, 1.0},
113  {0.3888, 0.39527, 0.43921, 0.62834, 0.8164, 0.79866, 0.82272, 0.90163, 1.0045, 1.055, 1.0849, 1.1005, 1.1153, 1.1253, 1.134, 1.1365, 1.1255, 1.0895, 1.0652, 1.0348, 1.0},
114  {0.32004, 0.31119, 0.30453, 0.30004, 0.31954, 0.40148, 0.5481, 0.74485, 0.99317, 1.1642, 1.2117, 1.2351, 1.2649, 1.3054, 1.375, 1.4992, 1.4098, 1.3191, 1.2232, 1.118, 1.0},
115  {0.10553, 0.14623, 0.20655, 0.26279, 0.19996, 0.40125, 0.5139, 0.71271, 0.89269, 1.0108, 1.1673, 1.3052, 1.4149, 1.429, 1.4521, 1.4886, 1.4006, 1.3116, 1.2177, 1.1151, 1.0},
116  {0.106, 0.14692, 0.20755, 0.26257, 0.20089, 0.40236, 0.51452, 0.71316, 0.89295, 1.0109, 1.1673, 1.3053, 1.4149, 1.429, 1.4521, 1.4886, 1.4006, 1.3116, 1.2177, 1.1151, 1.0},
117  {0.31991, 0.31111, 0.30445, 0.30004, 0.31995, 0.40221, 0.54884, 0.74534, 0.99364, 1.1644, 1.2117, 1.2351, 1.265, 1.3054, 1.375, 1.4992, 1.4098, 1.3191, 1.2232, 1.118, 1.0},
118  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
119  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
120  {0.37579, 0.39922, 0.37445, 0.32631, 0.39002, 0.42161, 0.54251, 0.69127, 0.90332, 1.0664, 1.1346, 1.1481, 1.1692, 1.2036, 1.2625, 1.3633, 1.2913, 1.2215, 1.1516, 1.0788, 1.0},
121  {0.31756, 0.33409, 0.25339, 0.35525, 0.52989, 0.63382, 0.7453, 0.93505, 1.1464, 1.2942, 1.3161, 1.328, 1.3393, 1.3525, 1.374, 1.4051, 1.3282, 1.2523, 1.1745, 1.0916, 1.0},
122  {0.38204, 0.39694, 0.36502, 0.33367, 0.39229, 0.43119, 0.54898, 0.70169, 0.91004, 1.0696, 1.1348, 1.1483, 1.1694, 1.2038, 1.2627, 1.3632, 1.2913, 1.2215, 1.1516, 1.0788, 1.0},
123  {0.38143, 0.39716, 0.36609, 0.33294, 0.39207, 0.43021, 0.54834, 0.70066, 0.90945, 1.0693, 1.1348, 1.1482, 1.1694, 1.2038, 1.2627, 1.3632, 1.2913, 1.2215, 1.1516, 1.0788, 1.0},
124  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
125  {0.29564, 0.32645, 0.29986, 0.30611, 0.48808, 0.59902, 0.71207, 0.8832, 1.1164, 1.2817, 1.3154, 1.3273, 1.3389, 1.3521, 1.3736, 1.4056, 1.3285, 1.2524, 1.1746, 1.0916, 1.0},
126  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
127  {0.3265, 0.3591, 0.39232, 0.42635, 0.46259, 0.50365, 0.55244, 0.61014, 0.67446, 0.74026, 0.80252, 0.85858, 0.90765, 0.94928, 0.9827, 1.0071, 1.0221, 1.0279, 1.0253, 1.0155, 1.0},
128  {0.13808, 0.15585, 0.17798, 0.2045, 0.22596, 0.25427, 0.33214, 0.44821, 0.5856, 0.74959, 0.89334, 1.0081, 1.0964, 1.1248, 1.173, 1.2548, 1.1952, 1.1406, 1.0903, 1.0437, 1.0},
129  {0.20585, 0.23253, 0.26371, 0.28117, 0.30433, 0.35417, 0.44902, 0.58211, 0.73486, 0.90579, 1.0395, 1.1488, 1.2211, 1.2341, 1.2553, 1.2877, 1.2245, 1.1654, 1.1093, 1.0547, 1.0},
130  {0.2852, 0.32363, 0.31419, 0.35164, 0.45463, 0.54331, 0.66908, 0.81735, 0.98253, 1.1557, 1.2557, 1.3702, 1.4186, 1.401, 1.374, 1.3325, 1.2644, 1.1991, 1.1348, 1.0694, 1.0},
131  {0.20928, 0.23671, 0.2664, 0.28392, 0.30584, 0.35929, 0.45725, 0.5893, 0.74047, 0.9101, 1.0407, 1.1503, 1.2212, 1.2342, 1.2554, 1.2876, 1.2245, 1.1654, 1.1093, 1.0547, 1.0},
132  {0.11897, 0.13611, 0.15796, 0.1797, 0.21335, 0.26367, 0.34705, 0.46115, 0.6016, 0.7759, 0.91619, 1.0523, 1.1484, 1.1714, 1.2098, 1.2721, 1.2106, 1.1537, 1.1004, 1.0496, 1.0},
133  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
134  {0.26663, 0.30469, 0.32886, 0.33487, 0.41692, 0.51616, 0.63323, 0.78162, 0.95551, 1.1372, 1.2502, 1.3634, 1.4189, 1.4013, 1.3743, 1.3329, 1.2646, 1.1992, 1.1349, 1.0694, 1.0},
135  {0.16553, 0.19066, 0.21468, 0.23609, 0.30416, 0.38821, 0.49644, 0.63386, 0.80299, 0.99907, 1.1304, 1.2724, 1.3535, 1.3475, 1.3381, 1.3219, 1.2549, 1.191, 1.1287, 1.0659, 1.0},
136  {0.37736, 0.41414, 0.45135, 0.48843, 0.52473, 0.55973, 0.59348, 0.62696, 0.66202, 0.70042, 0.74241, 0.786, 0.82819, 0.86688, 0.90128, 0.93107, 0.95589, 0.97532, 0.98908, 0.99719, 1.0},
137  {0.34354, 0.37692, 0.4109, 0.44492, 0.47873, 0.51296, 0.54937, 0.59047, 0.63799, 0.69117, 0.74652, 0.7998, 0.84832, 0.89111, 0.92783, 0.95798, 0.98095, 0.99635, 1.0043, 1.0052, 1.0},
138  {0.36364, 0.39792, 0.43277, 0.4676, 0.50186, 0.53538, 0.56884, 0.604, 0.64308, 0.68729, 0.73544, 0.7842, 0.83019, 0.87156, 0.90777, 0.93854, 0.96346, 0.98209, 0.99421, 1, 1.0},
139  {0.36362, 0.39791, 0.43276, 0.46759, 0.50185, 0.53537, 0.56883, 0.604, 0.64307, 0.68728, 0.73544, 0.7842, 0.83019, 0.87156, 0.90777, 0.93854, 0.96346, 0.98209, 0.99421, 1, 1.0},
140 };
141 
142 // inelastic interaction length in Silicon at 1 TeV per particle
143 const double nuclInelLength[numHadrons] = {
144 4.5606, 4.4916, 5.7511, 5.7856, 6.797, 6.8373, 6.8171, 6.8171, 0, 0, 4.6926, 4.6926, 4.6926, 4.6926, 0, 4.6926, 4.6926, 4.3171, 4.3171, 4.6926, 4.6926, 4.6926, 4.6926, 0, 4.6926, 4.6926, 2.509, 2.9048, 2.5479, 2.5479
145 };
146 
147 // elastic interaction length in Silicon at 1 TeV per particle
148 const double nuclElLength[numHadrons] = {
149 9.248, 9.451, 11.545, 11.671, 32.081, 34.373, 34.373, 32.081, 0, 0, 15.739, 20.348, 15.739, 15.739, 0, 20.349, 0, 9.7514, 12.864, 15.836, 20.516, 15.836, 15.744, 0, 20.517, 20.44, 4.129, 6.0904, 4.5204, 4.5204
150 };
151 
152 
154  unsigned int distAlgo, double distCut, double elimit, double eth) :
155  curr4Mom(0.,0.,0.,0.),
156  vectProj(0.,0.,1.),
157  theBoost(0.,0.,0.),
158  theBertiniLimit(elimit),
159  theEnergyLimit(eth),
160  theDistCut(distCut),
161  distMin(1E99),
162  theDistAlgo(distAlgo)
163 {
164  // FTF model
165  theHadronicModel = new G4TheoFSGenerator("FTF");
166  theStringModel = new G4FTFModel();
167  G4GeneratorPrecompoundInterface* cascade
168  = new G4GeneratorPrecompoundInterface(new CMSDummyDeexcitation());
169  theLund = new G4LundStringFragmentation();
170  theStringDecay = new G4ExcitedStringDecay(theLund);
171  theStringModel->SetFragmentationModel(theStringDecay);
172 
173  theHadronicModel->SetTransport(cascade);
174  theHadronicModel->SetHighEnergyGenerator(theStringModel);
175  theHadronicModel->SetMinEnergy(theEnergyLimit);
176 
177  // Bertini Cascade
178  theBertiniCascade = new G4CascadeInterface();
179 
180  theDiffuseElastic = new G4DiffuseElastic();
181 
182  // Geant4 particles and cross sections
183  std::call_once(initializeOnce, [this] () {
184  theG4Hadron[0] = G4Proton::Proton();
185  theG4Hadron[1] = G4Neutron::Neutron();
186  theG4Hadron[2] = G4PionPlus::PionPlus();
187  theG4Hadron[3] = G4PionMinus::PionMinus();
188  theG4Hadron[4] = G4KaonPlus::KaonPlus();
189  theG4Hadron[5] = G4KaonMinus::KaonMinus();
190  theG4Hadron[6] = G4KaonZeroLong::KaonZeroLong();
191  theG4Hadron[7] = G4KaonZeroShort::KaonZeroShort();
192  theG4Hadron[8] = G4KaonZero::KaonZero();
193  theG4Hadron[9] = G4AntiKaonZero::AntiKaonZero();
195  theG4Hadron[11]= G4OmegaMinus::OmegaMinus();
196  theG4Hadron[12]= G4SigmaMinus::SigmaMinus();
197  theG4Hadron[13]= G4SigmaPlus::SigmaPlus();
199  theG4Hadron[15]= G4XiMinus::XiMinus();
200  theG4Hadron[16]= G4XiZero::XiZero();
201  theG4Hadron[17]= G4AntiProton::AntiProton();
202  theG4Hadron[18]= G4AntiNeutron::AntiNeutron();
203  theG4Hadron[19]= G4AntiLambda::AntiLambda();
204  theG4Hadron[20]= G4AntiOmegaMinus::AntiOmegaMinus();
205  theG4Hadron[21]= G4AntiSigmaMinus::AntiSigmaMinus();
206  theG4Hadron[22]= G4AntiSigmaPlus::AntiSigmaPlus();
207  theG4Hadron[23]= G4AntiSigmaZero::AntiSigmaZero();
208  theG4Hadron[24]= G4AntiXiMinus::AntiXiMinus();
209  theG4Hadron[25]= G4AntiXiZero::AntiXiZero();
210  theG4Hadron[26]= G4AntiAlpha::AntiAlpha();
211  theG4Hadron[27]= G4AntiDeuteron::AntiDeuteron();
212  theG4Hadron[28]= G4AntiTriton::AntiTriton();
213  theG4Hadron[29]= G4AntiHe3::AntiHe3();
214 
215  // other Geant4 particles
216  G4ParticleDefinition* ion = G4GenericIon::GenericIon();
217  ion->SetProcessManager(new G4ProcessManager(ion));
218  G4DecayPhysics decays;
219  decays.ConstructParticle();
220  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
221  partTable->SetVerboseLevel(0);
222  partTable->SetReadiness();
223 
224  for(int i=0; i<numHadrons; ++i) {
225  theId[i] = theG4Hadron[i]->GetPDGEncoding();
226  }
227  });
228 
229  // local objects
230  vect = new G4PhysicsLogVector(npoints-1,100*MeV,TeV);
232  currIdx = 0;
233  index = 0;
234  currTrack = nullptr;
236 
237  // fill projectile particle definitions
238  dummyStep = new G4Step();
239  dummyStep->SetPreStepPoint(new G4StepPoint());
240 
241  // target is always Silicon
242  targetNucleus.SetParameters(28, 14);
243 }
244 
246 
247  delete theStringDecay;
248  delete theStringModel;
249  delete theLund;
250  delete vect;
251 }
252 
255 {
256  //std::cout << "#### Primary " << Particle.pid() << " E(GeV)= "
257  // << Particle.momentum().e() << std::endl;
258 
259  int thePid = Particle.pid();
260  if(thePid != theId[currIdx]) {
261  currParticle = nullptr;
262  currIdx = 0;
263  for(; currIdx<numHadrons; ++currIdx) {
264  if(theId[currIdx] == thePid) {
266  // neutral kaons
267  if(7 == currIdx || 8 == currIdx) {
269  if(random->flatShoot() > 0.5) { currParticle = theG4Hadron[10]; }
270  }
271  break;
272  }
273  }
274  }
275  if(!currParticle) { return; }
276 
277  // fill projectile for Geant4
278  double mass = currParticle->GetPDGMass();
279  double ekin = CLHEP::GeV*Particle.momentum().e() - mass;
280 
281  // check interaction length
284  if(0.0 == intLengthInelastic) { return; }
285 
286  // apply corrections
287  if(ekin <= vect->Energy(0)) {
290  } else if(ekin < vect->Energy(npoints-1)) {
291  index = vect->FindBin(ekin, index);
292  double e1 = vect->Energy(index);
293  double e2 = vect->Energy(index+1);
294  intLengthElastic *= ((corrfactors_el[currIdx][index]*(e2 - ekin) +
295  corrfactors_el[currIdx][index+1]*(ekin - e1))/(e2 - e1));
296  intLengthInelastic *= ((corrfactors_inel[currIdx][index]*(e2 - ekin) +
297  corrfactors_inel[currIdx][index+1]*(ekin - e1))/(e2 - e1));
298  }
299  /*
300  std::cout << " Primary " << currParticle->GetParticleName()
301  << " E(GeV)= " << e*fact << std::endl;
302  */
303 
304  double currInteractionLength = -G4Log(random->flatShoot())*intLengthElastic*intLengthInelastic
306  /*
307  std::cout << "*NuclearInteractionFTFSimulator::compute: R(X0)= " << radLengths
308  << " Rnuc(X0)= " << theNuclIntLength[currIdx] << " IntLength(X0)= "
309  << currInteractionLength << std::endl;
310  */
311  // Check position of nuclear interaction
312  if (currInteractionLength > radLengths) { return; }
313 
314  // fill projectile for Geant4
315  double px = Particle.momentum().px();
316  double py = Particle.momentum().py();
317  double pz = Particle.momentum().pz();
318  double ptot = sqrt(px*px + py*py + pz*pz);
319  double norm = 1./ptot;
320  G4ThreeVector dir(px*norm, py*norm, pz*norm);
321  /*
322  std::cout << " Primary " << currParticle->GetParticleName()
323  << " E(GeV)= " << e*fact << " P(GeV/c)= ("
324  << px << " " << py << " " << pz << ")" << std::endl;
325  */
326 
327  G4DynamicParticle* dynParticle = new G4DynamicParticle(theG4Hadron[currIdx],dir,ekin);
328  currTrack = new G4Track(dynParticle, 0.0, vectProj);
329  currTrack->SetStep(dummyStep);
330 
331  theProjectile.Initialise(*currTrack);
332  delete currTrack;
333 
334  G4HadFinalState* result;
335 
336  // elastic interaction
338 
339  result = theDiffuseElastic->ApplyYourself(theProjectile, targetNucleus);
340  G4ThreeVector dirnew = result->GetMomentumChange().unit();
341  double cost = (dir*dirnew);
342  double sint = std::sqrt((1. - cost)*(1. + cost));
343 
344  curr4Mom.set(ptot*dirnew.x(),ptot*dirnew.y(),ptot*dirnew.z(),Particle.momentum().e());
345 
346  // Always create a daughter if the kink is large engough
347  if (sint > theDistCut) {
348  saveDaughter(Particle, curr4Mom, thePid);
349  } else {
350  Particle.SetXYZT(curr4Mom.px(), curr4Mom.py(), curr4Mom.pz(), curr4Mom.e());
351  }
352 
353  // inelastic interaction
354  } else {
355 
356  // Bertini cascade for low-energy hadrons (except light anti-nuclei)
357  // FTFP is applied above energy limit and for all anti-hyperons and anti-ions
358  if(ekin <= theBertiniLimit && currIdx < 17) {
359  result = theBertiniCascade->ApplyYourself(theProjectile, targetNucleus);
360  } else {
361  result = theHadronicModel->ApplyYourself(theProjectile, targetNucleus);
362  }
363  if(result) {
364 
365  int nsec = result->GetNumberOfSecondaries();
366  if(0 < nsec) {
367 
368  result->SetTrafoToLab(theProjectile.GetTrafoToLab());
369  _theUpdatedState.clear();
370 
371  //std::cout << " " << nsec << " secondaries" << std::endl;
372  // Generate angle
373  double phi = random->flatShoot()*CLHEP::twopi;
375  distMin = 1e99;
376 
377  // rotate and store secondaries
378  for (int j=0; j<nsec; ++j) {
379 
380  const G4DynamicParticle* dp = result->GetSecondary(j)->GetParticle();
381  int thePid = dp->GetParticleDefinition()->GetPDGEncoding();
382 
383  // rotate around primary direction
384  curr4Mom = dp->Get4Momentum();
385  curr4Mom.rotate(phi, vectProj);
386  curr4Mom *= result->GetTrafoToLab();
387  /*
388  std::cout << j << ". " << dp->GetParticleDefinition()->GetParticleName()
389  << " " << thePid
390  << " " << curr4Mom*fact << std::endl;
391  */
392  // prompt 2-gamma decay for pi0, eta, eta'p
393  if(111 == thePid || 221 == thePid || 331 == thePid) {
394  theBoost = curr4Mom.boostVector();
395  double e = 0.5*dp->GetParticleDefinition()->GetPDGMass();
396  double fi = random->flatShoot()*CLHEP::twopi;
397  double cth = 2*random->flatShoot() - 1.0;
398  double sth = sqrt((1.0 - cth)*(1.0 + cth));
399  G4LorentzVector lv(e*sth*cos(fi),e*sth*sin(fi),e*cth,e);
400  lv.boost(theBoost);
401  saveDaughter(Particle, lv, 22);
402  curr4Mom -= lv;
403  if(curr4Mom.e() > theEnergyLimit) {
404  saveDaughter(Particle, curr4Mom, 22);
405  }
406  } else {
407  if(curr4Mom.e() > theEnergyLimit + dp->GetParticleDefinition()->GetPDGMass()) {
408  saveDaughter(Particle, curr4Mom, thePid);
409  }
410  }
411  }
412  result->Clear();
413  }
414  }
415  }
416 }
417 
419  const G4LorentzVector& lv, int pdgid)
420 {
421  unsigned int idx = _theUpdatedState.size();
422  _theUpdatedState.push_back(Particle);
423  _theUpdatedState[idx].SetXYZT(lv.px()*fact,lv.py()*fact,lv.pz()*fact,lv.e()*fact);
424  _theUpdatedState[idx].setID(pdgid);
425 
426  // Store the closest daughter index (for later tracking purposes, so charged particles only)
427  double distance = distanceToPrimary(Particle,_theUpdatedState[idx]);
428  // Find the closest daughter, if closer than a given upper limit.
429  if ( distance < distMin && distance < theDistCut ) {
430  distMin = distance;
432  }
433  // std::cout << _theUpdatedState[idx] << std::endl;
434 }
435 
436 double
438  const RawParticle& aDaughter) const
439 {
440  double distance = 2E99;
441  // Compute the distance only for charged primaries
442  if ( fabs(Particle.charge()) > 1E-6 ) {
443 
444  // The secondary must have the same charge
445  double chargeDiff = fabs(aDaughter.charge()-Particle.charge());
446  if ( fabs(chargeDiff) < 1E-6 ) {
447 
448  // Here are two distance definitions * to be tuned *
449  switch ( theDistAlgo ) {
450 
451  case 1:
452  // sin(theta12)
453  distance = (aDaughter.Vect().Unit().Cross(Particle.Vect().Unit())).R();
454  break;
455 
456  case 2:
457  // sin(theta12) * p1/p2
458  distance = (aDaughter.Vect().Cross(Particle.Vect())).R()
459  /aDaughter.Vect().Mag2();
460  break;
461 
462  default:
463  // Should not happen
464  break;
465  }
466  }
467  }
468  return distance;
469 }
const double nuclInelLength[numHadrons]
const double GeV
Definition: MathUtil.h:16
double flatShoot(double xmin=0.0, double xmax=1.0) const
void compute(ParticlePropagator &Particle, RandomEngineAndDistribution const *) override
Generate a nuclear interaction according to the probability that it happens.
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
int pid() const
get the HEP particle ID number
Definition: RawParticle.h:265
TRandom random
Definition: MVATrainer.cc:138
const double corrfactors_inel[numHadrons][npoints]
~NuclearInteractionFTFSimulator() override
Default Destructor.
const double MeV
const XYZTLorentzVector & momentum() const
the momentum fourvector
Definition: RawParticle.h:286
void saveDaughter(ParticlePropagator &Particle, const G4LorentzVector &lv, int pdgid)
T sqrt(T t)
Definition: SSEVec.h:18
static std::once_flag initializeOnce
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
const double nuclElLength[numHadrons]
static const int npoints
double charge() const
get the MEASURED charge
Definition: RawParticle.h:282
double distanceToPrimary(const RawParticle &Particle, const RawParticle &aDaughter) const
static const int numHadrons
const double corrfactors_el[numHadrons][npoints]
NuclearInteractionFTFSimulator(unsigned int distAlgo, double distCut, double elimit, double eth)
Constructor.
const G4ParticleDefinition * currParticle
Float e1
Definition: deltaR.h:20
const double fact
auto dp
Definition: deltaR.h:22
Float e2
Definition: deltaR.h:21
static const G4ParticleDefinition * theG4Hadron[numHadrons]
std::vector< RawParticle > _theUpdatedState
dbl *** dir
Definition: mlp_gen.cc:35