[1228] | 1 | // |
---|
| 2 | // ******************************************************************** |
---|
[968] | 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 | // |
---|
[1347] | 27 | // $Id: G4FTFParameters.cc,v 1.15 2010/11/15 10:02:38 vuzhinsk Exp $ |
---|
| 28 | // GEANT4 tag $Name: geant4-09-04-ref-00 $ |
---|
[968] | 29 | // |
---|
| 30 | |
---|
| 31 | #include "G4FTFParameters.hh" |
---|
| 32 | |
---|
[1196] | 33 | #include "G4ios.hh" |
---|
[1228] | 34 | #include <utility> |
---|
[1196] | 35 | |
---|
[968] | 36 | G4FTFParameters::G4FTFParameters() |
---|
| 37 | {;} |
---|
| 38 | |
---|
| 39 | |
---|
| 40 | G4FTFParameters::~G4FTFParameters() |
---|
| 41 | {;} |
---|
| 42 | //********************************************************************************************** |
---|
| 43 | |
---|
[1347] | 44 | //G4FTFParameters::G4FTFParameters(const G4ParticleDefinition * particle, |
---|
| 45 | // G4double theA, |
---|
| 46 | // G4double theZ, |
---|
| 47 | // G4double s) |
---|
[968] | 48 | G4FTFParameters::G4FTFParameters(const G4ParticleDefinition * particle, |
---|
[1347] | 49 | G4int theA, |
---|
| 50 | G4int theZ, |
---|
[968] | 51 | G4double s) |
---|
[1196] | 52 | { |
---|
[968] | 53 | G4int PDGcode = particle->GetPDGEncoding(); |
---|
| 54 | G4int absPDGcode = std::abs(PDGcode); |
---|
[1196] | 55 | G4double ProjectileMass = particle->GetPDGMass(); |
---|
| 56 | G4double TargetMass = G4Proton::Proton()->GetPDGMass(); |
---|
[968] | 57 | |
---|
[1196] | 58 | G4double Elab = (s - ProjectileMass*ProjectileMass - TargetMass*TargetMass)/ |
---|
| 59 | (2*TargetMass); |
---|
| 60 | G4double Plab = std::sqrt(Elab * Elab - ProjectileMass*ProjectileMass); |
---|
| 61 | |
---|
[1340] | 62 | G4double Ylab=0.5*std::log((Elab+Plab)/(Elab-Plab)); |
---|
| 63 | |
---|
| 64 | Plab/=GeV; // Uzhi 8.07.10 |
---|
[968] | 65 | G4double LogPlab = std::log( Plab ); |
---|
| 66 | G4double sqrLogPlab = LogPlab * LogPlab; |
---|
| 67 | |
---|
[1347] | 68 | G4int NumberOfTargetProtons = theZ; |
---|
| 69 | G4int NumberOfTargetNeutrons = theA-theZ; |
---|
| 70 | // G4int NumberOfTargetProtons = (G4int) theZ; |
---|
| 71 | // G4int NumberOfTargetNeutrons = (G4int) theA- (G4int) theZ; |
---|
[968] | 72 | G4int NumberOfTargetNucleons = NumberOfTargetProtons + NumberOfTargetNeutrons; |
---|
| 73 | |
---|
| 74 | G4double Xtotal, Xelastic; |
---|
| 75 | |
---|
[1347] | 76 | if( PDGcode > 1000 ) //------Projectile is baryon -------- |
---|
[968] | 77 | { |
---|
| 78 | G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) + 0.522*sqrLogPlab - 4.51*LogPlab; |
---|
| 79 | G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) + 0.513*sqrLogPlab - 4.27*LogPlab; |
---|
| 80 | |
---|
| 81 | G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab; |
---|
| 82 | G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab; |
---|
| 83 | |
---|
| 84 | Xtotal = ( NumberOfTargetProtons * XtotPP + |
---|
| 85 | NumberOfTargetNeutrons * XtotPN ) / NumberOfTargetNucleons; |
---|
| 86 | Xelastic = ( NumberOfTargetProtons * XelPP + |
---|
| 87 | NumberOfTargetNeutrons * XelPN ) / NumberOfTargetNucleons; |
---|
| 88 | } |
---|
[1347] | 89 | else if( PDGcode < -1000 ) //------Projectile is anti_baryon -------- |
---|
| 90 | { |
---|
| 91 | G4double XtotPP = 38.4 + 77.6*std::pow(Plab,-0.64) + 0.26*sqrLogPlab - 1.2*LogPlab; |
---|
| 92 | G4double XtotPN = 0. + 133.6*std::pow(Plab,-0.70) + 1.22*sqrLogPlab +13.7*LogPlab; |
---|
| 93 | |
---|
| 94 | G4double XelPP = 10.2 + 52.7*std::pow(Plab,-1.16) + 0.125*sqrLogPlab - 1.28*LogPlab; |
---|
| 95 | G4double XelPN = 36.5 + 0. *std::pow(Plab, 0. ) + 0. *sqrLogPlab -11.9 *LogPlab; |
---|
| 96 | |
---|
| 97 | Xtotal = ( NumberOfTargetProtons * XtotPP + |
---|
| 98 | NumberOfTargetNeutrons * XtotPN ) / NumberOfTargetNucleons; |
---|
| 99 | Xelastic = ( NumberOfTargetProtons * XelPP + |
---|
| 100 | NumberOfTargetNeutrons * XelPN ) / NumberOfTargetNucleons; |
---|
| 101 | } |
---|
[968] | 102 | else if( PDGcode == 211 ) //------Projectile is PionPlus ------- |
---|
| 103 | { |
---|
| 104 | G4double XtotPiP = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab; |
---|
| 105 | G4double XtotPiN = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab; |
---|
| 106 | |
---|
| 107 | G4double XelPiP = 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab; |
---|
| 108 | G4double XelPiN = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab; |
---|
| 109 | |
---|
| 110 | Xtotal = ( NumberOfTargetProtons * XtotPiP + |
---|
| 111 | NumberOfTargetNeutrons * XtotPiN ) / NumberOfTargetNucleons; |
---|
| 112 | Xelastic = ( NumberOfTargetProtons * XelPiP + |
---|
| 113 | NumberOfTargetNeutrons * XelPiN ) / NumberOfTargetNucleons; |
---|
| 114 | } |
---|
| 115 | else if( PDGcode == -211 ) //------Projectile is PionMinus ------- |
---|
| 116 | { |
---|
| 117 | G4double XtotPiP = 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab; |
---|
| 118 | G4double XtotPiN = 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab; |
---|
| 119 | |
---|
| 120 | G4double XelPiP = 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab; |
---|
| 121 | G4double XelPiN = 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab; |
---|
| 122 | |
---|
| 123 | Xtotal = ( NumberOfTargetProtons * XtotPiP + |
---|
| 124 | NumberOfTargetNeutrons * XtotPiN ) / NumberOfTargetNucleons; |
---|
| 125 | Xelastic = ( NumberOfTargetProtons * XelPiP + |
---|
| 126 | NumberOfTargetNeutrons * XelPiN ) / NumberOfTargetNucleons; |
---|
| 127 | } |
---|
| 128 | |
---|
| 129 | else if( PDGcode == 111 ) //------Projectile is PionZero ------- |
---|
| 130 | { |
---|
| 131 | G4double XtotPiP =(16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab + //Pi+ |
---|
| 132 | 33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab)/2; //Pi- |
---|
| 133 | |
---|
| 134 | G4double XtotPiN =(33.0 + 14.0 *std::pow(Plab,-1.36) + 0.456*sqrLogPlab - 4.03*LogPlab + //Pi+ |
---|
| 135 | 16.4 + 19.3 *std::pow(Plab,-0.42) + 0.19 *sqrLogPlab - 0.0 *LogPlab)/2; //Pi- |
---|
| 136 | |
---|
| 137 | G4double XelPiP =( 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab + //Pi+ |
---|
| 138 | 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab)/2; //Pi- |
---|
| 139 | G4double XelPiN =( 1.76 + 11.2*std::pow(Plab,-0.64) + 0.043*sqrLogPlab - 0.0 *LogPlab + //Pi+ |
---|
| 140 | 0.0 + 11.4*std::pow(Plab,-0.40) + 0.079*sqrLogPlab - 0.0 *LogPlab)/2; //Pi- |
---|
| 141 | |
---|
| 142 | Xtotal = ( NumberOfTargetProtons * XtotPiP + |
---|
| 143 | NumberOfTargetNeutrons * XtotPiN ) / NumberOfTargetNucleons; |
---|
| 144 | Xelastic = ( NumberOfTargetProtons * XelPiP + |
---|
| 145 | NumberOfTargetNeutrons * XelPiN ) / NumberOfTargetNucleons; |
---|
| 146 | } |
---|
| 147 | else if( PDGcode == 321 ) //------Projectile is KaonPlus ------- |
---|
| 148 | { |
---|
| 149 | G4double XtotKP = 18.1 + 0. *std::pow(Plab, 0. ) + 0.26 *sqrLogPlab - 1.0 *LogPlab; |
---|
| 150 | G4double XtotKN = 18.7 + 0. *std::pow(Plab, 0. ) + 0.21 *sqrLogPlab - 0.89*LogPlab; |
---|
| 151 | |
---|
| 152 | G4double XelKP = 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab; |
---|
| 153 | G4double XelKN = 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab; |
---|
| 154 | |
---|
| 155 | Xtotal = ( NumberOfTargetProtons * XtotKP + |
---|
| 156 | NumberOfTargetNeutrons * XtotKN ) / NumberOfTargetNucleons; |
---|
| 157 | Xelastic = ( NumberOfTargetProtons * XelKP + |
---|
| 158 | NumberOfTargetNeutrons * XelKN ) / NumberOfTargetNucleons; |
---|
| 159 | } |
---|
| 160 | else if( PDGcode ==-321 ) //------Projectile is KaonMinus ------ |
---|
| 161 | { |
---|
| 162 | G4double XtotKP = 32.1 + 0. *std::pow(Plab, 0. ) + 0.66 *sqrLogPlab - 5.6 *LogPlab; |
---|
| 163 | G4double XtotKN = 25.2 + 0. *std::pow(Plab, 0. ) + 0.38 *sqrLogPlab - 2.9 *LogPlab; |
---|
| 164 | |
---|
| 165 | G4double XelKP = 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab; |
---|
| 166 | G4double XelKN = 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab; |
---|
| 167 | |
---|
| 168 | Xtotal = ( NumberOfTargetProtons * XtotKP + |
---|
| 169 | NumberOfTargetNeutrons * XtotKN ) / NumberOfTargetNucleons; |
---|
| 170 | Xelastic = ( NumberOfTargetProtons * XelKP + |
---|
| 171 | NumberOfTargetNeutrons * XelKN ) / NumberOfTargetNucleons; |
---|
| 172 | } |
---|
[1196] | 173 | else if((PDGcode == 311) || (PDGcode == 130) || (PDGcode == 310))//Projectile is KaonZero |
---|
[968] | 174 | { |
---|
| 175 | G4double XtotKP =( 18.1 + 0. *std::pow(Plab, 0. ) + 0.26 *sqrLogPlab - 1.0 *LogPlab + //K+ |
---|
| 176 | 32.1 + 0. *std::pow(Plab, 0. ) + 0.66 *sqrLogPlab - 5.6 *LogPlab)/2; //K- |
---|
| 177 | G4double XtotKN =( 18.7 + 0. *std::pow(Plab, 0. ) + 0.21 *sqrLogPlab - 0.89*LogPlab + //K+ |
---|
| 178 | 25.2 + 0. *std::pow(Plab, 0. ) + 0.38 *sqrLogPlab - 2.9 *LogPlab)/2; //K- |
---|
| 179 | |
---|
| 180 | G4double XelKP =( 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab + //K+ |
---|
| 181 | 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab)/2; //K- |
---|
| 182 | G4double XelKN =( 7.3 + 0. *std::pow(Plab,-0. ) + 0.29 *sqrLogPlab - 2.4 *LogPlab + //K+ |
---|
| 183 | 5.0 + 8.1*std::pow(Plab,-1.8 ) + 0.16 *sqrLogPlab - 1.3 *LogPlab)/2; //K- |
---|
| 184 | Xtotal = ( NumberOfTargetProtons * XtotKP + |
---|
| 185 | NumberOfTargetNeutrons * XtotKN ) / NumberOfTargetNucleons; |
---|
| 186 | Xelastic = ( NumberOfTargetProtons * XelKP + |
---|
| 187 | NumberOfTargetNeutrons * XelKN ) / NumberOfTargetNucleons; |
---|
| 188 | } |
---|
| 189 | else //------Projectile is undefined, Nucleon assumed |
---|
| 190 | { |
---|
| 191 | G4double XtotPP = 48.0 + 0. *std::pow(Plab, 0. ) + 0.522*sqrLogPlab - 4.51*LogPlab; |
---|
| 192 | G4double XtotPN = 47.3 + 0. *std::pow(Plab, 0. ) + 0.513*sqrLogPlab - 4.27*LogPlab; |
---|
| 193 | |
---|
| 194 | G4double XelPP = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab; |
---|
| 195 | G4double XelPN = 11.9 + 26.9*std::pow(Plab,-1.21) + 0.169*sqrLogPlab - 1.85*LogPlab; |
---|
| 196 | |
---|
| 197 | Xtotal = ( NumberOfTargetProtons * XtotPP + |
---|
| 198 | NumberOfTargetNeutrons * XtotPN ) / NumberOfTargetNucleons; |
---|
| 199 | Xelastic = ( NumberOfTargetProtons * XelPP + |
---|
| 200 | NumberOfTargetNeutrons * XelPN ) / NumberOfTargetNucleons; |
---|
| 201 | }; |
---|
| 202 | |
---|
| 203 | // Xtotal and Xelastic in mb |
---|
| 204 | |
---|
[1347] | 205 | // For Pi- P interactions only! |
---|
| 206 | if(std::abs(Plab-1.4) < 0.05) {Xtotal=3.500599e+01; Xelastic= 1.150032e+01;} |
---|
| 207 | if(std::abs(Plab-1.5) < 0.05) {Xtotal=3.450591e+01; Xelastic= 1.050038e+01;} |
---|
| 208 | if(std::abs(Plab-1.6) < 0.05) {Xtotal=3.430576e+01; Xelastic= 9.800433e+00;} |
---|
| 209 | if(std::abs(Plab-1.7) < 0.05) {Xtotal=3.455560e+01; Xelastic= 9.300436e+00;} |
---|
| 210 | if(std::abs(Plab-1.8) < 0.05) {Xtotal=3.480545e+01; Xelastic= 8.800438e+00;} |
---|
| 211 | if(std::abs(Plab-2.0) < 0.05) {Xtotal=3.570503e+01; Xelastic= 8.200370e+00;} |
---|
| 212 | if(std::abs(Plab-2.2) < 0.05) {Xtotal=3.530495e+01; Xelastic= 7.800362e+00;} |
---|
| 213 | if(std::abs(Plab-2.5) < 0.05) {Xtotal=3.410484e+01; Xelastic= 7.350320e+00;} |
---|
| 214 | if(std::abs(Plab-2.75) < 0.05){Xtotal=3.280479e+01; Xelastic= 7.050273e+00;} |
---|
| 215 | if(std::abs(Plab-3.0) < 0.05) {Xtotal=3.180473e+01; Xelastic= 6.800258e+00;} |
---|
| 216 | if(std::abs(Plab-4.0) < 0.05) {Xtotal=2.910441e+01; Xelastic= 6.100229e+00;} |
---|
| 217 | if(std::abs(Plab-5.0) < 0.05) {Xtotal=2.820372e+01; Xelastic= 5.700275e+00;} |
---|
| 218 | if(std::abs(Plab-6.0) < 0.05) {Xtotal=2.760367e+01; Xelastic= 5.400255e+00;} |
---|
| 219 | if(std::abs(Plab-7.0) < 0.05) {Xtotal=2.725366e+01; Xelastic= 5.150256e+00;} |
---|
| 220 | if(std::abs(Plab-8.0) < 0.05) {Xtotal=2.690365e+01; Xelastic= 4.900258e+00;} |
---|
| 221 | if(std::abs(Plab-10.0) < 0.05){Xtotal=2.660342e+01; Xelastic= 4.600237e+00;} |
---|
| 222 | if(std::abs(Plab-12.0) < 0.05){Xtotal=2.632341e+01; Xelastic= 4.480229e+00;} |
---|
| 223 | if(std::abs(Plab-14.0) < 0.05){Xtotal=2.604340e+01; Xelastic= 4.360221e+00;} |
---|
| 224 | if(std::abs(Plab-20.0) < 0.05){Xtotal=2.520337e+01; Xelastic= 4.000197e+00;} |
---|
| 225 | if(std::abs(Plab-30.0) < 0.05){Xtotal=2.505334e+01; Xelastic= 3.912679e+00;} |
---|
| 226 | // |
---|
[968] | 227 | //----------- Geometrical parameters ------------------------------------------------ |
---|
| 228 | SetTotalCrossSection(Xtotal); |
---|
| 229 | SetElastisCrossSection(Xelastic); |
---|
| 230 | SetInelasticCrossSection(Xtotal-Xelastic); |
---|
| 231 | |
---|
[1347] | 232 | //G4cout<<"Plab Xtotal, Xelastic Xinel "<<Plab<<" "<<Xtotal<<" "<<Xelastic<<Xtotal-Xelastic)<<G4endl; |
---|
[1196] | 233 | // // Interactions with elastic and inelastic collisions |
---|
[968] | 234 | SetProbabilityOfElasticScatt(Xtotal, Xelastic); |
---|
| 235 | SetRadiusOfHNinteractions2(Xtotal/pi/10.); |
---|
| 236 | // |
---|
| 237 | /* //==== No elastic scattering ============================ |
---|
| 238 | SetProbabilityOfElasticScatt(Xtotal, 0.); |
---|
| 239 | SetRadiusOfHNinteractions2((Xtotal-Xelastic)/pi/10.); |
---|
| 240 | */ //======================================================= |
---|
| 241 | |
---|
| 242 | //----------------------------------------------------------------------------------- |
---|
[1340] | 243 | |
---|
[968] | 244 | SetSlope( Xtotal*Xtotal/16./pi/Xelastic/0.3894 ); // Slope parameter of elastic scattering |
---|
| 245 | // (GeV/c)^(-2)) |
---|
| 246 | //----------------------------------------------------------------------------------- |
---|
| 247 | SetGamma0( GetSlope()*Xtotal/10./2./pi ); |
---|
| 248 | |
---|
| 249 | //----------- Parameters of elastic scattering -------------------------------------- |
---|
| 250 | // Gaussian parametrization of |
---|
| 251 | // elastic scattering amplitude assumed |
---|
| 252 | SetAvaragePt2ofElasticScattering(1./(Xtotal*Xtotal/16./pi/Xelastic/0.3894)*GeV*GeV); |
---|
| 253 | |
---|
| 254 | //----------- Parameters of excitations --------------------------------------------- |
---|
[1347] | 255 | if( PDGcode > 1000 ) //------Projectile is baryon -------- |
---|
[968] | 256 | { |
---|
[1340] | 257 | SetMagQuarkExchange(1.84);//(3.63); |
---|
| 258 | SetSlopeQuarkExchange(0.7);//(1.2); |
---|
| 259 | SetDeltaProbAtQuarkExchange(0.); |
---|
[1347] | 260 | if(NumberOfTargetNucleons > 26) {SetProbOfSameQuarkExchange(1.);} |
---|
| 261 | else {SetProbOfSameQuarkExchange(0.);} |
---|
[1196] | 262 | |
---|
[1340] | 263 | SetProjMinDiffMass(1.16); // GeV |
---|
[1347] | 264 | SetProjMinNonDiffMass(1.16); // GeV |
---|
| 265 | SetProbabilityOfProjDiff(0.805*std::exp(-0.35*Ylab));// 0.5 |
---|
| 266 | |
---|
| 267 | SetTarMinDiffMass(1.16); // GeV |
---|
| 268 | SetTarMinNonDiffMass(1.16); // GeV |
---|
| 269 | SetProbabilityOfTarDiff(0.805*std::exp(-0.35*Ylab));// 0.5 |
---|
| 270 | |
---|
| 271 | SetAveragePt2(0.15); // 0.15 GeV^2 |
---|
| 272 | } |
---|
| 273 | if( PDGcode < -1000 ) //------Projectile is anti_baryon -------- |
---|
| 274 | { |
---|
| 275 | SetMagQuarkExchange(0.); |
---|
| 276 | SetSlopeQuarkExchange(0.); |
---|
| 277 | SetDeltaProbAtQuarkExchange(0.); |
---|
| 278 | SetProbOfSameQuarkExchange(0.); |
---|
| 279 | |
---|
| 280 | SetProjMinDiffMass(1.16); // GeV |
---|
[1340] | 281 | SetProjMinNonDiffMass(1.16); // GeV |
---|
[1347] | 282 | SetProbabilityOfProjDiff(0.805*std::exp(-0.35*Ylab));// 0.5 |
---|
[968] | 283 | |
---|
[1340] | 284 | SetTarMinDiffMass(1.16); // GeV |
---|
| 285 | SetTarMinNonDiffMass(1.16); // GeV |
---|
[1347] | 286 | SetProbabilityOfTarDiff(0.805*std::exp(-0.35*Ylab));// 0.5 |
---|
[968] | 287 | |
---|
[1340] | 288 | SetAveragePt2(0.15); // 0.15 GeV^2 |
---|
[968] | 289 | } |
---|
| 290 | else if( absPDGcode == 211 || PDGcode == 111) //------Projectile is Pion ----------- |
---|
| 291 | { |
---|
[1347] | 292 | SetMagQuarkExchange(240.); |
---|
| 293 | SetSlopeQuarkExchange(2.); |
---|
| 294 | SetDeltaProbAtQuarkExchange(0.56); //(0.35); |
---|
[1196] | 295 | |
---|
[968] | 296 | SetProjMinDiffMass(0.5); // GeV |
---|
[1347] | 297 | SetProjMinNonDiffMass(0.5); // GeV 0.3 |
---|
| 298 | SetProbabilityOfProjDiff(0.);//(0.*0.62*std::pow(s/GeV/GeV,-0.51)); // 40/32 X-dif/X-inel |
---|
[1228] | 299 | |
---|
[1347] | 300 | SetTarMinDiffMass(1.16); // GeV |
---|
| 301 | SetTarMinNonDiffMass(1.16); // GeV |
---|
| 302 | // SetProbabilityOfTarDiff(1.);//(2.*0.62*std::pow(s/GeV/GeV,-0.51)); |
---|
| 303 | // SetProbabilityOfTarDiff(2.6*std::exp(-0.46*Ylab)); |
---|
| 304 | SetProbabilityOfTarDiff(0.8*std::exp(-0.6*(Ylab-3.))); |
---|
[1228] | 305 | |
---|
[968] | 306 | SetAveragePt2(0.3); // GeV^2 |
---|
| 307 | } |
---|
[1196] | 308 | else if( (absPDGcode == 321) || (PDGcode == 311) || |
---|
| 309 | (PDGcode == 130) || (PDGcode == 310)) //Projectile is Kaon |
---|
[968] | 310 | { |
---|
[1196] | 311 | // Must be corrected, taken from PiN |
---|
| 312 | SetMagQuarkExchange(120.); |
---|
| 313 | SetSlopeQuarkExchange(2.0); |
---|
| 314 | SetDeltaProbAtQuarkExchange(0.6); |
---|
| 315 | |
---|
[968] | 316 | SetProjMinDiffMass(0.7); // GeV 1.1 |
---|
| 317 | SetProjMinNonDiffMass(0.7); // GeV |
---|
| 318 | SetProbabilityOfProjDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel |
---|
| 319 | |
---|
| 320 | SetTarMinDiffMass(1.1); // GeV |
---|
| 321 | SetTarMinNonDiffMass(1.1); // GeV |
---|
| 322 | SetProbabilityOfTarDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel |
---|
| 323 | |
---|
| 324 | SetAveragePt2(0.3); // GeV^2 |
---|
| 325 | } |
---|
| 326 | else //------Projectile is undefined, |
---|
| 327 | //------Nucleon assumed |
---|
| 328 | { |
---|
[1196] | 329 | SetMagQuarkExchange(3.5); |
---|
| 330 | SetSlopeQuarkExchange(1.0); |
---|
| 331 | SetDeltaProbAtQuarkExchange(0.1); |
---|
| 332 | |
---|
[968] | 333 | SetProjMinDiffMass((particle->GetPDGMass()+160.*MeV)/GeV); |
---|
| 334 | SetProjMinNonDiffMass((particle->GetPDGMass()+160.*MeV)/GeV); |
---|
| 335 | SetProbabilityOfProjDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel |
---|
| 336 | |
---|
| 337 | SetTarMinDiffMass(1.1); // GeV |
---|
| 338 | SetTarMinNonDiffMass(1.1); // GeV |
---|
| 339 | SetProbabilityOfTarDiff(0.95*std::pow(s/GeV/GeV,-0.35)); // 40/32 X-dif/X-inel |
---|
| 340 | |
---|
| 341 | SetAveragePt2(0.3); // GeV^2 |
---|
[1196] | 342 | } |
---|
[968] | 343 | |
---|
[1196] | 344 | // ---------- Set parameters of a string kink ------------------------------- |
---|
| 345 | SetPt2Kink(6.*GeV*GeV); |
---|
| 346 | G4double Puubar(1./3.), Pddbar(1./3.), Pssbar(1./3.); // SU(3) symmetry |
---|
| 347 | // G4double Puubar(0.41 ), Pddbar(0.41 ), Pssbar(0.18 ); // Broken SU(3) symmetry |
---|
| 348 | SetQuarkProbabilitiesAtGluonSplitUp(Puubar, Pddbar, Pssbar); |
---|
[968] | 349 | |
---|
[1196] | 350 | // --------- Set parameters of nuclear destruction-------------------- |
---|
[968] | 351 | |
---|
[1196] | 352 | if( absPDGcode < 1000 ) |
---|
| 353 | { |
---|
[1340] | 354 | SetMaxNumberOfCollisions(1000.,1.); //(Plab,2.); //3.); ############################## |
---|
| 355 | |
---|
[1347] | 356 | // SetCofNuclearDestruction(0.); //1.0); // for meson projectile |
---|
| 357 | // SetCofNuclearDestruction(1.*std::exp(4.*(Ylab-2.1))/(1.+std::exp(4.*(Ylab-2.1)))); |
---|
| 358 | //G4cout<<Ylab<<" "<<0.62*std::exp(4.*(Ylab-4.5))/(1.+std::exp(4.*(Ylab-4.5)))<<G4endl; |
---|
| 359 | //G4int Uzhi; G4cin>>Uzhi; |
---|
[1340] | 360 | |
---|
[1347] | 361 | // SetMaxNumberOfCollisions(Plab,2.); //4.); // ############################## |
---|
| 362 | |
---|
| 363 | SetCofNuclearDestruction(1.*std::exp(4.*(Ylab-2.1))/(1.+std::exp(4.*(Ylab-2.1)))); //0.62 1.0 |
---|
| 364 | //------------------------------------------ |
---|
| 365 | // SetDofNuclearDestruction(0.4); |
---|
| 366 | // SetPt2ofNuclearDestruction(0.17*GeV*GeV); |
---|
| 367 | // SetMaxPt2ofNuclearDestruction(1.0*GeV*GeV); |
---|
| 368 | |
---|
| 369 | // SetExcitationEnergyPerWoundedNucleon(100*MeV); |
---|
[1340] | 370 | SetDofNuclearDestruction(0.4); |
---|
[1347] | 371 | SetPt2ofNuclearDestruction((0.035+ |
---|
| 372 | 0.04*std::exp(4.*(Ylab-2.5))/(1.+std::exp(4.*(Ylab-2.5))))*GeV*GeV); //0.09 |
---|
[1340] | 373 | SetMaxPt2ofNuclearDestruction(1.0*GeV*GeV); |
---|
| 374 | |
---|
[1347] | 375 | SetExcitationEnergyPerWoundedNucleon(75.*MeV); |
---|
[1340] | 376 | } else // for baryon projectile |
---|
[1196] | 377 | { |
---|
[1347] | 378 | SetMaxNumberOfCollisions(Plab,2.); //4.); // ############################## |
---|
[1340] | 379 | |
---|
[1347] | 380 | SetCofNuclearDestruction(1.*std::exp(4.*(Ylab-2.1))/(1.+std::exp(4.*(Ylab-2.1)))); //0.62 1.0 |
---|
| 381 | //G4cout<<Ylab<<" "<<0.62*std::exp(4.*(Ylab-2.1))/(1.+std::exp(4.*(Ylab-2.1)))<<G4endl; |
---|
| 382 | //G4int Uzhi; G4cin>>Uzhi; |
---|
[1340] | 383 | |
---|
| 384 | SetDofNuclearDestruction(0.4); |
---|
| 385 | SetPt2ofNuclearDestruction((0.035+ |
---|
| 386 | 0.04*std::exp(4.*(Ylab-2.5))/(1.+std::exp(4.*(Ylab-2.5))))*GeV*GeV); //0.09 |
---|
| 387 | SetMaxPt2ofNuclearDestruction(1.0*GeV*GeV); |
---|
| 388 | |
---|
| 389 | SetExcitationEnergyPerWoundedNucleon(75.*MeV); |
---|
[1196] | 390 | } |
---|
| 391 | |
---|
| 392 | SetR2ofNuclearDestruction(1.5*fermi*fermi); |
---|
| 393 | |
---|
[1340] | 394 | //SetCofNuclearDestruction(0.47*std::exp(2.*(Ylab-2.5))/(1.+std::exp(2.*(Ylab-2.5)))); |
---|
| 395 | //SetPt2ofNuclearDestruction((0.035+0.1*std::exp(4.*(Ylab-3.))/(1.+std::exp(4.*(Ylab-3.))))*GeV*GeV); |
---|
[1196] | 396 | |
---|
[1347] | 397 | //SetMagQuarkExchange(120.); // 210. PipP |
---|
| 398 | //SetSlopeQuarkExchange(2.0); |
---|
| 399 | //SetDeltaProbAtQuarkExchange(0.6); |
---|
| 400 | //SetProjMinDiffMass(0.7); // GeV 1.1 |
---|
| 401 | //SetProjMinNonDiffMass(0.7); // GeV |
---|
| 402 | //SetProbabilityOfProjDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel |
---|
| 403 | //SetTarMinDiffMass(1.1); // GeV |
---|
| 404 | //SetTarMinNonDiffMass(1.1); // GeV |
---|
| 405 | //SetProbabilityOfTarDiff(0.85*std::pow(s/GeV/GeV,-0.5)); // 40/32 X-dif/X-inel |
---|
| 406 | // |
---|
| 407 | //SetAveragePt2(0.3); // GeV^2 |
---|
| 408 | //------------------------------------ |
---|
[1340] | 409 | //SetProbabilityOfElasticScatt(1.,1.); //(Xtotal, Xelastic); |
---|
| 410 | //SetProbabilityOfProjDiff(1.*0.62*std::pow(s/GeV/GeV,-0.51)); // 0->1 |
---|
| 411 | //SetProbabilityOfTarDiff(4.*0.62*std::pow(s/GeV/GeV,-0.51)); // 2->4 |
---|
| 412 | //SetAveragePt2(0.3); //(0.15); |
---|
| 413 | //SetAvaragePt2ofElasticScattering(0.); |
---|
| 414 | |
---|
[1347] | 415 | //SetMaxNumberOfCollisions(4.*(Plab+0.01),Plab); //6.); // ############################## |
---|
| 416 | //SetCofNuclearDestruction(0.2); //(0.4); |
---|
| 417 | //SetExcitationEnergyPerWoundedNucleon(0.*MeV); //(75.*MeV); |
---|
| 418 | //SetDofNuclearDestruction(0.4); //(0.4); |
---|
| 419 | //SetPt2ofNuclearDestruction(0.1*GeV*GeV); //(0.168*GeV*GeV); |
---|
[1340] | 420 | |
---|
[1196] | 421 | } |
---|
[968] | 422 | //********************************************************************************************** |
---|