CMS 3D CMS Logo

NuclearInteractionFTF.cc
Go to the documentation of this file.
1 // system headers
2 #include <mutex>
3 #include <cmath>
4 #include <memory>
5 
6 // Framework Headers
9 
10 
11 
12 // Fast Sim headers
20 
21 // Math
23 
24 // Geant4 headers
25 #include "G4Nucleus.hh"
26 #include "G4HadProjectile.hh"
27 #include "G4LorentzVector.hh"
28 #include "G4ThreeVector.hh"
29 
30 #include "G4ParticleDefinition.hh"
31 #include "G4DynamicParticle.hh"
32 #include "G4Track.hh"
33 #include "G4Step.hh"
34 #include "G4StepPoint.hh"
35 #include "G4TheoFSGenerator.hh"
36 #include "G4FTFModel.hh"
37 #include "G4ExcitedStringDecay.hh"
38 #include "G4LundStringFragmentation.hh"
39 #include "G4GeneratorPrecompoundInterface.hh"
40 #include "G4CascadeInterface.hh"
41 #include "G4DiffuseElastic.hh"
42 
43 #include "G4Proton.hh"
44 #include "G4Neutron.hh"
45 #include "G4PionPlus.hh"
46 #include "G4PionMinus.hh"
47 #include "G4AntiProton.hh"
48 #include "G4AntiNeutron.hh"
49 #include "G4KaonPlus.hh"
50 #include "G4KaonMinus.hh"
51 #include "G4KaonZeroLong.hh"
52 #include "G4KaonZeroShort.hh"
53 #include "G4KaonZero.hh"
54 #include "G4AntiKaonZero.hh"
55 #include "G4GenericIon.hh"
56 
57 #include "G4Lambda.hh"
58 #include "G4OmegaMinus.hh"
59 #include "G4SigmaMinus.hh"
60 #include "G4SigmaPlus.hh"
61 #include "G4SigmaZero.hh"
62 #include "G4XiMinus.hh"
63 #include "G4XiZero.hh"
64 #include "G4AntiLambda.hh"
65 #include "G4AntiOmegaMinus.hh"
66 #include "G4AntiSigmaMinus.hh"
67 #include "G4AntiSigmaPlus.hh"
68 #include "G4AntiSigmaZero.hh"
69 #include "G4AntiXiMinus.hh"
70 #include "G4AntiXiZero.hh"
71 #include "G4AntiAlpha.hh"
72 #include "G4AntiDeuteron.hh"
73 #include "G4AntiTriton.hh"
74 #include "G4AntiHe3.hh"
75 
76 #include "G4Material.hh"
77 #include "G4DecayPhysics.hh"
78 #include "G4ParticleTable.hh"
79 #include "G4IonTable.hh"
80 #include "G4ProcessManager.hh"
81 #include "G4PhysicsLogVector.hh"
82 #include "G4SystemOfUnits.hh"
83 
84 
86 // Author: Vladimir Ivanchenko
87 // Date: 20-Jan-2015
88 //
89 // Revision: Class structure modified to match SimplifiedGeometryPropagator
90 // Fixed a bug in which units where not properly converted from G4 to FastSim standard
91 // S. Kurz, 29 May 2017
93 
94 
96 
97 namespace fastsim
98 {
100 
105  {
106  public:
109 
111  ~NuclearInteractionFTF() override;
112 
114 
120  void interact(fastsim::Particle & particle, const SimplifiedGeometry & layer,std::vector<std::unique_ptr<fastsim::Particle> > & secondaries,const RandomEngineAndDistribution & random) override;
121 
122  private:
123  static const int numHadrons = 30;
124  static const int npoints = 21;
125 
126  static const G4ParticleDefinition* theG4Hadron[numHadrons];
127  static int theId[numHadrons];
128 
129  G4PhysicsLogVector* vect;
130  G4TheoFSGenerator* theHadronicModel;
131  G4FTFModel* theStringModel;
132  G4ExcitedStringDecay* theStringDecay;
133  G4LundStringFragmentation* theLund;
134  G4GeneratorPrecompoundInterface* theCascade;
135 
136  G4CascadeInterface* theBertiniCascade;
137  G4DiffuseElastic* theDiffuseElastic;
138 
139  G4Step* dummyStep;
140  G4Track* currTrack;
141  const G4ParticleDefinition* currParticle;
142 
143  G4Nucleus targetNucleus;
144  G4HadProjectile theProjectile;
145  G4LorentzVector curr4Mom;
146  G4ThreeVector vectProj;
147  G4ThreeVector theBoost;
148 
150  double theEnergyLimit;
151 
152  double theDistCut;
153 
156 
157  int currIdx;
158  size_t index;
159 
162  {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},
163  {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},
164  {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},
165  {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},
166  {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},
167  {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},
168  {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},
169  {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},
170  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
171  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
172  {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},
173  {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},
174  {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},
175  {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},
176  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
177  {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},
178  {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},
179  {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},
180  {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},
181  {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},
182  {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},
183  {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},
184  {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},
185  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
186  {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},
187  {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},
188  {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},
189  {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},
190  {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},
191  {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},
192  };
193 
195  const double corrfactors_el[numHadrons][npoints] = {
196  {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},
197  {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},
198  {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},
199  {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},
200  {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},
201  {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},
202  {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},
203  {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},
204  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
205  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
206  {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},
207  {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},
208  {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},
209  {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},
210  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
211  {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},
212  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
213  {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},
214  {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},
215  {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},
216  {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},
217  {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},
218  {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},
219  {0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 0, 1.0},
220  {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},
221  {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},
222  {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},
223  {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},
224  {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},
225  {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},
226  };
227 
229  const double nuclInelLength[numHadrons] = {
230  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
231  };
232 
234  const double nuclElLength[numHadrons] = {
235  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
236  };
237  };
238 
239  // Is this correct?
240  // Thread safety
241  static std::once_flag initializeOnce;
242  [[cms::thread_guard("initializeOnce")]] const G4ParticleDefinition* NuclearInteractionFTF::theG4Hadron[] = {nullptr};
243  [[cms::thread_guard("initializeOnce")]] int NuclearInteractionFTF::theId[] = {0};
244 }
245 
247  : fastsim::InteractionModel(name)
248  , curr4Mom(0.,0.,0.,0.)
249  , vectProj(0.,0.,1.)
250  , theBoost(0.,0.,0.)
251 {
252  // Angular distance of particle before/after interaction. If it too large, a daughter is created instead
253  theDistCut = cfg.getParameter<double>("distCut");
254  // Set energy limits for processes
255  theBertiniLimit = cfg.getParameter<double>("bertiniLimit") * CLHEP::GeV;
256  theEnergyLimit = cfg.getParameter<double>("energyLimit") * CLHEP::GeV;
257 
258  // FTF model
259  theHadronicModel = new G4TheoFSGenerator("FTF");
260  theStringModel = new G4FTFModel();
261  G4GeneratorPrecompoundInterface* cascade
262  = new G4GeneratorPrecompoundInterface(new CMSDummyDeexcitation());
263  theLund = new G4LundStringFragmentation();
264  theStringDecay = new G4ExcitedStringDecay(theLund);
265  theStringModel->SetFragmentationModel(theStringDecay);
266 
267  theHadronicModel->SetTransport(cascade);
268  theHadronicModel->SetHighEnergyGenerator(theStringModel);
269  theHadronicModel->SetMinEnergy(theEnergyLimit);
270 
271  // Bertini Cascade
272  theBertiniCascade = new G4CascadeInterface();
273 
274  theDiffuseElastic = new G4DiffuseElastic();
275 
276  // Geant4 particles and cross sections
277  std::call_once(initializeOnce, [this] () {
278  theG4Hadron[0] = G4Proton::Proton();
279  theG4Hadron[1] = G4Neutron::Neutron();
280  theG4Hadron[2] = G4PionPlus::PionPlus();
281  theG4Hadron[3] = G4PionMinus::PionMinus();
282  theG4Hadron[4] = G4KaonPlus::KaonPlus();
283  theG4Hadron[5] = G4KaonMinus::KaonMinus();
284  theG4Hadron[6] = G4KaonZeroLong::KaonZeroLong();
285  theG4Hadron[7] = G4KaonZeroShort::KaonZeroShort();
286  theG4Hadron[8] = G4KaonZero::KaonZero();
287  theG4Hadron[9] = G4AntiKaonZero::AntiKaonZero();
289  theG4Hadron[11]= G4OmegaMinus::OmegaMinus();
290  theG4Hadron[12]= G4SigmaMinus::SigmaMinus();
291  theG4Hadron[13]= G4SigmaPlus::SigmaPlus();
293  theG4Hadron[15]= G4XiMinus::XiMinus();
294  theG4Hadron[16]= G4XiZero::XiZero();
295  theG4Hadron[17]= G4AntiProton::AntiProton();
296  theG4Hadron[18]= G4AntiNeutron::AntiNeutron();
297  theG4Hadron[19]= G4AntiLambda::AntiLambda();
298  theG4Hadron[20]= G4AntiOmegaMinus::AntiOmegaMinus();
299  theG4Hadron[21]= G4AntiSigmaMinus::AntiSigmaMinus();
300  theG4Hadron[22]= G4AntiSigmaPlus::AntiSigmaPlus();
301  theG4Hadron[23]= G4AntiSigmaZero::AntiSigmaZero();
302  theG4Hadron[24]= G4AntiXiMinus::AntiXiMinus();
303  theG4Hadron[25]= G4AntiXiZero::AntiXiZero();
304  theG4Hadron[26]= G4AntiAlpha::AntiAlpha();
305  theG4Hadron[27]= G4AntiDeuteron::AntiDeuteron();
306  theG4Hadron[28]= G4AntiTriton::AntiTriton();
307  theG4Hadron[29]= G4AntiHe3::AntiHe3();
308 
309  // other Geant4 particles
310  G4ParticleDefinition* ion = G4GenericIon::GenericIon();
311  ion->SetProcessManager(new G4ProcessManager(ion));
312  G4DecayPhysics decays;
313  decays.ConstructParticle();
314  G4ParticleTable* partTable = G4ParticleTable::GetParticleTable();
315  partTable->SetVerboseLevel(0);
316  partTable->SetReadiness();
317 
318  for(int i=0; i<numHadrons; ++i) {
319  theId[i] = theG4Hadron[i]->GetPDGEncoding();
320  }
321  });
322 
323  // local objects
324  vect = new G4PhysicsLogVector(npoints-1,100*MeV,TeV);
326  currIdx = 0;
327  index = 0;
328  currTrack = nullptr;
330 
331  // fill projectile particle definitions
332  dummyStep = new G4Step();
333  dummyStep->SetPreStepPoint(new G4StepPoint());
334 
335  // target is always Silicon
336  targetNucleus.SetParameters(28, 14);
337 }
338 
340  delete theStringDecay;
341  delete theStringModel;
342  delete theLund;
343  delete vect;
344 }
345 
346 
347 void fastsim::NuclearInteractionFTF::interact(fastsim::Particle & particle, const SimplifiedGeometry & layer,std::vector<std::unique_ptr<fastsim::Particle> > & secondaries,const RandomEngineAndDistribution & random)
348 {
349  int thePid = particle.pdgId();
350  //
351  // no valid PDGid
352  //
353  if(abs(thePid) <= 100 || abs(thePid) >= 1000000)
354  {
355  return;
356  }
357 
358  // particle lost all its enrgy in previous interaction
359  if (particle.momentum().E() < 1E-10)
360  {
361  return;
362  }
363 
364  double radLengths = layer.getThickness(particle.position(),particle.momentum());
365  // TEC layers have some fudge factor for nuclear interactions
366  radLengths *= layer.getNuclearInteractionThicknessFactor();
367  //
368  // no material
369  //
370  if(radLengths < 1E-10)
371  {
372  return;
373  }
374 
375  // get the G4 hadron
376  if(thePid != theId[currIdx]) {
377  currParticle = nullptr;
378  currIdx = 0;
379  for(; currIdx<numHadrons; ++currIdx) {
380  if(theId[currIdx] == thePid) {
382  // neutral kaons
383  if(7 == currIdx || 8 == currIdx) {
385  if(random.flatShoot() > 0.5) { currParticle = theG4Hadron[10]; }
386  }
387  break;
388  }
389  }
390  }
391 
392  // no valid G4 hadron found
393  if(!currParticle){
394  return;
395  }
396 
397  // fill projectile for Geant4
398  double mass = currParticle->GetPDGMass();
399  double ekin = CLHEP::GeV * particle.momentum().e() - mass;
400 
401  // check interaction length
404  if(0.0 == intLengthInelastic){
405  return;
406  }
407 
408  // apply corrections
409  if(ekin <= vect->Energy(0)){
412  }else if(ekin < vect->Energy(npoints-1)){
413  index = vect->FindBin(ekin, index);
414  double e1 = vect->Energy(index);
415  double e2 = vect->Energy(index+1);
416  intLengthElastic *= ((corrfactors_el[currIdx][index] * (e2 - ekin) +
417  corrfactors_el[currIdx][index+1] * (ekin - e1)) / (e2 - e1));
418  intLengthInelastic *= ((corrfactors_inel[currIdx][index] * (e2 - ekin) +
419  corrfactors_inel[currIdx][index+1] * (ekin - e1)) / (e2 - e1));
420  }
421 
422  double currInteractionLength = -G4Log(random.flatShoot()) *
424 
425  // Check position of nuclear interaction
426  if(currInteractionLength > radLengths)
427  {
428  return;
429  }
430 
431  // fill projectile for Geant4
432  double px = particle.momentum().px();
433  double py = particle.momentum().py();
434  double pz = particle.momentum().pz();
435  double ptot = sqrt(px*px + py*py + pz*pz);
436  double norm = 1./ptot;
437  G4ThreeVector dir(px*norm, py*norm, pz*norm);
438 
439  G4DynamicParticle* dynParticle = new G4DynamicParticle(theG4Hadron[currIdx],dir,ekin);
440  currTrack = new G4Track(dynParticle, 0.0, vectProj);
441  currTrack->SetStep(dummyStep);
442 
443  theProjectile.Initialise(*currTrack);
444  delete currTrack;
445 
446  G4HadFinalState* result;
447 
448  // elastic interaction
450  result = theDiffuseElastic->ApplyYourself(theProjectile, targetNucleus);
451  G4ThreeVector dirnew = result->GetMomentumChange().unit();
452  double cost = (dir*dirnew);
453  double sint = std::sqrt((1. - cost) * (1. + cost));
454 
455  // This vector is already in units of GeV (FastSim standard even if it is a G4Vector)
456  curr4Mom.set(ptot*dirnew.x(),ptot*dirnew.y(),ptot*dirnew.z(),particle.momentum().e());
457 
458  // Always create a daughter if the kink is large enough
459  if(sint > theDistCut){
460  // Create secondary
461  secondaries.emplace_back(new fastsim::Particle(thePid
462  ,particle.position()
463  ,XYZTLorentzVector(curr4Mom.px(), curr4Mom.py(), curr4Mom.pz(), curr4Mom.e())));
464 
465  // ClosestChargedDaughter thing for FastSim tracking
466  if(particle.charge() != 0){
467  secondaries.back()->setMotherDeltaR(particle.momentum());
468  secondaries.back()->setMotherPdgId(thePid);
469  secondaries.back()->setMotherSimTrackIndex(particle.simTrackIndex());
470  }
471 
472  // The particle is destroyed
473  particle.momentum().SetXYZT(0., 0., 0., 0.);
474  }else{
475  // ... or just rotate particle
476  particle.momentum().SetXYZT(curr4Mom.px(), curr4Mom.py(), curr4Mom.pz(), curr4Mom.e());
477  }
478 
479  // inelastic interaction
480  }else{
481  // Bertini cascade for low-energy hadrons (except light anti-nuclei)
482  // FTFP is applied above energy limit and for all anti-hyperons and anti-ions
483  if(ekin <= theBertiniLimit && currIdx < 17){
484  result = theBertiniCascade->ApplyYourself(theProjectile, targetNucleus);
485  }else{
486  result = theHadronicModel->ApplyYourself(theProjectile, targetNucleus);
487  }
488 
489  if(result){
490  int nsec = result->GetNumberOfSecondaries();
491  if(nsec > 0){
492  result->SetTrafoToLab(theProjectile.GetTrafoToLab());
493 
494  // Generate angle
495  double phi = random.flatShoot()*CLHEP::twopi;
496 
497  // rotate and store secondaries
498  for (int j=0; j<nsec; ++j) {
499  const G4DynamicParticle* dp = result->GetSecondary(j)->GetParticle();
500  int daughterID = dp->GetParticleDefinition()->GetPDGEncoding();
501 
502  // rotate around primary direction
503  curr4Mom = dp->Get4Momentum();
504  curr4Mom.rotate(phi, vectProj);
505  curr4Mom *= result->GetTrafoToLab();
506 
507  // prompt 2-gamma decay for pi0, eta, eta'p
508  // don't have a charge so the closest daughter thing does not have to be done
509  if(111 == daughterID || 221 == daughterID || 331 == daughterID){
510  theBoost = curr4Mom.boostVector();
511  double e = 0.5*dp->GetParticleDefinition()->GetPDGMass();
512  double fi = random.flatShoot()*CLHEP::twopi;
513  double cth = 2.*random.flatShoot() - 1.0;
514  double sth = sqrt((1.0 - cth)*(1.0 + cth));
515  G4LorentzVector lv(e*sth*cos(fi),e*sth*sin(fi),e*cth,e);
516  lv.boost(theBoost);
517 
518  // create secondaries
519  secondaries.emplace_back(new fastsim::Particle(22
520  ,particle.position()
521  ,XYZTLorentzVector(lv.px()/CLHEP::GeV, lv.py()/CLHEP::GeV, lv.pz()/CLHEP::GeV, lv.e()/CLHEP::GeV)));
522 
523  curr4Mom -= lv;
524  if(curr4Mom.e() > theEnergyLimit) {
525  secondaries.emplace_back(new fastsim::Particle(22
526  ,particle.position()
528  }
529 
530  // The mother particle is destroyed
531  particle.momentum().SetXYZT(0., 0., 0., 0.);
532  }else{
533  // Copied from the original code, however I am not sure if all that is correct!
534  // Energy of particles increases during the interaction!
535  // std::cout << particle.momentum().e() << "(" << dynParticle->GetParticleDefinition()->GetPDGMass()/CLHEP::GeV << ") -> " << curr4Mom.e()/CLHEP::GeV << "(" << dp->GetParticleDefinition()->GetPDGMass()/CLHEP::GeV << ")" << std::endl;
536 
537  if(curr4Mom.e() > theEnergyLimit + dp->GetParticleDefinition()->GetPDGMass()){
538  // Create secondary
539  secondaries.emplace_back(new fastsim::Particle(daughterID
540  ,particle.position()
542 
543  // ClosestChargedDaughter thing for FastSim tracking
544  if(particle.charge() != 0 && std::abs(particle.charge()-dp->GetParticleDefinition()->GetPDGCharge()) < 1E-10){
545  secondaries.back()->setMotherDeltaR(particle.momentum());
546  secondaries.back()->setMotherPdgId(particle.getMotherDeltaR() == -1 ? particle.pdgId() : particle.getMotherPdgId());
547  secondaries.back()->setMotherSimTrackIndex(particle.simTrackIndex());
548  }
549 
550  // The mother particle is destroyed
551  particle.momentum().SetXYZT(0., 0., 0., 0.);
552  }
553  }
554  }
555  }
556 
557  // Clear the final state
558  result->Clear();
559  }
560  }
561 }
562 
566  "fastsim::NuclearInteractionFTF"
567  );
const double corrfactors_el[numHadrons][npoints]
elastic interaction length corrections per particle and energy
T getParameter(std::string const &) const
G4LundStringFragmentation * theLund
const double GeV
Definition: MathUtil.h:16
const G4ParticleDefinition * currParticle
Implementation of a generic detector layer (base class for forward/barrel layers).
const math::XYZTLorentzVector & position() const
Return position of the particle.
Definition: Particle.h:142
G4GeneratorPrecompoundInterface * theCascade
double flatShoot(double xmin=0.0, double xmax=1.0) const
double theDistCut
Cut on deltaR for ClosestChargedDaughter algorithm (FastSim tracking)
int getMotherPdgId() const
Get pdgIdto mother particle (only makes sense if mother and daughter charged).
Definition: Particle.h:213
size_t index
Index for energy of particle (vectors parametrized as a function of energy)
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const double getNuclearInteractionThicknessFactor() const
Some layers have a different thickness for nuclear interactions.
Needed as a dummy interface to Geant4 nuclear de-excitation module.
double intLengthElastic
Elastic interaction length.
int pdgId() const
Return pdgId of the particle.
Definition: Particle.h:136
TRandom random
Definition: MVATrainer.cc:138
const double MeV
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
Base class for any interaction model between a particle and a tracker layer.
double theEnergyLimit
Minimum energy of interacting particle.
NuclearInteractionFTF(const std::string &name, const edm::ParameterSet &cfg)
Constructor.
T sqrt(T t)
Definition: SSEVec.h:18
double theBertiniLimit
Bertini cascade for low-energy hadrons.
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
static const G4ParticleDefinition * theG4Hadron[numHadrons]
The Geant4 particles.
static const int npoints
Number of steps to cover range of particle energies.
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
int currIdx
Index of particle in vector of Geant4 particles.
double intLengthInelastic
Inelastic interaction length.
int simTrackIndex() const
Return index of the SimTrack.
Definition: Particle.h:155
Float e1
Definition: deltaR.h:20
const double corrfactors_inel[numHadrons][npoints]
inelastic interaction length corrections per particle and energy
const double nuclInelLength[numHadrons]
inelastic interaction length in Silicon at 1 TeV per particle
double getMotherDeltaR() const
Get delta R to mother particle (only makes sense if mother and daughter charged). ...
Definition: Particle.h:206
auto dp
Definition: deltaR.h:22
void interact(fastsim::Particle &particle, const SimplifiedGeometry &layer, std::vector< std::unique_ptr< fastsim::Particle > > &secondaries, const RandomEngineAndDistribution &random) override
Perform the interaction.
double charge() const
Return charge of the particle.
Definition: Particle.h:139
virtual const double getThickness(const math::XYZTLorentzVector &position) const =0
Return thickness of the layer at a given position.
math::XYZTLorentzVector XYZTLorentzVector
Float e2
Definition: deltaR.h:21
const double nuclElLength[numHadrons]
elastic interaction length in Silicon at 1 TeV per particle
static std::once_flag initializeOnce
static int theId[numHadrons]
The pdgIDs of the Geant4 particles.
const math::XYZTLorentzVector & momentum() const
Return momentum of the particle.
Definition: Particle.h:145
Implementation of nuclear interactions of hadrons in the tracker layers (based on FTF model of Geant4...
#define DEFINE_EDM_PLUGIN(factory, type, name)
dbl *** dir
Definition: mlp_gen.cc:35
Definition of a generic FastSim Particle which can be propagated through the detector (formerly Parti...
Definition: Particle.h:19
~NuclearInteractionFTF() override
Default destructor.
static const int numHadrons
Number of G4 hadrons.