Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FTFParameters.cc
Go to the documentation of this file.
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//
28
29#include <utility>
30
31#include "G4FTFParameters.hh"
32
33#include "G4ios.hh"
35#include "G4SystemOfUnits.hh"
36
38#include "G4Proton.hh"
39#include "G4Neutron.hh"
40#include "G4PionPlus.hh"
41#include "G4PionMinus.hh"
42#include "G4KaonPlus.hh"
43#include "G4KaonMinus.hh"
44
49
50#include "G4Exp.hh"
51#include "G4Log.hh"
52#include "G4Pow.hh"
53
55
56//============================================================================
57
58//#define debugFTFparams
59
60//============================================================================
61
63{
64 StringMass = new G4LundStringFragmentation; // for estimation of min. mass of diffr. states
65 Reset();
66 csGGinstance =
68 if (!csGGinstance) {
69 csGGinstance = new G4ComponentGGHadronNucleusXsc();
70 }
71
72 // Set parameters of a string kink
73 SetPt2Kink( 0.0*GeV*GeV ); // To switch off kinky strings (bad results obtained with 6.0*GeV*GeV)
74 G4double Puubar( 1.0/3.0 ), Pddbar( 1.0/3.0 ), Pssbar( 1.0/3.0 ); // SU(3) symmetry
75 //G4double Puubar( 0.41 ), Pddbar( 0.41 ), Pssbar( 0.18 ); // Broken SU(3) symmetry
76 SetQuarkProbabilitiesAtGluonSplitUp( Puubar, Pddbar, Pssbar );
77}
78
79//============================================================================
80
82 G4int theA, G4int theZ, G4double PlabPerParticle )
83{
84 Reset();
85
86 G4int ProjectilePDGcode = particle->GetPDGEncoding();
87 G4int ProjectileabsPDGcode = std::abs( ProjectilePDGcode );
88 G4double ProjectileMass = particle->GetPDGMass();
89 G4double ProjectileMass2 = ProjectileMass * ProjectileMass;
90
91 G4int ProjectileBaryonNumber( 0 ), AbsProjectileBaryonNumber( 0 ), AbsProjectileCharge( 0 );
92 G4bool ProjectileIsNucleus = false;
93
94 if ( std::abs( particle->GetBaryonNumber() ) > 1 ) { // The projectile is a nucleus
95 ProjectileIsNucleus = true;
96 ProjectileBaryonNumber = particle->GetBaryonNumber();
97 AbsProjectileBaryonNumber = std::abs( ProjectileBaryonNumber );
98 AbsProjectileCharge = std::abs( G4int( particle->GetPDGCharge() ) );
99 if ( ProjectileBaryonNumber > 1 ) {
100 ProjectilePDGcode = 2212; ProjectileabsPDGcode = 2212; // Proton
101 } else {
102 ProjectilePDGcode = -2212; ProjectileabsPDGcode = 2212; // Anti-Proton
103 }
104 ProjectileMass = G4Proton::Proton()->GetPDGMass();
105 ProjectileMass2 = sqr( ProjectileMass );
106 }
107
108 G4double TargetMass = G4Proton::Proton()->GetPDGMass();
109 G4double TargetMass2 = TargetMass * TargetMass;
110
111 G4double Plab = PlabPerParticle;
112 G4double Elab = std::sqrt( Plab*Plab + ProjectileMass2 );
113 G4double KineticEnergy = Elab - ProjectileMass;
114
115 G4double S = ProjectileMass2 + TargetMass2 + 2.0*TargetMass*Elab;
116
117 #ifdef debugFTFparams
118 G4cout << "--------- FTF Parameters --------------" << G4endl << "Proj Plab "
119 << ProjectilePDGcode << " " << Plab << G4endl << "Mass KinE " << ProjectileMass
120 << " " << KineticEnergy << G4endl << " A Z " << theA << " " << theZ << G4endl;
121 #endif
122
123 G4double Ylab, Xtotal( 0.0 ), Xelastic( 0.0 ), Xannihilation( 0.0 );
124 G4int NumberOfTargetNucleons;
125
126 Ylab = 0.5 * G4Log( (Elab + Plab)/(Elab - Plab) );
127
128 G4double ECMSsqr = S/GeV/GeV;
129 G4double SqrtS = std::sqrt( S )/GeV;
130
131 #ifdef debugFTFparams
132 G4cout << "Sqrt(s) " << SqrtS << G4endl;
133 #endif
134
135 TargetMass /= GeV; TargetMass2 /= (GeV*GeV);
136 ProjectileMass /= GeV; ProjectileMass2 /= (GeV*GeV);
137
138 Plab /= GeV;
139 G4double Xftf = 0.0;
140
141 G4int NumberOfTargetProtons = theZ;
142 G4int NumberOfTargetNeutrons = theA - theZ;
143 NumberOfTargetNucleons = NumberOfTargetProtons + NumberOfTargetNeutrons;
144
145 // ---------- hadron projectile ----------------
146 if ( AbsProjectileBaryonNumber <= 1 ) { // Projectile is hadron or baryon
147
148 // Interaction on P
149 G4double xTtP = csGGinstance->GetTotalIsotopeCrossSection( particle, KineticEnergy, 1, 1);
150 G4double xElP = csGGinstance->GetElasticIsotopeCrossSection(particle, KineticEnergy, 1, 1);
151
152 // Interaction on N
153 G4double xTtN = csGGinstance->GetTotalIsotopeCrossSection( particle, KineticEnergy, 0, 1);
154 G4double xElN = csGGinstance->GetElasticIsotopeCrossSection(particle, KineticEnergy, 0, 1);
155
156 // Average properties of h+N interactions
157 Xtotal = ( NumberOfTargetProtons * xTtP + NumberOfTargetNeutrons * xTtN ) / NumberOfTargetNucleons;
158 Xelastic = ( NumberOfTargetProtons * xElP + NumberOfTargetNeutrons * xElN ) / NumberOfTargetNucleons;
159 Xannihilation = 0.0;
160
161 Xtotal /= millibarn;
162 Xelastic /= millibarn;
163
164 #ifdef debugFTFparams
165 G4cout<<"Estimated cross sections (total and elastic) of h+N interactions "<<Xtotal<<" "<<Xelastic<<" (mb)"<<G4endl;
166 #endif
167 }
168
169 // ---------- nucleus projectile ----------------
170 if ( ProjectileIsNucleus && ProjectileBaryonNumber > 1 ) {
171
172 #ifdef debugFTFparams
173 G4cout<<"Projectile is a nucleus: A and Z - "<<ProjectileBaryonNumber<<" "<<ProjectileCharge<<G4endl;
174 #endif
175
177 // Interaction on P
178 G4double XtotPP = csGGinstance->GetTotalIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
179 G4double XelPP = csGGinstance->GetElasticIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
180
182 // Interaction on N
183 G4double XtotPN = csGGinstance->GetTotalIsotopeCrossSection(Neutron, KineticEnergy, 0, 1);
184 G4double XelPN = csGGinstance->GetElasticIsotopeCrossSection(Neutron, KineticEnergy, 0, 1);
185
186 #ifdef debugFTFparams
187 G4cout << "XsPP (total and elastic) " << XtotPP/millibarn << " " << XelPP/millibarn <<" (mb)"<< G4endl
188 << "XsPN (total and elastic) " << XtotPN/millibarn << " " << XelPN/millibarn <<" (mb)"<< G4endl;
189 #endif
190
191 Xtotal = (
192 AbsProjectileCharge * NumberOfTargetProtons * XtotPP +
193 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
194 NumberOfTargetNeutrons * XtotPP
195 +
196 ( AbsProjectileCharge * NumberOfTargetNeutrons +
197 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
198 NumberOfTargetProtons ) * XtotPN
199 ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
200 Xelastic= (
201 AbsProjectileCharge * NumberOfTargetProtons * XelPP +
202 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
203 NumberOfTargetNeutrons * XelPP
204 +
205 ( AbsProjectileCharge * NumberOfTargetNeutrons +
206 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
207 NumberOfTargetProtons ) * XelPN
208 ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
209
210 Xannihilation = 0.0;
211 Xtotal /= millibarn;
212 Xelastic /= millibarn;
213 }
214
215 // ---------- The projectile is anti-baryon or anti-nucleus ----------------
216 // anti Sigma^0_c anti Delta^-
217 if ( ProjectilePDGcode >= -4112 && ProjectilePDGcode <= -1114 ) {
218 // Only non-strange and strange baryons are considered
219
220 #ifdef debugFTFparams
221 G4cout<<"Projectile is a anti-baryon or anti-nucleus - "<<ProjectileBaryonNumber<<" "<<ProjectileCharge<<G4endl;
222 G4cout<<"(Only non-strange and strange baryons are considered)"<<G4endl;
223 #endif
224
225 G4double X_a( 0.0 ), X_b( 0.0 ), X_c( 0.0 ), X_d( 0.0 );
226 G4double MesonProdThreshold = ProjectileMass + TargetMass +
227 ( 2.0 * 0.14 + 0.016 ); // 2 Mpi + DeltaE;
228
229 if ( PlabPerParticle < 40.0*MeV ) { // Low energy limits. Projectile at rest.
230 Xtotal = 1512.9; // mb
231 Xelastic = 473.2; // mb
232 X_a = 625.1; // mb
233 X_b = 9.780; // mb
234 X_c = 49.989; // mb
235 X_d = 6.614; // mb
236 } else { // Total and elastic cross section of PbarP interactions a'la Arkhipov
237 G4double LogS = G4Log( ECMSsqr / 33.0625 );
238 G4double Xasmpt = 36.04 + 0.304*LogS*LogS; // mb
239 LogS = G4Log( SqrtS / 20.74 );
240 G4double Basmpt = 11.92 + 0.3036*LogS*LogS; // GeV^(-2)
241 G4double R0 = std::sqrt( 0.40874044*Xasmpt - Basmpt ); // GeV^(-1)
242
243 G4double FlowF = SqrtS / std::sqrt( ECMSsqr*ECMSsqr + ProjectileMass2*ProjectileMass2 +
244 TargetMass2*TargetMass2 - 2.0*ECMSsqr*ProjectileMass2
245 - 2.0*ECMSsqr*TargetMass2
246 - 2.0*ProjectileMass2*TargetMass2 );
247
248 Xtotal = Xasmpt * ( 1.0 + 13.55*FlowF/R0/R0/R0*
249 (1.0 - 4.47/SqrtS + 12.38/ECMSsqr - 12.43/SqrtS/ECMSsqr) ); // mb
250
251 Xasmpt = 4.4 + 0.101*LogS*LogS; // mb
252 Xelastic = Xasmpt * ( 1.0 + 59.27*FlowF/R0/R0/R0*
253 (1.0 - 6.95/SqrtS + 23.54/ECMSsqr - 25.34/SqrtS/ECMSsqr ) ); // mb
254
255 X_a = 25.0*FlowF; // mb, 3-shirts diagram
256
257 if ( SqrtS < MesonProdThreshold ) {
258 X_b = 3.13 + 140.0*G4Pow::GetInstance()->powA( MesonProdThreshold - SqrtS, 2.5 ); // mb anti-quark-quark annihilation
259 Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
260 } else {
261 X_b = 6.8/SqrtS; // mb anti-quark-quark annihilation
262 Xelastic -= 3.0*X_b; // Xel-X(PbarP->NNbar)
263 }
264
265 X_c = 2.0*FlowF*sqr( ProjectileMass + TargetMass )/ECMSsqr; // mb rearrangement
266
267 X_d = 23.3/ECMSsqr; // mb anti-quark-quark string creation
268 }
269
270 G4double Xann_on_P( 0.0), Xann_on_N( 0.0 );
271
272 if ( ProjectilePDGcode == -2212 ) { // Pbar+P/N
273 Xann_on_P = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
274 Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
275 } else if ( ProjectilePDGcode == -2112 ) { // NeutrBar+P/N
276 Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*4.0;
277 Xann_on_N = X_a + X_b*5.0 + X_c*5.0 + X_d*6.0;
278 } else if ( ProjectilePDGcode == -3122 ) { // LambdaBar+P/N
279 Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
280 Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
281 } else if ( ProjectilePDGcode == -3112 ) { // Sigma-Bar+P/N
282 Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
283 Xann_on_N = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
284 } else if ( ProjectilePDGcode == -3212 ) { // Sigma0Bar+P/N
285 Xann_on_P = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
286 Xann_on_N = X_a + X_b*3.0 + X_c*3.0 + X_d*2.0;
287 } else if ( ProjectilePDGcode == -3222 ) { // Sigma+Bar+P/N
288 Xann_on_P = X_a + X_b*4.0 + X_c*4.0 + X_d*2.0;
289 Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
290 } else if ( ProjectilePDGcode == -3312 ) { // Xi-Bar+P/N
291 Xann_on_P = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
292 Xann_on_N = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
293 } else if ( ProjectilePDGcode == -3322 ) { // Xi0Bar+P/N
294 Xann_on_P = X_a + X_b*2.0 + X_c*2.0 + X_d*0.0;
295 Xann_on_N = X_a + X_b*1.0 + X_c*1.0 + X_d*0.0;
296 } else if ( ProjectilePDGcode == -3334 ) { // Omega-Bar+P/N
297 Xann_on_P = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
298 Xann_on_N = X_a + X_b*0.0 + X_c*0.0 + X_d*0.0;
299 } else {
300 G4cout << "Unknown anti-baryon for FTF annihilation" << G4endl;
301 }
302
303 //G4cout << "Sum " << Xann_on_P << G4endl;
304
305 if ( ! ProjectileIsNucleus ) { // Projectile is anti-baryon
306 Xannihilation = ( NumberOfTargetProtons * Xann_on_P + NumberOfTargetNeutrons * Xann_on_N )
307 / NumberOfTargetNucleons;
308 } else { // Projectile is a nucleus
309 Xannihilation = (
310 ( AbsProjectileCharge * NumberOfTargetProtons +
311 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
312 NumberOfTargetNeutrons ) * Xann_on_P
313 +
314 ( AbsProjectileCharge * NumberOfTargetNeutrons +
315 ( AbsProjectileBaryonNumber - AbsProjectileCharge ) *
316 NumberOfTargetProtons ) * Xann_on_N
317 ) / ( AbsProjectileBaryonNumber * NumberOfTargetNucleons );
318 }
319
320 //G4double Xftf = 0.0;
321 MesonProdThreshold = ProjectileMass + TargetMass + (0.14 + 0.08); // Mpi + DeltaE
322 if ( SqrtS > MesonProdThreshold ) {
323 Xftf = 36.0 * ( 1.0 - MesonProdThreshold/SqrtS );
324 }
325
326 Xtotal = Xelastic + Xannihilation + Xftf;
327
328 #ifdef debugFTFparams
329 G4cout << "Plab Xtotal, Xelastic Xinel Xftf " << Plab << " " << Xtotal << " " << Xelastic
330 << " " << Xtotal - Xelastic << " " << Xtotal - Xelastic - Xannihilation << " (mb)"<< G4endl
331 << "Plab Xelastic/Xtotal, Xann/Xin " << Plab << " " << Xelastic/Xtotal << " "
332 << Xannihilation/(Xtotal - Xelastic) << G4endl;
333 #endif
334
335 }
336
337 if ( Xtotal == 0.0 ) { // Projectile is undefined, Nucleon assumed
338
340 // Interaction on P
341 G4double XtotPP = csGGinstance->GetTotalIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
342 G4double XelPP = csGGinstance->GetElasticIsotopeCrossSection(Proton, KineticEnergy, 1, 1);
343
344 // Interaction on N
345 G4double XtotPN = csGGinstance->GetTotalIsotopeCrossSection(Proton, KineticEnergy, 0, 1);
346 G4double XelPN = csGGinstance->GetElasticIsotopeCrossSection(Proton, KineticEnergy, 0, 1);
347
348 Xtotal = ( NumberOfTargetProtons * XtotPP + NumberOfTargetNeutrons * XtotPN )
349 / NumberOfTargetNucleons;
350 Xelastic = ( NumberOfTargetProtons * XelPP + NumberOfTargetNeutrons * XelPN )
351 / NumberOfTargetNucleons;
352 Xannihilation = 0.0;
353 Xtotal /= millibarn;
354 Xelastic /= millibarn;
355 };
356
357 // Geometrical parameters
358 SetTotalCrossSection( Xtotal );
359 SetElastisCrossSection( Xelastic );
360 SetInelasticCrossSection( Xtotal - Xelastic );
361
362 // Interactions with elastic and inelastic collisions
363 SetProbabilityOfElasticScatt( Xtotal, Xelastic );
364
365 SetRadiusOfHNinteractions2( Xtotal/pi/10.0 );
366
367 if ( ( Xtotal - Xelastic ) == 0.0 ) {
369 } else {
370 SetProbabilityOfAnnihilation( Xannihilation / (Xtotal - Xelastic) );
371 }
372
373 if(Xelastic > 0.0) {
374 SetSlope( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 );// Slope parameter of elastic scattering
375 // (GeV/c)^(-2))
376 // Parameters of elastic scattering
377 // Gaussian parametrization of elastic scattering amplitude assumed
378 SetAvaragePt2ofElasticScattering( 1.0/( Xtotal*Xtotal/16.0/pi/Xelastic/0.3894 )*GeV*GeV );
379 } else {
380 SetSlope(1.0);
382 }
383 SetGamma0( GetSlope()*Xtotal/10.0/2.0/pi );
384
385 G4double Xinel = Xtotal - Xelastic;
386
387 #ifdef debugFTFparams
388 G4cout<< "Slope of hN elastic scattering" << GetSlope() << G4endl;
389 G4cout << "AvaragePt2ofElasticScattering " << GetAvaragePt2ofElasticScattering() << G4endl;
390 G4cout<<"Parameters of excitation for projectile "<<ProjectilePDGcode<< G4endl;
391 #endif
392
393 if ( (ProjectilePDGcode == 2212) || (ProjectilePDGcode == 2112) ) { // Projectile is proton or neutron
394
395 // A process probability is parameterized as Prob = A_1*exp(-A_2*y) + A_3*exp(-A_4*y) + A_top
396 // y is a rapidity of a partcle in the target nucleus. Ymin is a minimal rapidity below it X=0
397
398 // Proc# A1 B1 A2 B2 A3 Atop Ymin
399 /* original hadr-string-diff-V10-03-07 (similar to 10.3.x)
400 SetParams( 0, 13.71, 1.75, -214.5, 4.25, 0.0, 0.5 , 1.1 ); // Qexchange without Exc.
401 SetParams( 1, 25.0, 1.0, -50.34, 1.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc.
402 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Projectile diffraction
403 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Target diffraction
404 SetParams( 4, 1.0, 0.0 , -2.01 , 0.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc. Additional multiplier
405 */
406
407 // Proc#
408 SetParams( 0, fParCollBaryonProj.GetProc0A1(), fParCollBaryonProj.GetProc0B1(),
409 fParCollBaryonProj.GetProc0A2(), fParCollBaryonProj.GetProc0B2(),
410 fParCollBaryonProj.GetProc0A3(), fParCollBaryonProj.GetProc0Atop(),
411 fParCollBaryonProj.GetProc0Ymin() ); // Qexchange without Exc.
412 SetParams( 1, fParCollBaryonProj.GetProc1A1(), fParCollBaryonProj.GetProc1B1(),
413 fParCollBaryonProj.GetProc1A2(), fParCollBaryonProj.GetProc1B2(),
414 fParCollBaryonProj.GetProc1A3(), fParCollBaryonProj.GetProc1Atop(),
415 fParCollBaryonProj.GetProc1Ymin() ); // Qexchange with Exc.
416 if ( Xinel > 0.0 ) {
417 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Projectile diffraction
418 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Target diffraction
419
420 SetParams( 4, fParCollBaryonProj.GetProc4A1(), fParCollBaryonProj.GetProc4B1(),
421 fParCollBaryonProj.GetProc4A2(), fParCollBaryonProj.GetProc4B2(),
422 fParCollBaryonProj.GetProc4A3(), fParCollBaryonProj.GetProc4Atop(),
423 fParCollBaryonProj.GetProc4Ymin() ); // Qexchange with Exc. Additional multiplier
424 } else { // if Xinel=0., zero everything out (obviously)
425 SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
426 SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
427 SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
428 }
429
430 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
431 //
432 // It is not decided what to do with diffraction dissociation in Had-Nucl and Nucl-Nucl interactions
433 // For the moment both ProjDiffDisso & TgtDiffDisso for A > 10 are set to false,
434 // so both projectile and target diffraction are turned OFF
435 //
436 if ( ! fParCollBaryonProj.IsProjDiffDissociation() )
437 SetParams( 2, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
438 if ( ! fParCollBaryonProj.IsTgtDiffDissociation() )
439 SetParams( 3, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
440 }
441
443
444 if ( NumberOfTargetNucleons > 26 ) {
446 } else {
448 }
449
450 SetProjMinDiffMass( fParCollBaryonProj.GetProjMinDiffMass() ); // GeV
451 SetProjMinNonDiffMass( fParCollBaryonProj.GetProjMinNonDiffMass() ); // GeV
452
453 SetTarMinDiffMass( fParCollBaryonProj.GetTgtMinDiffMass() ); // GeV
454 SetTarMinNonDiffMass( fParCollBaryonProj.GetTgtMinNonDiffMass() ); // GeV
455
456 SetAveragePt2( fParCollBaryonProj.GetAveragePt2() ); // GeV^2
457 SetProbLogDistrPrD( fParCollBaryonProj.GetProbLogDistrPrD() );
458 SetProbLogDistr( fParCollBaryonProj.GetProbLogDistr() );
459
460 } else if ( ProjectilePDGcode == -2212 || ProjectilePDGcode == -2112 ) { // Projectile is anti_proton or anti_neutron
461
462 // Proc# A1 B1 A2 B2 A3 Atop Ymin
463 SetParams( 0, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange without Exc.
464 SetParams( 1, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 1000.0 ); // Qexchange with Exc.
465 if ( Xinel > 0.) {
466 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Projectile diffraction
467 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93 ); // Target diffraction
468 SetParams( 4, 1.0, 0.0 , 0.0, 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
469 } else {
470 SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0 );
471 SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0 );
472 SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0 );
473 }
474
475 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
476 //
477 // It is not decided what to do with diffraction dissociation in Had-Nucl and Nucl-Nucl interactions
478 // For the moment both ProjDiffDisso & TgtDiffDisso are set to false,
479 // so both projectile and target diffraction are turned OFF
480 //
481 if ( ! fParCollBaryonProj.IsProjDiffDissociation() )
482 SetParams( 2, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
483 if ( ! fParCollBaryonProj.IsTgtDiffDissociation() )
484 SetParams( 3, 0.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
485 }
486
489 SetProjMinDiffMass( ProjectileMass + 0.22 ); // GeV
490 SetProjMinNonDiffMass( ProjectileMass + 0.22 ); // GeV
491 SetTarMinDiffMass( TargetMass + 0.22 ); // GeV
492 SetTarMinNonDiffMass( TargetMass + 0.22 ); // GeV
493 SetAveragePt2( 0.3 ); // GeV^2
494 SetProbLogDistrPrD( 0.55 );
495 SetProbLogDistr( 0.55 );
496
497 } else if ( ProjectileabsPDGcode == 211 || ProjectilePDGcode == 111 ) { // Projectile is Pion
498
499 // Proc# A1 B1 A2 B2 A3 Atop Ymin
500 /* --> original code
501 SetParams( 0, 150.0, 1.8 , -247.3, 2.3, 0., 1. , 2.3 );
502 SetParams( 1, 5.77, 0.6 , -5.77, 0.8, 0., 0. , 0.0 );
503 SetParams( 2, 2.27, 0.5 , -98052.0, 4.0, 0., 0. , 3.0 );
504 SetParams( 3, 7.0, 0.9, -85.28, 1.9, 0.08, 0. , 2.2 );
505 SetParams( 4, 1.0, 0.0 , -11.02, 1.0, 0.0, 0. , 2.4 ); // Qexchange with Exc. Additional multiply
506 */
507 // Proc#
508 SetParams( 0, fParCollPionProj.GetProc0A1(), fParCollPionProj.GetProc0B1(),
509 fParCollPionProj.GetProc0A2(), fParCollPionProj.GetProc0B2(),
510 fParCollPionProj.GetProc0A3(), fParCollPionProj.GetProc0Atop(),
511 fParCollPionProj.GetProc0Ymin() ); // Qexchange without Exc.
512 SetParams( 1, fParCollPionProj.GetProc1A1(), fParCollPionProj.GetProc1B1(),
513 fParCollPionProj.GetProc1A2(), fParCollPionProj.GetProc1B2(),
514 fParCollPionProj.GetProc1A3(), fParCollPionProj.GetProc1Atop(),
515 fParCollPionProj.GetProc1Ymin() ); // Qexchange with Exc.
516 SetParams( 2, fParCollPionProj.GetProc2A1(), fParCollPionProj.GetProc2B1(),
517 fParCollPionProj.GetProc2A2(), fParCollPionProj.GetProc2B2(),
518 fParCollPionProj.GetProc2A3(), fParCollPionProj.GetProc2Atop(),
519 fParCollPionProj.GetProc2Ymin() ); // Projectile diffraction
520 SetParams( 3, fParCollPionProj.GetProc3A1(), fParCollPionProj.GetProc3B1(),
521 fParCollPionProj.GetProc3A2(), fParCollPionProj.GetProc3B2(),
522 fParCollPionProj.GetProc3A3(), fParCollPionProj.GetProc3Atop(),
523 fParCollPionProj.GetProc3Ymin() ); // Target diffraction
524 SetParams( 4, fParCollPionProj.GetProc4A1(), fParCollPionProj.GetProc4B1(),
525 fParCollPionProj.GetProc4A2(), fParCollPionProj.GetProc4B2(),
526 fParCollPionProj.GetProc4A3(), fParCollPionProj.GetProc4Atop(),
527 fParCollPionProj.GetProc4Ymin() ); // Qexchange with Exc. Additional multiply
528
529 // NOTE: how can it be |ProjectileBaryonNumber| > 10 if projectile is a pion ???
530 //
531 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
532 if ( ! fParCollPionProj.IsProjDiffDissociation() )
533 SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
534 if ( ! fParCollPionProj.IsTgtDiffDissociation() )
535 SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
536 }
537
538 /* original code -->
539 SetDeltaProbAtQuarkExchange( 0.56 );
540 SetProjMinDiffMass( 1.0 ); // GeV
541 SetProjMinNonDiffMass( 1.0 ); // GeV
542 SetTarMinDiffMass( 1.16 ); // GeV
543 SetTarMinNonDiffMass( 1.16 ); // GeV
544 SetAveragePt2( 0.3 ); // GeV^2
545 SetProbLogDistrPrD( 0.55 );
546 SetProbLogDistr( 0.55 );
547 */
548
549 // JVY update, Aug.8, 2018 --> Feb.14, 2019
550 //
552 SetProjMinDiffMass( fParCollPionProj.GetProjMinDiffMass() ); // GeV
553 SetProjMinNonDiffMass( fParCollPionProj.GetProjMinNonDiffMass() ); // GeV
554 SetTarMinDiffMass( fParCollPionProj.GetTgtMinDiffMass() ); // GeV
555 SetTarMinNonDiffMass( fParCollPionProj.GetTgtMinNonDiffMass() ); // GeV
556 SetAveragePt2( fParCollPionProj.GetAveragePt2() ); // GeV^2
557 SetProbLogDistrPrD( fParCollPionProj.GetProbLogDistrPrD() );
558 SetProbLogDistr( fParCollPionProj.GetProbLogDistr() );
559
560 // ---> end update
561
562 } else if ( ProjectileabsPDGcode == 321 || ProjectileabsPDGcode == 311 ||
563 ProjectilePDGcode == 130 || ProjectilePDGcode == 310 ) { // Projectile is Kaon
564
565 // Proc# A1 B1 A2 B2 A3 Atop Ymin
566 SetParams( 0, 60.0 , 2.5 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Qexchange without Exc.
567 SetParams( 1, 6.0 , 1.0 , -24.33 , 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
568 SetParams( 2, 2.76, 1.2 , -22.5 , 2.7 ,0.04, 0.0 , 1.40 ); // Projectile diffraction
569 SetParams( 3, 1.09, 0.5 , -8.88 , 2. ,0.05, 0.0 , 1.40 ); // Target diffraction
570 SetParams( 4, 1.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
571 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
572 SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
573 SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
574 }
575
577 SetProjMinDiffMass( 0.7 ); // GeV
578 SetProjMinNonDiffMass( 0.7 ); // GeV
579 SetTarMinDiffMass( 1.16 ); // GeV
580 SetTarMinNonDiffMass( 1.16 ); // GeV
581 SetAveragePt2( 0.3 ); // GeV^2
582 SetProbLogDistrPrD( 0.55 );
583 SetProbLogDistr( 0.55 );
584
585 } else { // Projectile is not p, n, Pi0, Pi+, Pi-, K+, K-, K0 or their anti-particles
586
587 if ( ProjectileabsPDGcode > 1000 ) { // The projectile is a baryon as P or N
588 // Proc# A1 B1 A2 B2 A3 Atop Ymin
589 SetParams( 0, 13.71, 1.75, -30.69, 3.0 , 0.0, 1.0 , 0.93 ); // Qexchange without Exc.
590 SetParams( 1, 25.0, 1.0, -50.34, 1.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc.
591 if ( Xinel > 0.) {
592 SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Projectile diffraction
593 SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93); // Target diffraction
594 SetParams( 4, 1.0, 0.0 , -2.01 , 0.5 , 0.0, 0.0 , 1.4 ); // Qexchange with Exc. Additional multiply
595 } else {
596 SetParams( 2, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
597 SetParams( 3, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
598 SetParams( 4, 0.0, 0.0 ,0.0, 0.0 , 0.0, 0.0 , 0.0);
599 }
600
601 } else { // The projectile is a meson as K+-0
602 // Proc# A1 B1 A2 B2 A3 Atop Ymin
603 SetParams( 0, 60.0 , 2.5 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Qexchange without Exc.
604 SetParams( 1, 6.0 , 1.0 , -24.33 , 2.0 , 0.0, 0.0 , 1.40 ); // Qexchange with Exc.
605 SetParams( 2, 2.76, 1.2 , -22.5 , 2.7 ,0.04, 0.0 , 1.40 ); // Projectile diffraction
606 SetParams( 3, 1.09, 0.5 , -8.88 , 2. ,0.05, 0.0 , 1.40 ); // Target diffraction
607 SetParams( 4, 1.0, 0.0 , 0.0 , 0.0 , 0.0, 0.0 , 0.93 ); // Qexchange with Exc. Additional multiply
608 }
609
610 if ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 ) {
611 SetParams( 2, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Projectile diffraction
612 SetParams( 3, 0.0 , 0.0 , 0.0 , 0.0 , 0.0, 0.0 , -100.0 ); // Target diffraction
613 }
614
617
618 SetProjMinDiffMass( GetMinMass(particle)/GeV );
619 SetProjMinNonDiffMass( GetMinMass(particle)/GeV );
620
622 SetTarMinDiffMass( GetMinMass(Neutron)/GeV );
623 SetTarMinNonDiffMass( GetMinMass(Neutron)/GeV );
624
625 SetAveragePt2( 0.3 ); // GeV^2
626 SetProbLogDistrPrD( 0.55 );
627 SetProbLogDistr( 0.55 );
628
629 }
630
631 #ifdef debugFTFparams
632 G4cout<<"DeltaProbAtQuarkExchange "<< GetDeltaProbAtQuarkExchange() << G4endl;
633 G4cout<<"ProbOfSameQuarkExchange "<< GetProbOfSameQuarkExchange() << G4endl;
634 G4cout<<"ProjMinDiffMass "<< GetProjMinDiffMass()/GeV <<" GeV"<< G4endl;
635 G4cout<<"ProjMinNonDiffMass "<< GetProjMinNonDiffMass() <<" GeV"<< G4endl;
636 G4cout<<"TarMinDiffMass "<< GetTarMinDiffMass() <<" GeV"<< G4endl;
637 G4cout<<"TarMinNonDiffMass "<< GetTarMinNonDiffMass() <<" GeV"<< G4endl;
638 G4cout<<"AveragePt2 "<< GetAveragePt2() <<" GeV^2"<< G4endl;
639 G4cout<<"ProbLogDistrPrD "<< GetProbLogDistrPrD() << G4endl;
640 G4cout<<"ProbLogDistrTrD "<< GetProbLogDistr() << G4endl;
641 #endif
642
643 // Set parameters of nuclear destruction
644
645 if ( ProjectileabsPDGcode < 1000 ) { // Meson projectile
646
647 SetMaxNumberOfCollisions( Plab, 2.0 ); // 3.0 )
648 //
649 // target destruction
650 //
651 /* original code --->
652 SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)*
653 G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
654
655 SetR2ofNuclearDestruction( 1.5*fermi*fermi );
656 SetDofNuclearDestruction( 0.3 );
657 SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
658 ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
659 SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
660 SetExcitationEnergyPerWoundedNucleon( 40.0*MeV );
661 */
662 double coeff = fParCollMesonProj.GetNuclearTgtDestructP1();
663 //
664 // NOTE (JVY): Set this switch to false/true on line 138
665 //
666 if ( fParCollMesonProj.IsNuclearTgtDestructP1_ADEP() )
667 {
668 coeff *= G4double(NumberOfTargetNucleons);
669 }
670 double exfactor = G4Exp( fParCollMesonProj.GetNuclearTgtDestructP2()
671 * (Ylab-fParCollMesonProj.GetNuclearTgtDestructP3()) );
672 coeff *= exfactor;
673 coeff /= ( 1.+ exfactor );
674
676
678 SetDofNuclearDestruction( fParCollMesonProj.GetDofNuclearDestruct() );
679 coeff = fParCollMesonProj.GetPt2NuclearDestructP2();
680 exfactor = G4Exp( fParCollMesonProj.GetPt2NuclearDestructP3()
681 * (Ylab-fParCollMesonProj.GetPt2NuclearDestructP4()) );
682 coeff *= exfactor;
683 coeff /= ( 1. + exfactor );
684 SetPt2ofNuclearDestruction( (fParCollMesonProj.GetPt2NuclearDestructP1()+coeff)*CLHEP::GeV*CLHEP::GeV );
685
688
689 } else if ( ProjectilePDGcode == -2212 || ProjectilePDGcode == -2112 ) { // for anti-baryon projectile
690
691 SetMaxNumberOfCollisions( Plab, 2.0 );
692
693 SetCofNuclearDestruction( 0.00481*G4double(NumberOfTargetNucleons)*
694 G4Exp( 4.0*(Ylab - 2.1) )/( 1.0 + G4Exp( 4.0*(Ylab - 2.1) ) ) );
695 SetR2ofNuclearDestruction( 1.5*fermi*fermi );
697 SetPt2ofNuclearDestruction( ( 0.035 + 0.04*G4Exp( 4.0*(Ylab - 2.5) )/
698 ( 1.0 + G4Exp( 4.0*(Ylab - 2.5) ) ) )*GeV*GeV );
699 SetMaxPt2ofNuclearDestruction( 1.0*GeV*GeV );
701 if ( Plab < 2.0 ) { // 2 GeV/c
702 // For slow anti-baryon we have to garanty putting on mass-shell
704 SetR2ofNuclearDestruction( 1.5*fermi*fermi ); // this is equivalent to setting a few line above
705 // is it even necessary ?
707 SetPt2ofNuclearDestruction( 0.035*GeV*GeV );
708 SetMaxPt2ofNuclearDestruction( 0.04*GeV*GeV );
709 }
710
711 } else { // Projectile baryon assumed
712
713 // NOTE (JVY) FIXME !!! Will decide later how/if to make this one configurable...
714 //
715 SetMaxNumberOfCollisions( Plab, 2.0 );
716
717 // projectile destruction - does NOT really matter for particle projectile, only for a nucleus projectile
718 //
719 double coeff = 0.;
720 coeff = fParCollBaryonProj.GetNuclearProjDestructP1();
721 //
722 // NOTE (JVY): Set this switch to false/true on line 136
723 //
724 if ( fParCollBaryonProj.IsNuclearProjDestructP1_NBRNDEP() )
725 {
726 coeff *= G4double(AbsProjectileBaryonNumber);
727 }
728 double exfactor = G4Exp( fParCollBaryonProj.GetNuclearProjDestructP2()*(Ylab-fParCollBaryonProj.GetNuclearProjDestructP3()) );
729 coeff *= exfactor;
730 coeff /= ( 1.+ exfactor );
732
733 // target desctruction
734 //
735 coeff = fParCollBaryonProj.GetNuclearTgtDestructP1();
736 //
737 // NOTE (JVY): Set this switch to false/true on line 138
738 //
739 if ( fParCollBaryonProj.IsNuclearTgtDestructP1_ADEP() )
740 {
741 coeff *= G4double(NumberOfTargetNucleons);
742 }
743 exfactor = G4Exp( fParCollBaryonProj.GetNuclearTgtDestructP2()*(Ylab-fParCollBaryonProj.GetNuclearTgtDestructP3()) );
744 coeff *= exfactor;
745 coeff /= ( 1.+ exfactor );
747
749 SetDofNuclearDestruction( fParCollBaryonProj.GetDofNuclearDestruct() );
750
751 coeff = fParCollBaryonProj.GetPt2NuclearDestructP2();
752 exfactor = G4Exp( fParCollBaryonProj.GetPt2NuclearDestructP3()*(Ylab-fParCollBaryonProj.GetPt2NuclearDestructP4()) );
753 coeff *= exfactor;
754 coeff /= ( 1. + exfactor );
755 SetPt2ofNuclearDestruction( (fParCollBaryonProj.GetPt2NuclearDestructP1()+coeff)*CLHEP::GeV*CLHEP::GeV );
756
759
760 }
761
762 #ifdef debugFTFparams
763 G4cout<<"CofNuclearDestructionPr "<< GetCofNuclearDestructionPr() << G4endl;
764 G4cout<<"CofNuclearDestructionTr "<< GetCofNuclearDestruction() << G4endl;
765 G4cout<<"R2ofNuclearDestruction "<< GetR2ofNuclearDestruction()/fermi/fermi <<" fermi^2"<< G4endl;
766 G4cout<<"DofNuclearDestruction "<< GetDofNuclearDestruction() << G4endl;
767 G4cout<<"Pt2ofNuclearDestruction "<< GetPt2ofNuclearDestruction()/GeV/GeV <<" GeV^2"<< G4endl;
768 G4cout<<"ExcitationEnergyPerWoundedNucleon "<< GetExcitationEnergyPerWoundedNucleon() <<" MeV"<< G4endl;
769 #endif
770
771 //SetCofNuclearDestruction( 0.47*G4Exp( 2.0*(Ylab - 2.5) )/( 1.0 + G4Exp( 2.0*(Ylab - 2.5) ) ) );
772 //SetPt2ofNuclearDestruction( ( 0.035 + 0.1*G4Exp( 4.0*(Ylab - 3.0) )/( 1.0 + G4Exp( 4.0*(Ylab - 3.0) ) ) )*GeV*GeV );
773
774 //SetMagQuarkExchange( 120.0 ); // 210.0 PipP
775 //SetSlopeQuarkExchange( 2.0 );
776 //SetDeltaProbAtQuarkExchange( 0.6 );
777 //SetProjMinDiffMass( 0.7 ); // GeV 1.1
778 //SetProjMinNonDiffMass( 0.7 ); // GeV
779 //SetProbabilityOfProjDiff( 0.0); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
780 //SetTarMinDiffMass( 1.1 ); // GeV
781 //SetTarMinNonDiffMass( 1.1 ); // GeV
782 //SetProbabilityOfTarDiff( 0.0 ); // 0.85*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.5 ) ); // 40/32 X-dif/X-inel
783
784 //SetAveragePt2( 0.0 ); // GeV^2 0.3
785 //------------------------------------
786 //SetProbabilityOfElasticScatt( 1.0, 1.0); //(Xtotal, Xelastic);
787 //SetProbabilityOfProjDiff( 1.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 0->1
788 //SetProbabilityOfTarDiff( 4.0*0.62*G4Pow::GetInstance()->powA( s/GeV/GeV, -0.51 ) ); // 2->4
789 //SetAveragePt2( 0.3 ); // (0.15)
790 //SetAvaragePt2ofElasticScattering( 0.0 );
791
792 //SetMaxNumberOfCollisions( Plab, 6.0 ); //(4.0*(Plab + 0.01), Plab); // 6.0 );
793 //SetAveragePt2( 0.15 );
794 //SetCofNuclearDestruction(-1.);//( 0.75 ); // (0.25)
795 //SetExcitationEnergyPerWoundedNucleon(0.);//( 30.0*MeV ); // (75.0*MeV)
796 //SetDofNuclearDestruction(0.);//( 0.2 ); //0.4 // 0.3 0.5
797
798 /*
799 SetAveragePt2(0.3);
800 SetCofNuclearDestructionPr(0.);
801 SetCofNuclearDestruction(0.); //( 0.5 ); (0.25)
802 SetExcitationEnergyPerWoundedNucleon(0.); // 30.0*MeV; (75.0*MeV)
803 SetDofNuclearDestruction(0.); // 0.2; 0.4; 0.3; 0.5
804 SetPt2ofNuclearDestruction(0.); //(2.*0.075*GeV*GeV); ( 0.3*GeV*GeV ); (0.168*GeV*GeV)
805 */
806
807 //SetExcitationEnergyPerWoundedNucleon(0.001);
808 //SetPt2Kink( 0.0*GeV*GeV );
809
810 //SetRadiusOfHNinteractions2( Xtotal/pi/10.0 /2.);
811 //SetRadiusOfHNinteractions2( (Xtotal - Xelastic)/pi/10.0 );
812 //SetProbabilityOfElasticScatt( 1.0, 0.0);
813 /*
814 G4cout << "Pt2 " << GetAveragePt2()<<" "<<GetAveragePt2()/GeV/GeV<<G4endl;
815 G4cout << "Cnd " << GetCofNuclearDestruction() << G4endl;
816 G4cout << "Dnd " << GetDofNuclearDestruction() << G4endl;
817 G4cout << "Pt2 " << GetPt2ofNuclearDestruction()/GeV/GeV << G4endl;
818 */
819
820}
821
822//============================================================================
823
824G4double G4FTFParameters::GetMinMass( const G4ParticleDefinition* aParticle ) {
825 // The code is used for estimating the minimal string mass produced in diffraction dissociation.
826 // The indices used for minMassQDiQStr must be between 1 and 5, corresponding to the 5 considered
827 // quarks: d, u, s, c and b; enforcing this explicitly avoids compilation errors.
828 G4double EstimatedMass = 0.0;
829 G4int partID = std::abs(aParticle->GetPDGEncoding());
830 G4int Qleft = std::max( partID/100, 1 );
831 G4int Qright = std::max( (partID/ 10)%10, 1 );
832 if ( Qleft < 6 && Qright < 6 ) { // Q-Qbar string
833 EstimatedMass = StringMass->minMassQQbarStr[Qleft-1][Qright-1];
834 } else if ( Qleft < 6 && Qright > 6 ) { // Q - DiQ string
835 G4int q1 = std::max( std::min( Qright/10, 5 ), 1 );
836 G4int q2 = std::max( std::min( Qright%10, 5 ), 1 );
837 EstimatedMass = StringMass->minMassQDiQStr[Qleft-1][q1-1][q2-1];
838 } else if ( Qleft > 6 && Qright < 6 ) { // DiQ - Q string
839 G4int q1 = std::max( std::min( Qleft/10, 5 ), 1 );
840 G4int q2 = std::max( std::min( Qleft%10, 5 ), 1 );
841 EstimatedMass = StringMass->minMassQDiQStr[Qright-1][q1-1][q2-1];
842 }
843 return EstimatedMass;
844}
845
846//============================================================================
847
849 G4double Prob( 0.0 );
850 if ( y < ProcParams[ProcN][6] ) {
851 Prob = ProcParams[ProcN][5];
852 if (Prob < 0.) Prob=0.;
853 return Prob;
854 }
855 Prob = ProcParams[ProcN][0] * G4Exp( -ProcParams[ProcN][1]*y ) +
856 ProcParams[ProcN][2] * G4Exp( -ProcParams[ProcN][3]*y ) +
857 ProcParams[ProcN][4];
858 if (Prob < 0.) Prob=0.;
859 return Prob;
860}
861
862//============================================================================
863
865
866
868{
869 public:
870
871 // ctor
873
874 //==================================================================
875 //
876 // Cross sections for elementary processes
877 //
878 // these are for Inelastic interactions, i.e. Xinelastic=(Xtotal-Xelastix)>0.
879 // for elastic, all the A's & B's, Atop & Ymin are zeros
880 // general formula: Pp = A1*exp(B1*Y) + A2*exp(B2*Y) + A3
881 // but if Y<Ymin, then Pp=max(0.,Atop)
882 // for details, see also G4FTFParameters::GetProcProb( ProcN, y )
883 //
884 // baryons
885 //
886 /* JVY, Oct. 31, 2017: Per Alberto R. & Vladimir U., keep this group of parameters FIXED */
887 /* JVY, June 11, 2020: try to open up... */
888 // Process=0 --> Qexchg w/o excitation
889 //
890 HDP.SetDefault( "FTF_BARYON_PROC0_A1", 13.71 );
891 HDP.SetDefault( "FTF_BARYON_PROC0_B1", 1.75 );
892 HDP.SetDefault( "FTF_BARYON_PROC0_A2", -30.69 );
893 HDP.SetDefault( "FTF_BARYON_PROC0_B2", 3.0 );
894 HDP.SetDefault( "FTF_BARYON_PROC0_A3", 0.0 );
895 HDP.SetDefault( "FTF_BARYON_PROC0_ATOP", 1.0 );
896 HDP.SetDefault( "FTF_BARYON_PROC0_YMIN", 0.93 );
897 /* */
898 //
899 // Process=1 --> Qexchg w/excitation
900 //
901 HDP.SetDefault( "FTF_BARYON_PROC1_A1", 25. );
902 HDP.SetDefault( "FTF_BARYON_PROC1_B1", 1. );
903 HDP.SetDefault( "FTF_BARYON_PROC1_A2", -50.34 );
904 HDP.SetDefault( "FTF_BARYON_PROC1_B2", 1.5 );
905 HDP.SetDefault( "FTF_BARYON_PROC1_A3", 0. );
906 HDP.SetDefault( "FTF_BARYON_PROC1_ATOP", 0. );
907 HDP.SetDefault( "FTF_BARYON_PROC1_YMIN", 1.4 );
908 /* */
909 //
910 // NOTE: Process #2 & 3 are projectile & target diffraction
911 // they have more complex definition of A1 & A2
912 // (see example below)
913 // SetParams( 2, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Projectile diffraction
914 // SetParams( 3, 6.0/Xinel, 0.0 ,-6.0/Xinel*16.28, 3.0 , 0.0, 0.0 , 0.93);// Target diffraction
915 //
916 // Also, for ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 )
917 // projectile and/or target diffraction (dissociation) may be switched ON/OFF
918 //
919 HDP.SetDefault( "FTF_BARYON_DIFF_DISSO_PROJ", false );
920 HDP.SetDefault( "FTF_BARYON_DIFF_DISSO_TGT", false ); // as in hadr-string-diff-V10-03-07
921 //
922 /* JVY, Oct. 31, 2017: Per Alberto R. & Vladimir U., keep this group of parameters FIXED */
923 /* JVY, June 11, 2020: try to open up... */
924 // Process=4 --> Qexchg w/additional multiplier in excitation
925 //
926 HDP.SetDefault( "FTF_BARYON_PROC4_A1", 0.6 );
927 HDP.SetDefault( "FTF_BARYON_PROC4_B1", 0. );
928 HDP.SetDefault( "FTF_BARYON_PROC4_A2", -1.2 );
929 HDP.SetDefault( "FTF_BARYON_PROC4_B2", 0.5 );
930 HDP.SetDefault( "FTF_BARYON_PROC4_A3", 0. );
931 HDP.SetDefault( "FTF_BARYON_PROC4_ATOP",0. );
932 HDP.SetDefault( "FTF_BARYON_PROC4_YMIN",1.4 );
933 /* */
934 //
935 // Parameters of participating hadron (baryon) excitation
936 //
937 HDP.SetDefault( "FTF_BARYON_DELTA_PROB_QEXCHG", 0. );
938 HDP.SetDefault( "FTF_BARYON_PROB_SAME_QEXCHG", 0. );
939 HDP.SetDefault( "FTF_BARYON_DIFF_M_PROJ", 1.16, 1.16, 3. ); // it's supposed to be in GeV but do NOT do (*CLHEP::GeV)
940 // because it'll be done in the G4FTFParameters::SetProjMinDiffMass
941 HDP.SetDefault( "FTF_BARYON_NONDIFF_M_PROJ", 1.16, 1.16, 3. ); // do NOT (*CLHEP::GeV) - same as above
942 HDP.SetDefault( "FTF_BARYON_DIFF_M_TGT", 1.16, 1.16, 3. ); // do NOT (*CLHEP::GeV) - same as above
943 HDP.SetDefault( "FTF_BARYON_NONDIFF_M_TGT", 1.16, 1.16, 3. ); // do NOT (*CLHEP::GeV) - same as above
944 HDP.SetDefault( "FTF_BARYON_AVRG_PT2", 0.3, 0.08, 1. ); // do NOT (*CLHEP::GeV*CLHEP::GeV)
945 //
946 // JVY, Oct. 6, 2017: Per Alberto R., keep these two settings fixed (for now)
947 //
948 // HDP.SetDefault( "FTF_BARYON_PROB_DISTR_PROJ", 0.3 );
949 // HDP.SetDefault( "FTF_BARYON_PROB_DISTR_TGT", 0.3 );
950
951
952 // pions
953 //
954 // JVY, Aug.8, 2018 --> Feb.14, 2019 --> June 25, 2019:
955 // Parameters of participating hadron (pions) excitation
956 //
957 /* JVY, June 25, 2019: For now, keep this group of parameters FIXED */
958 // Process=0 --> Qexchg w/o excitation
959 //
960 HDP.SetDefault( "FTF_PION_PROC0_A1", 150.0 );
961 HDP.SetDefault( "FTF_PION_PROC0_B1", 1.8 );
962 HDP.SetDefault( "FTF_PION_PROC0_A2",-247.3 );
963 HDP.SetDefault( "FTF_PION_PROC0_B2", 2.3 );
964 HDP.SetDefault( "FTF_PION_PROC0_A3", 0.0 );
965 HDP.SetDefault( "FTF_PION_PROC0_ATOP", 1.0 );
966 HDP.SetDefault( "FTF_PION_PROC0_YMIN", 2.3 );
967 //
968 // Process=1 --> Qexchg w/excitation
969 //
970 HDP.SetDefault( "FTF_PION_PROC1_A1", 5.77 );
971 HDP.SetDefault( "FTF_PION_PROC1_B1", 0.6 );
972 HDP.SetDefault( "FTF_PION_PROC1_A2", -5.77 );
973 HDP.SetDefault( "FTF_PION_PROC1_B2", 0.8 );
974 HDP.SetDefault( "FTF_PION_PROC1_A3", 0. );
975 HDP.SetDefault( "FTF_PION_PROC1_ATOP", 0. );
976 HDP.SetDefault( "FTF_PION_PROC1_YMIN", 0.0 );
977 /*
978 //
979 // NOTE: Process #2 & 3 are projectile & target diffraction
980 //
981 // Process=2 --> Projectile diffraction
982 //
983 // Q: Would it even make sense to make these configurable ?
984 // The following is hadrcoded:
985 // Projectile Baryon Number > 10 (AbsProjectileBaryonNumber > 10)
986 // ... which is "strange" because projectile is a pion !!!... so it's always OFF
987 // (see also lines 1007-1016)
988 //
989 HDP.SetDefault( "FTF_PION_PROC2_A1", 2.27 );
990 HDP.SetDefault( "FTF_PION_PROC2_B1", 0.5 );
991 HDP.SetDefault( "FTF_PION_PROC2_A2", -98052.0);
992 HDP.SetDefault( "FTF_PION_PROC2_B2", 4.0 );
993 HDP.SetDefault( "FTF_PION_PROC2_A3", 0. );
994 HDP.SetDefault( "FTF_PION_PROC2_ATOP", 0. );
995 HDP.SetDefault( "FTF_PION_PROC2_YMIN", 3.0 );
996 */
997 //
998 // Process=3 --> Target diffraction
999 //
1000 HDP.SetDefault( "FTF_PION_PROC3_A1", 7.0 );
1001 HDP.SetDefault( "FTF_PION_PROC3_B1", 0.9 );
1002 HDP.SetDefault( "FTF_PION_PROC3_A2", -85.28);
1003 HDP.SetDefault( "FTF_PION_PROC3_B2", 1.9 );
1004 HDP.SetDefault( "FTF_PION_PROC3_A3", 0.08);
1005 HDP.SetDefault( "FTF_PION_PROC3_ATOP", 0. );
1006 HDP.SetDefault( "FTF_PION_PROC3_YMIN", 2.2 );
1007 //
1008 // projectile and/or target diffraction (dissociation) may be switched ON/OFF
1009 //
1010 // NOTE: Both projectile and target diffraction are turned OFF if
1011 // a) Number of Target Nucleons > 10 (NumberOfTargetNucleons > 10)
1012 // OR
1013 // b) Projectile Baryon Number > 10 (AbsProjectileBaryonNumber > 10)
1014 // ... which is "strange" because projectile is a pion !!!... so it's always OFF
1015 //
1016 HDP.SetDefault( "FTF_PION_DIFF_DISSO_PROJ", false );
1017 HDP.SetDefault( "FTF_PION_DIFF_DISSO_TGT", false );
1018 //
1019 /* JVY, June 25, 2019: For now keep this group of parameters FIXED */
1020 /* JVY, June 11, 2020: try to open up... */
1021 // Process=4 --> Qexchg w/additional multiplier in excitation
1022 //
1023 HDP.SetDefault( "FTF_PION_PROC4_A1", 1.0 );
1024 HDP.SetDefault( "FTF_PION_PROC4_B1", 0. );
1025 HDP.SetDefault( "FTF_PION_PROC4_A2",-11.02);
1026 HDP.SetDefault( "FTF_PION_PROC4_B2", 1.0 );
1027 HDP.SetDefault( "FTF_PION_PROC4_A3", 0. );
1028 HDP.SetDefault( "FTF_PION_PROC4_ATOP",0. );
1029 HDP.SetDefault( "FTF_PION_PROC4_YMIN",2.4 );
1030 /* */
1031 //
1032 // NOTE; As of geant4-10-05, all these settings beloe are correct
1033 // (and are the same as they were in 10.4.ref06)
1034
1035 HDP.SetDefault( "FTF_PION_DELTA_PROB_QEXCHG", 0.56 ); // in the past used to be 0.35
1036 HDP.SetDefault( "FTF_PION_DIFF_M_PROJ", 1.0, 0.5, 3. ); // in the past used to be 0.5...
1037 HDP.SetDefault( "FTF_PION_NONDIFF_M_PROJ", 1.0, 0.5, 3. ); // so why not set 0.5 as a low limit ?...
1038 // ... or perhaps even lower ?
1039 HDP.SetDefault( "FTF_PION_DIFF_M_TGT", 1.16, 1.16, 3. ); // All (NON)DIFF_M's are supposed to be in GeV but do NOT do (*CLHEP::GeV)
1040 // because it'll be done in the G4FTFParameters::SetProjMinDiffMass
1041 HDP.SetDefault( "FTF_PION_NONDIFF_M_TGT", 1.16, 1.16, 3. );
1042 HDP.SetDefault( "FTF_PION_AVRG_PT2", 0.3, 0.08, 1. ); // do NOT (*CLHEP::GeV*CLHEP::GeV)
1043
1044 //==================================================================
1045 //
1046 // nuclear destruction
1047 //
1048 // NOTE: Settings of most of these parameters are the same
1049 // for different types of projectile hadron
1050 // However, we decided to introduce separate variables
1051 // and configuration cards for each type of projectile
1052 //
1053 // baryons
1054 //
1055 // projectile destruction
1056 //
1057 HDP.SetDefault( "FTF_BARYON_NUCDESTR_P1_PROJ", 1., 0., 1. ); // in principle, it should be 1./NBRN - FIXME later !
1058 HDP.SetDefault( "FTF_BARYON_NUCDESTR_P1_NBRN_PROJ", false );
1059 //
1060 // for now, keep fixed p2 & p3 for the proj destruction
1061 // they're defined explicitly in G4FTFParamCollection ctor
1062 //
1063 // target destruction
1064 //
1065 HDP.SetDefault( "FTF_BARYON_NUCDESTR_P1_TGT", 1., 0., 1. );
1066 HDP.SetDefault( "FTF_BARYON_NUCDESTR_P1_ADEP_TGT", false );
1067 HDP.SetDefault( "FTF_BARYON_NUCDESTR_P2_TGT", 4.0, 2., 16. );
1068 HDP.SetDefault( "FTF_BARYON_NUCDESTR_P3_TGT", 2.1, 0., 4. );
1069 //
1070 HDP.SetDefault( "FTF_BARYON_PT2_NUCDESTR_P1", 0.035, 0., 0.25 );
1071 HDP.SetDefault( "FTF_BARYON_PT2_NUCDESTR_P2", 0.04, 0., 0.25 );
1072 HDP.SetDefault( "FTF_BARYON_PT2_NUCDESTR_P3", 4.0, 2., 16. );
1073 HDP.SetDefault( "FTF_BARYON_PT2_NUCDESTR_P4", 2.5, 0., 4. );
1074 //
1075 HDP.SetDefault( "FTF_BARYON_NUCDESTR_R2", 1.5*CLHEP::fermi*CLHEP::fermi, 0.5*CLHEP::fermi*CLHEP::fermi, 2.*CLHEP::fermi*CLHEP::fermi );
1076 HDP.SetDefault( "FTF_BARYON_EXCI_E_PER_WNDNUCLN", 40.*CLHEP::MeV, 0., 100.*CLHEP::MeV );
1077 HDP.SetDefault( "FTF_BARYON_NUCDESTR_DISP", 0.3, 0.1, 0.4 );
1078 //
1079 // JVY, Oct. 6, 2017: Per Alberto R., this is just a technical parameter,
1080 // and it should NOT be changed
1081 //
1082 // HDP.SetDefault( "FTF_BARYON_NUCDESTR_MAXPT2", 1. * CLHEP::GeV*CLHEP::GeV );
1083
1084
1085 // mesons - these parameters are common for pions, kaons, etc. (per original code)
1086 //
1087 // NOTE: *NO* projectile destruction for mesons !!!
1088 //
1089 // target destruction
1090 //
1091 HDP.SetDefault( "FTF_MESON_NUCDESTR_P1_TGT", 0.00481, 0., 1. );
1092 HDP.SetDefault( "FTF_MESON_NUCDESTR_P1_ADEP_TGT", true );
1093 HDP.SetDefault( "FTF_MESON_NUCDESTR_P2_TGT", 4.0, 2., 16. );
1094 HDP.SetDefault( "FTF_MESON_NUCDESTR_P3_TGT", 2.1, 0., 4. );
1095 //
1096 HDP.SetDefault( "FTF_MESON_PT2_NUCDESTR_P1", 0.035, 0., 0.25 );
1097 HDP.SetDefault( "FTF_MESON_PT2_NUCDESTR_P2", 0.04, 0., 0.25 );
1098 HDP.SetDefault( "FTF_MESON_PT2_NUCDESTR_P3", 4.0, 2., 16. );
1099 HDP.SetDefault( "FTF_MESON_PT2_NUCDESTR_P4", 2.5, 0., 4. );
1100 //
1101 HDP.SetDefault( "FTF_MESON_NUCDESTR_R2", 1.5*CLHEP::fermi*CLHEP::fermi,
1102 0.5*CLHEP::fermi*CLHEP::fermi,
1103 2.*CLHEP::fermi*CLHEP::fermi );
1104 HDP.SetDefault( "FTF_MESON_EXCI_E_PER_WNDNUCLN", 40.*CLHEP::MeV, 0., 100.*CLHEP::MeV );
1105 HDP.SetDefault( "FTF_MESON_NUCDESTR_DISP", 0.3, 0.1, 0.4 );
1106
1107 }
1108};
1109
1110
1112
1113
1115{
1116
1117 Reset(); // zero out everything
1118
1119 //
1120 // keep the 2 parameters below fixed for now (i.e. do not take them from HDP)
1121 //
1124
1125}
1126
1127
1128void G4FTFParamCollection::Reset()
1129{
1130 // parameters of excitation
1131
1132 // Proc=0 --> Qexchg w/o excitation
1133 fProc0A1 = 0.;
1134 fProc0B1 = 0.;
1135 fProc0A2 = 0.;
1136 fProc0B2 = 0.;
1137 fProc0A3 = 0.;
1138 fProc0Atop = 0.;
1139 fProc0Ymin = 0.;
1140
1141 // Proc=1 --> Qexchg w/excitation
1142 fProc1A1 = 0.;
1143 fProc1B1 = 0.;
1144 fProc1A2 = 0.;
1145 fProc1B2 = 0.;
1146 fProc1A3 = 0.;
1147 fProc1Atop = 0.;
1148 fProc1Ymin = 0.;
1149
1150 // Proc=2 & Proc=3 for ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 )
1151 // Do NOT do anything as it's set once and for all !!!
1152
1153 fProjDiffDissociation = false;
1154 fTgtDiffDissociation = false;
1155
1156 // Proc=4 --> Qexchg w/additional multiplier in excitation
1157 fProc4A1 = 0.;
1158 fProc4B1 = 0.;
1159 fProc4A2 = 0.;
1160 fProc4B2 = 0.;
1161 fProc4A3 = 0.;
1162 fProc4Atop = 0.;
1163 fProc4Ymin = 0.;
1164
1165 // parameters of participating baryon excitation
1166
1169 fProjMinDiffMass = 0.;
1171 fTgtMinDiffMass = 0.;
1172 fTgtMinNonDiffMass = 0.;
1173 fAveragePt2 = 0.;
1174 fProbLogDistrPrD = 0.;
1175 fProbLogDistr = 0.;
1176
1177 // parameters of nuclear distruction
1178
1179 // COMMONs
1190
1191 // baryons
1196
1197 return;
1198}
1199
1200//============================================================================
1201
1204{
1205
1206 // parameters of participating hadron (baryon) excitation
1207 //
1208 // baryons projectile
1209 //
1210 // Proc=0 --> Qexchg w/o excitation
1211 //
1212 /* As of Oct. 31, 2017 keep these fixed */
1213 HDP.DeveloperGet( "FTF_BARYON_PROC0_A1", fProc0A1 );
1214 HDP.DeveloperGet( "FTF_BARYON_PROC0_B1", fProc0B1 );
1215 HDP.DeveloperGet( "FTF_BARYON_PROC0_A2", fProc0A2 );
1216 HDP.DeveloperGet( "FTF_BARYON_PROC0_B2", fProc0B2 );
1217 HDP.DeveloperGet( "FTF_BARYON_PROC0_A3", fProc0A3 );
1218 HDP.DeveloperGet( "FTF_BARYON_PROC0_ATOP", fProc0Atop );
1219 HDP.DeveloperGet( "FTF_BARYON_PROC0_YMIN", fProc0Ymin );
1220 /* */
1221 //
1222 /* JVY, June 11, 2020: make configurable
1223 fProc0A1 = 13.71;
1224 fProc0B1 = 1.75;
1225 fProc0A2 = -30.69; // (or -214.5 as in Doc ?)
1226 fProc0B2 = 3.; // ( or 4. as in Doc ?)
1227 fProc0A3 = 0.;
1228 fProc0Atop = 1.; // ( or 0.5 as in Doc ?)
1229 fProc0Ymin = 0.93; // (or 1.1 as in Doc ?)
1230 */
1231 //
1232 // Proc=1 --> Qexchg w/excitation
1233 //
1234 /* As of Oct. 31, 2017 keep these fixed */
1235 HDP.DeveloperGet( "FTF_BARYON_PROC1_A1", fProc1A1 );
1236 HDP.DeveloperGet( "FTF_BARYON_PROC1_B1", fProc1B1 );
1237 HDP.DeveloperGet( "FTF_BARYON_PROC1_A2", fProc1A2 );
1238 HDP.DeveloperGet( "FTF_BARYON_PROC1_B2", fProc1B2 );
1239 HDP.DeveloperGet( "FTF_BARYON_PROC1_A3", fProc1A3 );
1240 HDP.DeveloperGet( "FTF_BARYON_PROC1_ATOP", fProc1Atop );
1241 HDP.DeveloperGet( "FTF_BARYON_PROC1_YMIN", fProc1Ymin );
1242 /* */
1243 //
1244 /* JVY, June 11, 2020: make configurable
1245 fProc1A1 = 25.;
1246 fProc1B1 = 1.;
1247 fProc1A2 = -50.34;
1248 fProc1B2 = 1.5;
1249 fProc1A3 = 0.;
1250 fProc1Atop = 0.;
1251 fProc1Ymin = 1.4;
1252 */
1253 //
1254 // Proc=2 & Proc=3 for the case ( AbsProjectileBaryonNumber > 10 || NumberOfTargetNucleons > 10 )
1255 // (diffraction dissociation)
1256 // NOTE-1: used to be ( AbsProjectileBaryonNumber > 1 || NumberOfTargetNucleons > 1 )...
1257 // NOTE-2: As of 10.5, both are set to false (via HDP)
1258 //
1259 HDP.DeveloperGet( "FTF_BARYON_DIFF_DISSO_PROJ", fProjDiffDissociation );
1260 HDP.DeveloperGet( "FTF_BARYON_DIFF_DISSO_TGT", fTgtDiffDissociation );
1261 //
1262 //
1263 // Proc=4 --> Qexchg "w/additional multiplier" in excitation
1264 //
1265 /* As of Oct. 31, 2017 keep these fixed */
1266 /* */
1267 HDP.DeveloperGet( "FTF_BARYON_PROC4_A1", fProc4A1 );
1268 HDP.DeveloperGet( "FTF_BARYON_PROC4_B1", fProc4B1 );
1269 HDP.DeveloperGet( "FTF_BARYON_PROC4_A2", fProc4A2 );
1270 HDP.DeveloperGet( "FTF_BARYON_PROC4_B2", fProc4B2 );
1271 HDP.DeveloperGet( "FTF_BARYON_PROC4_A3", fProc4A3 );
1272 HDP.DeveloperGet( "FTF_BARYON_PROC4_ATOP", fProc4Atop );
1273 HDP.DeveloperGet( "FTF_BARYON_PROC4_YMIN", fProc4Ymin );
1274 /* */
1275 //
1276 /* JVY, June 11, 2020: make configurable
1277 fProc4A1 = 0.6; // (or 1. as in Doc ?)
1278 fProc4B1 = 0.;
1279 fProc4A2 = -1.2; // (or -2.01 as in Doc ?)
1280 fProc4B2 = 0.5;
1281 fProc4A3 = 0.;
1282 fProc4Atop = 0.;
1283 fProc4Ymin = 1.4;
1284 */
1285 //
1286 //
1287 HDP.DeveloperGet( "FTF_BARYON_DELTA_PROB_QEXCHG", fDeltaProbAtQuarkExchange );
1288 HDP.DeveloperGet( "FTF_BARYON_PROB_SAME_QEXCHG", fProbOfSameQuarkExchange );
1289 HDP.DeveloperGet( "FTF_BARYON_DIFF_M_PROJ", fProjMinDiffMass );
1290 HDP.DeveloperGet( "FTF_BARYON_NONDIFF_M_PROJ", fProjMinNonDiffMass );
1291 HDP.DeveloperGet( "FTF_BARYON_DIFF_M_TGT", fTgtMinDiffMass );
1292 HDP.DeveloperGet( "FTF_BARYON_NONDIFF_M_TGT", fTgtMinNonDiffMass );
1293 HDP.DeveloperGet( "FTF_BARYON_AVRG_PT2", fAveragePt2 );
1294 //
1295 // JVY - Per Alberto R., we're curretly keeping these two settings fixed,
1296 // thus they're defined here explicitly, rather than via HDP
1297 //
1298 // HDP.DeveloperGet( "FTF_BARYON_PROB_DISTR_PROJ", fProbLogDistrPrD );
1299 // HDP.DeveloperGet( "FTF_BARYON_PROB_DISTR_TGT", fProbLogDistr );
1300 fProbLogDistrPrD = 0.55;
1301 fProbLogDistr = 0.55;
1302
1303 // nuclear destruction
1304 //
1305 // ---> LATER !!! ---> fBaryonMaxNumberOfCollisions = 2.;
1306 //
1307
1308 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_P1_PROJ", fNuclearProjDestructP1 );
1309 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_P1_NBRN_PROJ",fNuclearProjDestructP1_NBRNDEP );
1310 //
1311 // keep the 2 parameters below fixed for now (i.e. do not take them from HDP)
1312 //
1315
1316 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_P1_TGT", fNuclearTgtDestructP1 );
1317 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_P1_ADEP_TGT", fNuclearTgtDestructP1_ADEP );
1318 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_P2_TGT", fNuclearTgtDestructP2 );
1319 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_P3_TGT", fNuclearTgtDestructP3 );
1320 //
1321 HDP.DeveloperGet( "FTF_BARYON_PT2_NUCDESTR_P1", fPt2NuclearDestructP1 );
1322 HDP.DeveloperGet( "FTF_BARYON_PT2_NUCDESTR_P2", fPt2NuclearDestructP2 );
1323 HDP.DeveloperGet( "FTF_BARYON_PT2_NUCDESTR_P3", fPt2NuclearDestructP3 );
1324 HDP.DeveloperGet( "FTF_BARYON_PT2_NUCDESTR_P4", fPt2NuclearDestructP4 );
1325 //
1326 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_R2", fR2ofNuclearDestruct );
1327 HDP.DeveloperGet( "FTF_BARYON_EXCI_E_PER_WNDNUCLN", fExciEnergyPerWoundedNucleon );
1328 //
1329 HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_DISP", fDofNuclearDestruct ); // NOTE: "Dof" means "Dispersion of..."
1330 //
1331 // NOTE-1: this parameter has changed from 1. to 9. between 10.2 and 10.3.ref07 !!!
1332 // ... then it went back to 1. for the 10.4-candidate...
1333 // NOTE-2: this is a "technical" parameter, it should not be changed; this is why
1334 // it is defined explicitly rather than via HDP
1335 // --> HDP.DeveloperGet( "FTF_BARYON_NUCDESTR_MAXPT2", fMaxPt2ofNuclearDestruct );
1336 fMaxPt2ofNuclearDestruct = 9. * CLHEP::GeV*CLHEP::GeV;
1337}
1338
1339
1342{
1343 // nuclear destruction
1344 //
1345 // JVY, Aug.8, 2018 --> Feb.14, 2018 --> June 25, 2019:
1346 // These parameters are common for all mesons
1347 //
1348
1349 HDP.DeveloperGet( "FTF_MESON_NUCDESTR_P1_TGT", fNuclearTgtDestructP1 );
1350 HDP.DeveloperGet( "FTF_MESON_NUCDESTR_P1_ADEP_TGT", fNuclearTgtDestructP1_ADEP );
1351 HDP.DeveloperGet( "FTF_MESON_NUCDESTR_P2_TGT", fNuclearTgtDestructP2 );
1352 HDP.DeveloperGet( "FTF_MESON_NUCDESTR_P3_TGT", fNuclearTgtDestructP3 );
1353 //
1354 HDP.DeveloperGet( "FTF_MESON_PT2_NUCDESTR_P1", fPt2NuclearDestructP1 );
1355 HDP.DeveloperGet( "FTF_MESON_PT2_NUCDESTR_P2", fPt2NuclearDestructP2 );
1356 HDP.DeveloperGet( "FTF_MESON_PT2_NUCDESTR_P3", fPt2NuclearDestructP3 );
1357 HDP.DeveloperGet( "FTF_MESON_PT2_NUCDESTR_P4", fPt2NuclearDestructP4 );
1358 //
1359 HDP.DeveloperGet( "FTF_MESON_NUCDESTR_R2", fR2ofNuclearDestruct );
1360 HDP.DeveloperGet( "FTF_MESON_EXCI_E_PER_WNDNUCLN", fExciEnergyPerWoundedNucleon );
1361 HDP.DeveloperGet( "FTF_MESON_NUCDESTR_DISP", fDofNuclearDestruct ); // NOTE: "Dof" means "Dispersion of..."
1362 //
1363 // NOTE: it is a "technical" parameter, it should not be changed;
1364 // this is why it is defined explicitly rather than via HDP
1365 //
1366 fMaxPt2ofNuclearDestruct = 1. * CLHEP::GeV*CLHEP::GeV;
1367}
1368
1369
1372{
1373 // parameters of participating pion excitation (pi+/- or pi0)
1374 //
1375 // Proc=0 --> Qexchg w/o excitation
1376 //
1377 /* As of June 25, 2019 keep these fixed */
1378 HDP.DeveloperGet( "FTF_PION_PROC0_A1", fProc0A1 );
1379 HDP.DeveloperGet( "FTF_PION_PROC0_B1", fProc0B1 );
1380 HDP.DeveloperGet( "FTF_PION_PROC0_A2", fProc0A2 );
1381 HDP.DeveloperGet( "FTF_PION_PROC0_B2", fProc0B2 );
1382 HDP.DeveloperGet( "FTF_PION_PROC0_A3", fProc0A3 );
1383 HDP.DeveloperGet( "FTF_PION_PROC0_ATOP", fProc0Atop );
1384 HDP.DeveloperGet( "FTF_PION_PROC0_YMIN", fProc0Ymin );
1385 /* */
1386 //
1387 /* JVY, June 11, 2020: make configurable
1388 fProc0A1 = 150.0;
1389 fProc0B1 = 1.8;
1390 fProc0A2 =-247.3;
1391 fProc0B2 = 2.3; // ( or 4. as in Doc ?)
1392 fProc0A3 = 0.;
1393 fProc0Atop = 1.; // ( or 0.5 as in Doc ?)
1394 fProc0Ymin = 2.3; // (or 1.1 as in Doc ?)
1395 */
1396 //
1397 // Proc=1 --> Qexchg w/excitation
1398 //
1399 /* As of Oct. 31, 2017 keep these fixed */
1400 HDP.DeveloperGet( "FTF_PION_PROC1_A1", fProc1A1 );
1401 HDP.DeveloperGet( "FTF_PION_PROC1_B1", fProc1B1 );
1402 HDP.DeveloperGet( "FTF_PION_PROC1_A2", fProc1A2 );
1403 HDP.DeveloperGet( "FTF_PION_PROC1_B2", fProc1B2 );
1404 HDP.DeveloperGet( "FTF_PION_PROC1_A3", fProc1A3 );
1405 HDP.DeveloperGet( "FTF_PION_PROC1_ATOP", fProc1Atop );
1406 HDP.DeveloperGet( "FTF_PION_PROC1_YMIN", fProc1Ymin );
1407 /* */
1408 //
1409 /* JVY, June 11, 2020: make configurable
1410 fProc1A1 = 5.77;
1411 fProc1B1 = 0.6;
1412 fProc1A2 = -5.77;
1413 fProc1B2 = 0.8;
1414 fProc1A3 = 0.;
1415 fProc1Atop = 0.;
1416 fProc1Ymin = 0.0;
1417 */
1418 //
1419 // Proc=2 --> Projectile diffraction
1420 //
1421 // Q: Would it even make sense to make these configurable ?
1422 // The following is hadrcoded:
1423 // Projectile Baryon Number > 10 (AbsProjectileBaryonNumber > 10)
1424 // ... which is "strange" because projectile is a pion !!!... so it's always OFF
1425 // (see also lines 1007-1016)
1426 //
1427 /* As of Oct. 31, 2017 keep these fixed
1428 HDP.DeveloperGet( "FTF_PION_PROC2_A1", fProc2A1 );
1429 HDP.DeveloperGet( "FTF_PION_PROC2_B1", fProc2B1 );
1430 HDP.DeveloperGet( "FTF_PION_PROC2_A2", fProc2A2 );
1431 HDP.DeveloperGet( "FTF_PION_PROC2_B2", fProc2B2 );
1432 HDP.DeveloperGet( "FTF_PION_PROC2_A3", fProc2A3 );
1433 HDP.DeveloperGet( "FTF_PION_PROC2_ATOP", fProc2Atop );
1434 HDP.DeveloperGet( "FTF_PION_PROC2_YMIN", fProc2Ymin );
1435 */
1436 // keep fixed so far; see note above
1437 fProc2A1 = 2.27;
1438 fProc2B1 = 0.5;
1439 fProc2A2 =-98052.0;
1440 fProc2B2 = 4.0;
1441 fProc2A3 = 0.;
1442 fProc2Atop = 0.;
1443 fProc2Ymin = 3.0;
1444 //
1445 // Proc=3 --> Target diffraction
1446 //
1447 /* As of Oct. 31, 2017 keep these fixed */
1448 HDP.DeveloperGet( "FTF_PION_PROC3_A1", fProc3A1 );
1449 HDP.DeveloperGet( "FTF_PION_PROC3_B1", fProc3B1 );
1450 HDP.DeveloperGet( "FTF_PION_PROC3_A2", fProc3A2 );
1451 HDP.DeveloperGet( "FTF_PION_PROC3_B2", fProc3B2 );
1452 HDP.DeveloperGet( "FTF_PION_PROC3_A3", fProc3A3 );
1453 HDP.DeveloperGet( "FTF_PION_PROC3_ATOP", fProc3Atop );
1454 HDP.DeveloperGet( "FTF_PION_PROC3_YMIN", fProc3Ymin );
1455 /* */
1456 //
1457 /* JVY, June 11, 2020: make configurable
1458 fProc3A1 = 7.0;
1459 fProc3B1 = 0.9;
1460 fProc3A2 = -85.28;
1461 fProc3B2 = 1.9;
1462 fProc3A3 = 0.08;
1463 fProc3Atop = 0.;
1464 fProc3Ymin = 2.2;
1465 */
1466 //
1467 // for Proc2 & Proc3, pprojectile or target diffraction can be turned ON/OFF
1468 // if num.baryons >10 (which is strange for projectile which is pion !!!)
1469 //
1470 HDP.DeveloperGet( "FTF_PION_DIFF_DISSO_PROJ", fProjDiffDissociation );
1471 HDP.DeveloperGet( "FTF_PION_DIFF_DISSO_TGT", fTgtDiffDissociation );
1472 //
1473 // Proc=4 --> Qexchg "w/additional multiplier" in excitation
1474 //
1475 /* As of Oct. 31, 2017 keep these fixed */
1476 HDP.DeveloperGet( "FTF_PION_PROC4_A1", fProc4A1 );
1477 HDP.DeveloperGet( "FTF_PION_PROC4_B1", fProc4B1 );
1478 HDP.DeveloperGet( "FTF_PION_PROC4_A2", fProc4A2 );
1479 HDP.DeveloperGet( "FTF_PION_PROC4_B2", fProc4B2 );
1480 HDP.DeveloperGet( "FTF_PION_PROC4_A3", fProc4A3 );
1481 HDP.DeveloperGet( "FTF_PION_PROC4_ATOP", fProc4Atop );
1482 HDP.DeveloperGet( "FTF_PION_PROC4_YMIN", fProc4Ymin );
1483 /* */
1484 //
1485 /* JVY, June 11, 2020: make configurable
1486 fProc4A1 = 1.0;
1487 fProc4B1 = 0.;
1488 fProc4A2 = -11.02;
1489 fProc4B2 = 1.0;
1490 fProc4A3 = 0.;
1491 fProc4Atop = 0.;
1492 fProc4Ymin = 2.4;
1493 */
1494 //
1495 //
1496 HDP.DeveloperGet( "FTF_PION_DELTA_PROB_QEXCHG", fDeltaProbAtQuarkExchange );
1497 HDP.DeveloperGet( "FTF_PION_DIFF_M_PROJ", fProjMinDiffMass );
1498 HDP.DeveloperGet( "FTF_PION_NONDIFF_M_PROJ", fProjMinNonDiffMass );
1499 HDP.DeveloperGet( "FTF_PION_DIFF_M_TGT", fTgtMinDiffMass );
1500 HDP.DeveloperGet( "FTF_PION_NONDIFF_M_TGT", fTgtMinNonDiffMass );
1501 HDP.DeveloperGet( "FTF_PION_AVRG_PT2", fAveragePt2 );
1502 //
1503 fProbOfSameQuarkExchange = 0.; // This does NOT seem to apply to the pion case
1504 //
1505 // currently keep these two parameters fixed
1506 // thus they're defined here explicitly, rather than via HDP
1507 //
1508 fProbLogDistrPrD = 0.55;
1509 fProbLogDistr = 0.55;
1510}
1511
1512
1513//============================================================================
1514
1516 if ( StringMass ) delete StringMass;
1517}
1518
1519//============================================================================
1520
1521void G4FTFParameters::Reset()
1522{
1523 FTFhNcmsEnergy = 0.0;
1524 FTFXtotal = 0.0;
1525 FTFXelastic = 0.0;
1526 FTFXinelastic = 0.0;
1527 FTFXannihilation = 0.0;
1531 FTFSlope = 0.0;
1533 FTFGamma0 = 0.0;
1536 ProjMinDiffMass = 0.0;
1537 ProjMinNonDiffMass = 0.0;
1538 ProbLogDistrPrD = 0.0;
1539 TarMinDiffMass = 0.0;
1540 TarMinNonDiffMass = 0.0;
1541 AveragePt2 = 0.0;
1542 ProbLogDistr = 0.0;
1543 Pt2kink = 0.0;
1544 MaxNumberOfCollisions = 0.0;
1545 ProbOfInelInteraction = 0.0;
1550 DofNuclearDestruction = 0.0;
1553
1554 for ( G4int i = 0; i < 4; i++ ) {
1555 for ( G4int j = 0; j < 7; j++ ) {
1556 ProcParams[i][j] = 0.0;
1557 }
1558 }
1559
1560 return;
1561}
1562
double S(double temp)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4HadronicDeveloperParameters & HDP
G4FTFSettingDefaultHDP FTFDefaultsHDP
G4double G4Log(G4double x)
Definition: G4Log.hh:226
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4VComponentCrossSection * GetComponentCrossSection(const G4String &name)
static G4CrossSectionDataSetRegistry * Instance()
double GetProc3B2() const
double GetProc4A3() const
double GetNuclearTgtDestructP3() const
double GetProc3B1() const
bool IsNuclearProjDestructP1_NBRNDEP() const
double GetProc0A2() const
double GetExciEnergyPerWoundedNucleon() const
double GetProc1Atop() const
double GetProc2A3() const
double GetNuclearProjDestructP2() const
double GetProbLogDistrPrD() const
double GetPt2NuclearDestructP4() const
double GetPt2NuclearDestructP3() const
double GetNuclearProjDestructP3() const
double GetProc4A2() const
double GetPt2NuclearDestructP1() const
double GetProc2Ymin() const
double GetProjMinNonDiffMass() const
double GetProc3Ymin() const
double GetProc2A2() const
double GetProc4B1() const
double GetProc1A1() const
double GetProc0B2() const
double GetProbLogDistr() const
double GetProc1B2() const
double GetProc3Atop() const
double GetProc1A3() const
double GetProc1B1() const
double GetProc1Ymin() const
double GetProc2B1() const
double GetNuclearTgtDestructP2() const
double GetPt2NuclearDestructP2() const
double GetR2ofNuclearDestruct() const
double GetDeltaProbAtQuarkExchange() const
double GetDofNuclearDestruct() const
double GetProc0A3() const
bool IsProjDiffDissociation() const
double GetProc0Atop() const
double GetProc2B2() const
double GetProc2A1() const
double GetProc2Atop() const
double GetProc1A2() const
double GetProc4Atop() const
double GetProc3A1() const
double GetAveragePt2() const
double GetProc4Ymin() const
double GetProc0B1() const
double GetTgtMinNonDiffMass() const
bool IsTgtDiffDissociation() const
double GetProbOfSameQuarkExchange() const
double GetProc4B2() const
double GetProc0Ymin() const
double GetProc3A2() const
double GetProc3A3() const
bool IsNuclearTgtDestructP1_ADEP() const
double GetMaxPt2ofNuclearDestruct() const
double GetProc4A1() const
double GetTgtMinDiffMass() const
double GetNuclearTgtDestructP1() const
double GetProjMinDiffMass() const
double GetProc0A1() const
double GetNuclearProjDestructP1() const
void SetTarMinNonDiffMass(const G4double aValue)
G4double R2ofNuclearDestruction
void SetProjMinDiffMass(const G4double aValue)
G4double GetProbLogDistrPrD()
G4double GetCofNuclearDestructionPr()
void SetTotalCrossSection(const G4double Xtotal)
void SetExcitationEnergyPerWoundedNucleon(const G4double aValue)
G4double ExcitationEnergyPerWoundedNucleon
G4double GetAveragePt2()
void SetProbabilityOfElasticScatt(const G4double Xtotal, const G4double Xelastic)
void SetElastisCrossSection(const G4double Xelastic)
void SetQuarkProbabilitiesAtGluonSplitUp(const G4double Puubar, const G4double Pddbar, const G4double Pssbar)
G4double GetPt2ofNuclearDestruction()
G4double GetProbLogDistr()
void SetMaxPt2ofNuclearDestruction(const G4double aValue)
G4double GetProjMinNonDiffMass()
G4double ProcParams[5][7]
G4double ProbOfSameQuarkExchange
G4double GetAvaragePt2ofElasticScattering()
G4double TarMinNonDiffMass
G4double ProbOfInelInteraction
void SetPt2Kink(const G4double aValue)
void SetSlope(const G4double Slope)
G4double GetTarMinDiffMass()
G4double RadiusOfHNinteractions2
void SetDeltaProbAtQuarkExchange(const G4double aValue)
G4double MaxPt2ofNuclearDestruction
G4double GetTarMinNonDiffMass()
void SetAvaragePt2ofElasticScattering(const G4double aPt2)
G4double DeltaProbAtQuarkExchange
G4double GetProcProb(const G4int ProcN, const G4double y)
G4double CofNuclearDestruction
void SetTarMinDiffMass(const G4double aValue)
void SetProjMinNonDiffMass(const G4double aValue)
void SetCofNuclearDestructionPr(const G4double aValue)
void SetR2ofNuclearDestruction(const G4double aValue)
void SetProbLogDistrPrD(const G4double aValue)
G4double CofNuclearDestructionPr
void SetGamma0(const G4double Gamma0)
void SetAveragePt2(const G4double aValue)
G4double GetDeltaProbAtQuarkExchange()
G4double DofNuclearDestruction
G4double GetExcitationEnergyPerWoundedNucleon()
void SetProbabilityOfAnnihilation(const G4double aValue)
void InitForInteraction(const G4ParticleDefinition *, G4int theA, G4int theZ, G4double s)
void SetDofNuclearDestruction(const G4double aValue)
G4double ProbabilityOfAnnihilation
G4double GetProjMinDiffMass()
G4double GetProbOfSameQuarkExchange()
G4double ProbabilityOfElasticScatt
G4double GetDofNuclearDestruction()
void SetRadiusOfHNinteractions2(const G4double Radius2)
void SetParams(const G4int ProcN, const G4double A1, const G4double B1, const G4double A2, const G4double B2, const G4double A3, const G4double Atop, const G4double Ymin)
G4double Pt2ofNuclearDestruction
G4double ProbLogDistrPrD
void SetPt2ofNuclearDestruction(const G4double aValue)
G4double GetR2ofNuclearDestruction()
G4double GetCofNuclearDestruction()
void SetProbLogDistr(const G4double aValue)
void SetCofNuclearDestruction(const G4double aValue)
G4double ProjMinNonDiffMass
void SetMaxNumberOfCollisions(const G4double aValue, const G4double bValue)
G4double ProjMinDiffMass
G4double AvaragePt2ofElasticScattering
void SetProbOfSameQuarkExchange(const G4double aValue)
G4double MaxNumberOfCollisions
G4double FTFXannihilation
void SetInelasticCrossSection(const G4double Xinelastic)
G4bool DeveloperGet(const std::string name, G4bool &value)
static G4HadronicDeveloperParameters & GetInstance()
G4bool SetDefault(const std::string name, const G4bool value)
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4double GetPDGCharge() const
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
static G4Proton * Proton()
Definition: G4Proton.cc:92
virtual G4double GetTotalIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)=0
virtual G4double GetElasticIsotopeCrossSection(const G4ParticleDefinition *, G4double kinEnergy, G4int, G4int)=0
T sqr(const T &x)
Definition: templates.hh:128