1 | // |
---|
2 | // ******************************************************************** |
---|
3 | // * License and Disclaimer * |
---|
4 | // * * |
---|
5 | // * The Geant4 software is copyright of the Copyright Holders of * |
---|
6 | // * the Geant4 Collaboration. It is provided under the terms and * |
---|
7 | // * conditions of the Geant4 Software License, included in the file * |
---|
8 | // * LICENSE and available at http://cern.ch/geant4/license . These * |
---|
9 | // * include a list of copyright holders. * |
---|
10 | // * * |
---|
11 | // * Neither the authors of this software system, nor their employing * |
---|
12 | // * institutes,nor the agencies providing financial support for this * |
---|
13 | // * work make any representation or warranty, express or implied, * |
---|
14 | // * regarding this software system or assume any liability for its * |
---|
15 | // * use. Please see the license in the file LICENSE and URL above * |
---|
16 | // * for the full disclaimer and the limitation of liability. * |
---|
17 | // * * |
---|
18 | // * This code implementation is the result of the scientific and * |
---|
19 | // * technical work of the GEANT4 collaboration. * |
---|
20 | // * By using, copying, modifying or distributing the software (or * |
---|
21 | // * any work based on the software) you agree to acknowledge its * |
---|
22 | // * use in resulting scientific publications, and indicate your * |
---|
23 | // * acceptance of all terms of the Geant4 Software license. * |
---|
24 | // ******************************************************************** |
---|
25 | // |
---|
26 | // |
---|
27 | // $Id: G4VLongitudinalStringDecay.cc,v 1.14 2009/05/22 16:35:47 gunter Exp $ |
---|
28 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
29 | // |
---|
30 | // ----------------------------------------------------------------------------- |
---|
31 | // GEANT 4 class implementation file |
---|
32 | // |
---|
33 | // History: first implementation, Maxim Komogorov, 1-Jul-1998 |
---|
34 | // redesign Gunter Folger, August/September 2001 |
---|
35 | // ----------------------------------------------------------------------------- |
---|
36 | #include "G4ios.hh" |
---|
37 | #include "Randomize.hh" |
---|
38 | #include "G4VLongitudinalStringDecay.hh" |
---|
39 | #include "G4FragmentingString.hh" |
---|
40 | |
---|
41 | #include "G4ParticleDefinition.hh" |
---|
42 | #include "G4ParticleTypes.hh" |
---|
43 | #include "G4ParticleChange.hh" |
---|
44 | #include "G4VShortLivedParticle.hh" |
---|
45 | #include "G4ShortLivedConstructor.hh" |
---|
46 | #include "G4ParticleTable.hh" |
---|
47 | #include "G4ShortLivedTable.hh" |
---|
48 | #include "G4PhaseSpaceDecayChannel.hh" |
---|
49 | #include "G4VDecayChannel.hh" |
---|
50 | #include "G4DecayTable.hh" |
---|
51 | |
---|
52 | #include "G4DiQuarks.hh" |
---|
53 | #include "G4Quarks.hh" |
---|
54 | #include "G4Gluons.hh" |
---|
55 | |
---|
56 | //------------------------debug switches |
---|
57 | //#define DEBUG_LightFragmentationTest 1 |
---|
58 | |
---|
59 | |
---|
60 | //******************************************************************************** |
---|
61 | // Constructors |
---|
62 | |
---|
63 | G4VLongitudinalStringDecay::G4VLongitudinalStringDecay() |
---|
64 | { |
---|
65 | MassCut = 0.35*GeV; |
---|
66 | ClusterMass = 0.15*GeV; |
---|
67 | |
---|
68 | SmoothParam = 0.9; |
---|
69 | StringLoopInterrupt = 1000; |
---|
70 | ClusterLoopInterrupt = 500; |
---|
71 | |
---|
72 | // Changable Parameters below. |
---|
73 | |
---|
74 | SigmaQT = 0.5 * GeV; |
---|
75 | |
---|
76 | StrangeSuppress = 0.44; // 27 % strange quarks produced, ie. u:d:s=1:1:0.27 |
---|
77 | DiquarkSuppress = 0.07; |
---|
78 | DiquarkBreakProb = 0.1; |
---|
79 | |
---|
80 | //... pspin_meson is probability to create vector meson |
---|
81 | pspin_meson = 0.5; |
---|
82 | |
---|
83 | //... pspin_barion is probability to create 3/2 barion |
---|
84 | pspin_barion = 0.5; |
---|
85 | |
---|
86 | //... vectorMesonMix[] is quark mixing parameters for vector mesons (Variable spin = 3) |
---|
87 | vectorMesonMix.resize(6); |
---|
88 | vectorMesonMix[0] = 0.5; |
---|
89 | vectorMesonMix[1] = 0.0; |
---|
90 | vectorMesonMix[2] = 0.5; |
---|
91 | vectorMesonMix[3] = 0.0; |
---|
92 | vectorMesonMix[4] = 1.0; |
---|
93 | vectorMesonMix[5] = 1.0; |
---|
94 | |
---|
95 | //... scalarMesonMix[] is quark mixing parameters for scalar mesons (Variable spin=1) |
---|
96 | scalarMesonMix.resize(6); |
---|
97 | scalarMesonMix[0] = 0.5; |
---|
98 | scalarMesonMix[1] = 0.25; |
---|
99 | scalarMesonMix[2] = 0.5; |
---|
100 | scalarMesonMix[3] = 0.25; |
---|
101 | scalarMesonMix[4] = 1.0; |
---|
102 | scalarMesonMix[5] = 0.5; |
---|
103 | |
---|
104 | // Parameters may be changed until the first fragmentation starts |
---|
105 | PastInitPhase=false; |
---|
106 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion, |
---|
107 | scalarMesonMix,vectorMesonMix); |
---|
108 | Kappa = 1.0 * GeV/fermi; |
---|
109 | |
---|
110 | |
---|
111 | } |
---|
112 | |
---|
113 | |
---|
114 | G4VLongitudinalStringDecay::~G4VLongitudinalStringDecay() |
---|
115 | { |
---|
116 | delete hadronizer; |
---|
117 | } |
---|
118 | |
---|
119 | //============================================================================= |
---|
120 | |
---|
121 | // Operators |
---|
122 | |
---|
123 | //const & G4VLongitudinalStringDecay::operator=(const G4VLongitudinalStringDecay &) |
---|
124 | // { |
---|
125 | // } |
---|
126 | |
---|
127 | //----------------------------------------------------------------------------- |
---|
128 | |
---|
129 | int G4VLongitudinalStringDecay::operator==(const G4VLongitudinalStringDecay &) const |
---|
130 | { |
---|
131 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator== forbidden"); |
---|
132 | return false; |
---|
133 | } |
---|
134 | |
---|
135 | //------------------------------------------------------------------------------------- |
---|
136 | |
---|
137 | int G4VLongitudinalStringDecay::operator!=(const G4VLongitudinalStringDecay &) const |
---|
138 | { |
---|
139 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator!= forbidden"); |
---|
140 | return true; |
---|
141 | } |
---|
142 | |
---|
143 | //*********************************************************************************** |
---|
144 | |
---|
145 | // For changing Mass Cut used for selection of very small mass strings |
---|
146 | void G4VLongitudinalStringDecay::SetMassCut(G4double aValue){MassCut=aValue;} |
---|
147 | |
---|
148 | //----------------------------------------------------------------------------- |
---|
149 | |
---|
150 | // For handling a string with very low mass |
---|
151 | |
---|
152 | G4KineticTrackVector* G4VLongitudinalStringDecay::LightFragmentationTest(const |
---|
153 | G4ExcitedString * const string) |
---|
154 | { |
---|
155 | // Check string decay threshold |
---|
156 | |
---|
157 | G4KineticTrackVector * result=0; // return 0 when string exceeds the mass cut |
---|
158 | |
---|
159 | pDefPair hadrons((G4ParticleDefinition *)0,(G4ParticleDefinition *)0); |
---|
160 | |
---|
161 | G4FragmentingString aString(*string); |
---|
162 | |
---|
163 | if ( sqr(FragmentationMass(&aString,0,&hadrons)+MassCut) < aString.Mass2()) { |
---|
164 | return 0; |
---|
165 | } |
---|
166 | |
---|
167 | // The string mass is very low --------------------------- |
---|
168 | |
---|
169 | result=new G4KineticTrackVector; |
---|
170 | |
---|
171 | if ( hadrons.second ==0 ) |
---|
172 | { |
---|
173 | // Substitute string by light hadron, Note that Energy is not conserved here! |
---|
174 | |
---|
175 | /* |
---|
176 | #ifdef DEBUG_LightFragmentationTest |
---|
177 | G4cout << "VlongSF Warning replacing string by single hadron " <<G4endl; |
---|
178 | G4cout << hadrons.first->GetParticleName() |
---|
179 | << "string .. " << string->Get4Momentum() << " " |
---|
180 | << string->Get4Momentum().m() << G4endl; |
---|
181 | #endif |
---|
182 | G4cout << "VlongSF Warning replacing string by single hadron " <<G4endl; |
---|
183 | G4cout << hadrons.first->GetParticleName() << " string .. " <<G4endl; |
---|
184 | G4cout << string->Get4Momentum() << " " << string->Get4Momentum().m() << G4endl; |
---|
185 | */ |
---|
186 | G4ThreeVector Mom3 = string->Get4Momentum().vect(); |
---|
187 | G4LorentzVector Mom(Mom3, |
---|
188 | std::sqrt(Mom3.mag2() + |
---|
189 | sqr(hadrons.first->GetPDGMass()))); |
---|
190 | result->push_back(new G4KineticTrack(hadrons.first, 0, |
---|
191 | string->GetPosition(), |
---|
192 | Mom)); |
---|
193 | } else |
---|
194 | { |
---|
195 | //... string was qq--qqbar type: Build two stable hadrons, |
---|
196 | |
---|
197 | #ifdef DEBUG_LightFragmentationTest |
---|
198 | G4cout << "VlongSF Warning replacing qq-qqbar string by TWO hadrons " |
---|
199 | << hadrons.first->GetParticleName() << " / " |
---|
200 | << hadrons.second->GetParticleName() |
---|
201 | << "string .. " << string->Get4Momentum() << " " |
---|
202 | << string->Get4Momentum().m() << G4endl; |
---|
203 | #endif |
---|
204 | |
---|
205 | // Uzhi Formation time in the case??? |
---|
206 | |
---|
207 | G4LorentzVector Mom1, Mom2; |
---|
208 | Sample4Momentum(&Mom1, hadrons.first->GetPDGMass(), |
---|
209 | &Mom2,hadrons.second->GetPDGMass(), |
---|
210 | string->Get4Momentum().mag()); |
---|
211 | |
---|
212 | result->push_back(new G4KineticTrack(hadrons.first, 0, |
---|
213 | string->GetPosition(), |
---|
214 | Mom1)); |
---|
215 | result->push_back(new G4KineticTrack(hadrons.second, 0, |
---|
216 | string->GetPosition(), |
---|
217 | Mom2)); |
---|
218 | |
---|
219 | G4ThreeVector Velocity = string->Get4Momentum().boostVector(); |
---|
220 | result->Boost(Velocity); |
---|
221 | } |
---|
222 | |
---|
223 | return result; |
---|
224 | |
---|
225 | } |
---|
226 | |
---|
227 | //---------------------------------------------------------------------------------------- |
---|
228 | |
---|
229 | G4double G4VLongitudinalStringDecay::FragmentationMass( |
---|
230 | const G4FragmentingString * const string, |
---|
231 | Pcreate build, pDefPair * pdefs ) |
---|
232 | { |
---|
233 | |
---|
234 | G4double mass; |
---|
235 | static G4bool NeedInit(true); |
---|
236 | static std::vector<double> nomix; |
---|
237 | static G4HadronBuilder * minMassHadronizer; |
---|
238 | if ( NeedInit ) |
---|
239 | { |
---|
240 | NeedInit = false; |
---|
241 | nomix.resize(6); |
---|
242 | for ( G4int i=0; i<6 ; i++ ) nomix[i]=0; |
---|
243 | |
---|
244 | // minMassHadronizer=new G4HadronBuilder(pspin_meson,pspin_barion,nomix,nomix); |
---|
245 | minMassHadronizer=hadronizer; |
---|
246 | } |
---|
247 | |
---|
248 | if ( build==0 ) build=&G4HadronBuilder::BuildLowSpin; |
---|
249 | |
---|
250 | G4ParticleDefinition *Hadron1, *Hadron2=0; |
---|
251 | |
---|
252 | if (!string->FourQuarkString() ) |
---|
253 | { |
---|
254 | // spin 0 meson or spin 1/2 barion will be built |
---|
255 | |
---|
256 | Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(), |
---|
257 | string->GetRightParton()); |
---|
258 | mass= (Hadron1)->GetPDGMass(); |
---|
259 | } else |
---|
260 | { |
---|
261 | //... string is qq--qqbar: Build two stable hadrons, |
---|
262 | //... with extra uubar or ddbar quark pair |
---|
263 | G4int iflc = (G4UniformRand() < 0.5)? 1 : 2; |
---|
264 | if (string->GetLeftParton()->GetPDGEncoding() < 0) iflc = -iflc; |
---|
265 | |
---|
266 | //... theSpin = 4; spin 3/2 baryons will be built |
---|
267 | Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(), |
---|
268 | FindParticle(iflc) ); |
---|
269 | Hadron2 = (minMassHadronizer->*build)(string->GetRightParton(), |
---|
270 | FindParticle(-iflc) ); |
---|
271 | mass = (Hadron1)->GetPDGMass() + (Hadron2)->GetPDGMass(); |
---|
272 | } |
---|
273 | |
---|
274 | if ( pdefs != 0 ) |
---|
275 | { // need to return hadrons as well.... |
---|
276 | pdefs->first = Hadron1; |
---|
277 | pdefs->second = Hadron2; |
---|
278 | } |
---|
279 | |
---|
280 | return mass; |
---|
281 | } |
---|
282 | |
---|
283 | //---------------------------------------------------------------------------- |
---|
284 | |
---|
285 | G4ParticleDefinition* G4VLongitudinalStringDecay::FindParticle(G4int Encoding) |
---|
286 | { |
---|
287 | G4ParticleDefinition* ptr = G4ParticleTable::GetParticleTable()->FindParticle(Encoding); |
---|
288 | if (ptr == NULL) |
---|
289 | { |
---|
290 | G4cout << "Particle with encoding "<<Encoding<<" does not exist!!!"<<G4endl; |
---|
291 | throw G4HadronicException(__FILE__, __LINE__, "Check your particle table"); |
---|
292 | } |
---|
293 | return ptr; |
---|
294 | } |
---|
295 | |
---|
296 | //----------------------------------------------------------------------------- |
---|
297 | // virtual void Sample4Momentum(G4LorentzVector* Mom, G4double Mass, |
---|
298 | // G4LorentzVector* AntiMom, G4double AntiMass, |
---|
299 | // G4double InitialMass)=0; |
---|
300 | //----------------------------------------------------------------------------- |
---|
301 | |
---|
302 | //********************************************************************************* |
---|
303 | // For decision on continue or stop string fragmentation |
---|
304 | // virtual G4bool StopFragmenting(const G4FragmentingString * const string)=0; |
---|
305 | // virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0; |
---|
306 | |
---|
307 | // If a string can not fragment, make last break into 2 hadrons |
---|
308 | // virtual G4bool SplitLast(G4FragmentingString * string, |
---|
309 | // G4KineticTrackVector * LeftVector, |
---|
310 | // G4KineticTrackVector * RightVector)=0; |
---|
311 | //----------------------------------------------------------------------------- |
---|
312 | // |
---|
313 | // If a string fragments, do the following |
---|
314 | // |
---|
315 | // For transver of a string to its CMS frame |
---|
316 | //----------------------------------------------------------------------------- |
---|
317 | |
---|
318 | G4ExcitedString *G4VLongitudinalStringDecay::CPExcited(const G4ExcitedString & in) |
---|
319 | { |
---|
320 | G4Parton *Left=new G4Parton(*in.GetLeftParton()); |
---|
321 | G4Parton *Right=new G4Parton(*in.GetRightParton()); |
---|
322 | return new G4ExcitedString(Left,Right,in.GetDirection()); |
---|
323 | } |
---|
324 | |
---|
325 | //----------------------------------------------------------------------------- |
---|
326 | |
---|
327 | G4KineticTrack * G4VLongitudinalStringDecay::Splitup( |
---|
328 | G4FragmentingString *string, |
---|
329 | G4FragmentingString *&newString) |
---|
330 | { |
---|
331 | //G4cout<<"In G4VLong String Dec######################"<<G4endl; |
---|
332 | |
---|
333 | //... random choice of string end to use for creating the hadron (decay) |
---|
334 | SideOfDecay = (G4UniformRand() < 0.5)? 1: -1; |
---|
335 | if (SideOfDecay < 0) |
---|
336 | { |
---|
337 | string->SetLeftPartonStable(); |
---|
338 | } else |
---|
339 | { |
---|
340 | string->SetRightPartonStable(); |
---|
341 | } |
---|
342 | |
---|
343 | G4ParticleDefinition *newStringEnd; |
---|
344 | G4ParticleDefinition * HadronDefinition; |
---|
345 | if (string->DecayIsQuark()) |
---|
346 | { |
---|
347 | HadronDefinition= QuarkSplitup(string->GetDecayParton(), newStringEnd); |
---|
348 | } else { |
---|
349 | HadronDefinition= DiQuarkSplitup(string->GetDecayParton(), newStringEnd); |
---|
350 | } |
---|
351 | // create new String from old, ie. keep Left and Right order, but replace decay |
---|
352 | |
---|
353 | newString=new G4FragmentingString(*string,newStringEnd); // To store possible |
---|
354 | // quark containt of new string |
---|
355 | G4LorentzVector* HadronMomentum=SplitEandP(HadronDefinition, string, newString); |
---|
356 | |
---|
357 | delete newString; newString=0; // Uzhi 20.06.08 |
---|
358 | |
---|
359 | G4KineticTrack * Hadron =0; |
---|
360 | if ( HadronMomentum != 0 ) { |
---|
361 | |
---|
362 | G4ThreeVector Pos; |
---|
363 | Hadron = new G4KineticTrack(HadronDefinition, 0,Pos, *HadronMomentum); |
---|
364 | |
---|
365 | newString=new G4FragmentingString(*string,newStringEnd, |
---|
366 | HadronMomentum); |
---|
367 | |
---|
368 | //G4cout<<"Out G4VLong String Dec######################"<<G4endl; |
---|
369 | //G4cout<<"newString 4Mom"<<newString->Get4Momentum()<<G4endl; |
---|
370 | //G4cout<<"newString Pl "<<newString->LightConePlus()<<G4endl; |
---|
371 | //G4cout<<"newString Mi "<<newString->LightConeMinus()<<G4endl; |
---|
372 | //G4cout<<"newString Pts "<<newString->StablePt()<<G4endl; |
---|
373 | //G4cout<<"newString Ptd "<<newString->DecayPt()<<G4endl; |
---|
374 | //G4cout<<"newString M2 "<<newString->Mass2()<<G4endl; |
---|
375 | //G4cout<<"newString M2 "<<newString->Mass()<<G4endl; |
---|
376 | delete HadronMomentum; |
---|
377 | } |
---|
378 | return Hadron; |
---|
379 | } |
---|
380 | |
---|
381 | //-------------------------------------------------------------------------------------- |
---|
382 | |
---|
383 | G4ParticleDefinition * |
---|
384 | G4VLongitudinalStringDecay::QuarkSplitup(G4ParticleDefinition* |
---|
385 | decay, G4ParticleDefinition *&created) |
---|
386 | { |
---|
387 | G4int IsParticle=(decay->GetPDGEncoding()>0) ? -1 : +1; // if we have a quark, |
---|
388 | // we need antiquark |
---|
389 | // (or diquark) |
---|
390 | pDefPair QuarkPair = CreatePartonPair(IsParticle); |
---|
391 | created = QuarkPair.second; |
---|
392 | return hadronizer->Build(QuarkPair.first, decay); |
---|
393 | |
---|
394 | } |
---|
395 | |
---|
396 | //----------------------------------------------------------------------------- |
---|
397 | |
---|
398 | G4ParticleDefinition *G4VLongitudinalStringDecay::DiQuarkSplitup( |
---|
399 | G4ParticleDefinition* decay, |
---|
400 | G4ParticleDefinition *&created) |
---|
401 | { |
---|
402 | //... can Diquark break or not? |
---|
403 | if (G4UniformRand() < DiquarkBreakProb ){ |
---|
404 | //... Diquark break |
---|
405 | |
---|
406 | G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000; |
---|
407 | G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10; |
---|
408 | if (G4UniformRand() < 0.5) |
---|
409 | { |
---|
410 | G4int Swap = stableQuarkEncoding; |
---|
411 | stableQuarkEncoding = decayQuarkEncoding; |
---|
412 | decayQuarkEncoding = Swap; |
---|
413 | } |
---|
414 | |
---|
415 | G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1; |
---|
416 | // if we have a quark, we need antiquark) |
---|
417 | pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted |
---|
418 | //... Build new Diquark |
---|
419 | G4int QuarkEncoding=QuarkPair.second->GetPDGEncoding(); |
---|
420 | G4int i10 = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding)); |
---|
421 | G4int i20 = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding)); |
---|
422 | G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3; |
---|
423 | G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin); |
---|
424 | created = FindParticle(NewDecayEncoding); |
---|
425 | G4ParticleDefinition * decayQuark=FindParticle(decayQuarkEncoding); |
---|
426 | G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decayQuark); |
---|
427 | return had; |
---|
428 | // return hadronizer->Build(QuarkPair.first, decayQuark); |
---|
429 | |
---|
430 | } else { |
---|
431 | //... Diquark does not break |
---|
432 | |
---|
433 | G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1; |
---|
434 | // if we have a diquark, we need quark) |
---|
435 | pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted |
---|
436 | created = QuarkPair.second; |
---|
437 | |
---|
438 | G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay); |
---|
439 | return had; |
---|
440 | // return G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay); |
---|
441 | } |
---|
442 | } |
---|
443 | |
---|
444 | //----------------------------------------------------------------------------- |
---|
445 | |
---|
446 | G4int G4VLongitudinalStringDecay::SampleQuarkFlavor(void) |
---|
447 | { |
---|
448 | return (1 + (int)(G4UniformRand()/StrangeSuppress)); |
---|
449 | } |
---|
450 | |
---|
451 | //----------------------------------------------------------------------------- |
---|
452 | |
---|
453 | G4VLongitudinalStringDecay::pDefPair G4VLongitudinalStringDecay::CreatePartonPair(G4int NeedParticle,G4bool AllowDiquarks) |
---|
454 | { |
---|
455 | // NeedParticle = +1 for Particle, -1 for Antiparticle |
---|
456 | |
---|
457 | if ( AllowDiquarks && G4UniformRand() < DiquarkSuppress ) |
---|
458 | { |
---|
459 | // Create a Diquark - AntiDiquark pair , first in pair is anti to IsParticle |
---|
460 | G4int q1 = SampleQuarkFlavor(); |
---|
461 | G4int q2 = SampleQuarkFlavor(); |
---|
462 | G4int spin = (q1 != q2 && G4UniformRand() <= 0.5)? 1 : 3; |
---|
463 | // convention: quark with higher PDG number is first |
---|
464 | G4int PDGcode = (std::max(q1,q2) * 1000 + std::min(q1,q2) * 100 + spin) * NeedParticle; |
---|
465 | return pDefPair (FindParticle(-PDGcode),FindParticle(PDGcode)); |
---|
466 | |
---|
467 | |
---|
468 | } else { |
---|
469 | // Create a Quark - AntiQuark pair, first in pair IsParticle |
---|
470 | G4int PDGcode=SampleQuarkFlavor()*NeedParticle; |
---|
471 | return pDefPair (FindParticle(PDGcode),FindParticle(-PDGcode)); |
---|
472 | } |
---|
473 | |
---|
474 | } |
---|
475 | |
---|
476 | //----------------------------------------------------------------------------- |
---|
477 | |
---|
478 | // G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt() |
---|
479 | // { |
---|
480 | // G4double width_param= 2.0 * GeV*GeV; |
---|
481 | // G4double R = G4UniformRand(); |
---|
482 | // G4double Pt = std::sqrt(width_param*R/(1-R)); |
---|
483 | // G4double phi = 2.*pi*G4UniformRand(); |
---|
484 | // return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0); |
---|
485 | // } |
---|
486 | |
---|
487 | G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt(G4double ptMax) |
---|
488 | { |
---|
489 | G4double Pt; |
---|
490 | if ( ptMax < 0 ) { |
---|
491 | // sample full gaussian |
---|
492 | Pt = -std::log(G4UniformRand()); |
---|
493 | } else { |
---|
494 | // sample in limited range |
---|
495 | Pt = -std::log(CLHEP::RandFlat::shoot(exp(-sqr(ptMax)/sqr(SigmaQT)), 1.)); |
---|
496 | } |
---|
497 | Pt = SigmaQT * std::sqrt(Pt); |
---|
498 | G4double phi = 2.*pi*G4UniformRand(); |
---|
499 | return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0); |
---|
500 | } |
---|
501 | |
---|
502 | //****************************************************************************** |
---|
503 | |
---|
504 | void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons) |
---|
505 | { |
---|
506 | |
---|
507 | // `yo-yo` formation time |
---|
508 | // const G4double kappa = 1.0 * GeV/fermi/4.; // Uzhi String tension 1.06.08 |
---|
509 | G4double kappa = GetStringTensionParameter(); |
---|
510 | //G4cout<<"Kappa "<<kappa<<G4endl; // Uzhi 20.06.08 |
---|
511 | //G4int Uzhi; G4cin>>Uzhi; // Uzhi 20.06.08 |
---|
512 | for(size_t c1 = 0; c1 < Hadrons->size(); c1++) |
---|
513 | { |
---|
514 | G4double SumPz = 0; |
---|
515 | G4double SumE = 0; |
---|
516 | for(size_t c2 = 0; c2 < c1; c2++) |
---|
517 | { |
---|
518 | SumPz += Hadrons->operator[](c2)->Get4Momentum().pz(); |
---|
519 | SumE += Hadrons->operator[](c2)->Get4Momentum().e(); |
---|
520 | } |
---|
521 | G4double HadronE = Hadrons->operator[](c1)->Get4Momentum().e(); |
---|
522 | G4double HadronPz = Hadrons->operator[](c1)->Get4Momentum().pz(); |
---|
523 | Hadrons->operator[](c1)->SetFormationTime((theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa)); |
---|
524 | G4ThreeVector aPosition(0, 0, (theInitialStringMass - 2.*SumE - HadronE + HadronPz)/(2.*kappa)); |
---|
525 | Hadrons->operator[](c1)->SetPosition(aPosition); |
---|
526 | } |
---|
527 | } |
---|
528 | |
---|
529 | //----------------------------------------------------------------------------- |
---|
530 | |
---|
531 | /* |
---|
532 | void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons) |
---|
533 | { |
---|
534 | // 'constituent' formation time |
---|
535 | const G4double kappa = 1.0 * GeV/fermi; |
---|
536 | for(G4int c1 = 0; c1 < Hadrons->length(); c1++) |
---|
537 | { |
---|
538 | G4double SumPz = 0; |
---|
539 | G4double SumE = 0; |
---|
540 | for(G4int c2 = 0; c2 <= c1; c2++) |
---|
541 | { |
---|
542 | SumPz += Hadrons->at(c2)->Get4Momentum().pz(); |
---|
543 | SumE += Hadrons->at(c2)->Get4Momentum().e(); |
---|
544 | } |
---|
545 | Hadrons->at(c1)->SetFormationTime((theInitialStringMass - 2.*SumPz)/(2.*kappa)); |
---|
546 | G4ThreeVector aPosition(0, 0, (theInitialStringMass - 2.*SumE)/(2.*kappa)); |
---|
547 | Hadrons->at(c1)->SetPosition(aPosition); |
---|
548 | } |
---|
549 | c1 = Hadrons->length()-1; |
---|
550 | Hadrons->at(c1)->SetFormationTime(Hadrons->at(c1-1)->GetFormationTime()); |
---|
551 | Hadrons->at(c1)->SetPosition(Hadrons->at(c1-1)->GetPosition()); |
---|
552 | } |
---|
553 | */ |
---|
554 | |
---|
555 | |
---|
556 | //***************************************************************************** |
---|
557 | |
---|
558 | void G4VLongitudinalStringDecay::SetSigmaTransverseMomentum(G4double aValue) |
---|
559 | { |
---|
560 | if ( PastInitPhase ) { |
---|
561 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetSigmaTransverseMomentum after FragmentString() not allowed"); |
---|
562 | } else { |
---|
563 | SigmaQT = aValue; |
---|
564 | } |
---|
565 | } |
---|
566 | |
---|
567 | //---------------------------------------------------------------------------------------------------------- |
---|
568 | |
---|
569 | void G4VLongitudinalStringDecay::SetStrangenessSuppression(G4double aValue) |
---|
570 | { |
---|
571 | if ( PastInitPhase ) { |
---|
572 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetStrangenessSuppression after FragmentString() not allowed"); |
---|
573 | } else { |
---|
574 | StrangeSuppress = aValue; |
---|
575 | } |
---|
576 | } |
---|
577 | |
---|
578 | //---------------------------------------------------------------------------------------------------------- |
---|
579 | |
---|
580 | void G4VLongitudinalStringDecay::SetDiquarkSuppression(G4double aValue) |
---|
581 | { |
---|
582 | if ( PastInitPhase ) { |
---|
583 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkSuppression after FragmentString() not allowed"); |
---|
584 | } else { |
---|
585 | DiquarkSuppress = aValue; |
---|
586 | } |
---|
587 | } |
---|
588 | |
---|
589 | //---------------------------------------------------------------------------------------- |
---|
590 | |
---|
591 | void G4VLongitudinalStringDecay::SetDiquarkBreakProbability(G4double aValue) |
---|
592 | { |
---|
593 | if ( PastInitPhase ) { |
---|
594 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkBreakProbability after FragmentString() not allowed"); |
---|
595 | } else { |
---|
596 | DiquarkBreakProb = aValue; |
---|
597 | } |
---|
598 | } |
---|
599 | |
---|
600 | //---------------------------------------------------------------------------------------------------------- |
---|
601 | |
---|
602 | void G4VLongitudinalStringDecay::SetVectorMesonProbability(G4double aValue) |
---|
603 | { |
---|
604 | if ( PastInitPhase ) { |
---|
605 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonProbability after FragmentString() not allowed"); |
---|
606 | } else { |
---|
607 | pspin_meson = aValue; |
---|
608 | delete hadronizer; |
---|
609 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion, |
---|
610 | scalarMesonMix,vectorMesonMix); |
---|
611 | } |
---|
612 | } |
---|
613 | |
---|
614 | //---------------------------------------------------------------------------------------------------------- |
---|
615 | |
---|
616 | void G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability(G4double aValue) |
---|
617 | { |
---|
618 | if ( PastInitPhase ) { |
---|
619 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability after FragmentString() not allowed"); |
---|
620 | } else { |
---|
621 | pspin_barion = aValue; |
---|
622 | delete hadronizer; |
---|
623 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion, |
---|
624 | scalarMesonMix,vectorMesonMix); |
---|
625 | } |
---|
626 | } |
---|
627 | |
---|
628 | //---------------------------------------------------------------------------------------------------------- |
---|
629 | |
---|
630 | void G4VLongitudinalStringDecay::SetScalarMesonMixings(std::vector<G4double> aVector) |
---|
631 | { |
---|
632 | if ( PastInitPhase ) { |
---|
633 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings after FragmentString() not allowed"); |
---|
634 | } else { |
---|
635 | if ( aVector.size() < 6 ) |
---|
636 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings( argument Vector too small"); |
---|
637 | scalarMesonMix[0] = aVector[0]; |
---|
638 | scalarMesonMix[1] = aVector[1]; |
---|
639 | scalarMesonMix[2] = aVector[2]; |
---|
640 | scalarMesonMix[3] = aVector[3]; |
---|
641 | scalarMesonMix[4] = aVector[4]; |
---|
642 | scalarMesonMix[5] = aVector[5]; |
---|
643 | delete hadronizer; |
---|
644 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion, |
---|
645 | scalarMesonMix,vectorMesonMix); |
---|
646 | } |
---|
647 | } |
---|
648 | |
---|
649 | //---------------------------------------------------------------------------------------------------------- |
---|
650 | |
---|
651 | void G4VLongitudinalStringDecay::SetVectorMesonMixings(std::vector<G4double> aVector) |
---|
652 | { |
---|
653 | if ( PastInitPhase ) { |
---|
654 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings after FragmentString() not allowed"); |
---|
655 | } else { |
---|
656 | if ( aVector.size() < 6 ) |
---|
657 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings( argument Vector too small"); |
---|
658 | vectorMesonMix[0] = aVector[0]; |
---|
659 | vectorMesonMix[1] = aVector[1]; |
---|
660 | vectorMesonMix[2] = aVector[2]; |
---|
661 | vectorMesonMix[3] = aVector[3]; |
---|
662 | vectorMesonMix[4] = aVector[4]; |
---|
663 | vectorMesonMix[5] = aVector[5]; |
---|
664 | delete hadronizer; |
---|
665 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion, |
---|
666 | scalarMesonMix,vectorMesonMix); |
---|
667 | |
---|
668 | } |
---|
669 | } |
---|
670 | |
---|
671 | //------------------------------------------------------------------------------------------- |
---|
672 | void G4VLongitudinalStringDecay::SetStringTensionParameter(G4double aValue)// Uzhi 20 June 08 |
---|
673 | { |
---|
674 | Kappa = aValue * GeV/fermi; |
---|
675 | } |
---|
676 | //************************************************************************************** |
---|
677 | |
---|