1 | /*
|
---|
2 | ---- Delphes ----
|
---|
3 | A Fast Simulator for general purpose LHC detector
|
---|
4 | S. Ovyn ~~~~ severine.ovyn@uclouvain.be
|
---|
5 |
|
---|
6 | Center for Particle Physics and Phenomenology (CP3)
|
---|
7 | Universite Catholique de Louvain (UCL)
|
---|
8 | Louvain-la-Neuve, Belgium
|
---|
9 | */
|
---|
10 |
|
---|
11 | /// \file SmearUtil.cc
|
---|
12 | /// \brief RESOLution class, and some generic definitions
|
---|
13 |
|
---|
14 |
|
---|
15 | #include "interface/SmearUtil.h"
|
---|
16 | #include "TRandom.h"
|
---|
17 |
|
---|
18 | #include <iostream>
|
---|
19 | #include <sstream>
|
---|
20 | #include <fstream>
|
---|
21 | using namespace std;
|
---|
22 |
|
---|
23 | //------------------------------------------------------------------------------
|
---|
24 |
|
---|
25 | RESOLution::RESOLution() {
|
---|
26 |
|
---|
27 | MAX_TRACKER = 2.5; // tracker coverage
|
---|
28 | MAX_CALO_CEN = 3.0; // central calorimeter coverage
|
---|
29 | MAX_CALO_FWD = 5.0; // forward calorimeter pseudorapidity coverage
|
---|
30 | MAX_MU = 2.4; // muon chambers pseudorapidity coverage
|
---|
31 | MIN_CALO_VFWD= 5.2; // very forward calorimeter (if any), like CASTOR
|
---|
32 | MAX_CALO_VFWD= 6.6; // very forward calorimeter (if any), like CASTOR
|
---|
33 | MIN_ZDC = 8.3; // zero-degree calorimeter, coverage
|
---|
34 |
|
---|
35 | ZDC_S = 140.; // ZDC distance to IP
|
---|
36 | RP220_S = 220; // distance of the RP to the IP, in meters
|
---|
37 | RP220_X = 0.002;// distance of the RP to the beam, in meters
|
---|
38 | FP420_S = 420; // distance of the RP to the IP, in meters
|
---|
39 | FP420_X = 0.004;// distance of the RP to the beam, in meters
|
---|
40 |
|
---|
41 |
|
---|
42 | ELG_Scen = 0.028; // S term for central ECAL
|
---|
43 | ELG_Ncen = 0.124 ; // N term for central ECAL
|
---|
44 | ELG_Ccen = 0.0026 ; // C term for central ECAL
|
---|
45 | ELG_Cfwd = 0.107 ; // S term for forward ECAL
|
---|
46 | ELG_Sfwd = 2.084 ; // C term for forward ECAL
|
---|
47 | ELG_Nfwd = 0.0 ; // N term for central ECAL
|
---|
48 |
|
---|
49 | HAD_Shcal = 0.91 ; // S term for central HCAL // hadronic calorimeter
|
---|
50 | HAD_Nhcal = 0. ; // N term for central HCAL
|
---|
51 | HAD_Chcal = 0.038 ; // C term for central HCAL
|
---|
52 | HAD_Shf = 2.7 ; // S term for central HF // forward calorimeter
|
---|
53 | HAD_Nhf = 0. ; // N term for central HF
|
---|
54 | HAD_Chf = 0.13 ; // C term for central HF
|
---|
55 |
|
---|
56 | MU_SmearPt = 0.01 ;
|
---|
57 |
|
---|
58 | ELEC_pt = 10.0;
|
---|
59 | MUON_pt = 10.0;
|
---|
60 | JET_pt = 20.0;
|
---|
61 | TAUJET_pt = 10.0;
|
---|
62 |
|
---|
63 |
|
---|
64 | TAU_CONE_ENERGY = 0.15 ; // Delta R = radius of the cone // for "electromagnetic collimation"
|
---|
65 | TAU_EM_COLLIMATION = 0.95;
|
---|
66 | TAU_CONE_TRACKS= 0.4 ; //Delta R for tracker isolation for tau's
|
---|
67 | PT_TRACK_TAU = 2.0 ; // GeV // 6 GeV ????
|
---|
68 |
|
---|
69 |
|
---|
70 | PT_TRACKS_MIN = 0.9 ; // minimal pt needed to reach the calorimeter, in GeV
|
---|
71 | PT_QUARKS_MIN = 2.0 ; // minimal pt needed by quarks to reach the tracker, in GeV (??????)
|
---|
72 | TRACKING_EFF = 90;
|
---|
73 |
|
---|
74 |
|
---|
75 | TAGGING_B = 40;
|
---|
76 | MISTAGGING_C = 10;
|
---|
77 | MISTAGGING_L = 1;
|
---|
78 |
|
---|
79 |
|
---|
80 | CONERADIUS = 0.7; // generic jet radius ; not for tau's !!!
|
---|
81 | JETALGO = 1; // 1 for Cone algorithm, 2 for MidPoint algorithm, 3 for SIScone algorithm, 4 for kt algorithm
|
---|
82 |
|
---|
83 | //General jet parameters
|
---|
84 | SEEDTHRESHOLD = 1.0;
|
---|
85 | OVERLAPTHRESHOLD = 0.75;
|
---|
86 |
|
---|
87 | // Define Cone algorithm.
|
---|
88 | C_ADJACENCYCUT = 2;
|
---|
89 | C_MAXITERATIONS = 100;
|
---|
90 | C_IRATCH = 1;
|
---|
91 |
|
---|
92 | //Define MidPoint algorithm.
|
---|
93 | M_CONEAREAFRACTION = 0.25;
|
---|
94 | M_MAXPAIRSIZE = 2;
|
---|
95 | M_MAXITERATIONS = 100;
|
---|
96 |
|
---|
97 | }
|
---|
98 |
|
---|
99 | //------------------------------------------------------------------------------
|
---|
100 | void RESOLution::ReadDataCard(const string datacard) {
|
---|
101 |
|
---|
102 | string temp_string;
|
---|
103 | istringstream curstring;
|
---|
104 |
|
---|
105 | ifstream fichier_a_lire(datacard.c_str());
|
---|
106 | if(!fichier_a_lire.good()) {
|
---|
107 | cout << datacard << "Datadard " << datacard << " not found, use default values" << endl;
|
---|
108 | return;
|
---|
109 | }
|
---|
110 |
|
---|
111 | while (getline(fichier_a_lire,temp_string)) {
|
---|
112 | curstring.clear(); // needed when using several times istringstream::str(string)
|
---|
113 | curstring.str(temp_string);
|
---|
114 | string varname;
|
---|
115 | float value;
|
---|
116 |
|
---|
117 | if(strstr(temp_string.c_str(),"#")) { }
|
---|
118 | else if(strstr(temp_string.c_str(),"MAX_TRACKER")){curstring >> varname >> value; MAX_TRACKER = value;}
|
---|
119 | else if(strstr(temp_string.c_str(),"MAX_CALO_CEN")){curstring >> varname >> value; MAX_CALO_CEN = value;}
|
---|
120 | else if(strstr(temp_string.c_str(),"MAX_CALO_FWD")){curstring >> varname >> value; MAX_CALO_FWD = value;}
|
---|
121 | else if(strstr(temp_string.c_str(),"MAX_MU")){curstring >> varname >> value; MAX_MU = value;}
|
---|
122 | else if(strstr(temp_string.c_str(),"ELG_Scen")){curstring >> varname >> value; ELG_Scen = value;}
|
---|
123 | else if(strstr(temp_string.c_str(),"ELG_Ncen")){curstring >> varname >> value; ELG_Ncen = value;}
|
---|
124 | else if(strstr(temp_string.c_str(),"ELG_Ccen")){curstring >> varname >> value; ELG_Ccen = value;}
|
---|
125 | else if(strstr(temp_string.c_str(),"ELG_Sfwd")){curstring >> varname >> value; ELG_Sfwd = value;}
|
---|
126 | else if(strstr(temp_string.c_str(),"ELG_Cfwd")){curstring >> varname >> value; ELG_Cfwd = value;}
|
---|
127 | else if(strstr(temp_string.c_str(),"ELG_Nfwd")){curstring >> varname >> value; ELG_Nfwd = value;}
|
---|
128 | else if(strstr(temp_string.c_str(),"HAD_Shcal")){curstring >> varname >> value; HAD_Shcal = value;}
|
---|
129 | else if(strstr(temp_string.c_str(),"HAD_Nhcal")){curstring >> varname >> value; HAD_Nhcal = value;}
|
---|
130 | else if(strstr(temp_string.c_str(),"HAD_Chcal")){curstring >> varname >> value; HAD_Chcal = value;}
|
---|
131 | else if(strstr(temp_string.c_str(),"HAD_Shf")){curstring >> varname >> value; HAD_Shf = value;}
|
---|
132 | else if(strstr(temp_string.c_str(),"HAD_Nhf")){curstring >> varname >> value; HAD_Nhf = value;}
|
---|
133 | else if(strstr(temp_string.c_str(),"HAD_Chf")){curstring >> varname >> value; HAD_Chf = value;}
|
---|
134 | else if(strstr(temp_string.c_str(),"MU_SmearPt")){curstring >> varname >> value; MU_SmearPt = value;}
|
---|
135 | else if(strstr(temp_string.c_str(),"TAU_CONE_ENERGY")){curstring >> varname >> value; TAU_CONE_ENERGY = value;}
|
---|
136 | else if(strstr(temp_string.c_str(),"TAU_CONE_TRACKS")){curstring >> varname >> value; TAU_CONE_TRACKS = value;}
|
---|
137 | else if(strstr(temp_string.c_str(),"PT_TRACK_TAU")){curstring >> varname >> value; PT_TRACK_TAU = value;}
|
---|
138 | else if(strstr(temp_string.c_str(),"PT_TRACKS_MIN")){curstring >> varname >> value; PT_TRACKS_MIN = value;}
|
---|
139 | else if(strstr(temp_string.c_str(),"TAGGING_B")){curstring >> varname >> value; TAGGING_B = (int)value;}
|
---|
140 | else if(strstr(temp_string.c_str(),"MISTAGGING_C")){curstring >> varname >> value; MISTAGGING_C = (int)value;}
|
---|
141 | else if(strstr(temp_string.c_str(),"MISTAGGING_L")){curstring >> varname >> value; MISTAGGING_L = (int)value;}
|
---|
142 | else if(strstr(temp_string.c_str(),"CONERADIUS")){curstring >> varname >> value; CONERADIUS = value;}
|
---|
143 | else if(strstr(temp_string.c_str(),"JETALGO")){curstring >> varname >> value; JETALGO = (int)value;}
|
---|
144 | else if(strstr(temp_string.c_str(),"TRACKING_EFF")){curstring >> varname >> value; TRACKING_EFF = (int)value;}
|
---|
145 | else if(strstr(temp_string.c_str(),"ELEC_pt")){curstring >> varname >> value; ELEC_pt = value;}
|
---|
146 | else if(strstr(temp_string.c_str(),"MUON_pt")){curstring >> varname >> value; MUON_pt = value;}
|
---|
147 | else if(strstr(temp_string.c_str(),"JET_pt")){curstring >> varname >> value; JET_pt = value;}
|
---|
148 | else if(strstr(temp_string.c_str(),"TAUJET_pt")){curstring >> varname >> value; TAUJET_pt = value;}
|
---|
149 |
|
---|
150 | }
|
---|
151 |
|
---|
152 | // General jet variables
|
---|
153 | SEEDTHRESHOLD = 1.0;
|
---|
154 | OVERLAPTHRESHOLD = 0.75;
|
---|
155 |
|
---|
156 | // Define Cone algorithm.
|
---|
157 | C_ADJACENCYCUT = 2;
|
---|
158 | C_MAXITERATIONS = 100;
|
---|
159 | C_IRATCH = 1;
|
---|
160 |
|
---|
161 | //Define MidPoint algorithm.
|
---|
162 | M_CONEAREAFRACTION = 0.25;
|
---|
163 | M_MAXPAIRSIZE = 2;
|
---|
164 | M_MAXITERATIONS = 100;
|
---|
165 |
|
---|
166 | }
|
---|
167 |
|
---|
168 |
|
---|
169 | // **********Provides the smeared TLorentzVector for the electrons********
|
---|
170 | // Smears the electron energy, and changes the 4-momentum accordingly
|
---|
171 | // different smearing if the electron is central (eta < 2.5) or forward
|
---|
172 | void RESOLution::SmearElectron(TLorentzVector &electron) {
|
---|
173 | // the 'electron' variable will be changed by the function
|
---|
174 | float energy = electron.E(); // before smearing
|
---|
175 | float energyS = 0.0; // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
|
---|
176 |
|
---|
177 | if(fabs(electron.Eta()) < MAX_TRACKER) { // if the electron is inside the tracker
|
---|
178 | energyS = gRandom->Gaus(energy, sqrt(
|
---|
179 | pow(ELG_Ncen,2) +
|
---|
180 | pow(ELG_Ccen*energy,2) +
|
---|
181 | pow(ELG_Scen*sqrt(energy),2) ));
|
---|
182 | } else { // outside the tracker
|
---|
183 | energyS = gRandom->Gaus(energy, sqrt(
|
---|
184 | pow(ELG_Nfwd,2) +
|
---|
185 | pow(ELG_Cfwd*energy,2) +
|
---|
186 | pow(ELG_Sfwd*sqrt(energy),2) ) );
|
---|
187 | }
|
---|
188 | electron.SetPtEtaPhiE(energyS/cosh(electron.Eta()), electron.Eta(), electron.Phi(), energyS);
|
---|
189 | if(electron.E() < 0)electron.SetPxPyPzE(0,0,0,0); // no negative values after smearing !
|
---|
190 | }
|
---|
191 |
|
---|
192 |
|
---|
193 | // **********Provides the smeared TLorentzVector for the muons********
|
---|
194 | // Smears the muon pT and changes the 4-momentum accordingly
|
---|
195 | void RESOLution::SmearMu(TLorentzVector &muon) {
|
---|
196 | // the 'muon' variable will be changed by the function
|
---|
197 | float pt = muon.Pt(); // before smearing
|
---|
198 | float ptS = gRandom->Gaus(pt, MU_SmearPt*pt ); // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
|
---|
199 |
|
---|
200 | muon.SetPtEtaPhiE(ptS, muon.Eta(), muon.Phi(), ptS*cosh(muon.Eta()));
|
---|
201 |
|
---|
202 | if(muon.E() < 0)muon.SetPxPyPzE(0,0,0,0); // no negative values after smearing !
|
---|
203 | }
|
---|
204 |
|
---|
205 |
|
---|
206 | // **********Provides the smeared TLorentzVector for the hadrons********
|
---|
207 | // Smears the hadron 4-momentum
|
---|
208 | void RESOLution::SmearHadron(TLorentzVector &hadron, const float frac)
|
---|
209 | // the 'hadron' variable will be changed by the function
|
---|
210 | // the 'frac' variable describes the long-living particles. Should be 0.7 for K0S and Lambda, 1. otherwise
|
---|
211 | {
|
---|
212 | float energy = hadron.E(); // before smearing
|
---|
213 | float energyS = 0.0; // after smearing // \sigma/E = C + N/E + S/\sqrt{E}
|
---|
214 | float energy_ecal = (1.0 - frac)*energy; // electromagnetic calorimeter
|
---|
215 | float energy_hcal = frac*energy; // hadronic calorimeter
|
---|
216 | // frac takes into account the decay of long-living particles, that decay in the calorimeters
|
---|
217 | // some of the particles decay mostly in the ecal, some mostly in the hcal
|
---|
218 |
|
---|
219 | float energyS1,energyS2;
|
---|
220 | if(fabs(hadron.Eta()) < MAX_CALO_CEN) {
|
---|
221 | energyS1 = gRandom->Gaus(energy_hcal, sqrt(
|
---|
222 | pow(HAD_Nhcal,2) +
|
---|
223 | pow(HAD_Chcal*energy_hcal,2) +
|
---|
224 | pow(HAD_Shcal*sqrt(energy_hcal),2) )) ;
|
---|
225 |
|
---|
226 |
|
---|
227 | energyS2 = gRandom->Gaus(energy_ecal, sqrt(
|
---|
228 | pow(ELG_Ncen,2) +
|
---|
229 | pow(ELG_Ccen*energy_ecal,2) +
|
---|
230 | pow(ELG_Scen*sqrt(energy_ecal),2) ) );
|
---|
231 |
|
---|
232 | energyS = ((energyS1>0)?energyS1:0) + ((energyS2>0)?energyS2:0);
|
---|
233 | } else {
|
---|
234 | energyS = gRandom->Gaus(energy, sqrt(
|
---|
235 | pow(HAD_Nhf,2) +
|
---|
236 | pow(HAD_Chf*energy,2) +
|
---|
237 | pow(HAD_Shf*sqrt(energy),2) ));
|
---|
238 | }
|
---|
239 |
|
---|
240 |
|
---|
241 | hadron.SetPtEtaPhiE(energyS/cosh(hadron.Eta()),hadron.Eta(), hadron.Phi(), energyS);
|
---|
242 |
|
---|
243 | if(hadron.E() < 0)hadron.SetPxPyPzE(0,0,0,0);
|
---|
244 | }
|
---|
245 |
|
---|
246 | // **********Provides the energy in the cone of radius TAU_CONE_ENERGY for the tau identification********
|
---|
247 | // to be taken into account, a calo tower should
|
---|
248 | // 1) have a transverse energy \f$ E_T = \sqrt{E_X^2 + E_Y^2} \f$ above a given threshold
|
---|
249 | // 2) be inside a cone with a radius R and the axis defined by (eta,phi)
|
---|
250 | double RESOLution::EnergySmallCone(const vector<PhysicsTower> &towers, const float eta, const float phi) {
|
---|
251 | double Energie=0;
|
---|
252 | for(unsigned int i=0; i < towers.size(); i++) {
|
---|
253 | if(towers[i].fourVector.pt() < SEEDTHRESHOLD) continue;
|
---|
254 | if((DeltaR(phi,eta,towers[i].fourVector.phi(),towers[i].fourVector.eta()) < TAU_CONE_ENERGY)) {
|
---|
255 | Energie += towers[i].fourVector.E;
|
---|
256 | }
|
---|
257 | }
|
---|
258 | return Energie;
|
---|
259 | }
|
---|
260 |
|
---|
261 |
|
---|
262 | // **********Provides the number of tracks in the cone of radius TAU_CONE_TRACKS for the tau identification********
|
---|
263 | // to be taken into account, a track should
|
---|
264 | // 1) avec a transverse momentum \$f p_T \$ above a given threshold
|
---|
265 | // 2) be inside a cone with a radius R and the axis defined by (eta,phi)
|
---|
266 | // IMPORTANT REMARK !!!!!
|
---|
267 | // previously, the argument 'phi' was before the argument 'eta'
|
---|
268 | // this has been changed for consistency with the other functions
|
---|
269 | // double check your running code that uses NumTracks !
|
---|
270 | unsigned int RESOLution::NumTracks(const vector<TLorentzVector> &tracks, const float pt_track, const float eta, const float phi) {
|
---|
271 | unsigned int numtrack=0;
|
---|
272 | for(unsigned int i=0; i < tracks.size(); i++) {
|
---|
273 | if((tracks[i].Pt() < pt_track )||
|
---|
274 | (DeltaR(phi,eta,tracks[i].Phi(),tracks[i].Eta()) > TAU_CONE_TRACKS)
|
---|
275 | )continue;
|
---|
276 | numtrack++;
|
---|
277 | }
|
---|
278 | return numtrack;
|
---|
279 | }
|
---|
280 |
|
---|
281 |
|
---|
282 | //*** Returns the PID of the particle with the highest energy, in a cone with a radius CONERADIUS and an axis (eta,phi) *********
|
---|
283 | //used by Btaggedjet
|
---|
284 | ///// Attention : bug removed => CONERADIUS/2 -> CONERADIUS !!
|
---|
285 | int RESOLution::Bjets(const TSimpleArray<TRootGenParticle> &subarray, const float eta, const float phi) {
|
---|
286 | float emax=0;
|
---|
287 | int Ppid=0;
|
---|
288 | if(subarray.GetEntries()>0) {
|
---|
289 | for(int i=0; i < subarray.GetEntries();i++) { // should have pt>PT_JETMIN and a small cone radius (r<CONE_JET)
|
---|
290 | float genDeltaR = DeltaR(subarray[i]->Phi,subarray[i]->Eta,phi,eta);
|
---|
291 | if(genDeltaR < CONERADIUS && subarray[i]->E > emax) {
|
---|
292 | emax=subarray[i]->E;
|
---|
293 | Ppid=abs(subarray[i]->PID);
|
---|
294 | }
|
---|
295 | }
|
---|
296 | }
|
---|
297 | return Ppid;
|
---|
298 | }
|
---|
299 |
|
---|
300 |
|
---|
301 | //******************** Simulates the b-tagging efficiency for real bjet, or the misendentification for other jets****************
|
---|
302 | bool RESOLution::Btaggedjet(const TLorentzVector &JET, const TSimpleArray<TRootGenParticle> &subarray) {
|
---|
303 | if( rand()%100 < (TAGGING_B+1) && Bjets(subarray,JET.Eta(),JET.Phi())==pB ) return true; // b-tag of b-jets is 40%
|
---|
304 | else if( rand()%100 < (MISTAGGING_C+1) && Bjets(subarray,JET.Eta(),JET.Phi())==pC ) return true; // b-tag of c-jets is 10%
|
---|
305 | else if( rand()%100 < (MISTAGGING_L+1) && Bjets(subarray,JET.Eta(),JET.Phi())!=0) return true; // b-tag of light jets is 1%
|
---|
306 | return false;
|
---|
307 | }
|
---|
308 |
|
---|
309 | //***********************Isolation criteria***********************
|
---|
310 | //****************************************************************
|
---|
311 | bool RESOLution::Isolation(Float_t phi,Float_t eta,const vector<TLorentzVector> &tracks,float PT_TRACK2)
|
---|
312 | {
|
---|
313 | bool isolated = false;
|
---|
314 | Float_t deltar=5000.; // Initial value; should be high; no further repercussion
|
---|
315 | // loop on all final charged particles, with p_t >2, close enough from the electron
|
---|
316 | for(unsigned int i=0; i < tracks.size(); i++)
|
---|
317 | {
|
---|
318 | if(tracks[i].Pt() < PT_TRACK2)continue;
|
---|
319 | Float_t genDeltaR = DeltaR(phi,eta,tracks[i].Phi(),tracks[i].Eta()); // slower to evaluate
|
---|
320 | if(
|
---|
321 | (genDeltaR > deltar) ||
|
---|
322 | (genDeltaR==0)
|
---|
323 | ) continue ;
|
---|
324 | deltar=genDeltaR;
|
---|
325 | }
|
---|
326 | if(deltar > 0.5)isolated = true; // returns the closest distance
|
---|
327 | return isolated;
|
---|
328 | }
|
---|
329 |
|
---|
330 |
|
---|
331 | //**************************** Returns the delta Phi ****************************
|
---|
332 | float DeltaPhi(const float phi1, const float phi2) {
|
---|
333 | float deltaphi=phi1-phi2; // in here, -PI < phi < PI
|
---|
334 | if(fabs(deltaphi) > PI) deltaphi=2.*PI-fabs(deltaphi);// put deltaphi between 0 and PI
|
---|
335 | else deltaphi=fabs(deltaphi);
|
---|
336 |
|
---|
337 | return deltaphi;
|
---|
338 | }
|
---|
339 |
|
---|
340 | //**************************** Returns the delta R****************************
|
---|
341 | float DeltaR(const float phi1, const float eta1, const float phi2, const float eta2) {
|
---|
342 | return sqrt(pow(DeltaPhi(phi1,phi2),2) + pow(eta1-eta2,2));
|
---|
343 | }
|
---|
344 |
|
---|
345 | int sign(const int myint) {
|
---|
346 | if (myint >0) return 1;
|
---|
347 | else if (myint <0) return -1;
|
---|
348 | else return 0;
|
---|
349 | }
|
---|
350 |
|
---|
351 | int sign(const float myfloat) {
|
---|
352 | if (myfloat >0) return 1;
|
---|
353 | else if (myfloat <0) return -1;
|
---|
354 | else return 0;
|
---|
355 | }
|
---|
356 |
|
---|
357 |
|
---|
358 | float Charge(const long int pid) {
|
---|
359 | // source: RPP chap 34 Monte Carlo Particle Numbering Scheme
|
---|
360 | /* switch (abs(pid)) {
|
---|
361 | case 1: case 3: case 5: case 7: return (float) sign(pid)*(-1/3); break; // d, s, b, b'
|
---|
362 | case 2: case 4: case 6: case 8: return (float) sign(pid)*2/3; break; // u, c, t, t'
|
---|
363 |
|
---|
364 | case 11: case 13: case 15: return (float) sign(pid)*(-1); break; // e, mu, tau
|
---|
365 | case 12: case 14: case 16: return (float) 0; break; // nu_e, nu_mu, nu_tau
|
---|
366 |
|
---|
367 | case 9: case 21: case 22: case 23: case 25:
|
---|
368 | case 32: case 33: case 35: case 36: return (float) 0; break; // neutral gauge/higgs bosons
|
---|
369 | case 24: case 34: case 37: return (float) sign(pid); break; // charged gauge/higgs bosons
|
---|
370 | }
|
---|
371 | */
|
---|
372 | return 0;
|
---|
373 | }
|
---|