source: trunk/source/processes/hadronic/models/parton_string/diffraction/src/G4FTFParameters.cc@ 1201

Last change on this file since 1201 was 1196, checked in by garnier, 16 years ago

update CVS release candidate geant4.9.3.01

File size: 16.1 KB
Line 
1//********************************************************************
2// * License and Disclaimer *
3// * *
4// * The Geant4 software is copyright of the Copyright Holders of *
5// * the Geant4 Collaboration. It is provided under the terms and *
6// * conditions of the Geant4 Software License, included in the file *
7// * LICENSE and available at http://cern.ch/geant4/license . These *
8// * include a list of copyright holders. *
9// * *
10// * Neither the authors of this software system, nor their employing *
11// * institutes,nor the agencies providing financial support for this *
12// * work make any representation or warranty, express or implied, *
13// * regarding this software system or assume any liability for its *
14// * use. Please see the license in the file LICENSE and URL above *
15// * for the full disclaimer and the limitation of liability. *
16// * *
17// * This code implementation is the result of the scientific and *
18// * technical work of the GEANT4 collaboration. *
19// * By using, copying, modifying or distributing the software (or *
20// * any work based on the software) you agree to acknowledge its *
21// * use in resulting scientific publications, and indicate your *
22// * acceptance of all terms of the Geant4 Software license. *
23// ********************************************************************
24//
25//
26// $Id: G4FTFParameters.cc,v 1.11 2009/10/25 10:50:54 vuzhinsk Exp $
27// GEANT4 tag $Name: geant4-09-03-cand-01 $
28//
29
30#include "G4FTFParameters.hh"
31
32#include "G4ios.hh"
33#include <utility> // Uzhi 29.03.08
34
35G4FTFParameters::G4FTFParameters()
36{;}
37
38
39G4FTFParameters::~G4FTFParameters()
40{;}
41//**********************************************************************************************
42
43G4FTFParameters::G4FTFParameters(const G4ParticleDefinition * particle,
44 G4double theA,
45 G4double theZ,
46 G4double s)
47{
48 G4int PDGcode = particle->GetPDGEncoding();
49 G4int absPDGcode = std::abs(PDGcode);
50 G4double ProjectileMass = particle->GetPDGMass();
51 G4double TargetMass = G4Proton::Proton()->GetPDGMass();
52
53 G4double Elab = (s - ProjectileMass*ProjectileMass - TargetMass*TargetMass)/
54 (2*TargetMass);
55 G4double Plab = std::sqrt(Elab * Elab - ProjectileMass*ProjectileMass);
56
57 G4double LogPlab = std::log( Plab );
58 G4double sqrLogPlab = LogPlab * LogPlab;
59
60 G4int NumberOfTargetProtons = (G4int) theZ;
61 G4int NumberOfTargetNeutrons = (G4int) theA- (G4int) theZ;
62 G4int NumberOfTargetNucleons = NumberOfTargetProtons + NumberOfTargetNeutrons;
63
64 G4double Xtotal, Xelastic;
65
66 if( absPDGcode > 1000 ) //------Projectile is baryon --------
67 {
68 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) + 0.522*sqrLogPlab - 4.51*LogPlab;
69 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) + 0.513*sqrLogPlab - 4.27*LogPlab;
70
71 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
72 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
73
74 Xtotal = ( NumberOfTargetProtons * XtotPP +
75 NumberOfTargetNeutrons * XtotPN ) / NumberOfTargetNucleons;
76 Xelastic = ( NumberOfTargetProtons * XelPP +
77 NumberOfTargetNeutrons * XelPN ) / NumberOfTargetNucleons;
78 }
79 else if( PDGcode == 211 ) //------Projectile is PionPlus -------
80 {
81 G4double XtotPiP = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab;
82 G4double XtotPiN = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab;
83
84 G4double XelPiP = 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab;
85 G4double XelPiN = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab;
86
87 Xtotal = ( NumberOfTargetProtons * XtotPiP +
88 NumberOfTargetNeutrons * XtotPiN ) / NumberOfTargetNucleons;
89 Xelastic = ( NumberOfTargetProtons * XelPiP +
90 NumberOfTargetNeutrons * XelPiN ) / NumberOfTargetNucleons;
91 }
92 else if( PDGcode == -211 ) //------Projectile is PionMinus -------
93 {
94 G4double XtotPiP = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab;
95 G4double XtotPiN = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab;
96
97 G4double XelPiP = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab;
98 G4double XelPiN = 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab;
99
100 Xtotal = ( NumberOfTargetProtons * XtotPiP +
101 NumberOfTargetNeutrons * XtotPiN ) / NumberOfTargetNucleons;
102 Xelastic = ( NumberOfTargetProtons * XelPiP +
103 NumberOfTargetNeutrons * XelPiN ) / NumberOfTargetNucleons;
104 }
105
106 else if( PDGcode == 111 ) //------Projectile is PionZero -------
107 {
108 G4double XtotPiP =(16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+
109 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
110
111 G4double XtotPiN =(33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab + //Pi+
112 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
113
114 G4double XelPiP =( 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+
115 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
116 G4double XelPiN =( 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+
117 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
118
119 Xtotal = ( NumberOfTargetProtons * XtotPiP +
120 NumberOfTargetNeutrons * XtotPiN ) / NumberOfTargetNucleons;
121 Xelastic = ( NumberOfTargetProtons * XelPiP +
122 NumberOfTargetNeutrons * XelPiN ) / NumberOfTargetNucleons;
123 }
124 else if( PDGcode == 321 ) //------Projectile is KaonPlus -------
125 {
126 G4double XtotKP = 18.1 + 0. *std::pow(Plab, 0. ) + 0.26 *sqrLogPlab - 1.0 *LogPlab;
127 G4double XtotKN = 18.7 + 0. *std::pow(Plab, 0. ) + 0.21 *sqrLogPlab - 0.89*LogPlab;
128
129 G4double XelKP = 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab;
130 G4double XelKN = 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab;
131
132 Xtotal = ( NumberOfTargetProtons * XtotKP +
133 NumberOfTargetNeutrons * XtotKN ) / NumberOfTargetNucleons;
134 Xelastic = ( NumberOfTargetProtons * XelKP +
135 NumberOfTargetNeutrons * XelKN ) / NumberOfTargetNucleons;
136 }
137 else if( PDGcode ==-321 ) //------Projectile is KaonMinus ------
138 {
139 G4double XtotKP = 32.1 + 0. *std::pow(Plab, 0. ) + 0.66 *sqrLogPlab - 5.6 *LogPlab;
140 G4double XtotKN = 25.2 + 0. *std::pow(Plab, 0. ) + 0.38 *sqrLogPlab - 2.9 *LogPlab;
141
142 G4double XelKP = 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab;
143 G4double XelKN = 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab;
144
145 Xtotal = ( NumberOfTargetProtons * XtotKP +
146 NumberOfTargetNeutrons * XtotKN ) / NumberOfTargetNucleons;
147 Xelastic = ( NumberOfTargetProtons * XelKP +
148 NumberOfTargetNeutrons * XelKN ) / NumberOfTargetNucleons;
149 }
150 else if((PDGcode == 311) || (PDGcode == 130) || (PDGcode == 310))//Projectile is KaonZero
151 {
152 G4double XtotKP =( 18.1 + 0. *std::pow(Plab, 0. ) + 0.26 *sqrLogPlab - 1.0 *LogPlab + //K+
153 32.1 + 0. *std::pow(Plab, 0. ) + 0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
154 G4double XtotKN =( 18.7 + 0. *std::pow(Plab, 0. ) + 0.21 *sqrLogPlab - 0.89*LogPlab + //K+
155 25.2 + 0. *std::pow(Plab, 0. ) + 0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
156
157 G4double XelKP =( 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+
158 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
159 G4double XelKN =( 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab + //K+
160 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
161 Xtotal = ( NumberOfTargetProtons * XtotKP +
162 NumberOfTargetNeutrons * XtotKN ) / NumberOfTargetNucleons;
163 Xelastic = ( NumberOfTargetProtons * XelKP +
164 NumberOfTargetNeutrons * XelKN ) / NumberOfTargetNucleons;
165 }
166 else //------Projectile is undefined, Nucleon assumed
167 {
168 G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) + 0.522*sqrLogPlab - 4.51*LogPlab;
169 G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) + 0.513*sqrLogPlab - 4.27*LogPlab;
170
171 G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
172 G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
173
174 Xtotal = ( NumberOfTargetProtons * XtotPP +
175 NumberOfTargetNeutrons * XtotPN ) / NumberOfTargetNucleons;
176 Xelastic = ( NumberOfTargetProtons * XelPP +
177 NumberOfTargetNeutrons * XelPN ) / NumberOfTargetNucleons;
178 };
179
180// Xtotal and Xelastic in mb
181
182//----------- Geometrical parameters ------------------------------------------------
183 SetTotalCrossSection(Xtotal);
184 SetElastisCrossSection(Xelastic);
185 SetInelasticCrossSection(Xtotal-Xelastic);
186
187// // Interactions with elastic and inelastic collisions
188 SetProbabilityOfElasticScatt(Xtotal, Xelastic);
189 SetRadiusOfHNinteractions2(Xtotal/pi/10.);
190//
191/* //==== No elastic scattering ============================
192 SetProbabilityOfElasticScatt(Xtotal, 0.);
193 SetRadiusOfHNinteractions2((Xtotal-Xelastic)/pi/10.);
194*/ //=======================================================
195
196//-----------------------------------------------------------------------------------
197 SetSlope( Xtotal*Xtotal/16./pi/Xelastic/0.3894 ); // Slope parameter of elastic scattering
198 // (GeV/c)^(-2))
199//-----------------------------------------------------------------------------------
200 SetGamma0( GetSlope()*Xtotal/10./2./pi );
201
202//----------- Parameters of elastic scattering --------------------------------------
203 // Gaussian parametrization of
204 // elastic scattering amplitude assumed
205 SetAvaragePt2ofElasticScattering(1./(Xtotal*Xtotal/16./pi/Xelastic/0.3894)*GeV*GeV);
206
207//----------- Parameters of excitations ---------------------------------------------
208 if( absPDGcode > 1000 ) //------Projectile is baryon --------
209 {
210 SetMagQuarkExchange(3.4); //3.8);
211 SetSlopeQuarkExchange(1.2);
212 SetDeltaProbAtQuarkExchange(0.1); //(0.1*4.);
213
214 SetProjMinDiffMass(1.1); // GeV
215 SetProjMinNonDiffMass(1.1); // GeV
216 SetProbabilityOfProjDiff(0.76*std::pow(s/GeV/GeV,-0.35));
217
218 SetTarMinDiffMass(1.1); // GeV
219 SetTarMinNonDiffMass(1.1); // GeV
220 SetProbabilityOfTarDiff(0.76*std::pow(s/GeV/GeV,-0.35));
221
222 SetAveragePt2(0.3); // GeV^2
223 }
224 else if( absPDGcode == 211 || PDGcode == 111) //------Projectile is Pion -----------
225 {
226 SetMagQuarkExchange(120.); // 210.
227 SetSlopeQuarkExchange(2.0);
228 SetDeltaProbAtQuarkExchange(0.6);
229
230 SetProjMinDiffMass(0.5); // GeV
231 SetProjMinNonDiffMass(0.3); // GeV
232 SetProbabilityOfProjDiff(0.*0.62*std::pow(s/GeV/GeV,-0.51)); // 40/32 X-dif/X-inel
233//G4cout<<"Params "<<0.6*0.62*std::pow(s/GeV/GeV,-0.51)<<G4endl;
234 SetTarMinDiffMass(1.1); // GeV
235 SetTarMinNonDiffMass(1.1); // GeV
236//G4cout<<" "<<2.7*0.62*std::pow(s/GeV/GeV,-0.51)<<G4endl; // was 2
237//G4int Uzhi; G4cin>>Uzhi;
238 SetProbabilityOfTarDiff(2.*0.62*std::pow(s/GeV/GeV,-0.51)); // 40/32 X-dif/X-inel
239
240/*
241SetProjMinDiffMass(0.5);
242SetProjMinNonDiffMass(0.3); // Uzhi 12.06.08
243SetProbabilityOfProjDiff(0.05);
244SetProbabilityOfTarDiff(0.05);
245*/
246 SetAveragePt2(0.3); // GeV^2
247 }
248 else if( (absPDGcode == 321) || (PDGcode == 311) ||
249 (PDGcode == 130) || (PDGcode == 310)) //Projectile is Kaon
250 {
251// Must be corrected, taken from PiN
252 SetMagQuarkExchange(120.);
253 SetSlopeQuarkExchange(2.0);
254 SetDeltaProbAtQuarkExchange(0.6);
255
256 SetProjMinDiffMass(0.7); // GeV 1.1
257 SetProjMinNonDiffMass(0.7); // GeV
258 SetProbabilityOfProjDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel
259
260 SetTarMinDiffMass(1.1); // GeV
261 SetTarMinNonDiffMass(1.1); // GeV
262 SetProbabilityOfTarDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel
263
264 SetAveragePt2(0.3); // GeV^2
265 }
266 else //------Projectile is undefined,
267 //------Nucleon assumed
268 {
269 SetMagQuarkExchange(3.5);
270 SetSlopeQuarkExchange(1.0);
271 SetDeltaProbAtQuarkExchange(0.1);
272
273 SetProjMinDiffMass((particle->GetPDGMass()+160.*MeV)/GeV);
274 SetProjMinNonDiffMass((particle->GetPDGMass()+160.*MeV)/GeV);
275 SetProbabilityOfProjDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel
276
277 SetTarMinDiffMass(1.1); // GeV
278 SetTarMinNonDiffMass(1.1); // GeV
279 SetProbabilityOfTarDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel
280
281 SetAveragePt2(0.3); // GeV^2
282 }
283
284// ---------- Set parameters of a string kink -------------------------------
285 SetPt2Kink(6.*GeV*GeV);
286 G4double Puubar(1./3.), Pddbar(1./3.), Pssbar(1./3.); // SU(3) symmetry
287// G4double Puubar(0.41 ), Pddbar(0.41 ), Pssbar(0.18 ); // Broken SU(3) symmetry
288 SetQuarkProbabilitiesAtGluonSplitUp(Puubar, Pddbar, Pssbar);
289
290// --------- Set parameters of nuclear destruction--------------------
291
292 if( absPDGcode < 1000 )
293 {
294 SetCofNuclearDestruction(1.); //1.0); // for meson projectile
295 } else if( theA > 20. )
296 {
297 SetCofNuclearDestruction(0.2); //2); // for baryon projectile and heavy target
298 } else
299 {
300 SetCofNuclearDestruction(0.2); //1.0); // for baryon projectile and light target
301 }
302
303 SetR2ofNuclearDestruction(1.5*fermi*fermi);
304
305 SetExcitationEnergyPerWoundedNucleon(100*MeV);
306
307 SetDofNuclearDestruction(0.4);
308 SetPt2ofNuclearDestruction(0.17*GeV*GeV);
309 SetMaxPt2ofNuclearDestruction(1.0*GeV*GeV);
310}
311//**********************************************************************************************
Note: See TracBrowser for help on using the repository browser.