source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4QuasiFreeRatios.cc @ 968

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

update processes

File size: 39.8 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//
27// $Id: G4QuasiFreeRatios.cc,v 1.19 2008/03/21 21:44:39 dennis Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30//
31// G4 Physics class: G4QuasiFreeRatios for N+A elastic cross sections
32// Created: M.V. Kossov, CERN/ITEP(Moscow), 10-OCT-01
33// The last update: M.V. Kossov, CERN/ITEP (Moscow) 15-Oct-06
34//
35//================================================================================
36
37//#define debug
38//#define pdebug
39//#define ppdebug
40//#define nandebug
41
42#include "G4QuasiFreeRatios.hh"
43
44// initialisation of statics
45std::vector<G4double*> G4QuasiFreeRatios::vT; // Vector of pointers to LinTable in C++ heap
46std::vector<G4double*> G4QuasiFreeRatios::vL; // Vector of pointers to LogTable in C++ heap
47std::vector<std::pair<G4double,G4double>*> G4QuasiFreeRatios::vX; // ETPointers to LogTable
48
49G4QuasiFreeRatios::G4QuasiFreeRatios()
50{
51#ifdef pdebug
52                G4cout<<"***^^^*** G4QuasiFreeRatios singletone is created ***^^^***"<<G4endl;
53#endif
54}
55
56G4QuasiFreeRatios::~G4QuasiFreeRatios()
57{
58  std::vector<G4double*>::iterator pos;
59  for(pos=vT.begin(); pos<vT.end(); pos++)
60  { delete [] *pos; }
61  vT.clear();
62  for(pos=vL.begin(); pos<vL.end(); pos++)
63  { delete [] *pos; }
64  vL.clear();
65
66  std::vector<std::pair<G4double,G4double>*>::iterator pos2;
67  for(pos2=vX.begin(); pos2<vX.end(); pos2++)
68  { delete [] *pos2; }
69  vX.clear();
70}
71
72// Returns Pointer to the G4VQCrossSection class
73G4QuasiFreeRatios* G4QuasiFreeRatios::GetPointer()
74{
75  static G4QuasiFreeRatios theRatios;   // *** Static body of the QEl Cross Section ***
76  return &theRatios;
77}
78
79// Calculation of pair(QuasiFree/Inelastic,QuasiElastic/QuasiFree)
80std::pair<G4double,G4double> G4QuasiFreeRatios::GetRatios(G4double pIU, G4int pPDG,
81                                                           G4int tgZ,    G4int tgN)
82{
83                std::pair<G4double,G4double> ElTot=GetElTot(pIU, pPDG, tgZ, tgN); // mean hN El&Tot(IU)
84  G4double R=0.;
85  G4double QF2In=1.;               // Prototype of QuasiFree/Inelastic ratio for hN_tot
86  if(ElTot.second>0.)
87  {
88    R=ElTot.first/ElTot.second;    // Elastic/Total ratio (does not depend on units
89    QF2In=GetQF2IN_Ratio(ElTot.second/millibarn, tgZ+tgN); // QuasiFree/Inelastic ratio
90  }
91  return std::make_pair(QF2In,R);
92}
93
94// Calculatio QasiFree/Inelastic Ratio as a function of total hN cross-section (mb) and A
95G4double G4QuasiFreeRatios::GetQF2IN_Ratio(G4double s, G4int A)
96{
97  static const G4int    nps=150;        // Number of steps in the R(s) LinTable
98  static const G4int    mps=nps+1;      // Number of elements in the R(s) LinTable
99  static const G4double sma=150.;       // The first LinTabEl(s=0)=1., s>sma -> logTab
100  static const G4double ds=sma/nps;     // Step of the linear Table
101  static const G4int    nls=100;        // Number of steps in the R(lns) logTable
102  static const G4int    mls=nls+1;      // Number of elements in the R(lns) logTable
103  static const G4double lsi=5.;         // The min ln(s) logTabEl(s=148.4 < sma=150.)
104  static const G4double lsa=9.;         // The max ln(s) logTabEl(s=148.4 - 8103. mb)
105  static const G4double mi=std::exp(lsi);// The min s of logTabEl(~ 148.4 mb)
106  static const G4double ms=std::exp(lsa);// The max s of logTabEl(~ 8103. mb)
107  static const G4double dl=(lsa-lsi)/nls;// Step of the logarithmic Table
108  static const G4double edl=std::exp(dl);// Multiplication step of the logarithmic Table
109  static const G4double toler=.01;      // The tolarence mb defining the same cross-section
110  static G4double lastS=0.;             // The last sigma value for which R was calculated
111  static G4double lastR=0.;             // The last ratio R which was calculated
112  // Local Associative Data Base:
113  static std::vector<G4int>     vA;     // Vector of calculated A
114  static std::vector<G4double>  vH;     // Vector of max s initialized in the LinTable
115  static std::vector<G4int>     vN;     // Vector of topBin number initialized in LinTable
116  static std::vector<G4double>  vM;     // Vector of rel max ln(s) initialized in LogTable
117  static std::vector<G4int>     vK;     // Vector of topBin number initialized in LogTable
118  // Last values of the Associative Data Base:
119  static G4int     lastA=0;             // theLast of calculated A
120  static G4double  lastH=0.;            // theLast of max s initialized in the LinTable
121  static G4int     lastN=0;             // theLast of topBin number initialized in LinTable
122  static G4double  lastM=0.;            // theLast of rel max ln(s) initialized in LogTable
123  static G4int     lastK=0;             // theLast of topBin number initialized in LogTable
124  static G4double* lastT=0;             // theLast of pointer to LinTable in the C++ heap
125  static G4double* lastL=0;             // theLast of pointer to LogTable in the C++ heap
126  // LogTable is created only if necessary. The ratio R(s>8100 mb) = 0 for any nuclei
127  if(s<toler || A<2) return 1.;
128  if(s>ms) return 0.;
129  if(A>238)
130  {
131    G4cout<<"-Warning-G4QuasiFreeRatio::GetQF2IN_Ratio:A="<<A<<">238, return zero"<<G4endl;
132    return 0.;
133  }
134  G4int nDB=vA.size();                  // A number of nuclei already initialized in AMDB
135  //  if(nDB && lastA==A && std::fabs(s-lastS)<toler) return lastR;
136  if(nDB && lastA==A && s==lastS) return lastR;  // VI do not use tolerance
137  G4bool found=false;
138  G4int i=-1;
139                if(nDB) for (i=0; i<nDB; i++) if(A==vA[i]) // Sirch for this A in AMDB
140  {
141    found=true;                         // The A value is found
142    break;
143  }
144  if(!nDB || !found)                    // Create new line in the AMDB
145         {
146    lastA = A;
147#ifdef pdebug
148                                G4cout<<"G4QuasiFreeRatios::GetQF2IN_Ratio: NewT, A="<<A<<", nDB="<<nDB<<G4endl;
149#endif
150    lastT = new G4double[mps];          // Create the linear Table
151    lastN = static_cast<int>(s/ds)+1;   // MaxBin to be initialized
152    if(lastN>nps)
153    {
154      lastN=nps;
155      lastH=sma;
156    }
157    else lastH = lastN*ds;              // Calculate max initialized s for LinTab
158    G4double sv=0;
159    lastT[0]=1.;
160    for(G4int j=1; j<=lastN; j++)       // Calculate LogTab values
161    {
162      sv+=ds;
163      lastT[j]=CalcQF2IN_Ratio(sv,A);
164    }
165    if(s>sma)                           // Initialize the logarithmic Table
166    {
167#ifdef pdebug
168                                G4cout<<"G4QuasiFreeRatios::GetQF2IN_Ratio: NewL, A="<<A<<", nDB="<<nDB<<G4endl;
169#endif
170      lastL=new G4double[mls];          // Create the logarithmic Table
171      G4double ls=std::log(s);
172      lastK = static_cast<int>((ls-lsi)/dl)+1; // MaxBin to be initialized in LogTaB
173      if(lastK>nls)
174      {
175        lastK=nls;
176        lastM=lsa-lsi;
177      }
178      else lastM = lastK*dl;            // Calculate max initialized ln(s)-lsi for LogTab
179      sv=mi;
180      for(G4int j=0; j<=lastK; j++)     // Calculate LogTab values
181      {
182        lastL[j]=CalcQF2IN_Ratio(sv,A);
183               if(j!=lastK) sv*=edl;
184      }
185    }
186    else                                // LogTab is not initialized
187    {
188      lastL = 0;
189      lastK = 0;
190      lastM = 0.;
191    }
192    i++;                                // Make a new record to AMDB and position on it
193    vA.push_back(lastA);
194    vH.push_back(lastH);
195    vN.push_back(lastN);
196    vM.push_back(lastM);
197    vK.push_back(lastK);
198    vT.push_back(lastT);
199    vL.push_back(lastL);
200  }
201  else                                  // The A value was found in AMDB
202         {
203    lastA=vA[i];
204    lastH=vH[i];
205    lastN=vN[i];
206    lastM=vM[i];
207    lastK=vK[i];
208    lastT=vT[i];
209    lastL=vL[i];
210    if(s>lastM)                          // At least LinTab must be updated
211    {
212      G4int nextN=lastN+1;               // The next bin to be initialized
213      if(lastN<nps)
214      {
215        lastN = static_cast<int>(s/ds)+1;// MaxBin to be initialized
216        if(lastN>nps)
217        {
218          lastN=nps;
219          lastH=sma;
220        }
221        else lastH = lastN*ds;           // Calculate max initialized s for LinTab
222        G4double sv=lastM;
223        for(G4int j=nextN; j<=lastN; j++)// Calculate LogTab values
224        {
225          sv+=ds;
226          lastT[j]=CalcQF2IN_Ratio(sv,A);
227        }
228      } // End of LinTab update
229      if(lastN>=nextN)
230                                                {
231        vH[i]=lastH;
232        vN[i]=lastN;
233      }
234      G4int nextK=lastK+1;
235      if(s>sma && lastK<nls)             // LogTab must be updated
236                                                {
237        G4double sv=std::exp(lastM+lsi); // Define starting poit (lastM will be changed)
238        G4double ls=std::log(s);
239        lastK = static_cast<int>((ls-lsi)/dl)+1; // MaxBin to be initialized in LogTaB
240        if(lastK>nls)
241        {
242          lastK=nls;
243          lastM=lsa-lsi;
244        }
245        else lastM = lastK*dl;           // Calculate max initialized ln(s)-lsi for LogTab
246        for(G4int j=nextK; j<=lastK; j++)// Calculate LogTab values
247        {
248                 sv*=edl;
249          lastL[j]=CalcQF2IN_Ratio(sv,A);
250        }
251      } // End of LogTab update
252      if(lastK>=nextK)
253                                                {
254        vM[i]=lastM;
255        vK[i]=lastK;
256      }
257    }
258  }
259  // Now one can use tabeles to calculate the value
260  if(s<sma)                             // Use linear table
261                {
262    G4int n=static_cast<int>(s/ds);     // Low edge number of the bin
263    G4double d=s-n*ds;                  // Linear shift
264    G4double v=lastT[n];                // Base
265    lastR=v+d*(lastT[n+1]-v)/ds;        // Result
266  }
267                else                                  // Use log table
268                {
269    G4double ls=std::log(s)-lsi;        // ln(s)-l_min
270    G4int n=static_cast<int>(ls/dl);    // Low edge number of the bin
271    G4double d=ls-n*dl;                 // Log shift
272    G4double v=lastL[n];                // Base
273    lastR=v+d*(lastL[n+1]-v)/dl;        // Result
274  }
275  if(lastR<0.) lastR=0.;
276  if(lastR>1.) lastR=1.;
277  return lastR;
278} // End of CalcQF2IN_Ratio
279
280// Calculatio QasiFree/Inelastic Ratio as a function of total hN cross-section and A
281G4double G4QuasiFreeRatios::CalcQF2IN_Ratio(G4double s, G4int A)
282{
283  static const G4double C=1.246;
284                G4double s2=s*s;
285  G4double s4=s2*s2;
286                G4double ss=std::sqrt(std::sqrt(s));
287  G4double P=7.48e-5*s2/(1.+8.77e12/s4/s4/s2);
288  G4double E=.2644+.016/(1.+std::exp((29.54-s)/2.49));
289  G4double F=ss*.1526*std::exp(-s2*ss*.0000859);
290         return C*std::exp(-E*std::pow(G4double(A-1.),F))/std::pow(G4double(A),P);
291} // End of CalcQF2IN_Ratio
292
293// Calculatio pair(hN_el,hN_tot) (mb): p in GeV/c, index(PDG,F) (see FetchElTot)
294std::pair<G4double,G4double> G4QuasiFreeRatios::CalcElTot(G4double p, G4int I)
295{
296  // ---------> Each parameter set can have not more than nPoints=128 parameters
297  static const G4double lmi=3.5;       // min of (lnP-lmi)^2 parabola
298  static const G4double pbe=.0557;     // elastic (lnP-lmi)^2 parabola coefficient
299  static const G4double pbt=.3;        // total (lnP-lmi)^2 parabola coefficient
300  static const G4double pmi=.1;        // Below that fast LE calculation is made
301  static const G4double pma=1000.;     // Above that fast HE calculation is made
302  // ======================================================================================
303  G4double El=0.;                      // prototype of the elastic hN cross-section
304  G4double To=0.;                      // prototype of the total hN cross-section
305  if(p<=0.)
306  {
307    G4cout<<"-Warning-G4QuasiFreeRatios::CalcElTot: p="<<p<<" is zero or negative"<<G4endl;
308    return std::make_pair(El,To);
309  }
310  if     (!I)                          // pp/nn
311                {
312#ifdef debug
313                                G4cout<<"G4QuasiFreeR::CalcElTot:I=0, p="<<p<<", pmi="<<pmi<<", pma="<<pma<<G4endl;
314#endif
315    if(p<pmi)
316    {
317      G4double p2=p*p;
318      El=1./(.00012+p2*.2);
319      To=El;
320#ifdef debug
321                                  G4cout<<"G4QuasiFreeR::CalcElTot:I=0i, El="<<El<<", To="<<To<<", p2="<<p2<<G4endl;
322#endif
323    }
324    else if(p>pma)
325    {
326      G4double lp=std::log(p)-lmi;
327      G4double lp2=lp*lp;
328      El=pbe*lp2+6.72;
329      To=pbt*lp2+38.2;
330#ifdef debug
331                                  G4cout<<"G4QuasiFreeR::CalcElTot:I=0a, El="<<El<<", To="<<To<<", lp2="<<lp2<<G4endl;
332#endif
333    }
334    else
335    {
336      G4double p2=p*p;
337      G4double LE=1./(.00012+p2*.2);
338      G4double lp=std::log(p)-lmi;
339      G4double lp2=lp*lp;
340      G4double rp2=1./p2;
341      El=LE+(pbe*lp2+6.72+32.6/p)/(1.+rp2/p);
342      To=LE+(pbt*lp2+38.2+52.7*rp2)/(1.+2.72*rp2*rp2);
343#ifdef debug
344                                  G4cout<<"G4QuasiFreeR::CalcElTot:0,E="<<El<<",T="<<To<<",s="<<p2<<",l="<<lp2<<G4endl;
345#endif
346    }
347  }
348  else if(I==1)                        // np/pn
349                {
350    if(p<pmi)
351    {
352      G4double p2=p*p;
353      El=1./(.00012+p2*(.051+.1*p2));
354      To=El;
355    }
356    else if(p>pma)
357    {
358      G4double lp=std::log(p)-lmi;
359      G4double lp2=lp*lp;
360      El=pbe*lp2+6.72;
361      To=pbt*lp2+38.2;
362    }
363    else
364    {
365      G4double p2=p*p;
366      G4double LE=1./(.00012+p2*(.051+.1*p2));
367      G4double lp=std::log(p)-lmi;
368      G4double lp2=lp*lp;
369      G4double rp2=1./p2;
370      El=LE+(pbe*lp2+6.72+30./p)/(1.+.49*rp2/p);
371      To=LE+(pbt*lp2+38.2)/(1.+.54*rp2*rp2);
372    }
373  }
374  else if(I==2)                        // pimp/pipn
375                {
376    G4double lp=std::log(p);
377    if(p<pmi)
378    {
379      G4double lr=lp+1.27;
380      El=1.53/(lr*lr+.0676);
381      To=El*3;
382    }
383    else if(p>pma)
384    {
385      G4double ld=lp-lmi;
386      G4double ld2=ld*ld;
387      G4double sp=std::sqrt(p);
388      El=pbe*ld2+2.4+7./sp;
389      To=pbt*ld2+22.3+12./sp;
390    }
391    else
392    {
393      G4double lr=lp+1.27;
394      G4double LE=1.53/(lr*lr+.0676);
395      G4double ld=lp-lmi;
396      G4double ld2=ld*ld;
397      G4double p2=p*p;
398      G4double p4=p2*p2;
399      G4double sp=std::sqrt(p);
400      G4double lm=lp+.36;
401      G4double md=lm*lm+.04;
402      G4double lh=lp-.017;
403      G4double hd=lh*lh+.0025;
404      El=LE+(pbe*ld2+2.4+7./sp)/(1.+.7/p4)+.6/md+.05/hd;
405      To=LE*3+(pbt*ld2+22.3+12./sp)/(1.+.4/p4)+1./md+.06/hd;
406    }
407  }
408  else if(I==3)                        // pipp/pimn
409                {
410    G4double lp=std::log(p);
411    if(p<pmi)
412    {
413      G4double lr=lp+1.27;
414      G4double lr2=lr*lr;
415      El=13./(lr2+lr2*lr2+.0676);
416      To=El;
417    }
418    else if(p>pma)
419    {
420      G4double ld=lp-lmi;
421      G4double ld2=ld*ld;
422      G4double sp=std::sqrt(p);
423      El=pbe*ld2+2.4+6./sp;
424      To=pbt*ld2+22.3+5./sp;
425    }
426    else
427    {
428      G4double lr=lp+1.27;
429      G4double lr2=lr*lr;
430      G4double LE=13./(lr2+lr2*lr2+.0676);
431      G4double ld=lp-lmi;
432      G4double ld2=ld*ld;
433      G4double p2=p*p;
434      G4double p4=p2*p2;
435      G4double sp=std::sqrt(p);
436      G4double lm=lp-.32;
437      G4double md=lm*lm+.0576;
438      El=LE+(pbe*ld2+2.4+6./sp)/(1.+3./p4)+.7/md;
439      To=LE+(pbt*ld2+22.3+5./sp)/(1.+1./p4)+.8/md;
440    }
441  }
442                else if(I==4)                        // Kmp/Kmn/K0p/K0n
443                {
444
445    if(p<pmi)
446    {
447      G4double psp=p*std::sqrt(p);
448      El=5.2/psp;
449      To=14./psp;
450    }
451    else if(p>pma)
452    {
453      G4double ld=std::log(p)-lmi;
454      G4double ld2=ld*ld;
455      El=pbe*ld2+2.23;
456      To=pbt*ld2+19.5;
457    }
458    else
459    {
460      G4double ld=std::log(p)-lmi;
461      G4double ld2=ld*ld;
462      G4double sp=std::sqrt(p);
463      G4double psp=p*sp;
464      G4double p2=p*p;
465      G4double p4=p2*p2;
466      G4double lm=p-.39;
467      G4double md=lm*lm+.000156;
468      G4double lh=p-1.;
469      G4double hd=lh*lh+.0156;
470      El=5.2/psp+(pbe*ld2+2.23)/(1.-.7/sp+.075/p4)+.004/md+.15/hd;
471      To=14./psp+(pbt*ld2+19.5)/(1.-.21/sp+.52/p4)+.006/md+.30/hd;
472    }
473  }
474  else if(I==5)                        // Kpp/Kpn/aKp/aKn
475                {
476    if(p<pmi)
477    {
478      G4double lr=p-.38;
479      G4double lm=p-1.;
480      G4double md=lm*lm+.372;   
481      El=.7/(lr*lr+.0676)+2./md;
482      To=El+.6/md;
483    }
484    else if(p>pma)
485    {
486      G4double ld=std::log(p)-lmi;
487      G4double ld2=ld*ld;
488      El=pbe*ld2+2.23;
489      To=pbt*ld2+19.5;
490    }
491    else
492    {
493      G4double ld=std::log(p)-lmi;
494      G4double ld2=ld*ld;
495      G4double lr=p-.38;
496      G4double LE=.7/(lr*lr+.0676);
497      G4double sp=std::sqrt(p);
498      G4double p2=p*p;
499      G4double p4=p2*p2;
500      G4double lm=p-1.;
501      G4double md=lm*lm+.372;
502      El=LE+(pbe*ld2+2.23)/(1.-.7/sp+.1/p4)+2./md;
503      To=LE+(pbt*ld2+19.5)/(1.+.46/sp+1.6/p4)+2.6/md;
504    }
505  }
506  else if(I==6)                        // hyperon-N
507                {
508    if(p<pmi)
509    {
510      G4double p2=p*p;
511      El=1./(.002+p2*(.12+p2));
512      To=El;
513    }
514    else if(p>pma)
515    {
516      G4double lp=std::log(p)-lmi;
517      G4double lp2=lp*lp;
518      G4double sp=std::sqrt(p);
519      El=(pbe*lp2+6.72)/(1.+2./sp);
520      To=(pbt*lp2+38.2+900./sp)/(1.+27./sp);
521    }
522    else
523    {
524      G4double p2=p*p;
525      G4double LE=1./(.002+p2*(.12+p2));
526      G4double lp=std::log(p)-lmi;
527      G4double lp2=lp*lp;
528      G4double p4=p2*p2;
529      G4double sp=std::sqrt(p);
530      El=LE+(pbe*lp2+6.72+99./p2)/(1.+2./sp+2./p4);
531      To=LE+(pbt*lp2+38.2+900./sp)/(1.+27./sp+3./p4);
532    }
533  }
534  else if(I==7)                        // antibaryon-N
535                {
536    if(p>pma)
537    {
538      G4double lp=std::log(p)-lmi;
539      G4double lp2=lp*lp;
540      El=pbe*lp2+6.72;
541      To=pbt*lp2+38.2;
542    }
543    else
544    {
545      G4double ye=std::pow(p,1.25);
546      G4double yt=std::pow(p,.35);
547      G4double lp=std::log(p)-lmi;
548      G4double lp2=lp*lp;
549      El=80./(ye+1.)+pbe*lp2+6.72;
550      To=(80./yt+.3)/yt+pbt*lp2+38.2;
551    }
552  }
553  else
554  {
555    G4cout<<"*Error*G4QuasiFreeRatios::CalcElTot:ind="<<I<<" is not defined (0-7)"<<G4endl;
556    G4Exception("G4QuasiFreeRatios::CalcElTot:","23",FatalException,"CHIPScrash");
557  }
558  if(El>To) El=To;
559  return std::make_pair(El,To);
560} // End of CalcElTot
561
562// Calculatio pair(hN_el,hN_tot)(mb): p in GeV/c, F=true -> N=proton, F=false -> N=neutron
563std::pair<G4double,G4double> G4QuasiFreeRatios::FetchElTot(G4double p, G4int PDG, G4bool F)
564{
565  static const G4int    nlp=300;         // Number of steps in the S(lnp) logTable(5% step)
566  static const G4int    mlp=nlp+1;       // Number of elements in the S(lnp) logTable
567  static const G4double lpi=-5.;         // The min ln(p) logTabEl(p=6.7 MeV/c - 22. TeV/c)
568  static const G4double lpa=10.;         // The max ln(p) logTabEl(p=6.7 MeV/c - 22. TeV/c)
569  static const G4double mi=std::exp(lpi);// The min p of logTabEl(~ 6.7 MeV/c)
570  static const G4double ma=std::exp(lpa);// The max p of logTabEl(~ 22. TeV)
571  static const G4double dl=(lpa-lpi)/nlp;// Step of the logarithmic Table
572  static const G4double edl=std::exp(dl);// Multiplication step of the logarithmic Table
573  //  static const G4double toler=.001;      // Relative tolarence defining "the same momentum"
574  static G4double lastP=0.;              // The last momentum for which XS was calculated
575  static G4int    lastH=0;               // The last projPDG for which XS was calculated
576  static G4bool   lastF=true;            // The last nucleon for which XS was calculated
577  static std::pair<G4double,G4double> lastR=std::make_pair(0.,0.); // The last result
578  // Local Associative Data Base:
579  static std::vector<G4int>     vI;      // Vector of index for which XS was calculated
580  static std::vector<G4double>  vM;      // Vector of rel max ln(p) initialized in LogTable
581  static std::vector<G4int>     vK;      // Vector of topBin number initialized in LogTable
582  // Last values of the Associative Data Base:
583  static G4int     lastI=0;              // The Last index for which XS was calculated
584  static G4double  lastM=0.;             // The Last rel max ln(p) initialized in LogTable
585  static G4int     lastK=0;             // The Last topBin number initialized in LogTable
586  static std::pair<G4double,G4double>* lastX=0; // The Last ETPointers to LogTable in heap
587  // LogTable is created only if necessary. The ratio R(s>8100 mb) = 0 for any nuclei
588  G4int nDB=vI.size();                   // A number of hadrons already initialized in AMDB
589#ifdef pdebug
590  G4cout<<"G4QuasiFreeR::FetchElTot:p="<<p<<",PDG="<<PDG<<",F="<<F<<",nDB="<<nDB<<G4endl;
591#endif
592  if(nDB && lastH==PDG && lastF==F && p>0. && p==lastP) return lastR; // VI do not use toler
593  //  if(nDB && lastH==PDG && lastF==F && p>0. && std::fabs(p-lastP)/p<toler) return lastR;
594  lastH=PDG;
595  lastF=F;
596  G4int ind=-1;                          // Prototipe of the index of the PDG/F combination
597  // i=0: pp(nn), i=1: np(pn), i=2: pimp(pipn), i=3: pipp(pimn), i=4: Kmp(Kmn,K0n,K0p),
598                // i=5: Kpp(Kpn,aK0n,aK0p), i=6: Hp(Hn), i=7: app(apn,ann,anp)
599  G4bool kfl=true;                             // Flag of K0/aK0 oscillation
600  G4bool kf=false;
601  if(PDG==130||PDG==310)
602  {
603    kf=true;
604    if(G4UniformRand()>.5) kfl=false;
605  }
606  if ( (PDG == 2212 && F) || (PDG == 2112 && !F) ) {
607    ind=0; // pp/nn
608  } else if ( (PDG == 2112 && F) || (PDG == 2212 && !F) ) {
609    ind=1; // np/pn
610  } else if ( (PDG == -211 && F) || (PDG == 211 && !F) ) {
611    ind=2; // pimp/pipn
612  } else if ( (PDG == 211 && F) || (PDG == -211 && !F) ) {
613    ind=3; // pipp/pimn
614  } else if ( PDG == -321 || PDG == -311 || (kf && !kfl) ) {
615    ind=4; // KmN/K0N
616  } else if ( PDG == 321 || PDG == 311 || (kf && kfl) ) {
617    ind=5; // KpN/aK0N
618  } else if ( PDG > 3000 && PDG < 3335) {
619    ind=6;        // @@ for all hyperons - take Lambda
620  } else if ( PDG < -2000 && PDG > -3335 ) {
621    ind=7;        // @@ for all anti-baryons - anti-p/anti-n
622  } else {
623    G4cout<<"*Error*G4QuasiFreeRatios::FetchElTot: PDG="<<PDG
624          <<", while it is defined only for p,n,hyperons,anti-baryons,pi,K/antiK"<<G4endl;
625    G4Exception("G4QuasiFreeRatio::FetchElTot:","22",FatalException,"CHIPScrash");
626  }
627  if(nDB && lastI==ind && p>0. && p==lastP) return lastR;  // VI do not use toler
628  //  if(nDB && lastI==ind && p>0. && std::fabs(p-lastP)/p<toler) return lastR;
629  if(p<=mi || p>=ma) return CalcElTot(p,ind);   // @@ Slow calculation ! (Warning?)
630  G4bool found=false;
631  G4int i=-1;
632                if(nDB) for (i=0; i<nDB; i++) if(ind==vI[i])  // Sirch for this index in AMDB
633  {
634    found=true;                                 // The index is found
635    break;
636  }
637  G4double lp=std::log(p);
638#ifdef pdebug
639                G4cout<<"G4QuasiFreeR::FetchElTot:I="<<ind<<",i="<<i<<",fd="<<found<<",lp="<<lp<<G4endl;
640#endif
641  if(!nDB || !found)                            // Create new line in the AMDB
642         {
643#ifdef pdebug
644                                G4cout<<"G4QuasiFreeRatios::FetchElTot: NewX, ind="<<ind<<", nDB="<<nDB<<G4endl;
645#endif
646    lastX = new std::pair<G4double,G4double>[mlp]; // Create logarithmic Table for ElTot
647    lastI = ind;                                // Remember the initialized inex
648    lastK = static_cast<int>((lp-lpi)/dl)+1;    // MaxBin to be initialized in LogTaB
649    if(lastK>nlp)
650    {
651      lastK=nlp;
652      lastM=lpa-lpi;
653    }
654    else lastM = lastK*dl;               // Calculate max initialized ln(p)-lpi for LogTab
655    G4double pv=mi;
656    for(G4int j=0; j<=lastK; j++)        // Calculate LogTab values
657    {
658      lastX[j]=CalcElTot(pv,ind);
659#ifdef pdebug
660                    G4cout<<"G4QuasiFreeR::FetchElTot:I,j="<<j<<",pv="<<pv<<",E="<<lastX[j].first<<",T="
661            <<lastX[j].second<<G4endl;
662#endif
663      if(j!=lastK) pv*=edl;
664    }
665    i++;                                 // Make a new record to AMDB and position on it
666    vI.push_back(lastI);
667    vM.push_back(lastM);
668    vK.push_back(lastK);
669    vX.push_back(lastX);
670  }
671  else                                   // The A value was found in AMDB
672         {
673    lastI=vI[i];
674    lastM=vM[i];
675    lastK=vK[i];
676    lastX=vX[i];
677    G4int nextK=lastK+1;
678    G4double lpM=lastM+lpi;
679#ifdef pdebug
680                  G4cout<<"G4QuasiFreeR::FetchElTo:M="<<lpM<<",l="<<lp<<",K="<<lastK<<",n="<<nlp<<G4endl;
681#endif
682    if(lp>lpM && lastK<nlp)              // LogTab must be updated
683                                {
684      lastK = static_cast<int>((lp-lpi)/dl)+1; // MaxBin to be initialized in LogTab
685#ifdef pdebug
686                    G4cout<<"G4QuasiFreeR::FetET:K="<<lastK<<",lp="<<lp<<",li="<<lpi<<",dl="<<dl<<G4endl;
687#endif
688      if(lastK>nlp)
689      {
690        lastK=nlp;
691        lastM=lpa-lpi;
692      }
693      else lastM = lastK*dl;           // Calculate max initialized ln(p)-lpi for LogTab
694      G4double pv=std::exp(lpM);       // momentum of the last calculated beam
695      for(G4int j=nextK; j<=lastK; j++)// Calculate LogTab values
696      {
697        pv*=edl;
698        lastX[j]=CalcElTot(pv,ind);
699#ifdef pdebug
700                      G4cout<<"G4QuasiFreeR::FetchElTot:U:j="<<j<<",p="<<pv<<",E="<<lastX[j].first<<",T="
701              <<lastX[j].second<<G4endl;
702#endif
703      }
704    } // End of LogTab update
705    if(lastK>=nextK)                   // The AMDB was apdated
706                                {
707      vM[i]=lastM;
708      vK[i]=lastK;
709    }
710  }
711  // Now one can use tabeles to calculate the value
712  G4double dlp=lp-lpi;                       // Shifted log(p) value
713  G4int n=static_cast<int>(dlp/dl);          // Low edge number of the bin
714  G4double d=dlp-n*dl;                       // Log shift
715  G4double e=lastX[n].first;                 // E-Base
716  lastR.first=e+d*(lastX[n+1].first-e)/dl;   // E-Result
717  if(lastR.first<0.)  lastR.first = 0.;
718  G4double t=lastX[n].second;                // T-Base
719  lastR.second=t+d*(lastX[n+1].second-t)/dl; // T-Result
720  if(lastR.second<0.) lastR.second= 0.;
721  if(lastR.first>lastR.second) lastR.first = lastR.second;
722  return lastR;
723} // End of FetchElTot
724
725// (Mean Elastic and Mean Total) Cross-Sections (mb) for PDG+(Z,N) at P=p[GeV/c]
726std::pair<G4double,G4double> G4QuasiFreeRatios::GetElTot(G4double pIU, G4int hPDG,
727                                                         G4int Z,       G4int N)
728{
729  G4double pGeV=pIU/gigaelectronvolt;
730  if(Z<1 && N<1)
731  {
732    G4cout<<"-Warning-G4QuasiFreeRatio::GetElTot:Z="<<Z<<",N="<<N<<", return zero"<<G4endl;
733    return std::make_pair(0.,0.);
734  }
735  std::pair<G4double,G4double> hp=FetchElTot(pGeV, hPDG, true);
736  std::pair<G4double,G4double> hn=FetchElTot(pGeV, hPDG, false);
737  G4double A=(Z+N)/millibarn;                // To make the result in independent units(IU)
738  return std::make_pair((Z*hp.first+N*hn.first)/A,(Z*hp.second+N*hn.second)/A);
739} // End of GetElTot
740
741// (Mean Elastic and Mean Total) Cross-Sections (mb) for PDG+(Z,N) at P=p[GeV/c]
742std::pair<G4double,G4double> G4QuasiFreeRatios::GetChExFactor(G4double pIU, G4int hPDG,
743                                                              G4int Z, G4int N)
744{
745  G4double pGeV=pIU/gigaelectronvolt;
746  G4double resP=0.;
747  G4double resN=0.;
748  if(Z<1 && N<1)
749  {
750    G4cout<<"-Warning-G4QuasiFreeRatio::GetChExF:Z="<<Z<<",N="<<N<<", return zero"<<G4endl;
751    return std::make_pair(resP,resN);
752  }
753  G4double A=Z+N;
754  G4double pf=0.;                              // Possibility to interact with a proton
755  G4double nf=0.;                              // Possibility to interact with a neutron
756  if   (hPDG==-211||hPDG==-321||hPDG==3112||hPDG==3212||hPDG==3312) pf=Z/(A+N);
757  else if(hPDG==211||hPDG==321||hPDG==3222||hPDG==3212||hPDG==3322) nf=N/(A+Z);
758  else if(hPDG==-311||hPDG==311||hPDG==130||hPDG==310)
759  {
760    G4double dA=A+A;
761    pf=Z/(dA+N+N);
762    nf=N/(dA+Z+Z);
763  }
764  G4double mult=1.;  // Factor of increasing multiplicity ( ? @@)
765  if(pGeV>.5)
766  {
767    mult=1./(1.+std::log(pGeV+pGeV))/pGeV;
768    if(mult>1.) mult=1.;
769  }
770  if(pf)
771  {
772    std::pair<G4double,G4double> hp=FetchElTot(pGeV, hPDG, true);
773    resP=pf*(hp.second/hp.first-1.)*mult;
774  }
775  if(nf)
776  {
777    std::pair<G4double,G4double> hn=FetchElTot(pGeV, hPDG, false);
778    resN=nf*(hn.second/hn.first-1.)*mult;
779  }
780  return std::make_pair(resP,resN);
781} // End of GetChExFactor
782
783// scatter (pPDG,p4M) on a virtual nucleon (NPDG,N4M), result: final pair(newN4M,newp4M)
784// if(newN4M.e()==0.) - below threshold, XS=0, no scattering of the progectile happened
785std::pair<G4LorentzVector,G4LorentzVector> G4QuasiFreeRatios::Scatter(G4int NPDG,
786                                     G4LorentzVector N4M, G4int pPDG, G4LorentzVector p4M)
787{
788  static const G4double mNeut= G4QPDGCode(2112).GetMass();
789  static const G4double mProt= G4QPDGCode(2212).GetMass();
790  static const G4double mDeut= G4QPDGCode(2112).GetNuclMass(1,1,0);// Mass of deuteron
791  static const G4double mTrit= G4QPDGCode(2112).GetNuclMass(1,2,0);// Mass of tritium
792  static const G4double mHel3= G4QPDGCode(2112).GetNuclMass(2,1,0);// Mass of Helium3
793  static const G4double mAlph= G4QPDGCode(2112).GetNuclMass(2,2,0);// Mass of alpha
794  G4LorentzVector pr4M=p4M/megaelectronvolt;   // Convert 4-momenta in MeV (keep p4M)
795  N4M/=megaelectronvolt;
796  G4LorentzVector tot4M=N4M+p4M;
797  G4double mT=mNeut;
798  G4int Z=0;
799  G4int N=1;
800  if(NPDG==2212||NPDG==90001000)
801  {
802    mT=mProt;
803    Z=1;
804    N=0;
805  }
806  else if(NPDG==90001001)
807  {
808    mT=mDeut;
809    Z=1;
810    N=1;
811  }
812  else if(NPDG==90002001)
813  {
814    mT=mHel3;
815    Z=2;
816    N=1;
817  }
818  else if(NPDG==90001002)
819  {
820    mT=mTrit;
821    Z=1;
822    N=2;
823  }
824  else if(NPDG==90002002)
825  {
826    mT=mAlph;
827    Z=2;
828    N=2;
829  }
830  else if(NPDG!=2112&&NPDG!=90000001)
831  {
832    G4cout<<"Error:G4QuasiFreeRatios::Scatter:NPDG="<<NPDG<<" is not 2212 or 2112"<<G4endl;
833    G4Exception("G4QuasiFreeRatios::Scatter:","21",FatalException,"CHIPScomplain");
834    //return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M);// Use this if not exception
835  }
836  G4double mT2=mT*mT;
837  G4double mP2=pr4M.m2();
838  G4double E=(tot4M.m2()-mT2-mP2)/(mT+mT);
839#ifdef pdebug
840  G4cerr<<"G4QFR::Scat:qM="<<mT<<",qM2="<<mT2<<",pM2="<<mP2<<",totM2="<<tot4M.m2()<<G4endl;
841#endif
842  G4double E2=E*E;
843  if(E<0. || E2<mP2)
844  {
845#ifdef pdebug
846    G4cerr<<"-Warning-G4QFR::Scat:*Negative Energy*E="<<E<<",E2="<<E2<<"<M2="<<mP2<<G4endl;
847#endif
848    return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); // Do Nothing Action
849  }
850                G4double P=std::sqrt(E2-mP2);                   // Momentum in pseudo laboratory system
851  G4VQCrossSection* CSmanager=G4QElasticCrossSection::GetPointer();
852#ifdef debug
853  G4cout<<"G4QFR::Scatter: Before XS, P="<<P<<", Z="<<Z<<", N="<<N<<", PDG="<<pPDG<<G4endl;
854#endif
855  // @@ Temporary NN t-dependence for all hadrons
856  if(pPDG>3400 || pPDG<-3400) G4cout<<"-Warning-G4QElast::Scatter: pPDG="<<pPDG<<G4endl;
857  G4int PDG=2212;                                                // *TMP* instead of pPDG
858  if(pPDG==2112||pPDG==-211||pPDG==-321) PDG=2112;               // *TMP* instead of pPDG
859  G4double xSec=CSmanager->GetCrossSection(false, P, Z, N, PDG); // Rec.CrossSect *TMP*
860  //G4double xSec=CSmanager->GetCrossSection(false, P, Z, N, pPDG); // Rec.CrossSect
861#ifdef debug
862  G4cout<<"G4QElast::Scatter:pPDG="<<pPDG<<",P="<<P<<",CS="<<xSec/millibarn<<G4endl;
863#endif
864#ifdef nandebug
865  if(xSec>0. || xSec<0. || xSec==0);
866  else  G4cout<<"*Warning*G4QElast::Scatter: xSec="<<xSec/millibarn<<G4endl;
867#endif
868  // @@ check a possibility to separate p, n, or alpha (!)
869  if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing
870  {
871#ifdef pdebug
872    G4cerr<<"-Warning-G4QFR::Scat:**Zero XS**PDG="<<pPDG<<",NPDG="<<NPDG<<",P="<<P<<G4endl;
873#endif
874    return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); //Do Nothing Action
875  }
876  G4double mint=CSmanager->GetExchangeT(Z,N,PDG); // functional randomized -t (MeV^2) *TMP*
877  //G4double mint=CSmanager->GetExchangeT(Z,N,pPDG); // functional randomized -t in MeV^2
878  G4double maxt=CSmanager->GetHMaxT();            // max possible -t
879#ifdef ppdebug
880  G4cout<<"G4QFR::Scat:PDG="<<PDG<<",P="<<P<<",X="<<xSec<<",-t="<<mint<<"<"<<maxt<<", Z="<<Z<<",N="<<N<<G4endl;
881#endif
882#ifdef nandebug
883  if(mint>-.0000001);
884  else  G4cout<<"*Warning*G4QFR::Scat: -t="<<mint<<G4endl;
885#endif
886  G4double cost=1.-(mint+mint)/maxt; // cos(theta) in CMS
887#ifdef ppdebug
888  G4cout<<"G4QFR::Scat:-t="<<mint<<"<"<<maxt<<", cost="<<cost<<", Z="<<Z<<",N="<<N<<G4endl;
889#endif
890  if(cost>1. || cost<-1. || !(cost>-1. || cost<=1.))
891  {
892    if     (cost>1.)  cost=1.;
893    else if(cost<-1.) cost=-1.;
894    else
895                                {
896      G4cerr<<"G4QFR::S:*NAN*c="<<cost<<",t="<<mint<<",tm="<<CSmanager->GetHMaxT()<<G4endl;
897      return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); // Do Nothing Action
898    }
899  }
900  G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,mT);      // 4mom of the recoil nucleon
901  G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-mT)*.01);
902  if(!G4QHadron(tot4M).RelDecayIn2(pr4M, reco4M, dir4M, cost, cost))
903  {
904    G4cerr<<"G4QFR::Scat:t="<<tot4M<<tot4M.m()<<",mT="<<mT<<",mP="<<std::sqrt(mP2)<<G4endl;
905    //G4Exception("G4QFR::Scat:","009",FatalException,"Decay of ElasticComp");
906    return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); // Do Nothing Action
907  }
908#ifdef ppdebug
909                G4cout<<"G4QFR::Scat:p4M="<<pr4M<<"+r4M="<<reco4M<<",dr="<<dir4M<<",t4M="<<tot4M<<G4endl;
910#endif
911                return std::make_pair(reco4M*megaelectronvolt,pr4M*megaelectronvolt); // Result
912} // End of Scatter
913
914// scatter (pPDG,p4M) on a virtual nucleon (NPDG,N4M), result: final pair(newN4M,newp4M)
915// if(newN4M.e()==0.) - below threshold, XS=0, no scattering of the progectile happened
916// User should himself change the charge (PDG) (e.g. pn->np, pi+n->pi0p, pi-p->pi0n etc.)
917std::pair<G4LorentzVector,G4LorentzVector> G4QuasiFreeRatios::ChExer(G4int NPDG,
918                                     G4LorentzVector N4M, G4int pPDG, G4LorentzVector p4M)
919{
920  static const G4double mNeut= G4QPDGCode(2112).GetMass();
921  static const G4double mProt= G4QPDGCode(2212).GetMass();
922  G4LorentzVector pr4M=p4M/megaelectronvolt;          // Convert 4-momenta in MeV(keep p4M)
923  N4M/=megaelectronvolt;
924  G4LorentzVector tot4M=N4M+p4M;
925  G4int Z=0;
926  G4int N=1;
927  G4int RPDG=2212;                                     // PDG of the recoil nucleon
928  G4int sPDG=0;                                        // PDG code of the scattered hadron
929  G4double mS=0.;                                      // proto of mass of scattered hadron
930  G4double mT=mProt;                                   // mass of the recoil nucleon
931  if(NPDG==2212)
932  {
933    mT=mNeut;
934    Z=1;
935    N=0;
936    RPDG=2112;                                         // PDG of the recoil nucleon
937                if(pPDG==-211) sPDG=111;                           // pi+    -> pi0
938                else if(pPDG==-321)
939    {
940      sPDG=310;                                        // K+     -> K0S
941      if(G4UniformRand()>.5) sPDG=130;                 // K+     -> K0L
942    }
943    else if(pPDG==-311||pPDG==311||pPDG==130||pPDG==310) sPDG=321;  // K0     -> K+ (?)
944    else if(pPDG==3112) sPDG=3212;                     // Sigma- -> Sigma0
945    else if(pPDG==3212) sPDG=3222;                     // Sigma0 -> Sigma+
946    else if(pPDG==3312) sPDG=3322;                     // Xi-    -> Xi0
947  }
948  else if(NPDG==2112) // Default
949  {
950                if(pPDG==211)  sPDG=111;                           // pi+    -> pi0
951                else if(pPDG==321)
952    {
953      sPDG=310;                                        // K+     -> K0S
954      if(G4UniformRand()>.5) sPDG=130;                 // K+     -> K0L
955    }
956    else if(pPDG==-311||pPDG==311||pPDG==130||pPDG==310) sPDG=-321; // K0     -> K- (?)
957    else if(pPDG==3222) sPDG=3212;                     // Sigma+ -> Sigma0
958    else if(pPDG==3212) sPDG=3112;                     // Sigma0 -> Sigma-
959    else if(pPDG==3322) sPDG=3312;                     // Xi0    -> Xi-
960  }
961  else
962  {
963    G4cout<<"Error:G4QuasiFreeRatios::ChExer: NPDG="<<NPDG<<" is not 2212 or 2112"<<G4endl;
964    G4Exception("G4QuasiFreeRatios::ChExer:","21",FatalException,"CHIPS complain");
965    //return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M);// Use this if not exception
966  }
967  if(sPDG) mS=G4QPDGCode(2112).GetMass();
968  else
969  {
970    G4cout<<"Error:G4QuasiFreeRatios::ChExer: BAD pPDG="<<pPDG<<", NPDG="<<NPDG<<G4endl;
971    G4Exception("G4QuasiFreeRatios::ChExer:","21",FatalException,"CHIPS complain");
972    //return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M);// Use this if not exception
973  }
974  G4double mT2=mT*mT;
975  G4double mS2=mS*mS;
976  G4double E=(tot4M.m2()-mT2-mS2)/(mT+mT);
977  G4double E2=E*E;
978  if(E<0. || E2<mS2)
979  {
980#ifdef pdebug
981    G4cerr<<"-Warning-G4QFR::ChEx:*Negative Energy*E="<<E<<",E2="<<E2<<"<M2="<<mS2<<G4endl;
982#endif
983    return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); // Do Nothing Action
984  }
985                G4double P=std::sqrt(E2-mS2);                   // Momentum in pseudo laboratory system
986  G4VQCrossSection* CSmanager=G4QElasticCrossSection::GetPointer();
987#ifdef debug
988  G4cout<<"G4QFR::ChExer: Before XS, P="<<P<<", Z="<<Z<<", N="<<N<<", PDG="<<pPDG<<G4endl;
989#endif
990  // @@ Temporary NN t-dependence for all hadrons
991  G4int PDG=2212;                                                // *TMP* instead of pPDG
992  G4double xSec=CSmanager->GetCrossSection(false, P, Z, N, PDG); // Rec.CrossSect *TMP*
993  //G4double xSec=CSmanager->GetCrossSection(false, P, Z, N, pPDG); // Rec.CrossSect
994#ifdef debug
995  G4cout<<"G4QElast::ChExer:pPDG="<<pPDG<<",P="<<P<<",CS="<<xSec/millibarn<<G4endl;
996#endif
997#ifdef nandebug
998  if(xSec>0. || xSec<0. || xSec==0);
999  else  G4cout<<"*Warning*G4QElast::ChExer: xSec="<<xSec/millibarn<<G4endl;
1000#endif
1001  // @@ check a possibility to separate p, n, or alpha (!)
1002  if(xSec <= 0.) // The cross-section iz 0 -> Do Nothing
1003  {
1004#ifdef pdebug
1005    G4cerr<<"-Warning-G4QFR::ChEx:**Zero XS**PDG="<<pPDG<<",NPDG="<<NPDG<<",P="<<P<<G4endl;
1006#endif
1007    return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); //Do Nothing Action
1008  }
1009  G4double mint=CSmanager->GetExchangeT(Z,N,PDG); // functional randomized -t (MeV^2) *TMP*
1010  //G4double mint=CSmanager->GetExchangeT(Z,N,pPDG); // functional randomized -t in MeV^2
1011#ifdef pdebug
1012  G4cout<<"G4QFR::ChEx:PDG="<<pPDG<<", P="<<P<<", CS="<<xSec<<", -t="<<mint<<G4endl;
1013#endif
1014#ifdef nandebug
1015  if(mint>-.0000001);
1016  else  G4cout<<"*Warning*G4QFR::ChExer: -t="<<mint<<G4endl;
1017#endif
1018  G4double cost=1.-mint/CSmanager->GetHMaxT(); // cos(theta) in CMS
1019#ifdef pdebug
1020  G4cout<<"G4QFR::ChEx:-t="<<mint<<",dpc2="<<CSmanager->GetHMaxT()<<",cost="<<cost<<G4endl;
1021#endif
1022  if(cost>1. || cost<-1. || !(cost>-1. || cost<=1.))
1023  {
1024    if     (cost>1.)  cost=1.;
1025    else if(cost<-1.) cost=-1.;
1026    else
1027                                {
1028      G4cerr<<"G4QFR::C:*NAN*c="<<cost<<",t="<<mint<<",tm="<<CSmanager->GetHMaxT()<<G4endl;
1029      return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); // Do Nothing Action
1030    }
1031  }
1032  G4LorentzVector reco4M=G4LorentzVector(0.,0.,0.,mT);      // 4mom of the recoil nucleon
1033  pr4M=G4LorentzVector(0.,0.,0.,mS);                        // 4mom of the scattered hadron
1034  G4LorentzVector dir4M=tot4M-G4LorentzVector(0.,0.,0.,(tot4M.e()-mT)*.01);
1035  if(!G4QHadron(tot4M).RelDecayIn2(pr4M, reco4M, dir4M, cost, cost))
1036  {
1037    G4cerr<<"G4QFR::ChEx:t="<<tot4M<<tot4M.m()<<",mT="<<mT<<",mP="<<mS<<G4endl;
1038    //G4Exception("G4QFR::ChExer:","009",FatalException,"Decay of ElasticComp");
1039    return std::make_pair(G4LorentzVector(0.,0.,0.,0.),p4M); // Do Nothing Action
1040  }
1041#ifdef debug
1042                G4cout<<"G4QFR::ChEx:p4M="<<p4M<<"+r4M="<<reco4M<<"="<<p4M+reco4M<<"="<<tot4M<<G4endl;
1043#endif
1044                return std::make_pair(reco4M*megaelectronvolt,pr4M*megaelectronvolt); // Result
1045} // End of ChExer
1046
1047// Calculate ChEx/El ratio (p is in independent units, (Z,N) is target, pPDG is projectile)
1048G4double G4QuasiFreeRatios::ChExElCoef(G4double p, G4int Z, G4int N, G4int pPDG) 
1049{
1050  p/=MeV;                                // Converted from independent units
1051  G4double A=Z+N;
1052  if(A<1.5) return 0.;
1053  G4double C=0.;
1054  if     (pPDG==2212) C=N/(A+Z);
1055  else if(pPDG==2112) C=Z/(A+N);
1056  else G4cout<<"*Warning*G4QCohChrgExchange::ChExElCoef: wrong PDG="<<pPDG<<G4endl;
1057  C*=C;                         // Coherent processes squares the amplitude
1058  // @@ This is true only for nucleons: other projectiles must be treated differently
1059  G4double sp=std::sqrt(p);
1060  G4double p2=p*p;           
1061  G4double p4=p2*p2;
1062                G4double dl1=std::log(p)-5.;
1063  G4double T=(6.75+.14*dl1*dl1+13./p)/(1.+.14/p4)+.6/(p4+.00013);
1064  G4double U=(6.25+8.33e-5/p4/p)*(p*sp+.34)/p2/p; 
1065  G4double R=U/T;
1066  return C*R*R;
1067}
Note: See TracBrowser for help on using the repository browser.