source: trunk/source/processes/hadronic/models/parton_string/hadronization/src/G4LundStringFragmentation.cc@ 880

Last change on this file since 880 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 25.0 KB
Line 
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: G4LundStringFragmentation.cc,v 1.7 2007/04/24 14:55:23 gunter Exp $
28// GEANT4 tag $Name: geant4-09-01-patch-02 $
29//
30// -----------------------------------------------------------------------------
31// GEANT 4 class implementation file
32//
33// History: first implementation, Maxim Komogorov, 10-Jul-1998
34// -----------------------------------------------------------------------------
35#include "G4LundStringFragmentation.hh"
36#include "G4FragmentingString.hh"
37#include "G4DiQuarks.hh"
38#include "G4Quarks.hh"
39
40#include "Randomize.hh"
41
42// Class G4LundStringFragmentation
43//****************************************************************************************
44
45G4LundStringFragmentation::G4LundStringFragmentation()
46 {
47 MinimalStringMass = 0.; // Uzhi
48 MinimalStringMass2 = 0.; // Uzhi
49 WminLUND = 1.*GeV; // Uzhi
50 SmoothParam = 0.2; // Uzhi
51
52 }
53
54// G4LundStringFragmentation::G4LundStringFragmentation(G4double sigmaPt)
55// : G4VLongitudinalStringDecay(sigmaPt)
56// {
57// }
58
59G4LundStringFragmentation::G4LundStringFragmentation(const G4LundStringFragmentation &) : G4VLongitudinalStringDecay()
60 {
61 }
62
63
64G4LundStringFragmentation::~G4LundStringFragmentation()
65 {
66 }
67
68//****************************************************************************************
69
70const G4LundStringFragmentation & G4LundStringFragmentation::operator=(const G4LundStringFragmentation &)
71 {
72 throw G4HadronicException(__FILE__, __LINE__, "G4LundStringFragmentation::operator= meant to not be accessable");
73 return *this;
74 }
75
76int G4LundStringFragmentation::operator==(const G4LundStringFragmentation &right) const
77 {
78 return !memcmp(this, &right, sizeof(G4LundStringFragmentation));
79 }
80
81int G4LundStringFragmentation::operator!=(const G4LundStringFragmentation &right) const
82 {
83 return memcmp(this, &right, sizeof(G4LundStringFragmentation));
84 }
85
86//****************************************************************************************
87//----------------------------------------------------------------------------------------------------------
88
89G4KineticTrackVector* G4LundStringFragmentation::FragmentString(const G4ExcitedString& theString)
90{
91
92//G4cout<<"In FragmentString"<<G4endl;
93
94// Can no longer modify Parameters for Fragmentation.
95 PastInitPhase=true;
96
97// check if string has enough mass to fragment...
98 G4KineticTrackVector * LeftVector=LightFragmentationTest(&theString);
99 if ( LeftVector != 0 ) {
100//G4cout<<"Return single hadron from string"<<G4endl;
101 return LeftVector;}
102
103 LeftVector = new G4KineticTrackVector;
104 G4KineticTrackVector * RightVector=new G4KineticTrackVector;
105
106// this should work but its only a semi deep copy. %GF G4ExcitedString theStringInCMS(theString);
107 G4ExcitedString *theStringInCMS=CPExcited(theString);
108 G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
109
110 G4bool success=false, inner_sucess=true;
111 G4int attempt=0;
112 while ( !success && attempt++ < StringLoopInterrupt )
113 {
114 G4FragmentingString *currentString=new G4FragmentingString(*theStringInCMS);
115
116//G4cout<<"Main FragmentString cur M2 "<<std::sqrt(currentString->Mass2())<<G4endl;
117
118 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
119 LeftVector->clear();
120 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
121 RightVector->clear();
122
123 inner_sucess=true; // set false on failure..
124 while (! StopFragmenting(currentString) )
125 { // Split current string into hadron + new string
126
127// G4FragmentingString *PreviousString=currentString; // Uzhi
128
129 G4FragmentingString *newString=0; // used as output from SplitUp...
130
131//G4cout<<"FragmentString to Splitup ===================================="<<G4endl;
132//G4int Uzhi; G4cin>>Uzhi; // Uzhi
133 G4KineticTrack * Hadron=Splitup(currentString,newString);
134//G4cout<<" Hadron "<<Hadron<<G4endl;
135
136// if ( Hadron != 0 && IsFragmentable(newString)) // Uzhi
137 if ( Hadron != 0 ) // Uzhi
138 {
139 if ( currentString->GetDecayDirection() > 0 )
140 LeftVector->push_back(Hadron);
141 else
142 RightVector->push_back(Hadron);
143 delete currentString;
144 currentString=newString;
145 } /* else { // Uzhi
146 // abandon ... start from the beginning
147 if (newString) delete newString; // ??? Uzhi local?
148 if (Hadron) delete Hadron;
149// currentString = PreviousString; // Uzhi
150 inner_sucess=false;
151 break;
152 } */ // Uzhi
153// delete PreviousString; // ??? Uzhi local?
154 };
155 // Split current string into 2 final Hadrons
156//G4cout<<"FragmentString to SplitLast if inner_sucess#0"<<inner_sucess<<G4endl;
157 if ( inner_sucess && // Uzhi
158 SplitLast(currentString,LeftVector, RightVector) )
159 {
160 success=true;
161 }
162 delete currentString;
163 }
164
165 delete theStringInCMS;
166
167 if ( ! success )
168 {
169 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
170 LeftVector->clear();
171 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
172 delete RightVector;
173 return LeftVector;
174 }
175
176 // Join Left- and RightVector into LeftVector in correct order.
177 while(!RightVector->empty())
178 {
179 LeftVector->push_back(RightVector->back());
180 RightVector->erase(RightVector->end()-1);
181 }
182 delete RightVector;
183
184//G4cout<<"CalculateHadronTimePosition"<<G4endl;
185
186 CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
187
188 G4LorentzRotation toObserverFrame(toCms.inverse());
189
190 for(size_t C1 = 0; C1 < LeftVector->size(); C1++)
191 {
192 G4KineticTrack* Hadron = LeftVector->operator[](C1);
193 G4LorentzVector Momentum = Hadron->Get4Momentum();
194 Momentum = toObserverFrame*Momentum;
195 Hadron->Set4Momentum(Momentum);
196 G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
197 Momentum = toObserverFrame*Coordinate;
198 Hadron->SetFormationTime(Momentum.e());
199 G4ThreeVector aPosition(Momentum.vect());
200 Hadron->SetPosition(theString.GetPosition()+aPosition);
201 }
202
203//G4cout<<"Out FragmentString"<<G4endl;
204 return LeftVector;
205
206
207
208}
209
210//----------------------------------------------------------------------------------------------------------
211
212//G4double G4LundStringFragmentation::GetLightConeZ(G4double zmin, G4double zmax, // Uzhi
213// G4int , G4ParticleDefinition* pHadron, // Uzhi
214G4double G4LundStringFragmentation::GetLightConeZ(G4double zmin, G4double zmax,
215 G4int, G4ParticleDefinition* pHadron, // Uzhi
216G4double Px, G4double Py)
217 {
218 const G4double alund = 0.7/GeV/GeV;
219
220// If blund get restored, you MUST adapt the calculation of zOfMaxyf.
221// const G4double blund = 1;
222
223 G4double z, yf;
224 G4double Mass = pHadron->GetPDGMass();
225
226 G4double Mt2 = Px*Px + Py*Py + Mass*Mass;
227 G4double zOfMaxyf=alund*Mt2/(alund*Mt2 + 1.);
228 G4double maxYf=(1-zOfMaxyf)/zOfMaxyf * std::exp(-alund*Mt2/zOfMaxyf);
229
230// G4double N=1.; // Uzhi
231// G4double OverN=1./N; // Uzhi
232// G4double ZminN=std::pow(zmin,N); // Uzhi
233// G4double ZmaxN=std::pow(zmax,N); // Uzhi
234// G4double Brac=ZmaxN-ZminN; // Uzhi
235
236//G4cout<<" ZminN ZmaxN Brac Code "<<ZminN<<" "<< ZmaxN<<" "<<Brac<<" "<<PartonEncoding<<G4endl;
237
238// if(std::abs(PartonEncoding) < 1000) // Uzhi
239 { // Uzhi q or q-bar
240//G4cout<<" quark "<<G4endl; // Vova
241 do // Uzhi
242 {
243 z = zmin + G4UniformRand()*(zmax-zmin);
244// yf = std::pow(1. - z, blund)/z*std::exp(-alund*Mt2/z);
245 yf = (1-z)/z * std::exp(-alund*Mt2/z);
246 }
247 while (G4UniformRand()*maxYf > yf);
248 } // Uzhi
249// else // Uzhi
250// { // Uzhi qq or qq-bar
251// //G4cout<<"Di-quark"<<G4endl; // Vova
252// z = std::pow(Brac * G4UniformRand() + ZminN, OverN); // Uzhi
253// }; // Uzhi
254//
255//G4cout<<" test z "<<std::pow(2.,3.)<<" "<<z<<G4endl; // Vova
256 return z;
257 }
258//-----------------------------------------------------------------------------------------
259
260G4LorentzVector * G4LundStringFragmentation::SplitEandP(G4ParticleDefinition * pHadron,
261 G4FragmentingString * string)
262{
263 G4double HadronMass = pHadron->GetPDGMass();
264 SetMinimalStringMass(string); // Uzhi
265 G4double StringMass2 = string->Mass2(); // Uzhi
266
267//G4cout<<"SplitEandP string mass "<<string->Mass()<<" Hadron mass "<<HadronMass<<pHadron->GetParticleName()<<G4endl; // Uzhi
268//G4cout<<string->GetLeftParton()->GetPDGEncoding()<<" "<<G4endl;
269//G4cout<<string->GetRightParton()->GetPDGEncoding()<<" "<<G4endl;
270//G4cout<<" Min string mass "<<MinimalStringMass<<G4endl;
271
272 // calculate and assign hadron transverse momentum component HadronPx andHadronPy
273 G4ThreeVector thePt;
274 thePt=SampleQuarkPt();
275
276 G4ThreeVector HadronPt = thePt +string->DecayPt();
277 HadronPt.setZ(0);
278 //... sample z to define hadron longitudinal momentum and energy
279 //... but first check the available phase space
280
281// G4double DecayQuarkMass2 = sqr(string->GetDecayParton()->GetPDGMass());
282
283//G4cout<<" QuarkMass "<<string->GetDecayParton()->GetPDGMass()<<G4endl; // Uzhi
284
285 G4double HadronMass2T = sqr(HadronMass) + HadronPt.mag2();
286// G4double ResidualMass2T=sqr(MinimalStringMass + WminLUND) + HadronPt.mag2(); // Uzhi
287 G4double ResidualMass2T=sqr(MinimalStringMass + WminLUND) + HadronPt.mag2(); // Uzhi
288
289//G4cout<<" Mt h res str "<<std::sqrt(HadronMass2T)<<" "<<std::sqrt(ResidualMass2T)<<" srt mass"<<string->Mass()<<G4endl;
290
291// if (DecayQuarkMass2 + HadronMass2T >= SmoothParam*(string->Mass2()) ) // Uzhi
292
293 G4double Pz2 = (sqr(StringMass2 - HadronMass2T - ResidualMass2T) - // Uzhi
294 4*HadronMass2T * ResidualMass2T)/4./StringMass2; // Uzhi
295
296//G4cout<<" Pz**2 "<<Pz2<<G4endl;
297
298 if(Pz2 < 0 ) {return 0;} // have to start all over! // Uzhi
299
300 //... then compute allowed z region z_min <= z <= z_max
301
302 G4double Pz = std::sqrt(Pz2); // Uzhi
303 G4double zMin = (std::sqrt(HadronMass2T+Pz2) - Pz)/std::sqrt(StringMass2); // Uzhi
304 G4double zMax = (std::sqrt(HadronMass2T+Pz2) + Pz)/std::sqrt(StringMass2); // Uzhi
305
306//G4cout<<" Zmin max "<<zMin<<" "<<zMax<<G4endl; // Uzhi
307
308// G4double zMax = 1. - DecayQuarkMass2/(string->Mass2()); // Uzhi
309 if (zMin >= zMax) return 0; // have to start all over!
310
311 G4double z = GetLightConeZ(zMin, zMax,
312 string->GetDecayParton()->GetPDGEncoding(), pHadron,
313 HadronPt.x(), HadronPt.y());
314
315 //... now compute hadron longitudinal momentum and energy
316 // longitudinal hadron momentum component HadronPz
317
318 HadronPt.setZ(0.5* string->GetDecayDirection() *
319 (z * string->LightConeDecay() -
320 HadronMass2T/(z * string->LightConeDecay())));
321
322 G4double HadronE = 0.5* (z * string->LightConeDecay() +
323 HadronMass2T/(z * string->LightConeDecay()));
324
325 G4LorentzVector * a4Momentum= new G4LorentzVector(HadronPt,HadronE);
326
327//G4cout<<"Out of SplitEandP Pz E "<<HadronPt.getZ()<<" "<<0.5* (z * string->LightConeDecay() + HadronMass2T/(z * string->LightConeDecay()))<<G4endl;
328
329 return a4Momentum;
330}
331
332
333//-----------------------------------------------------------------------------------------
334
335G4bool G4LundStringFragmentation::SplitLast(G4FragmentingString * string,
336 G4KineticTrackVector * LeftVector,
337 G4KineticTrackVector * RightVector)
338{
339 //... perform last cluster decay
340//G4cout<<"SplitLast String mass "<<string->Mass()<<G4endl;
341//G4cout<<string->GetLeftParton()->GetPDGEncoding()<<" "<<G4endl;
342//G4cout<<string->GetRightParton()->GetPDGEncoding()<<" "<<G4endl;
343
344 G4ThreeVector ClusterVel =string->Get4Momentum().boostVector();
345
346 G4double ResidualMass = string->Mass();
347// G4double ClusterMassCut = ClusterMass;
348
349 G4int cClusterInterrupt = 0;
350
351 G4ParticleDefinition * LeftHadron, * RightHadron;
352 do
353 {
354//G4cout<<" Cicle "<<cClusterInterrupt<<" "<< ClusterLoopInterrupt<<G4endl;
355
356 if (cClusterInterrupt++ >= ClusterLoopInterrupt)
357 {
358 return false;
359 }
360 G4ParticleDefinition * quark = NULL;
361 string->SetLeftPartonStable(); // to query quark contents..
362
363 if (string->DecayIsQuark() && string->StableIsQuark() )
364 {
365 //... there are quarks on cluster ends
366 LeftHadron= QuarkSplitup(string->GetLeftParton(), quark);
367 } else {
368 //... there is a Diquark on cluster ends
369 G4int IsParticle;
370
371 if ( string->StableIsQuark() ) {
372 IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? -1 : +1;
373 } else {
374 IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? +1 : -1;
375 }
376
377 pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
378 quark = QuarkPair.second;
379
380 LeftHadron=hadronizer->Build(QuarkPair.first, string->GetLeftParton());
381 }
382
383 RightHadron = hadronizer->Build(string->GetRightParton(), quark);
384
385//G4cout<<"SplitLast Left Right hadrons "<<LeftHadron->GetPDGEncoding()<<" "<<RightHadron->GetPDGEncoding()<<G4endl;
386//G4cout<<"SplitLast Left Right hadrons "<<LeftHadron->GetPDGMass()<<" "<<RightHadron->GetPDGMass()<<G4endl;
387//G4cout<<" Sum H mass Str Mass "<<LeftHadron->GetPDGMass() + RightHadron->GetPDGMass()<<" "<<ResidualMass<<G4endl;
388
389 //... repeat procedure, if mass of cluster is too low to produce hadrons
390 //... ClusterMassCut = 0.15*GeV model parameter
391// if ( quark->GetParticleSubType()== "quark" ) {ClusterMassCut = 0.;} // Uzhi
392// else {ClusterMassCut = ClusterMass;} // Uzhi
393
394 }
395 while (ResidualMass <= LeftHadron->GetPDGMass() + RightHadron->GetPDGMass()); // Uzhi VOVA
396// while (ResidualMass <= LeftHadron->GetPDGMass() + RightHadron->GetPDGMass() + ClusterMassCut); // Uzhi
397 //... compute hadron momenta and energies
398
399 G4LorentzVector LeftMom, RightMom;
400 G4ThreeVector Pos;
401
402//G4cout<<"Sample4Momentum"<<G4endl;
403
404 Sample4Momentum(&LeftMom, LeftHadron->GetPDGMass(), &RightMom, RightHadron->GetPDGMass(), ResidualMass);
405
406 LeftMom.boost(ClusterVel);
407 RightMom.boost(ClusterVel);
408
409 LeftVector->push_back(new G4KineticTrack(LeftHadron, 0, Pos, LeftMom));
410 RightVector->push_back(new G4KineticTrack(RightHadron, 0, Pos, RightMom));
411
412 return true;
413
414}
415//----------------------------------------------------------------------------------------------------------
416
417G4bool G4LundStringFragmentation::IsFragmentable(const G4FragmentingString * const string)
418{
419//G4cout<<"In IsFragmentable"<<G4endl;
420 SetMinimalStringMass(string); // Uzhi
421//G4cout<<"Out IsFragmentable MinMass"<<MinimalStringMass<<" String Mass"<<std::sqrt(string->Get4Momentum().mag2())<<G4endl;
422 return sqr(MinimalStringMass + WminLUND) < string->Get4Momentum().mag2(); // Uzhi
423
424// return sqr(FragmentationMass(string)+MassCut) < // Uzhi
425// string->Mass2(); // Uzhi
426}
427
428//----------------------------------------------------------------------------------------------------------
429
430G4bool G4LundStringFragmentation::StopFragmenting(const G4FragmentingString * const string)
431{
432//G4cout<<"StopFragmenting"<<G4endl;
433
434 SetMinimalStringMass(string); // Uzhi
435//G4cout<<"StopFragm MinMass "<<MinimalStringMass<<" String Mass "<<std::sqrt(string->Get4Momentum().mag2())<<G4endl;
436 return sqr((MinimalStringMass + WminLUND)*(1 + SmoothParam * (1.-2*G4UniformRand()))) > // Uzhi
437 string->Get4Momentum().mag2(); // Uzhi
438// sqr(FragmentationMass(string,&G4HadronBuilder::BuildHighSpin)+MassCut) > // Uzhi
439// string->Get4Momentum().mag2(); // Uzhi
440}
441
442//----------------------------------------------------------------------------------------------------------
443
444void G4LundStringFragmentation::Sample4Momentum(G4LorentzVector* Mom, G4double Mass, G4LorentzVector* AntiMom, G4double AntiMass, G4double InitialMass)
445 {
446 G4ThreeVector Pt; // Uzhi
447 G4double MassMt2, AntiMassMt2; // Uzhi
448 G4double AvailablePz, AvailablePz2; // Uzhi
449
450//G4cout<<" Smpl4Mom "<<Mass<<" "<<AntiMass<<" "<<InitialMass<<G4endl;
451 // Uzhi
452 do // Uzhi
453 { // Uzhi
454 Pt=SampleQuarkPt(); Pt.setZ(0); G4double Pt2=Pt.mag2(); // Uzhi
455
456//G4cout<<"Sample4Momentum Pt x y "<<Pt.getX()<<" "<<Pt.getY()<<G4endl;
457
458 MassMt2 = Mass * Mass + Pt2; // Uzhi
459 AntiMassMt2= AntiMass * AntiMass + Pt2; // Uzhi
460
461//G4cout<<"Mts "<<MassMt2<<" "<<AntiMassMt2<<" "<<InitialMass*InitialMass<<G4endl;
462
463 AvailablePz2= sqr(InitialMass*InitialMass - MassMt2 - AntiMassMt2) -
464 4.*MassMt2*AntiMassMt2; // Uzhi
465 } // Uzhi
466 while(AvailablePz2 < 0.); // Uzhi
467 // Uzhi
468 AvailablePz2 /=(4.*InitialMass*InitialMass); // Uzhi
469 // Uzhi
470 AvailablePz = std::sqrt(AvailablePz2); // Uzhi
471
472//G4cout<<"AvailablePz "<<AvailablePz<<G4endl;
473
474
475 G4double Px=Pt.getX(); // Uzhi
476 G4double Py=Pt.getY(); // Uzhi
477 // Uzhi
478 Mom->setPx(Px); Mom->setPy(Py); Mom->setPz(AvailablePz); // Uzhi
479 Mom->setE(std::sqrt(MassMt2+AvailablePz2)); // Uzhi
480
481//G4cout<<" 1 part "<<Px<<" "<<Py<<" "<<AvailablePz<<" "<<std::sqrt(MassMt2+AvailablePz2)<<G4endl;
482
483 // Uzhi
484 AntiMom->setPx(-Px); AntiMom->setPy(-Py); AntiMom->setPz(-AvailablePz); // Uzhi
485 AntiMom->setE (std::sqrt(AntiMassMt2+AvailablePz2)); // Uzhi
486
487//G4cout<<" 2 part "<<-Px<<" "<<-Py<<" "<<-AvailablePz<<" "<<std::sqrt(AntiMassMt2+AvailablePz2)<<G4endl;
488
489// Maybe it must be inversed! // Uzhi
490/* // Uzhi
491 G4double r_val = sqr(InitialMass*InitialMass - Mass*Mass - AntiMass*AntiMass) - sqr(2.*Mass*AntiMass);
492 G4double Pabs = (r_val > 0.)? std::sqrt(r_val)/(2.*InitialMass) : 0;
493
494 //... sample unit vector
495 G4double pz = 1. - 2.*G4UniformRand();
496 G4double st = std::sqrt(1. - pz * pz)*Pabs;
497 G4double phi = 2.*pi*G4UniformRand();
498 G4double px = st*std::cos(phi);
499 G4double py = st*std::sin(phi);
500 pz *= Pabs;
501
502 Mom->setPx(px); Mom->setPy(py); Mom->setPz(pz);
503 Mom->setE(std::sqrt(Pabs*Pabs + Mass*Mass));
504
505 AntiMom->setPx(-px); AntiMom->setPy(-py); AntiMom->setPz(-pz);
506 AntiMom->setE (std::sqrt(Pabs*Pabs + AntiMass*AntiMass));
507*/ // Uzhi
508 }
509
510
511void G4LundStringFragmentation::SetMinimalStringMass(const G4FragmentingString * const string) // Uzhi
512{
513//G4cout<<"In SetMinMass -------------------"<<std::sqrt(string->Mass2())<<G4endl;
514//G4cout<<string->GetLeftParton()->GetPDGEncoding()<<G4endl;
515//G4cout<<string->GetRightParton()->GetPDGEncoding()<<G4endl;
516
517 G4double EstimatedMass=0.750* GeV; // 2*m_q
518
519 G4int Qleft =std::abs(string->GetLeftParton()->GetPDGEncoding());
520
521 if( Qleft > 1000)
522 {
523 G4int q1=Qleft/1000;
524 if( q1 < 3) {EstimatedMass += 0.325* GeV;}
525 if( q1 > 2) {EstimatedMass += 0.500* GeV;}
526
527 G4int q2=(Qleft/100)%10;
528 if( q2 < 3) {EstimatedMass += 0.325* GeV;}
529 if( q2 > 2) {EstimatedMass += 0.500* GeV;}
530 }
531 else
532 {
533 if( Qleft < 3) {EstimatedMass += 0.325* GeV;}
534 if( Qleft > 2) {EstimatedMass += 0.500* GeV;}
535 }
536
537 G4int Qright=std::abs(string->GetRightParton()->GetPDGEncoding());
538
539 if( Qright > 1000)
540 {
541 G4int q1=Qright/1000;
542 if( q1 < 3) {EstimatedMass += 0.325* GeV;}
543 if( q1 > 2) {EstimatedMass += 0.500* GeV;}
544
545 G4int q2=(Qright/100)%10;
546 if( q2 < 3) {EstimatedMass += 0.325* GeV;}
547 if( q2 > 2) {EstimatedMass += 0.500* GeV;}
548 }
549 else
550 {
551 if( Qright < 3) {EstimatedMass += 0.325* GeV;}
552 if( Qright > 2) {EstimatedMass += 0.500* GeV;}
553 }
554
555 MinimalStringMass=EstimatedMass;
556 SetMinimalStringMass2(EstimatedMass);
557
558/*
559 Pcreate build=&G4HadronBuilder::BuildLowSpin;
560
561 G4ParticleDefinition *Hadron1, *Hadron2=0;
562
563 G4int iflc = (G4UniformRand() < 0.5)? 1 : 2;
564
565 if (string->GetLeftParton()->GetParticleSubType() == "quark") iflc = -iflc;
566
567 if (string->GetLeftParton()->GetPDGEncoding() < 0) iflc = -iflc;
568
569 // 1/2 baryon (anti-baryon) and scalar meson (QQ-q or QbarQbar-Qbar),
570 // or 2 scalar mesons (Q-Qbar),
571 // or 2 1/2 baryons (anti-baryons) will be built (QQ-QbarQbar)
572
573//G4cout<<"In SetMinMass -------------------"<<std::sqrt(string->Mass2())<<G4endl;
574//G4cout<<string->GetLeftParton()->GetPDGEncoding()<<" "<<FindParticle(iflc)->GetPDGEncoding()<<G4endl;
575//G4cout<<string->GetRightParton()->GetPDGEncoding()<<" "<<FindParticle(-iflc)->GetPDGEncoding()<<G4endl;
576
577 Hadron1 = (hadronizer->*build)(string->GetLeftParton(),FindParticle(iflc));
578 Hadron2 =(hadronizer->*build)(string->GetRightParton(),FindParticle(-iflc));
579 MinimalStringMass = (Hadron1)->GetPDGMass() + (Hadron2)->GetPDGMass();
580
581//G4cout<<(Hadron1)->GetPDGEncoding()<<" "<<(Hadron2)->GetPDGEncoding()<<G4endl;
582//G4cout<<"Out SetMinMass "<<MinimalStringMass<<G4endl;
583*/
584// SetMinimalStringMass2(MinimalStringMass);
585}
586//*******************************************************************************************************
587
588void G4LundStringFragmentation::SetMinimalStringMass2(const G4double aValue) // Uzhi
589{
590 MinimalStringMass2=aValue * aValue;
591}
592//*******************************************************************************************************
593
594//****************************************************************************************
Note: See TracBrowser for help on using the repository browser.