Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4QGSMFragmentation.cc
Go to the documentation of this file.
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$
28//
29// -----------------------------------------------------------------------------
30// GEANT 4 class implementation file
31//
32// History: first implementation, Maxim Komogorov, 10-Jul-1998
33// -----------------------------------------------------------------------------
36#include "Randomize.hh"
37#include "G4ios.hh"
39#include "G4DiQuarks.hh"
40#include "G4Quarks.hh"
41
42// Class G4QGSMFragmentation
43//****************************************************************************************
44
46arho(0.5), aphi(0.), an(-0.5), ala(-0.75), aksi(-1.), alft(0.5)
47 {
48 }
49
51 {
52 }
53
54//----------------------------------------------------------------------------------------------------------
55
57{
58// Can no longer modify Parameters for Fragmentation.
59 PastInitPhase=true;
60
61// check if string has enough mass to fragment...
62 G4KineticTrackVector * LeftVector=LightFragmentationTest(&theString);
63 if ( LeftVector != 0 ) return LeftVector;
64
65 LeftVector = new G4KineticTrackVector;
67
68// this should work but its only a semi deep copy. %GF G4ExcitedString theStringInCMS(theString);
69 G4ExcitedString *theStringInCMS=CPExcited(theString);
70 G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms();
71
72 G4bool success=false, inner_sucess=true;
73 G4int attempt=0;
74 while ( !success && attempt++ < StringLoopInterrupt )
75 {
76 G4FragmentingString *currentString=new G4FragmentingString(*theStringInCMS);
77
78 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
79 LeftVector->clear();
80 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
81 RightVector->clear();
82
83 inner_sucess=true; // set false on failure..
84 while (! StopFragmenting(currentString) )
85 { // Split current string into hadron + new string
86 G4FragmentingString *newString=0; // used as output from SplitUp...
87 G4KineticTrack * Hadron=Splitup(currentString,newString);
88 if ( Hadron != 0 && IsFragmentable(newString))
89 {
90 if ( currentString->GetDecayDirection() > 0 )
91 LeftVector->push_back(Hadron);
92 else
93 RightVector->push_back(Hadron);
94 delete currentString;
95 currentString=newString;
96 } else {
97 // abandon ... start from the beginning
98 if (newString) delete newString; // Uzhi restore 20.06.08
99 if (Hadron) delete Hadron;
100 inner_sucess=false;
101 break;
102 }
103 }
104 // Split current string into 2 final Hadrons
105 if ( inner_sucess &&
106 SplitLast(currentString,LeftVector, RightVector) )
107 {
108 success=true;
109 }
110 delete currentString;
111 }
112
113 delete theStringInCMS;
114
115 if ( ! success )
116 {
117 std::for_each(LeftVector->begin(), LeftVector->end(), DeleteKineticTrack());
118 LeftVector->clear();
119 std::for_each(RightVector->begin(), RightVector->end(), DeleteKineticTrack());
120 delete RightVector;
121 return LeftVector;
122 }
123
124 // Join Left- and RightVector into LeftVector in correct order.
125 while(!RightVector->empty())
126 {
127 LeftVector->push_back(RightVector->back());
128 RightVector->erase(RightVector->end()-1);
129 }
130 delete RightVector;
131
132 CalculateHadronTimePosition(theString.Get4Momentum().mag(), LeftVector);
133
134 G4LorentzRotation toObserverFrame(toCms.inverse());
135
136 for(size_t C1 = 0; C1 < LeftVector->size(); C1++)
137 {
138 G4KineticTrack* Hadron = LeftVector->operator[](C1);
139 G4LorentzVector Momentum = Hadron->Get4Momentum();
140 Momentum = toObserverFrame*Momentum;
141 Hadron->Set4Momentum(Momentum);
142 G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime());
143 Momentum = toObserverFrame*Coordinate;
144 Hadron->SetFormationTime(Momentum.e());
145 G4ThreeVector aPosition(Momentum.vect());
146 Hadron->SetPosition(theString.GetPosition()+aPosition);
147 }
148 return LeftVector;
149
150
151
152}
153
154//----------------------------------------------------------------------------------------------------------
155
156G4double G4QGSMFragmentation::GetLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding, G4ParticleDefinition* , G4double , G4double )
157{
158 G4double z;
159 G4double theA(0), d1, d2, yf;
160 G4int absCode = std::abs( PartonEncoding );
161 if (absCode < 10)
162 {
163 if(absCode == 1 || absCode == 2) theA = arho;
164 else if(absCode == 3) theA = aphi;
165 else throw G4HadronicException(__FILE__, __LINE__, "Unknown PDGencoding in G4QGSMFragmentation::G4LightConeZ");
166
167 do
168 {
169 z = zmin + G4UniformRand() * (zmax - zmin);
170 d1 = (1. - z);
171 d2 = (alft - theA);
172 yf = std::pow(d1, d2);
173 }
174 while (G4UniformRand() > yf);
175 }
176 else
177 {
178 if(absCode == 1103 || absCode == 2101 ||
179 absCode == 2203 || absCode == 2103)
180 {
181 d2 = (alft - (2.*an - arho));
182 }
183 else if(absCode == 3101 || absCode == 3103 ||
184 absCode == 3201 || absCode == 3203)
185 {
186 d2 = (alft - (2.*ala - arho));
187 }
188 else
189 {
190 d2 = (alft - (2.*aksi - arho));
191 }
192
193 do
194 {
195 z = zmin + G4UniformRand() * (zmax - zmin);
196 d1 = (1. - z);
197 yf = std::pow(d1, d2);
198 }
199 while (G4UniformRand() > yf);
200 }
201 return z;
202}
203//-----------------------------------------------------------------------------------------
204
205G4LorentzVector * G4QGSMFragmentation::SplitEandP(G4ParticleDefinition * pHadron,
206 G4FragmentingString * string, // Uzhi
207 G4FragmentingString * ) // Uzhi
208{
209 G4double HadronMass = pHadron->GetPDGMass();
210
211 // calculate and assign hadron transverse momentum component HadronPx andHadronPy
212 G4ThreeVector thePt;
213 thePt=SampleQuarkPt();
214 G4ThreeVector HadronPt = thePt +string->DecayPt();
215 HadronPt.setZ(0);
216 //... sample z to define hadron longitudinal momentum and energy
217 //... but first check the available phase space
218 G4double DecayQuarkMass2 = sqr(string->GetDecayParton()->GetPDGMass());
219 G4double HadronMass2T = sqr(HadronMass) + HadronPt.mag2();
220 if (DecayQuarkMass2 + HadronMass2T >= SmoothParam*(string->Mass2()) )
221 return 0; // have to start all over!
222
223 //... then compute allowed z region z_min <= z <= z_max
224
225 G4double zMin = HadronMass2T/(string->Mass2());
226 G4double zMax = 1. - DecayQuarkMass2/(string->Mass2());
227 if (zMin >= zMax) return 0; // have to start all over!
228
229 G4double z = GetLightConeZ(zMin, zMax,
230 string->GetDecayParton()->GetPDGEncoding(), pHadron,
231 HadronPt.x(), HadronPt.y());
232
233 //... now compute hadron longitudinal momentum and energy
234 // longitudinal hadron momentum component HadronPz
235
236 HadronPt.setZ(0.5* string->GetDecayDirection() *
237 (z * string->LightConeDecay() -
238 HadronMass2T/(z * string->LightConeDecay())));
239 G4double HadronE = 0.5* (z * string->LightConeDecay() +
240 HadronMass2T/(z * string->LightConeDecay()));
241
242 G4LorentzVector * a4Momentum= new G4LorentzVector(HadronPt,HadronE);
243
244 return a4Momentum;
245}
246
247
248//-----------------------------------------------------------------------------------------
249
250G4bool G4QGSMFragmentation::SplitLast(G4FragmentingString * string,
251 G4KineticTrackVector * LeftVector,
252 G4KineticTrackVector * RightVector)
253{
254 //... perform last cluster decay
255 G4ThreeVector ClusterVel =string->Get4Momentum().boostVector();
256 G4double ResidualMass =string->Mass();
257 G4double ClusterMassCut = ClusterMass;
258 G4int cClusterInterrupt = 0;
259 G4ParticleDefinition * LeftHadron, * RightHadron;
260 do
261 {
262 if (cClusterInterrupt++ >= ClusterLoopInterrupt)
263 {
264 return false;
265 }
266 G4ParticleDefinition * quark = NULL;
267 string->SetLeftPartonStable(); // to query quark contents..
268 if (string->DecayIsQuark() && string->StableIsQuark() )
269 {
270 //... there are quarks on cluster ends
271 LeftHadron= QuarkSplitup(string->GetLeftParton(), quark);
272 } else {
273 //... there is a Diquark on cluster ends
274 G4int IsParticle;
275 if ( string->StableIsQuark() ) {
276 IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? -1 : +1;
277 } else {
278 IsParticle=(string->GetLeftParton()->GetPDGEncoding()>0) ? +1 : -1;
279 }
280 pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
281 quark = QuarkPair.second;
282 LeftHadron=hadronizer->Build(QuarkPair.first, string->GetLeftParton());
283 }
284 RightHadron = hadronizer->Build(string->GetRightParton(), quark);
285
286 //... repeat procedure, if mass of cluster is too low to produce hadrons
287 //... ClusterMassCut = 0.15*GeV model parameter
288 if ( quark->GetParticleSubType()== "quark" ) {ClusterMassCut = 0.;}
289 else {ClusterMassCut = ClusterMass;}
290 }
291 while (ResidualMass <= LeftHadron->GetPDGMass() + RightHadron->GetPDGMass() + ClusterMassCut);
292
293 //... compute hadron momenta and energies
294 G4LorentzVector LeftMom, RightMom;
296 Sample4Momentum(&LeftMom, LeftHadron->GetPDGMass(), &RightMom, RightHadron->GetPDGMass(), ResidualMass);
297 LeftMom.boost(ClusterVel);
298 RightMom.boost(ClusterVel);
299 LeftVector->push_back(new G4KineticTrack(LeftHadron, 0, Pos, LeftMom));
300 RightVector->push_back(new G4KineticTrack(RightHadron, 0, Pos, RightMom));
301
302 return true;
303
304}
305
306//----------------------------------------------------------------------------------------------------------
307
308G4bool G4QGSMFragmentation::IsFragmentable(const G4FragmentingString * const string)
309{
310 return sqr(FragmentationMass(string)+MassCut) <
311 string->Mass2();
312}
313
314//----------------------------------------------------------------------------------------------------------
315
316G4bool G4QGSMFragmentation::StopFragmenting(const G4FragmentingString * const string)
317{
318 return
320 string->Get4Momentum().mag2();
321}
322
323//----------------------------------------------------------------------------------------------------------
324
325//----------------------------------------------------------------------------------------------------------
326
327void G4QGSMFragmentation::Sample4Momentum(G4LorentzVector* Mom, G4double Mass, G4LorentzVector* AntiMom, G4double AntiMass, G4double InitialMass)
328 {
329 G4double r_val = sqr(InitialMass*InitialMass - Mass*Mass - AntiMass*AntiMass) - sqr(2.*Mass*AntiMass);
330 G4double Pabs = (r_val > 0.)? std::sqrt(r_val)/(2.*InitialMass) : 0;
331
332 //... sample unit vector
333 G4double pz = 1. - 2.*G4UniformRand();
334 G4double st = std::sqrt(1. - pz * pz)*Pabs;
335 G4double phi = 2.*pi*G4UniformRand();
336 G4double px = st*std::cos(phi);
337 G4double py = st*std::sin(phi);
338 pz *= Pabs;
339
340 Mom->setPx(px); Mom->setPy(py); Mom->setPz(pz);
341 Mom->setE(std::sqrt(Pabs*Pabs + Mass*Mass));
342
343 AntiMom->setPx(-px); AntiMom->setPy(-py); AntiMom->setPz(-pz);
344 AntiMom->setE (std::sqrt(Pabs*Pabs + AntiMass*AntiMass));
345 }
346
347
348//*********************************************************************************************
CLHEP::HepLorentzVector G4LorentzVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define C1
#define G4UniformRand()
Definition: Randomize.hh:53
double x() const
double mag2() const
double y() const
void setZ(double)
HepLorentzRotation inverse() const
HepLorentzVector & boost(double, double, double)
Hep3Vector vect() const
const G4ThreeVector & GetPosition() const
G4LorentzRotation TransformToAlignedCms()
G4LorentzVector Get4Momentum() const
G4ParticleDefinition * GetLeftParton(void) const
G4ParticleDefinition * GetRightParton(void) const
G4ParticleDefinition * GetDecayParton() const
G4int GetDecayDirection() const
G4ParticleDefinition * Build(G4ParticleDefinition *black, G4ParticleDefinition *white)
G4ParticleDefinition * BuildHighSpin(G4ParticleDefinition *black, G4ParticleDefinition *white)
G4double GetFormationTime() const
void SetPosition(const G4ThreeVector aPosition)
void Set4Momentum(const G4LorentzVector &a4Momentum)
const G4ThreeVector & GetPosition() const
void SetFormationTime(G4double aFormationTime)
const G4LorentzVector & Get4Momentum() const
const G4String & GetParticleSubType() const
virtual G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)
G4ThreeVector SampleQuarkPt(G4double ptMax=-1.)
G4KineticTrack * Splitup(G4FragmentingString *string, G4FragmentingString *&newString)
std::pair< G4ParticleDefinition *, G4ParticleDefinition * > pDefPair
G4ParticleDefinition * QuarkSplitup(G4ParticleDefinition *decay, G4ParticleDefinition *&created)
pDefPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true)
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
G4ExcitedString * CPExcited(const G4ExcitedString &string)
G4double FragmentationMass(const G4FragmentingString *const string, Pcreate build=0, pDefPair *pdefs=0)
G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString *const theString)
ush Pos
Definition: deflate.h:82
const G4double pi
T sqr(const T &x)
Definition: templates.hh:145