source: trunk/source/processes/hadronic/models/parton_string/hadronization/src/G4VKinkyStringDecay.cc@ 1036

Last change on this file since 1036 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 4.5 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: G4VKinkyStringDecay.cc,v 1.4 2008/04/25 14:20:14 vuzhinsk Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29// Maxim Komogorov
30//
31// -----------------------------------------------------------------------------
32// GEANT 4 class implementation file
33//
34// History: first implementation, Maxim Komogorov, 10-Oct-1998
35// -----------------------------------------------------------------------------
36
37#include "G4VKinkyStringDecay.hh"
38#include "G4KineticTrackVector.hh"
39#include "G4KineticTrack.hh"
40#include "Randomize.hh"
41
42//*****************************************************************************************************
43
44G4VKinkyStringDecay::G4VKinkyStringDecay(G4VLongitudinalStringDecay* theModal)
45 {
46 this->SetLongitudinalStringDecay(theModal);
47 }
48
49//*****************************************************************************************************
50
51G4double G4VKinkyStringDecay::GetLightConeGluonZ(G4double zmin, G4double zmax)
52 {
53 G4double z, yf;
54 do {
55 z = zmin + G4UniformRand()*(zmax-zmin);
56 yf = z*z +sqr(1 - z);
57 }
58 while (G4UniformRand() > yf);
59 return z;
60 }
61
62//*****************************************************************************************************
63
64G4KineticTrackVector* G4VKinkyStringDecay::FragmentString(const G4ExcitedString& String)
65 {
66 G4LorentzVector Mom = String.GetGluon()->Get4Momentum();
67 G4ThreeVector Pos = String.GetGluon()->GetPosition();
68 G4int QuarkEncoding = theLongitudinalStringDecay->SampleQuarkFlavor();
69 G4ThreeVector Pquark=theLongitudinalStringDecay->SampleQuarkPt();
70 G4double Pt2 = Pquark.mag2();
71 G4double z = GetLightConeGluonZ(0, 1);
72 G4double w = Mom.e() + Mom.pz();
73 //... now compute quark longitudinal momentum and energy
74
75 Pquark.setZ( (z*w - Pt2/(z*w))*0.5);
76 G4double E = (z*w + Pt2/(z*w))*0.5;
77
78 G4Parton* AntiColor = new G4Parton(-QuarkEncoding);
79 AntiColor->SetPosition(Pos);
80 G4LorentzVector AntiColorMom(-Pquark, E);
81 AntiColor->Set4Momentum(AntiColorMom);
82 G4Parton* Color = new G4Parton(*String.GetColorParton());
83 G4ExcitedString Str1(Color, AntiColor, String.GetDirection());
84 G4KineticTrackVector* KTV1 = theLongitudinalStringDecay->FragmentString(Str1);
85
86 Color = new G4Parton(QuarkEncoding);
87 Color->SetPosition(Pos);
88 G4LorentzVector ColorMom(Pquark, E);
89 Color->Set4Momentum(ColorMom);
90 AntiColor = new G4Parton(*String.GetAntiColorParton());
91 G4ExcitedString Str2(Color, AntiColor, String.GetDirection());
92 G4KineticTrackVector* KTV2 = theLongitudinalStringDecay->FragmentString(Str2);
93
94 if (KTV1 && KTV2)
95 while(!KTV2->empty())
96 {
97 KTV1->push_back(KTV2->back());
98 KTV1->erase(KTV1->end()-1);
99 }
100 return KTV1;
101 }
102
103//*****************************************************************************************************
104
105
106
107
108
109
110
111
112
113
114
115
116
117
Note: See TracBrowser for help on using the repository browser.