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

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

update ti head

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