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

Last change on this file since 1337 was 1337, checked in by garnier, 14 years ago

tag geant4.9.4 beta 1 + modifs locales

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