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

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

update to geant4.9.2

File size: 14.4 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.4 2008/12/18 13:02:00 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02 $
29//
30
31#include "G4FTFParameters.hh"
32
33G4FTFParameters::G4FTFParameters()
34{;}
35
36
37G4FTFParameters::~G4FTFParameters()
38{;}
39//**********************************************************************************************
40
41G4FTFParameters::G4FTFParameters(const G4ParticleDefinition * particle, 
42                                                   G4double   theA,
43                                                   G4double   theZ,
44                                                   G4double   s) 
45    {
46    G4int PDGcode = particle->GetPDGEncoding();
47    G4int absPDGcode = std::abs(PDGcode);
48    G4double Elab = (s - 2*0.88*GeV*GeV)/(2*0.939*GeV)/GeV;
49    G4double Plab = std::sqrt(Elab * Elab - 0.88);
50
51    G4double LogPlab = std::log( Plab );
52    G4double sqrLogPlab = LogPlab * LogPlab;
53
54//G4cout<<"G4FTFParameters Plab "<<Plab<<G4endl;
55
56    G4int NumberOfTargetProtons  = (G4int) theZ; 
57    G4int NumberOfTargetNeutrons = (G4int) theA- (G4int) theZ;
58    G4int NumberOfTargetNucleons = NumberOfTargetProtons + NumberOfTargetNeutrons;
59
60    G4double Xtotal, Xelastic;
61
62    if( absPDGcode > 1000 )                        //------Projectile is baryon --------
63      {       
64       G4double XtotPP = 48.0 +  0. *std::pow(Plab, 0.  ) + 0.522*sqrLogPlab - 4.51*LogPlab;
65       G4double XtotPN = 47.3 +  0. *std::pow(Plab, 0.  ) + 0.513*sqrLogPlab - 4.27*LogPlab;
66
67       G4double XelPP  = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
68       G4double XelPN  = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
69
70       Xtotal          = ( NumberOfTargetProtons  * XtotPP + 
71                           NumberOfTargetNeutrons * XtotPN  ) / NumberOfTargetNucleons;
72       Xelastic        = ( NumberOfTargetProtons  * XelPP  + 
73                           NumberOfTargetNeutrons * XelPN   ) / NumberOfTargetNucleons;
74      }
75    else if( PDGcode ==  211 )                     //------Projectile is PionPlus -------
76      {
77       G4double XtotPiP = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab;
78       G4double XtotPiN = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab;
79           
80       G4double XelPiP  =  0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab;
81       G4double XelPiN  = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab;
82
83       Xtotal           = ( NumberOfTargetProtons  * XtotPiP + 
84                            NumberOfTargetNeutrons * XtotPiN  ) / NumberOfTargetNucleons;
85       Xelastic         = ( NumberOfTargetProtons  * XelPiP  + 
86                            NumberOfTargetNeutrons * XelPiN   ) / NumberOfTargetNucleons; 
87      }
88    else if( PDGcode == -211 )                     //------Projectile is PionMinus -------
89      {
90       G4double XtotPiP = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab;
91       G4double XtotPiN = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab;
92           
93       G4double XelPiP  = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab;
94       G4double XelPiN  =  0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab;
95
96       Xtotal           = ( NumberOfTargetProtons  * XtotPiP + 
97                            NumberOfTargetNeutrons * XtotPiN  ) / NumberOfTargetNucleons;
98       Xelastic         = ( NumberOfTargetProtons  * XelPiP  + 
99                            NumberOfTargetNeutrons * XelPiN   ) / NumberOfTargetNucleons;
100      }
101
102    else if( PDGcode ==  111 )                     //------Projectile is PionZero  -------
103      {
104       G4double XtotPiP =(16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab +   //Pi+
105                          33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi-
106
107       G4double XtotPiN =(33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab +   //Pi+
108                          16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
109           
110       G4double XelPiP  =( 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab +    //Pi+
111                           1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
112       G4double XelPiN  =( 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab +   //Pi+
113                           0.0  + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi-
114
115       Xtotal           = ( NumberOfTargetProtons  * XtotPiP + 
116                            NumberOfTargetNeutrons * XtotPiN  ) / NumberOfTargetNucleons;
117       Xelastic         = ( NumberOfTargetProtons  * XelPiP  + 
118                            NumberOfTargetNeutrons * XelPiN   ) / NumberOfTargetNucleons; 
119      }
120    else if( PDGcode == 321 )                      //------Projectile is KaonPlus -------
121      {
122       G4double XtotKP = 18.1 +  0. *std::pow(Plab, 0.  ) + 0.26 *sqrLogPlab - 1.0 *LogPlab;
123       G4double XtotKN = 18.7 +  0. *std::pow(Plab, 0.  ) + 0.21 *sqrLogPlab - 0.89*LogPlab;
124
125       G4double XelKP  =  5.0 +  8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab;
126       G4double XelKN  =  7.3 +  0. *std::pow(Plab,-0.  ) + 0.29 *sqrLogPlab - 2.4 *LogPlab;
127
128       Xtotal          = ( NumberOfTargetProtons  * XtotKP + 
129                           NumberOfTargetNeutrons * XtotKN  ) / NumberOfTargetNucleons;
130       Xelastic        = ( NumberOfTargetProtons  * XelKP  + 
131                           NumberOfTargetNeutrons * XelKN   ) / NumberOfTargetNucleons;
132      }
133    else if( PDGcode ==-321 )                      //------Projectile is KaonMinus ------
134      {
135       G4double XtotKP = 32.1 +  0. *std::pow(Plab, 0.  ) + 0.66 *sqrLogPlab - 5.6 *LogPlab;
136       G4double XtotKN = 25.2 +  0. *std::pow(Plab, 0.  ) + 0.38 *sqrLogPlab - 2.9 *LogPlab;
137
138       G4double XelKP  =  7.3 +  0. *std::pow(Plab,-0.  ) + 0.29 *sqrLogPlab - 2.4 *LogPlab;
139       G4double XelKN  =  5.0 +  8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab;
140
141       Xtotal          = ( NumberOfTargetProtons  * XtotKP + 
142                           NumberOfTargetNeutrons * XtotKN  ) / NumberOfTargetNucleons;
143       Xelastic        = ( NumberOfTargetProtons  * XelKP  + 
144                           NumberOfTargetNeutrons * XelKN   ) / NumberOfTargetNucleons;
145      }
146    else if( PDGcode == 311 )                      //------Projectile is KaonZero ------
147      {
148       G4double XtotKP =( 18.1 +  0. *std::pow(Plab, 0.  ) + 0.26 *sqrLogPlab - 1.0 *LogPlab +   //K+
149                          32.1 +  0. *std::pow(Plab, 0.  ) + 0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K-
150       G4double XtotKN =( 18.7 +  0. *std::pow(Plab, 0.  ) + 0.21 *sqrLogPlab - 0.89*LogPlab +   //K+
151                          25.2 +  0. *std::pow(Plab, 0.  ) + 0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K-
152
153       G4double XelKP  =(  5.0 +  8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab +   //K+
154                           7.3 +  0. *std::pow(Plab,-0.  ) + 0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K-
155       G4double XelKN  =(  7.3 +  0. *std::pow(Plab,-0.  ) + 0.29 *sqrLogPlab - 2.4 *LogPlab +   //K+
156                           5.0 +  8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K-
157       Xtotal          = ( NumberOfTargetProtons  * XtotKP + 
158                           NumberOfTargetNeutrons * XtotKN  ) / NumberOfTargetNucleons;
159       Xelastic        = ( NumberOfTargetProtons  * XelKP  + 
160                           NumberOfTargetNeutrons * XelKN   ) / NumberOfTargetNucleons;
161      }
162    else                 //------Projectile is undefined, Nucleon assumed
163      {
164       G4double XtotPP = 48.0 +  0. *std::pow(Plab, 0.  ) + 0.522*sqrLogPlab - 4.51*LogPlab;
165       G4double XtotPN = 47.3 +  0. *std::pow(Plab, 0.  ) + 0.513*sqrLogPlab - 4.27*LogPlab;
166
167       G4double XelPP  = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
168       G4double XelPN  = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab;
169
170       Xtotal          = ( NumberOfTargetProtons  * XtotPP + 
171                           NumberOfTargetNeutrons * XtotPN  ) / NumberOfTargetNucleons;
172       Xelastic        = ( NumberOfTargetProtons  * XelPP  + 
173                           NumberOfTargetNeutrons * XelPN   ) / NumberOfTargetNucleons;
174      };
175
176//      Xtotal and Xelastic in mb
177
178//----------- Geometrical parameters ------------------------------------------------
179      SetTotalCrossSection(Xtotal);
180      SetElastisCrossSection(Xelastic);
181      SetInelasticCrossSection(Xtotal-Xelastic);
182
183//  // Interactions with elastic ans inelastic collisions
184      SetProbabilityOfElasticScatt(Xtotal, Xelastic);
185      SetRadiusOfHNinteractions2(Xtotal/pi/10.);
186//
187/* //==== No elastic scattering ============================
188      SetProbabilityOfElasticScatt(Xtotal, 0.);
189      SetRadiusOfHNinteractions2((Xtotal-Xelastic)/pi/10.);
190*/ //=======================================================
191
192//G4cout<<" Rnn "<<Xtotal/pi/10.<<" "<<Xtotal/pi/10.*fermi*fermi<<G4endl;
193//G4cout<<"G4FTFParameters Xt Xel MeV "<<Xtotal<<" "<<Xelastic<<" "<<GeV<<G4endl;
194
195//----------------------------------------------------------------------------------- 
196      SetSlope( Xtotal*Xtotal/16./pi/Xelastic/0.3894 ); // Slope parameter of elastic scattering
197                                                        //      (GeV/c)^(-2))
198//G4cout<<"G4FTFParameters Slope "<<GetSlope()<<G4endl;
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              SetProjMinDiffMass(1.1);                    // GeV
211              SetProjMinNonDiffMass(1.1);                 // GeV
212              SetProbabilityOfProjDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel
213
214              SetTarMinDiffMass(1.1);                     // GeV
215              SetTarMinNonDiffMass(1.1);                  // GeV
216              SetProbabilityOfTarDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel
217
218              SetAveragePt2(0.3);                         // GeV^2
219             }
220           else if( absPDGcode == 211 || PDGcode ==  111) //------Projectile is Pion -----------
221             {
222              SetProjMinDiffMass(0.5);                    // GeV
223              SetProjMinNonDiffMass(0.3);                 // GeV
224              SetProbabilityOfProjDiff(0.62*std::pow(s/GeV/GeV,-0.51)); // 40/32 X-dif/X-inel
225
226              SetTarMinDiffMass(1.1);                     // GeV
227              SetTarMinNonDiffMass(1.1);                  // GeV
228              SetProbabilityOfTarDiff(0.62*std::pow(s/GeV/GeV,-0.51)); // 40/32 X-dif/X-inel
229
230/*
231SetProjMinDiffMass(0.5);
232SetProjMinNonDiffMass(0.3);   // Uzhi 12.06.08
233SetProbabilityOfProjDiff(0.05);
234SetProbabilityOfTarDiff(0.05);
235*/
236              SetAveragePt2(0.3);                         // GeV^2
237             }
238           else if( absPDGcode == 321 || PDGcode == -311) //------Projectile is Kaon -----------
239             {
240              SetProjMinDiffMass(0.7);                    // GeV 1.1
241              SetProjMinNonDiffMass(0.7);                 // GeV
242              SetProbabilityOfProjDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel
243
244              SetTarMinDiffMass(1.1);                     // GeV
245              SetTarMinNonDiffMass(1.1);                  // GeV
246              SetProbabilityOfTarDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel
247
248              SetAveragePt2(0.3);                         // GeV^2
249             }
250           else                                           //------Projectile is undefined,
251                                                          //------Nucleon assumed
252             {
253              SetProjMinDiffMass((particle->GetPDGMass()+160.*MeV)/GeV);
254              SetProjMinNonDiffMass((particle->GetPDGMass()+160.*MeV)/GeV);
255              SetProbabilityOfProjDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel
256
257              SetTarMinDiffMass(1.1);                     // GeV
258              SetTarMinNonDiffMass(1.1);                  // GeV
259              SetProbabilityOfTarDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel
260
261              SetAveragePt2(0.3);                         // GeV^2
262             };
263
264
265//G4cout<<"G4FTFParameters Out"<<G4endl;
266
267    } 
268//**********************************************************************************************
Note: See TracBrowser for help on using the repository browser.