]>
Commit | Line | Data |
---|---|---|
b1c2e580 | 1 | /* |
2 | ||
3 | July 2008 BW mass is limited by "PYTHIA method", by I. Lokhtin and L. Malinina | |
4 | ||
5 | ||
6 | Nikolai Amelin, Ludmila Malinina, Timur Pocheptsov (C) JINR/Dubna | |
7 | amelin@sunhe.jinr.ru, malinina@sunhe.jinr.ru, pocheptsov@sunhe.jinr.ru | |
8 | November. 2, 2005 | |
9 | ||
10 | */ | |
11 | ||
12 | #include <functional> | |
13 | #include <algorithm> | |
14 | #include <vector> | |
15 | #include <iostream> | |
16 | #include <TRandom.h> | |
17 | #include <TError.h> | |
18 | #include <TMath.h> | |
19 | ||
3fa37a65 | 20 | #ifndef DATABASEPDG_H |
b1c2e580 | 21 | #include "DatabasePDG.h" |
22 | #endif | |
23 | #ifndef PARTICLE_PDG | |
24 | #include "ParticlePDG.h" | |
25 | #endif | |
26 | #ifndef DECAY_CHANNEL | |
27 | #include "DecayChannel.h" | |
28 | #endif | |
29 | #ifndef HADRONDECAYER_INCLUDED | |
30 | #include "HadronDecayer.h" | |
31 | #endif | |
32 | #ifndef UKUTILITY_INCLUDED | |
33 | #include "UKUtility.h" | |
34 | #endif | |
35 | #ifndef PARTICLE_INCLUDED | |
36 | #include "Particle.h" | |
37 | #endif | |
38 | #include "HYJET_COMMONS.h" | |
39 | ||
40 | //calculates decay time in fm/c | |
41 | //calculates 1,2 and 3 body decays | |
42 | ||
43 | using std::cout; | |
44 | using std::endl; | |
45 | ||
46 | Double_t GetDecayTime(const Particle &parent, Double_t weakDecayLimit) { | |
47 | ParticlePDG *pDef = parent.Def(); | |
48 | Double_t width = pDef->GetWidth(); //GeV | |
7b7936e9 | 49 | |
50 | // if particle is set to be stable then return 0 | |
51 | if(pDef->GetStableStatus()) | |
52 | return 0.; | |
53 | ||
54 | if(width > weakDecayLimit && weakDecayLimit>=0.0) { | |
b1c2e580 | 55 | const Double_t slope = parent.E() * 0.1973 / (pDef->GetMass() * width); |
b1c2e580 | 56 | return -slope * TMath::Log(gRandom->Rndm());//in fm/c |
57 | } | |
58 | ||
59 | return 0.; | |
60 | } | |
61 | ||
62 | ||
63 | extern "C" void mydelta_(); | |
64 | extern SERVICEEVCommon SERVICEEV; | |
65 | ||
66 | void Decay(List_t &output, Particle &parent, ParticleAllocator &allocator, DatabasePDG* database) { | |
b1c2e580 | 67 | |
68 | // Get the PDG properties of the particle | |
69 | ParticlePDG *pDef = parent.Def(); | |
70 | ||
71 | // Get the number of posible decay channels | |
72 | Int_t nDecayChannel = pDef->GetNDecayChannels(); | |
73 | ||
74 | // check the full branching of this specie | |
75 | Double_t fullBranching = pDef->GetFullBranching(); // Only 3 or less body decays | |
76 | ||
77 | // return if particle has no branching | |
78 | if(fullBranching < 0.00001) | |
79 | return; | |
80 | ||
81 | // get the PDG mass of the specie | |
82 | Double_t PDGmass = pDef->GetMass(); | |
3fa37a65 | 83 | Int_t ComprCodePyth=0; |
b1c2e580 | 84 | Float_t Delta =0; |
85 | ||
86 | Bool_t success = kFALSE; | |
87 | Int_t iterations = 0; | |
88 | // Try to decay the particle | |
89 | while(!success) { | |
7b7936e9 | 90 | // std::cout << "HadronDecayer::Decay() iteration #" << iterations << std::endl; |
91 | if(iterations>1000) { | |
92 | std::cout << "HadronDecayer::Decay() more than 1000 iterations to decay particle with code " | |
93 | << pDef->GetPDG() << std::endl; | |
94 | std::cout << " Will be left undecayed ... check it out!" << std::endl; | |
95 | return; | |
96 | } | |
3fa37a65 | 97 | |
b1c2e580 | 98 | // get a random mass using the Breit Wigner distribution |
3fa37a65 | 99 | Double_t BWmass = gRandom->BreitWigner(PDGmass, pDef->GetWidth()); |
b1c2e580 | 100 | // BWmass = PDGmass; |
101 | // Try to cut the Breit Wigner tail of the particle using the cuts from pythia | |
102 | // The Delta variable is obtained from pythia based on the specie | |
103 | int encoding =pDef->GetPDG(); | |
104 | SERVICEEV.ipdg = encoding; | |
105 | mydelta_(); | |
106 | ComprCodePyth=SERVICEEV.KC; | |
107 | Delta = SERVICEEV.delta;// PYDAT2.PMAS[KC][3]; | |
108 | ||
109 | //if there are no such particle in PYTHIA particle table, we take Delta=0.4 | |
110 | if(ComprCodePyth==0){ | |
111 | BWmass=PDGmass; | |
112 | Delta=0.0; | |
113 | } | |
114 | ||
115 | //bad delta - an exception | |
116 | if(ComprCodePyth==254){ | |
117 | BWmass=PDGmass; | |
118 | Delta=0.0; | |
119 | } | |
120 | ||
7b7936e9 | 121 | // K0 decay into K0s or K0l |
122 | if(TMath::Abs(encoding)==311) { | |
123 | BWmass=PDGmass; | |
124 | Delta=0.0; | |
125 | } | |
126 | ||
b1c2e580 | 127 | //for particles from PYTHIA table only, if the BW mass is outside the cut range then quit this iteration and generate another BW mass |
128 | if(ComprCodePyth!=0 && Delta>0 && (BWmass<PDGmass-Delta || BWmass>PDGmass+Delta)){ | |
7b7936e9 | 129 | iterations++; |
b1c2e580 | 130 | continue; |
131 | } | |
b1c2e580 | 132 | |
133 | // check how many decay channels are allowed with the generated mass | |
134 | Int_t nAllowedChannels = database->GetNAllowedChannels(pDef, BWmass); | |
135 | // if no decay channels are posible with this mass, then generate another BW mass | |
136 | if(nAllowedChannels==0) { | |
137 | iterations++; | |
3fa37a65 | 138 | std::cout << "HadronDecayer::Decay() no decays allowed at this BW mass" << std::endl; |
b1c2e580 | 139 | continue; |
140 | } | |
141 | ||
142 | std::vector<Particle> apDaughter; | |
143 | std::vector<Double_t> dMass; //daughters'mass | |
144 | std::vector<Double_t> dMom; | |
145 | std::vector<Double_t> sm; | |
146 | std::vector<Double_t> rd; | |
147 | ||
148 | // we need to choose an allowed decay channel | |
149 | Double_t randValue = gRandom->Rndm() * fullBranching; | |
150 | Int_t chosenChannel = 1000; | |
151 | Bool_t found = kFALSE; | |
152 | Int_t channelIterations = 0; | |
153 | while(!found) { | |
3fa37a65 | 154 | for(Int_t nChannel = 0; nChannel < nDecayChannel; nChannel++) { |
b1c2e580 | 155 | randValue -= pDef->GetDecayChannel(nChannel)->GetBranching(); |
156 | if(randValue <= 0. && database->IsChannelAllowed(pDef->GetDecayChannel(nChannel), BWmass)) { | |
157 | chosenChannel = nChannel; | |
158 | found = kTRUE; | |
159 | break; | |
160 | } | |
161 | } | |
162 | channelIterations++; | |
163 | } | |
164 | ||
165 | // get the PDG information for the chosen decay channel | |
166 | DecayChannel *dc = pDef->GetDecayChannel(chosenChannel); | |
167 | Int_t nSec = dc->GetNDaughters(); | |
168 | ||
169 | // Adjust the parent momentum four-vector for the MC generated Breit-Wigner mass | |
170 | Particle parentBW(database->GetPDGParticle(parent.Encoding())); | |
171 | parentBW.Pos(parent.Pos()); | |
172 | Double_t BWenergy = TMath::Sqrt(parent.Mom().X()*parent.Mom().X() + | |
173 | parent.Mom().Y()*parent.Mom().Y() + | |
174 | parent.Mom().Z()*parent.Mom().Z() + | |
175 | BWmass*BWmass); | |
176 | ||
177 | Int_t NB = (Int_t)parent.GetType(); //particle from jets | |
178 | ||
179 | TLorentzVector MomparentBW(parent.Mom().X(), parent.Mom().Y(), parent.Mom().Z(), BWenergy); | |
180 | parentBW.Mom(MomparentBW); | |
181 | // take into account BW when calculating boost velocity (for wide resonances it matters) | |
182 | TVector3 velocityBW(parentBW.Mom().BoostVector()); | |
183 | ||
184 | // now we have an allowed decay | |
185 | // first case: one daughter particle | |
186 | if(nSec == 1) { | |
187 | // initialize the daughter particle | |
188 | Particle p1(database->GetPDGParticle(dc->GetDaughterPDG(0))); | |
7b7936e9 | 189 | p1.Pos(parentBW.Pos()); |
190 | p1.Mom(parent.Mom()); | |
b1c2e580 | 191 | p1.SetLastMotherPdg(parentBW.Encoding()); |
192 | p1.SetLastMotherDecayCoor(parentBW.Pos()); | |
193 | p1.SetLastMotherDecayMom(parentBW.Mom()); | |
7b7936e9 | 194 | p1.SetType(NB); |
b1c2e580 | 195 | |
3fa37a65 | 196 | // add the daughter to the list of secondaries |
197 | Int_t parentIndex = parent.GetIndex(); | |
b1c2e580 | 198 | Int_t p1Index = p1.SetIndex(); // set the daughter index |
199 | p1.SetMother(parentIndex); // set the mother index for this daughter | |
3fa37a65 | 200 | parent.SetFirstDaughterIndex(p1Index); |
201 | parent.SetLastDaughterIndex(p1Index); | |
b1c2e580 | 202 | allocator.AddParticle(p1, output); |
203 | success = kTRUE; | |
204 | } | |
205 | // second case: two daughter particles | |
206 | else if(nSec == 2) { | |
207 | // initialize the daughter particles | |
208 | Particle p1(database->GetPDGParticle(dc->GetDaughterPDG(0))); | |
209 | p1.Pos(parentBW.Pos()); | |
210 | Particle p2(database->GetPDGParticle(dc->GetDaughterPDG(1))); | |
211 | p2.Pos(parentBW.Pos()); | |
212 | ||
213 | // calculate the momenta in rest frame of mother for the two particles (theta and phi are isotropic) | |
214 | MomAntiMom(p1.Mom(), p1.TableMass(), p2.Mom(), p2.TableMass(), BWmass); | |
215 | ||
216 | // boost to the laboratory system (to the mother velocity) | |
217 | p1.Mom().Boost(velocityBW); | |
218 | p2.Mom().Boost(velocityBW); | |
219 | ||
220 | //store information about mother | |
221 | p1.SetLastMotherPdg(parentBW.Encoding()); | |
222 | p1.SetLastMotherDecayCoor(parentBW.Pos()); | |
223 | p1.SetLastMotherDecayMom(parentBW.Mom()); | |
224 | p2.SetLastMotherPdg(parentBW.Encoding()); | |
225 | p2.SetLastMotherDecayCoor(parentBW.Pos()); | |
226 | p2.SetLastMotherDecayMom(parentBW.Mom()); | |
b1c2e580 | 227 | //set to daughters the same type as has mother |
228 | p1.SetType(NB); | |
229 | p2.SetType(NB); | |
230 | ||
231 | ||
232 | // check the kinematics in the lab system | |
233 | Double_t deltaS = TMath::Sqrt((parentBW.Mom().X()-p1.Mom().X()-p2.Mom().X())*(parentBW.Mom().X()-p1.Mom().X()-p2.Mom().X())+ | |
234 | (parentBW.Mom().Y()-p1.Mom().Y()-p2.Mom().Y())*(parentBW.Mom().Y()-p1.Mom().Y()-p2.Mom().Y())+ | |
235 | (parentBW.Mom().Z()-p1.Mom().Z()-p2.Mom().Z())*(parentBW.Mom().Z()-p1.Mom().Z()-p2.Mom().Z())+ | |
236 | (parentBW.Mom().E()-p1.Mom().E()-p2.Mom().E())*(parentBW.Mom().E()-p1.Mom().E()-p2.Mom().E())); | |
237 | // if deltaS is too big then repeat the kinematic procedure | |
238 | ||
b1c2e580 | 239 | if(deltaS>0.001) { |
b1c2e580 | 240 | iterations++; |
241 | continue; | |
242 | } | |
243 | // push particles to the list of secondaries | |
3fa37a65 | 244 | Int_t parentIndex = parent.GetIndex(); |
b1c2e580 | 245 | p1.SetIndex(); |
246 | p2.SetIndex(); | |
247 | p1.SetMother(parentIndex); | |
248 | p2.SetMother(parentIndex); | |
3fa37a65 | 249 | parent.SetFirstDaughterIndex(p1.GetIndex()); |
250 | parent.SetLastDaughterIndex(p2.GetIndex()); | |
b1c2e580 | 251 | allocator.AddParticle(p1, output); |
252 | allocator.AddParticle(p2, output); | |
253 | success = kTRUE; | |
3fa37a65 | 254 | |
b1c2e580 | 255 | } |
256 | ||
257 | // third case: three daughter particle | |
258 | else if(nSec == 3) { | |
259 | // initialize the daughter particle | |
260 | Particle p1(database->GetPDGParticle(dc->GetDaughterPDG(0))); | |
261 | p1.Pos(parentBW.Pos()); | |
262 | Particle p2(database->GetPDGParticle(dc->GetDaughterPDG(1))); | |
263 | p2.Pos(parentBW.Pos()); | |
264 | Particle p3(database->GetPDGParticle(dc->GetDaughterPDG(2))); | |
265 | p3.Pos(parentBW.Pos()); | |
266 | // calculate the momenta in the rest frame of the mother particle | |
267 | Double_t pAbs1 = 0., pAbs2 = 0., pAbs3 = 0., sumPabs = 0., maxPabs = 0.; | |
268 | Double_t mass1 = p1.TableMass(), mass2 = p2.TableMass(), mass3 = p3.TableMass(); | |
269 | TLorentzVector &mom1 = p1.Mom(), &mom2 = p2.Mom(), &mom3 = p3.Mom(); | |
270 | Double_t deltaMass = BWmass - mass1 - mass2 - mass3; | |
271 | ||
272 | do { | |
273 | Double_t rd1 = gRandom->Rndm(); | |
274 | Double_t rd2 = gRandom->Rndm(); | |
275 | if (rd2 > rd1) | |
276 | std::swap(rd1, rd2); | |
277 | // 1 | |
278 | Double_t e = rd2*deltaMass; | |
279 | pAbs1 = TMath::Sqrt(e*e + 2*e*mass1); | |
280 | sumPabs = pAbs1; | |
281 | maxPabs = sumPabs; | |
282 | // 2 | |
283 | e = (1-rd1)*deltaMass; | |
284 | pAbs2 = TMath::Sqrt(e*e + 2*e*mass2); | |
285 | ||
286 | if(pAbs2 > maxPabs) | |
287 | maxPabs = pAbs2; | |
288 | ||
289 | sumPabs += pAbs2; | |
290 | // 3 | |
291 | e = (rd1-rd2)*deltaMass; | |
292 | pAbs3 = TMath::Sqrt(e*e + 2*e*mass3); | |
293 | ||
294 | if (pAbs3 > maxPabs) | |
295 | maxPabs = pAbs3; | |
296 | sumPabs += pAbs3; | |
297 | } while(maxPabs > sumPabs - maxPabs); | |
298 | ||
299 | // isotropic sample first particle 3-momentum | |
300 | Double_t cosTheta = 2*(gRandom->Rndm()) - 1; | |
301 | Double_t sinTheta = TMath::Sqrt(1 - cosTheta*cosTheta); | |
302 | Double_t phi = TMath::TwoPi()*(gRandom->Rndm()); | |
303 | Double_t sinPhi = TMath::Sin(phi); | |
304 | Double_t cosPhi = TMath::Cos(phi); | |
305 | ||
306 | mom1.SetPxPyPzE(sinTheta*cosPhi, sinTheta*sinPhi, cosTheta, 0); | |
307 | mom1 *= pAbs1; | |
308 | // sample rest particle 3-momentum | |
309 | Double_t cosThetaN = (pAbs2*pAbs2 - pAbs3*pAbs3 - pAbs1*pAbs1)/(2*pAbs1*pAbs3); | |
310 | Double_t sinThetaN = TMath::Sqrt(1 - cosThetaN*cosThetaN); | |
311 | Double_t phiN = TMath::TwoPi()*(gRandom->Rndm()); | |
312 | Double_t sinPhiN = TMath::Sin(phiN); | |
313 | Double_t cosPhiN = TMath::Cos(phiN); | |
314 | ||
315 | mom3.SetPxPyPzE(sinThetaN*cosPhiN*cosTheta*cosPhi - sinThetaN*sinPhiN*sinPhi + cosThetaN*sinTheta*cosPhi, | |
316 | sinThetaN*cosPhiN*cosTheta*sinPhi + sinThetaN*sinPhiN*cosPhi + cosThetaN*sinTheta*sinPhi, | |
317 | -sinThetaN*cosPhiN*sinTheta + cosThetaN*cosTheta, | |
318 | 0.); | |
319 | ||
320 | mom3 *= pAbs3*mom3.P(); | |
321 | mom2 = mom1; | |
322 | mom2 += mom3; | |
323 | mom2 *= -1.; | |
324 | // calculate energy | |
325 | mom1.SetE(TMath::Sqrt(mom1.P()*mom1.P() + mass1*mass1)); | |
326 | mom2.SetE(TMath::Sqrt(mom2.P()*mom2.P() + mass2*mass2)); | |
327 | mom3.SetE(TMath::Sqrt(mom3.P()*mom3.P() + mass3*mass3)); | |
328 | ||
329 | // boost to Lab system | |
330 | mom1.Boost(velocityBW); | |
331 | mom2.Boost(velocityBW); | |
332 | mom3.Boost(velocityBW); | |
333 | ||
334 | p1.SetLastMotherPdg(parentBW.Encoding()); | |
335 | p1.SetLastMotherDecayCoor(parentBW.Pos()); | |
336 | p1.SetLastMotherDecayMom(parentBW.Mom()); | |
337 | p2.SetLastMotherPdg(parentBW.Encoding()); | |
338 | p2.SetLastMotherDecayCoor(parentBW.Pos()); | |
339 | p2.SetLastMotherDecayMom(parentBW.Mom()); | |
340 | p3.SetLastMotherPdg(parentBW.Encoding()); | |
341 | p3.SetLastMotherDecayCoor(parentBW.Pos()); | |
342 | p3.SetLastMotherDecayMom(parentBW.Mom()); | |
343 | ||
344 | //set to daughters the same type as has mother | |
345 | p1.SetType(NB); | |
346 | p2.SetType(NB); | |
347 | p3.SetType(NB); | |
b1c2e580 | 348 | |
b1c2e580 | 349 | // energy conservation check in the lab system |
350 | Double_t deltaS = TMath::Sqrt((parentBW.Mom().X()-p1.Mom().X()-p2.Mom().X()-p3.Mom().X())*(parentBW.Mom().X()-p1.Mom().X()-p2.Mom().X()-p3.Mom().X()) + | |
351 | (parentBW.Mom().Y()-p1.Mom().Y()-p2.Mom().Y()-p3.Mom().Y())*(parentBW.Mom().Y()-p1.Mom().Y()-p2.Mom().Y()-p3.Mom().Y()) + | |
352 | (parentBW.Mom().Z()-p1.Mom().Z()-p2.Mom().Z()-p3.Mom().Z())*(parentBW.Mom().Z()-p1.Mom().Z()-p2.Mom().Z()-p3.Mom().Z()) + | |
353 | (parentBW.Mom().E()-p1.Mom().E()-p2.Mom().E()-p3.Mom().E())*(parentBW.Mom().E()-p1.Mom().E()-p2.Mom().E()-p3.Mom().E())); | |
354 | // if deltaS is too big then repeat the kinematic procedure | |
355 | if(deltaS>0.001) { | |
7b7936e9 | 356 | // cout << "3-body decay kinematic check in lab system: " << pDef->GetPDG() << " --> " << p1.Encoding() << " + " << p2.Encoding() << " + " << p3.Encoding() << endl; |
357 | // cout << "Mother (e,px,py,pz): " << parentBW.Mom().E() << "\t" << parentBW.Mom().X() << "\t" << parentBW.Mom().Y() << "\t" << parentBW.Mom().Z() << endl; | |
358 | // cout << "Daughter1 (e,px,py,pz): " << p1.Mom().E() << "\t" << p1.Mom().X() << "\t" << p1.Mom().Y() << "\t" << p1.Mom().Z() << endl; | |
359 | // cout << "Daughter2 (e,px,py,pz): " << p2.Mom().E() << "\t" << p2.Mom().X() << "\t" << p2.Mom().Y() << "\t" << p2.Mom().Z() << endl; | |
360 | // cout << "Daughter3 (e,px,py,pz): " << p3.Mom().E() << "\t" << p3.Mom().X() << "\t" << p3.Mom().Y() << "\t" << p3.Mom().Z() << endl; | |
361 | // cout << "3-body decay delta(sqrtS) = " << deltaS << endl; | |
362 | // cout << "Repeating the decay algorithm..." << endl; | |
b1c2e580 | 363 | |
364 | iterations++; | |
365 | continue; | |
366 | } | |
367 | ||
3fa37a65 | 368 | Int_t parentIndex = parent.GetIndex(); |
b1c2e580 | 369 | p1.SetIndex(); |
370 | p2.SetIndex(); | |
371 | p3.SetIndex(); | |
372 | p1.SetMother(parentIndex); | |
373 | p2.SetMother(parentIndex); | |
374 | p3.SetMother(parentIndex); | |
3fa37a65 | 375 | parent.SetFirstDaughterIndex(p1.GetIndex()); |
376 | parent.SetLastDaughterIndex(p3.GetIndex()); | |
b1c2e580 | 377 | allocator.AddParticle(p1, output); |
378 | allocator.AddParticle(p2, output); | |
379 | allocator.AddParticle(p3, output); | |
380 | success = kTRUE; | |
381 | } | |
382 | } | |
383 | return; | |
384 | } |