source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/cross_sections/src/G4QElasticCrossSection.cc @ 1228

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

update geant4.9.3 tag

File size: 62.0 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: G4QElasticCrossSection.cc,v 1.1 2009/11/16 18:15:43 mkossov Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30//
31// G4 Physics class: G4QElasticCrossSection 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// Short description: Interaction cross-sections for the G4QElastic process
37// -------------------------------------------------------------------------------
38
39//#define debug
40//#define isodebug
41//#define pdebug
42//#define ppdebug
43//#define tdebug
44//#define sdebug
45
46#include "G4QElasticCrossSection.hh"
47
48// Initialization of the static parameters
49const G4int G4QElasticCrossSection::nPoints=128;//#of points in the AMDB tables(>anyPar)(D)
50const G4int G4QElasticCrossSection::nLast=nPoints-1; // the Last element in the table   (D)
51G4double  G4QElasticCrossSection::lPMin=-8.;  // Min tabulated logarithmic Momentum     (D)
52G4double  G4QElasticCrossSection::lPMax= 8.;  // Max tabulated logarithmic Momentum     (D)
53G4double  G4QElasticCrossSection::dlnP=(lPMax-lPMin)/nLast;// Log step in the table     (D)
54G4bool    G4QElasticCrossSection::onlyCS=true;// Flag to calculate only CS (not Si/Bi)  (L)
55G4double  G4QElasticCrossSection::lastSIG=0.; // Last calculated cross section          (L)
56G4double  G4QElasticCrossSection::lastLP=-10.;// Last log(mom_of_the_incident_hadron)   (L)
57G4double  G4QElasticCrossSection::lastTM=0.;  // Last t_maximum                         (L)
58G4double  G4QElasticCrossSection::theSS=0.;   // The Last sq.slope of first difruction  (L)
59G4double  G4QElasticCrossSection::theS1=0.;   // The Last mantissa of first difruction  (L)
60G4double  G4QElasticCrossSection::theB1=0.;   // The Last slope of first difruction     (L)
61G4double  G4QElasticCrossSection::theS2=0.;   // The Last mantissa of second difruction (L)
62G4double  G4QElasticCrossSection::theB2=0.;   // The Last slope of second difruction    (L)
63G4double  G4QElasticCrossSection::theS3=0.;   // The Last mantissa of third difruction  (L)
64G4double  G4QElasticCrossSection::theB3=0.;   // The Last slope of third difruction     (L)
65G4double  G4QElasticCrossSection::theS4=0.;   // The Last mantissa of 4-th difruction   (L)
66G4double  G4QElasticCrossSection::theB4=0.;   // The Last slope of 4-th difruction      (L)
67G4int     G4QElasticCrossSection::lastTZ=0;   // Last atomic number of the target   
68G4int     G4QElasticCrossSection::lastTN=0;   // Last number of neutrons of the target
69G4double  G4QElasticCrossSection::lastPIN=0.; // Last initialized max momentum
70G4double* G4QElasticCrossSection::lastCST=0;  // Elastic cross-section table               
71G4double* G4QElasticCrossSection::lastPAR=0;  // Parameters for functional calculation     
72G4double* G4QElasticCrossSection::lastSST=0;  // E-dep of sq.slope of the first difruction
73G4double* G4QElasticCrossSection::lastS1T=0;  // E-dep of mantissa of the first difruction
74G4double* G4QElasticCrossSection::lastB1T=0;  // E-dep of the slope of the first difruction
75G4double* G4QElasticCrossSection::lastS2T=0;  // E-dep of mantissa of the second difruction
76G4double* G4QElasticCrossSection::lastB2T=0;  // E-dep of the slope of theSecond difruction
77G4double* G4QElasticCrossSection::lastS3T=0;  // E-dep of mantissa of the third difruction
78G4double* G4QElasticCrossSection::lastB3T=0;  // E-dep of the slope of the third difruction
79G4double* G4QElasticCrossSection::lastS4T=0;  // E-dep of mantissa of the 4-th difruction
80G4double* G4QElasticCrossSection::lastB4T=0;  // E-dep of the slope of the 4-th difruction
81G4int     G4QElasticCrossSection::lastPDG=0;  // The last PDG code of the projectile
82G4int     G4QElasticCrossSection::lastN=0;    // The last N of calculated nucleus
83G4int     G4QElasticCrossSection::lastZ=0;    // The last Z of calculated nucleus
84G4double  G4QElasticCrossSection::lastP=0.;   // Last used in cross section Momentum
85G4double  G4QElasticCrossSection::lastTH=0.;  // Last threshold momentum
86G4double  G4QElasticCrossSection::lastCS=0.;  // Last value of the Cross Section
87G4int     G4QElasticCrossSection::lastI=0;    // The last position in the DAMDB
88
89std::vector<G4double*> G4QElasticCrossSection::PAR; // Vector of parameters for functCalcul
90std::vector<G4double*> G4QElasticCrossSection::CST; // Vector of cross-section table
91std::vector<G4double*> G4QElasticCrossSection::SST; // Vector of the first squared slope
92std::vector<G4double*> G4QElasticCrossSection::S1T; // Vector of the first mantissa
93std::vector<G4double*> G4QElasticCrossSection::B1T; // Vector of the first slope
94std::vector<G4double*> G4QElasticCrossSection::S2T; // Vector of the secon mantissa
95std::vector<G4double*> G4QElasticCrossSection::B2T; // Vector of the second slope
96std::vector<G4double*> G4QElasticCrossSection::S3T; // Vector of the third mantissa
97std::vector<G4double*> G4QElasticCrossSection::B3T; // Vector of the third slope
98std::vector<G4double*> G4QElasticCrossSection::S4T; // Vector of the 4-th mantissa (gloria)
99std::vector<G4double*> G4QElasticCrossSection::B4T; // Vector of the 4-th slope    (gloria)
100
101G4QElasticCrossSection::G4QElasticCrossSection()
102{
103}
104
105G4QElasticCrossSection::~G4QElasticCrossSection()
106{
107  std::vector<G4double*>::iterator pos;
108  for (pos=CST.begin(); pos<CST.end(); pos++)
109  { delete [] *pos; }
110  CST.clear();
111  for (pos=PAR.begin(); pos<PAR.end(); pos++)
112  { delete [] *pos; }
113  PAR.clear();
114  for (pos=SST.begin(); pos<SST.end(); pos++)
115  { delete [] *pos; }
116  SST.clear();
117  for (pos=S1T.begin(); pos<S1T.end(); pos++)
118  { delete [] *pos; }
119  S1T.clear();
120  for (pos=B1T.begin(); pos<B1T.end(); pos++)
121  { delete [] *pos; }
122  B1T.clear();
123  for (pos=S2T.begin(); pos<S2T.end(); pos++)
124  { delete [] *pos; }
125  S2T.clear();
126  for (pos=B2T.begin(); pos<B2T.end(); pos++)
127  { delete [] *pos; }
128  B2T.clear();
129  for (pos=S3T.begin(); pos<S3T.end(); pos++)
130  { delete [] *pos; }
131  S3T.clear();
132  for (pos=B3T.begin(); pos<B3T.end(); pos++)
133  { delete [] *pos; }
134  B3T.clear();
135  for (pos=S4T.begin(); pos<S4T.end(); pos++)
136  { delete [] *pos; }
137  S4T.clear();
138  for (pos=B4T.begin(); pos<B4T.end(); pos++)
139  { delete [] *pos; }
140  B4T.clear();
141}
142
143// Returns Pointer to the G4VQCrossSection class
144G4VQCrossSection* G4QElasticCrossSection::GetPointer()
145{
146  static G4QElasticCrossSection theCrossSection; //**Static body of the QEl Cross Section**
147  return &theCrossSection;
148}
149
150// The main member function giving the collision cross section (P is in IU, CS is in mb)
151// Make pMom in independent units ! (Now it is MeV)
152G4double G4QElasticCrossSection::GetCrossSection(G4bool fCS, G4double pMom, G4int tgZ,
153                                                 G4int tgN, G4int pPDG)
154{
155  static std::vector <G4int>    colPDG;// Vector of the projectile PDG code
156  static std::vector <G4int>    colN;  // Vector of N for calculated nuclei (isotops)
157  static std::vector <G4int>    colZ;  // Vector of Z for calculated nuclei (isotops)
158  static std::vector <G4double> colP;  // Vector of last momenta for the reaction
159  static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
160  static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
161  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
162  if(tgZ==0 && tgN==1)                 // Temporary change for Quasi-Elastic
163  {
164    tgZ=1;
165    tgN=1;
166    if     (pPDG==2212) pPDG=2112;
167    else if(pPDG==2112) pPDG=2212;
168  }
169  G4double pEn=pMom;
170  onlyCS=fCS;
171#ifdef pdebug
172  G4cout<<"G4QElCS::GetCS:>>> f="<<fCS<<", p="<<pMom<<", Z="<<tgZ<<"("<<lastZ<<") ,N="<<tgN
173        <<"("<<lastN<<"),PDG="<<pPDG<<"("<<lastPDG<<"), T="<<pEn<<"("<<lastTH<<")"<<",Sz="
174        <<colN.size()<<G4endl;
175  //CalculateCrossSection(fCS,-27,j,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
176#endif
177  if(!pPDG)
178  {
179#ifdef pdebug
180    G4cout<<"G4QElCS::GetCS: *** Found pPDG="<<pPDG<<" ====> CS=0"<<G4endl;
181    //CalculateCrossSection(fCS,-27,j,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
182#endif
183    return 0.;                         // projectile PDG=0 is a mistake (?!) @@
184  }
185  G4bool in=false;                     // By default the isotope must be found in the AMDB
186//GF the block will update parameters, needed for quasi-eleastic
187//GF if(tgN!=lastN || tgZ!=lastZ || pPDG!=lastPDG)// The nucleus isn't the LastUsedIsotope
188  {
189    in = false;                        // By default the isotope haven't be found in AMDB 
190    lastP   = 0.;                      // New momentum history (nothing to compare with)
191    lastPDG = pPDG;                    // The last PDG of the projectile
192    lastN   = tgN;                     // The last N of the calculated nucleus
193    lastZ   = tgZ;                     // The last Z of the calculated nucleus
194    lastI   = colN.size();             // Size of the Associative Memory DB in the heap
195    if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
196    {                                  // The nucleus with projPDG is found in AMDB
197      if(colPDG[i]==pPDG && colN[i]==tgN && colZ[i]==tgZ)
198      {
199        lastI=i;
200        lastTH =colTH[i];                // Last THreshold (A-dependent)
201#ifdef pdebug
202        G4cout<<"G4QElCS::GetCS:*Found* P="<<pMom<<",Threshold="<<lastTH<<",i="<<i<<G4endl;
203        //CalculateCrossSection(fCS,-27,i,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
204#endif
205        if(pEn<=lastTH)
206        {
207#ifdef pdebug
208          G4cout<<"G4QElCS::GetCS:Found T="<<pEn<<" < Threshold="<<lastTH<<",CS=0"<<G4endl;
209          //CalculateCrossSection(fCS,-27,i,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
210#endif
211          return 0.;                     // Energy is below the Threshold value
212        }
213        lastP  =colP [i];                // Last Momentum  (A-dependent)
214        lastCS =colCS[i];                // Last CrossSect (A-dependent)
215        //  if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
216        if(lastP == pMom)   // V.Ivanchenko safe solution
217        {
218#ifdef pdebug
219          G4cout<<"G4QElCS::GetCS:P="<<pMom<<",CS="<<lastCS*millibarn<<G4endl;
220#endif
221          CalculateCrossSection(fCS,-1,i,lastPDG,lastZ,lastN,pMom); // Update param's only
222          return lastCS*millibarn;     // Use theLastCS
223        }
224        in = true;                       // This is the case when the isotop is found in DB
225        // Momentum pMom is in IU ! @@ Units
226#ifdef pdebug
227        G4cout<<"G4QElCS::G:UpdateDB P="<<pMom<<",f="<<fCS<<",I="<<lastI<<",i="<<i<<G4endl;
228#endif
229        lastCS=CalculateCrossSection(fCS,-1,i,lastPDG,lastZ,lastN,pMom); // read & update
230#ifdef pdebug
231        G4cout<<"G4QElCS::GetCrosSec: *****> New (inDB) Calculated CS="<<lastCS<<G4endl;
232        //CalculateCrossSection(fCS,-27,i,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
233#endif
234        if(lastCS<=0. && pEn>lastTH)    // Correct the threshold
235        {
236#ifdef pdebug
237          G4cout<<"G4QElCS::GetCS: New T="<<pEn<<"(CS=0) > Threshold="<<lastTH<<G4endl;
238#endif
239          lastTH=pEn;
240        }
241        break;                           // Go out of the LOOP with found lastI
242      }
243#ifdef pdebug
244      G4cout<<"---G4QElCrossSec::GetCrosSec:pPDG="<<pPDG<<",i="<<i<<",N="<<colN[i]
245            <<",Z["<<i<<"]="<<colZ[i]<<",cPDG="<<colPDG[i]<<G4endl;
246      //CalculateCrossSection(fCS,-27,i,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
247#endif
248    }
249    if(!in)                            // This nucleus has not been calculated previously
250    {
251#ifdef pdebug
252      G4cout<<"G4QElCS::GetCrosSec:CalcNew P="<<pMom<<",f="<<fCS<<",lastI="<<lastI<<G4endl;
253#endif
254      //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
255      lastCS=CalculateCrossSection(fCS,0,lastI,lastPDG,lastZ,lastN,pMom);//calculate&create
256      if(lastCS<=0.)
257      {
258        lastTH = ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
259#ifdef pdebug
260        G4cout<<"G4QElCrossSection::GetCrossSect: NewThresh="<<lastTH<<",T="<<pEn<<G4endl;
261#endif
262        if(pEn>lastTH)
263        {
264#ifdef pdebug
265          G4cout<<"G4QElCS::GetCS: First T="<<pEn<<"(CS=0) > Threshold="<<lastTH<<G4endl;
266#endif
267          lastTH=pEn;
268        }
269      }
270#ifdef pdebug
271      G4cout<<"G4QElCS::GetCrosSec: New CS="<<lastCS<<",lZ="<<lastN<<",lN="<<lastZ<<G4endl;
272      //CalculateCrossSection(fCS,-27,lastI,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
273#endif
274      colN.push_back(tgN);
275      colZ.push_back(tgZ);
276      colPDG.push_back(pPDG);
277      colP.push_back(pMom);
278      colTH.push_back(lastTH);
279      colCS.push_back(lastCS);
280#ifdef pdebug
281      G4cout<<"G4QElCS::GetCS:1st,P="<<pMom<<"(MeV),CS="<<lastCS*millibarn<<"(mb)"<<G4endl;
282      //CalculateCrossSection(fCS,-27,lastI,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
283#endif
284      return lastCS*millibarn;
285    } // End of creation of the new set of parameters
286    else
287    {
288#ifdef pdebug
289      G4cout<<"G4QElCS::GetCS: Update lastI="<<lastI<<G4endl;
290#endif
291      colP[lastI]=pMom;
292      colPDG[lastI]=pPDG;
293      colCS[lastI]=lastCS;
294    }
295  } // End of parameters udate
296
297//GF
298// else if(pEn<=lastTH)
299// {
300//#ifdef pdebug
301//  G4cout<<"G4QElCS::GetCS: Current T="<<pEn<<" < Threshold="<<lastTH<<", CS=0"<<G4endl;
302//  //CalculateCrossSection(fCS,-27,lastI,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
303//#endif
304//  return 0.;                         // Momentum is below the Threshold Value -> CS=0
305// }
306// //  else if(std::fabs(lastP/pMom-1.)<tolerance)
307// else if(lastP == pMom) // V.Ivanchenko safe solution
308// {
309//#ifdef pdebug
310//  G4cout<<"G4QElCS::GetCS:OldCur P="<<pMom<<"="<<pMom<<", CS="<<lastCS*millibarn<<G4endl;
311//  //CalculateCrossSection(fCS,-27,lastI,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
312//  G4cout<<"G4QElCS::GetCrSec:***SAME***, onlyCS="<<onlyCS<<G4endl;
313//#endif
314//  return lastCS*millibarn;     // Use theLastCS
315// }
316// else
317// {
318//#ifdef pdebug
319//  G4cout<<"G4QElCS::GetCS:UpdatCur P="<<pMom<<",f="<<fCS<<",I="<<lastI<<",j="<<j<<G4endl;
320//#endif
321//  lastCS=CalculateCrossSection(fCS,1,lastI,lastPDG,lastZ,lastN,pMom); // Only UpdateDB
322//  lastP=pMom;
323// }
324//GF
325
326#ifdef pdebug
327  G4cout<<"G4QElCS::GetCrSec:End,P="<<pMom<<"(MeV),CS="<<lastCS*millibarn<<"(mb)"<<G4endl;
328  //CalculateCrossSection(fCS,-27,lastI,lastPDG,lastZ,lastN,pMom); // DUMMY TEST
329  G4cout<<"G4QElCS::GetCrSec:***End***, onlyCS="<<onlyCS<<G4endl;
330#endif
331  return lastCS*millibarn;
332}
333
334// Calculation of total elastic cross section (p in IU, CS in mb) @@ Units (?)
335// F=0 - create AMDB, F=-1 - read&update AMDB, F=1 - update AMDB (sinchro with higher AMDB)
336G4double G4QElasticCrossSection::CalculateCrossSection(G4bool CS,G4int F,G4int I,G4int PDG,
337                                                       G4int tgZ, G4int tgN, G4double pIU)
338{
339  // *** Begin of Associative Memory DB for acceleration of the cross section calculations
340  static std::vector <G4double>  PIN;   // Vector of max initialized log(P) in the table
341  // *** End of Static Definitions (Associative Memory Data Base) ***
342  if(tgZ==0 && tgN==1)                 // Temporary change for Quasi-Elastic
343  {
344    tgZ=1;
345    tgN=1;
346    if     (PDG==2212) PDG=2112;
347    else if(PDG==2112) PDG=2212;
348  }
349  G4double pMom=pIU/GeV;                // All calculations are in GeV
350  onlyCS=CS;                            // Flag to calculate only CS (not Si/Bi)
351#ifdef pdebug
352  G4cout<<"G4QElasticCrossSection::CalcCS:->onlyCS="<<onlyCS<<",F="<<F<<",p="<<pIU<<G4endl;
353#endif
354  lastLP=std::log(pMom);                // Make a logarithm of the momentum for calculation
355  if(F)                                 // This isotope was found in AMDB =>RETRIEVE/UPDATE
356  {
357    if(F<0)                             // the AMDB must be loded
358    {
359      lastPIN = PIN[I];                 // Max log(P) initialised for this table set
360      lastPAR = PAR[I];                 // Pointer to the parameter set
361      lastCST = CST[I];                 // Pointer to the total sross-section table
362      lastSST = SST[I];                 // Pointer to the first squared slope
363      lastS1T = S1T[I];                 // Pointer to the first mantissa
364      lastB1T = B1T[I];                 // Pointer to the first slope
365      lastS2T = S2T[I];                 // Pointer to the second mantissa
366      lastB2T = B2T[I];                 // Pointer to the second slope
367      lastS3T = S3T[I];                 // Pointer to the third mantissa
368      lastB3T = B3T[I];                 // Pointer to the rhird slope
369      lastS4T = S4T[I];                 // Pointer to the 4-th mantissa
370      lastB4T = B4T[I];                 // Pointer to the 4-th slope
371#ifdef pdebug
372      G4cout<<"G4QElasticCS::CalcCS: DB is updated for I="<<I<<",*,PIN4="<<PIN[4]<<G4endl;
373#endif
374    }
375#ifdef pdebug
376    G4cout<<"G4QElasticCrossSection::CalcCS:*read*, LP="<<lastLP<<",PIN="<<lastPIN<<G4endl;
377#endif
378    if(lastLP>lastPIN && lastLP<lPMax)
379    {
380      lastPIN=GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);// Can update upper logP-Limit in tabs
381#ifdef pdebug
382      G4cout<<"G4QElCS::CalcCS:*updated(I)*,LP="<<lastLP<<"<IN["<<I<<"]="<<lastPIN<<G4endl;
383#endif
384      PIN[I]=lastPIN;                   // Remember the new P-Limit of the tables
385    }
386  }
387  else                                  // This isotope wasn't initialized => CREATE
388  {
389    lastPAR = new G4double[nPoints];    // Allocate memory for parameters of CS function
390    lastPAR[nLast]=0;                   // Initialization for VALGRIND
391    lastCST = new G4double[nPoints];    // Allocate memory for Tabulated CS function   
392    lastSST = new G4double[nPoints];    // Allocate memory for Tabulated first sqaredSlope
393    lastS1T = new G4double[nPoints];    // Allocate memory for Tabulated first mantissa
394    lastB1T = new G4double[nPoints];    // Allocate memory for Tabulated first slope   
395    lastS2T = new G4double[nPoints];    // Allocate memory for Tabulated second mantissa
396    lastB2T = new G4double[nPoints];    // Allocate memory for Tabulated second slope   
397    lastS3T = new G4double[nPoints];    // Allocate memory for Tabulated third mantissa
398    lastB3T = new G4double[nPoints];    // Allocate memory for Tabulated third slope   
399    lastS4T = new G4double[nPoints];    // Allocate memory for Tabulated 4-th mantissa
400    lastB4T = new G4double[nPoints];    // Allocate memory for Tabulated 4-th slope   
401#ifdef pdebug
402    G4cout<<"G4QElasticCrossSection::CalcCS:*ini*,lastLP="<<lastLP<<",min="<<lPMin<<G4endl;
403#endif
404    lastPIN = GetPTables(lastLP,lPMin,PDG,tgZ,tgN); // Returns the new P-limit for tables
405#ifdef pdebug
406    G4cout<<"G4QElCS::CalcCS:*i*Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<",LP"<<lastPIN<<G4endl;
407#endif
408    PIN.push_back(lastPIN);             // Fill parameters of CS function to AMDB
409    PAR.push_back(lastPAR);             // Fill parameters of CS function to AMDB
410    CST.push_back(lastCST);             // Fill Tabulated CS function to AMDB   
411    SST.push_back(lastSST);             // Fill Tabulated first sq.slope to AMDB
412    S1T.push_back(lastS1T);             // Fill Tabulated first mantissa to AMDB
413    B1T.push_back(lastB1T);             // Fill Tabulated first slope to AMDB   
414    S2T.push_back(lastS2T);             // Fill Tabulated second mantissa to AMDB
415    B2T.push_back(lastB2T);             // Fill Tabulated second slope to AMDB   
416    S3T.push_back(lastS3T);             // Fill Tabulated third mantissa to AMDB
417    B3T.push_back(lastB3T);             // Fill Tabulated third slope to AMDB   
418    S4T.push_back(lastS4T);             // Fill Tabulated 4-th mantissa to AMDB
419    B4T.push_back(lastB4T);             // Fill Tabulated 4-th slope to AMDB   
420  } // End of creation/update of the new set of parameters and tables
421  // ============= NOW Update (if necessary) and Calculate the Cross Section ===========
422#ifdef pdebug
423  G4cout<<"G4QElCS::CalcCS:?update?,LP="<<lastLP<<",IN="<<lastPIN<<",ML="<<lPMax<<G4endl;
424#endif
425  if(lastLP>lastPIN && lastLP<lPMax)
426  {
427    lastPIN = GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);
428#ifdef pdebug
429    G4cout<<"G4QElCS::CalcCS: *updated(O)*, LP="<<lastLP<<" < IN="<<lastPIN<<G4endl;
430#endif
431  }
432#ifdef pdebug
433  G4cout<<"G4QElastCS::CalcCS: lastLP="<<lastLP<<",lPM="<<lPMin<<",lPIN="<<lastPIN<<G4endl;
434#endif
435  if(!onlyCS) lastTM=GetQ2max(PDG, tgZ, tgN, pMom); // Calculate (-t)_max=Q2_max (GeV2)
436#ifdef pdebug
437  G4cout<<"G4QElasticCrosSec::CalcCS:oCS="<<onlyCS<<",-t="<<lastTM<<", p="<<lastLP<<G4endl;
438#endif
439  if(lastLP>lPMin && lastLP<=lastPIN)   // Linear fit is made using precalculated tables
440  {
441    if(lastLP==lastPIN) // VI do not use tolerance
442      //    if(std::fabs(lastLP-lastPIN)<.0001) // Just take the highest tabulated value
443    {
444      G4double shift=(lastLP-lPMin)/dlnP+.000001; // Log distance from lPMin
445      G4int    blast=static_cast<int>(shift); // this is a bin number of the lower edge (0)
446      if(blast<0 || blast>=nLast) G4cout<<"G4QEleastCS::CCS:b="<<blast<<","<<nLast<<G4endl;
447      lastSIG = lastCST[blast];
448      if(!onlyCS)                       // Skip the differential cross-section parameters
449      {
450        theSS  = lastSST[blast];
451        theS1  = lastS1T[blast];
452        theB1  = lastB1T[blast];
453        theS2  = lastS2T[blast];
454        theB2  = lastB2T[blast];
455        theS3  = lastS3T[blast];
456        theB3  = lastB3T[blast];
457        theS4  = lastS4T[blast];
458        theB4  = lastB4T[blast];
459      }
460#ifdef pdebug
461      G4cout<<"G4QElasticCrossSection::CalculateCS:(E) S1="<<theS1<<", B1="<<theB1<<G4endl;
462#endif
463    }
464    else
465    {
466      G4double shift=(lastLP-lPMin)/dlnP;        // a shift from the beginning of the table
467      G4int    blast=static_cast<int>(shift);    // the lower bin number
468      if(blast<0)   blast=0;
469      if(blast>=nLast) blast=nLast-1;            // low edge of the last bin
470      shift-=blast;                              // step inside the unit bin
471      G4int lastL=blast+1;                       // the upper bin number
472      G4double SIGL=lastCST[blast];              // the basic value of the cross-section
473      lastSIG= SIGL+shift*(lastCST[lastL]-SIGL); // calculated total elastic cross-section
474#ifdef pdebug
475      G4cout<<"G4QElCS::CalcCrossSection: Sig="<<lastSIG<<", P="<<pMom<<", Z="<<tgZ<<", N="
476            <<tgN<<", PDG="<<PDG<<", onlyCS="<<onlyCS<<G4endl;
477#endif
478      if(!onlyCS)                       // Skip the differential cross-section parameters
479      {
480        G4double SSTL=lastSST[blast];           // the low bin of the first squared slope
481        theSS=SSTL+shift*(lastSST[lastL]-SSTL); // the basic value of the first sq.slope
482        G4double S1TL=lastS1T[blast];           // the low bin of the first mantissa
483        theS1=S1TL+shift*(lastS1T[lastL]-S1TL); // the basic value of the first mantissa
484        G4double B1TL=lastB1T[blast];           // the low bin of the first slope
485#ifdef pdebug
486        G4cout<<"G4QElCS::CalcCrossSection:bl="<<blast<<",ls="<<lastL<<",SL="<<S1TL<<",SU="
487              <<lastS1T[lastL]<<",BL="<<B1TL<<",BU="<<lastB1T[lastL]<<G4endl;
488#endif
489        theB1=B1TL+shift*(lastB1T[lastL]-B1TL); // the basic value of the first slope
490        G4double S2TL=lastS2T[blast];           // the low bin of the second mantissa
491        theS2=S2TL+shift*(lastS2T[lastL]-S2TL); // the basic value of the second mantissa
492        G4double B2TL=lastB2T[blast];           // the low bin of the second slope
493        theB2=B2TL+shift*(lastB2T[lastL]-B2TL); // the basic value of the second slope
494        G4double S3TL=lastS3T[blast];           // the low bin of the third mantissa
495        theS3=S3TL+shift*(lastS3T[lastL]-S3TL); // the basic value of the third mantissa
496#ifdef pdebug
497        G4cout<<"G4QElCS::CCS: s3l="<<S3TL<<",sh3="<<shift<<",s3h="<<lastS3T[lastL]<<",b="
498              <<blast<<",l="<<lastL<<G4endl;
499#endif
500        G4double B3TL=lastB3T[blast];           // the low bin of the third slope
501        theB3=B3TL+shift*(lastB3T[lastL]-B3TL); // the basic value of the third slope
502        G4double S4TL=lastS4T[blast];           // the low bin of the 4-th mantissa
503        theS4=S4TL+shift*(lastS4T[lastL]-S4TL); // the basic value of the 4-th mantissa
504#ifdef pdebug
505        G4cout<<"G4QElCS::CCS: s4l="<<S4TL<<",sh4="<<shift<<",s4h="<<lastS4T[lastL]<<",b="
506              <<blast<<",l="<<lastL<<G4endl;
507#endif
508        G4double B4TL=lastB4T[blast];           // the low bin of the 4-th slope
509        theB4=B4TL+shift*(lastB4T[lastL]-B4TL); // the basic value of the 4-th slope
510      }
511#ifdef pdebug
512      G4cout<<"G4QElasticCrossSection::CalculateCS:(I) S1="<<theS1<<", B1="<<theB1<<G4endl;
513#endif
514    }
515  }
516  else lastSIG=GetTabValues(lastLP, PDG, tgZ, tgN); // Direct calculation beyond the table
517  if(lastSIG<0.) lastSIG = 0.;                   // @@ a Warning print can be added
518#ifdef pdebug
519  G4cout<<"G4QElasticCrossSection::CalculateCS: END, onlyCS="<<onlyCS<<G4endl;
520#endif
521  return lastSIG;
522}
523
524// It has parameter sets for all tZ/tN/PDG, using them the tables can be created/updated
525G4double G4QElasticCrossSection::GetPTables(G4double LP,G4double ILP, G4int PDG, G4int tgZ,
526                                                                                 G4int tgN)
527{
528  // @@ At present all nA==pA ---------> Each neucleus can have not more than 51 parameters
529  static const G4double pwd=2727;
530  const G4int n_npel=24;                // #of parameters for np-elastic (<nPoints=128)
531  const G4int n_ppel=32;                // #of parameters for pp-elastic (<nPoints=128)
532  //                      -0- -1-  -2- -3- -4-  -5- -6- -7- -8- -9--10--11--12--13- -14-
533  G4double np_el[n_npel]={12.,.05,.0001,5.,.35,6.75,.14,19.,.6,6.75,.14,13.,.14,.6,.00013,
534                          75.,.001,7.2,4.32,.012,2.5,0.0,12.,.34};
535  //                      -15--16--17- -18- -19--20--21--22--23-
536  //                       -0-   -1-  -2- -3- -4- -5-  -6-  -7-  -8--9--10--11--12--13-
537  G4double pp_el[n_ppel]={2.865,18.9,.6461,3.,9.,.425,.4276,.0022,5.,74.,3.,3.4,.2,.17,
538                          .001,8.,.055,3.64,5.e-5,4000.,1500.,.46,1.2e6,3.5e6,5.e-5,1.e10,
539                          8.5e8,1.e10,1.1,3.4e6,6.8e6,0.};
540  //                      -14--15- -16- -17- -18-  -19- -20- -21- -22-  -23-   -24-  -25-
541  //                       -26- -27-  -28- -29- -30- -31-
542  //                      -0- -1--2- -3-  -4-  -5-  -6-    -7-  -8- -9- -10--11--12--13-
543  if(tgZ==0 && tgN==1)                 // Temporary change for Quasi-Elastic
544  {
545    tgZ=1;
546    tgN=1;
547    if     (PDG==2212) PDG=2112;
548    else if(PDG==2112) PDG=2212;
549  }
550  if(PDG==2212 || PDG==2112)
551  {
552    // --- Total np elastic cross section cs & s1/b1 (t), s2/b2 (u) --- NotTuned for highE
553    //p2=p*p;p3=p2*p;sp=sqrt(p);p2s=p2*sp;lp=log(p);dl1=lp-(5.=par(3));p4=p2*p2; p=|3-mom|
554    //CS=12./(p2s+.05*p+.0001/sqrt(sp))+.35/p+(6.75+.14*dl1*dl1+19./p)/(1.+.6/p4);
555    //  par(0)   par(1) par(2)        par(4) par(5) par(6)     par(7)     par(8)
556    //s1=(6.75+.14*dl2*dl2+13./p)/(1.+.14/p4)+.6/(p4+.00013), s2=(75.+.001/p4/p)/p3
557    //  par(9) par(10)   par(11)     par(12) par(13) par(14)  par(15) par(16)
558    //b1=(7.2+4.32/(p4*p4+.012*p3))/(1.+2.5/p4), ss=0., b2=12./(p*sp+.34)
559    //par(17) par(18)    par(19)       par(20)  par(21)   par(22)   par(23)
560    //
561    // -- Total pp elastic cross section cs & s1/b1 (main), s2/b2 (tail1), s3/b3 (tail2) --
562    //p2=p*p;p3=p2*p;sp=sqrt(p);p2s=p2*sp;lp=log(p);dl1=lp-(3.=par(3));p4=p2*p2; p=|3-mom|
563    //CS=2.865/p2s/(1+.0022/p2s)+(18.9+.6461*dl1*dl1+9./p)/(1.+.425*lp)/(1.+.4276/p4);
564    //   par(0)       par(7)     par(1) par(2)      par(4)      par(5)         par(6)
565    //dl2=lp-5., s1=(74.+3.*dl2*dl2)/(1+3.4/p4/p)+(.2/p2+17.*p)/(p4+.001*sp),
566    //     par(8) par(9) par(10)        par(11)   par(12)par(13)    par(14)
567    // b1=8.*p**.055/(1.+3.64/p3); s2=5.e-5+4000./(p4+1500.*p); b2=.46+1.2e6/(p4+3.5e6/sp);
568    // par(15) par(16)  par(17)     par(18) par(19)  par(20)   par(21) par(22)  par(23)
569    // s3=5.e-5+1.e10/(p4*p4+8.5e8*p2+1.e10); b3=1.1+3.4e6/(p4+6.8e6); ss=0.
570    //  par(24) par(25)     par(26)  par(27) par(28) par(29)  par(30)   par(31)
571    //
572    if(lastPAR[nLast]!=pwd) // A unique flag to avoid the repeatable definition
573    {
574      if ( (PDG == 2112 && tgZ == 1 && tgN == 0) || 
575           (PDG == 2212 && tgZ == 0 && tgN == 1) ) {
576        for (G4int ip=0; ip<n_npel; ip++) lastPAR[ip]=np_el[ip]; //np/pn
577
578      } else if ( (PDG == 2212 && tgZ == 1 && tgN == 0) || 
579                  (PDG == 2112 && tgZ == 0 && tgN == 1) ) {
580        for (G4int ip=0; ip<n_ppel; ip++) lastPAR[ip]=pp_el[ip]; //pp/nn
581
582      } else {
583        G4double a=tgZ+tgN;
584        G4double sa=std::sqrt(a);
585        G4double ssa=std::sqrt(sa);
586        G4double asa=a*sa;
587        G4double a2=a*a;
588        G4double a3=a2*a;
589        G4double a4=a3*a;
590        G4double a5=a4*a;
591        G4double a6=a4*a2;
592        G4double a7=a6*a;
593        G4double a8=a7*a;
594        G4double a9=a8*a;
595        G4double a10=a5*a5;
596        G4double a12=a6*a6;
597        G4double a14=a7*a7;
598        G4double a16=a8*a8;
599        G4double a17=a16*a;
600        G4double a32=a16*a16;
601        G4double r1a16=3.e16/a16;
602        G4double r2a16=2.e13/a16;
603        G4double r3a16=6.e13/a16;
604        G4double pa10=5.e-27*a10;
605        // Reaction cross-section parameters (pel=peh_fit.f)
606        lastPAR[0]=a/3.;                                                       // p1(p/n)
607        //lastPAR[0]=.28*a;                                                      // p1(p/n)
608        lastPAR[1]=4.5*a*ssa/(1.+6./a3)/(1.+1.e-2*a);                          // p2(p/n)
609        //lastPAR[1]=4.8*std::pow(a,1.14)/(1.+3.6/a3);                           // p2(p/n)
610        if(PDG==2112)                           // Neutron projectile
611        {
612          lastPAR[2]=3.3/a+.5*ssa/(1.+2.e7/a8);                                // p3 (old)
613          if(a<6.5)
614          {
615            lastPAR[3]=1./(1.+.00123*a4);                                      // p4
616            lastPAR[4]=1.5e-4/a2/(a+1.2e-6*a12)+1.5e-6;                        // p5
617            lastPAR[5]=.0062/(a+5.e-11*a16);                                   // p6
618            lastPAR[6]=3.2e-14/a6/(1.+4.e-10*a17);                             // p7
619            lastPAR[7]=.847/a4/(a+.0005*a8)+.00045;                            // p8
620            lastPAR[8]=.413/a16+7.e-7;                                         // p9
621          }
622          else
623          {
624            lastPAR[3]=a*(.5+1.e-10*a4)/(1.+7.e5/a6);                          // p4
625            lastPAR[4]=a*3.e-5/(1.+2.e12/a12);                                 // p5
626            lastPAR[5]=(.0006+1.e-14*a5)/(1.+4.e5/a3)/(1.+3.e23/a16/a8);       // p6
627            lastPAR[6]=1.e-22*a4/(1.+5.e16*(1.+5e31/a32)/a16);                 // p7
628            lastPAR[7]=(8.+.00016*a2)/(1.+1.4e14/a16);                         // p8
629            lastPAR[8]=.0013;                                                  // p9
630          }
631        }
632        else                                   // Proton projectile
633        {
634          lastPAR[2]=3.3/a+.8*ssa/(1.+2.e7/a8);                                // p3
635          //lastPAR[2]=3.3/a+.7*ssa/(1.+2.e7/a8);                                // p3(old)
636          lastPAR[3]=1./(1.+6.e-9*a12)+.6*a/(1.+1.6e15/a16);                   // p4
637          lastPAR[4]=6.e-4/(a4+3.e-6*a12)+.16/(a+9.e5/a3+r1a16*r1a16);         // p5
638          lastPAR[5]=3.e-4/a2+(1.5e-3+3.e-14*a7)/(1.+r2a16*r2a16+3.e-12*a6);   // p6
639          lastPAR[6]=4.e-30+(1.2e-28*a3+pa10*pa10)/(1.+r3a16*r3a16+4.e-26*a14);// p7
640          lastPAR[7]=.1/(1.+1.7e-8*a16)+1./(1.+1.7e18/a16+2.e-6*a4);         // p8
641          //lastPAR[7]=.07/(1.+1.7e-8*a16)+.5/(1.+1.7e18/a16+4.5e-7*a4);         // p8(old)
642          lastPAR[8]=(1.5e-10+2.e-18*a8)/(1.+3.e-25*a16);                      // p9
643        }
644        // @@ the differential cross-section is parameterized separately for A>6 & A<7
645        if(a<6.5)
646        {
647          G4double a28=a16*a12;
648          // The main pre-exponent      (pel_sg)
649          lastPAR[ 9]=4000*a;                                // p1
650          lastPAR[10]=1.2e7*a8+380*a17;                      // p2
651          lastPAR[11]=.7/(1.+4.e-12*a16);                    // p3
652          lastPAR[12]=2.5/a8/(a4+1.e-16*a32);                // p4
653          lastPAR[13]=.28*a;                                 // p5
654          lastPAR[14]=1.2*a2+2.3;                            // p6
655          lastPAR[15]=3.8/a;                                 // p7
656          // The main slope             (pel_sl)
657          lastPAR[16]=.01/(1.+.0024*a5);                     // p1
658          lastPAR[17]=.2*a;                                  // p2
659          lastPAR[18]=9.e-7/(1.+.035*a5);                    // p3
660          lastPAR[19]=(42.+2.7e-11*a16)/(1.+.14*a);          // p4
661          // The main quadratic         (pel_sh)
662          lastPAR[20]=2.25*a3;                               // p1
663          lastPAR[21]=18.;                                   // p2
664          lastPAR[22]=2.4e-3*a8/(1.+2.6e-4*a7);              // p3
665          lastPAR[23]=3.5e-36*a32*a8/(1.+5.e-15*a32/a);      // p4
666          // The 1st max pre-exponent   (pel_qq)
667          lastPAR[24]=1.e5/(a8+2.5e12/a16);                  // p1
668          //lastPAR[24]=8.e4/(a8+2.5e12/a16);                  // p1 (old)
669          lastPAR[25]=8.e7/(a12+1.e-27*a28*a28);             // p2
670          lastPAR[26]=.0006*a3;                              // p3
671          //lastPAR[26]=.0011*a3;                              // p3 (old)
672          // The 1st max slope          (pel_qs)
673          lastPAR[27]=10.+4.e-8*a12*a;                       // p1
674          lastPAR[28]=.114;                                  // p2
675          lastPAR[29]=.003;                                  // p3
676          lastPAR[30]=2.e-23;                                // p4
677          // The effective pre-exponent (pel_ss)
678          lastPAR[31]=1./(1.+.0001*a8);                      // p1
679          lastPAR[32]=1.5e-4/(1.+5.e-6*a12);                 // p2
680          lastPAR[33]=.03;                                   // p3
681          // The effective slope        (pel_sb)
682          lastPAR[34]=a/2;                                   // p1
683          lastPAR[35]=2.e-7*a4;                              // p2
684          lastPAR[36]=4.;                                    // p3
685          lastPAR[37]=64./a3;                                // p4
686          // The gloria pre-exponent    (pel_us)
687          lastPAR[38]=1.e8*std::exp(.32*asa);                // p1
688          //lastPAR[38]=1.05e8*std::exp(.32*asa);              // p1 (old)
689          lastPAR[39]=20.*std::exp(.45*asa);                // p2
690          //lastPAR[39]=19.5*std::exp(.45*asa);                // p2 (old)
691          lastPAR[40]=7.e3+2.4e6/a5;                         // p3
692          lastPAR[41]=2.5e5*std::exp(.085*a3);               // p4
693          lastPAR[42]=2.5*a;                                 // p5
694          // The gloria slope           (pel_ub)
695          lastPAR[43]=920.+.03*a8*a3;                        // p1
696          lastPAR[44]=93.+.0023*a12;                         // p2
697#ifdef pdebug
698         G4cout<<"G4QElCS::CalcCS:la "<<lastPAR[38]<<", "<<lastPAR[39]<<", "<<lastPAR[40]
699               <<", "<<lastPAR[42]<<", "<<lastPAR[43]<<", "<<lastPAR[44]<<G4endl;
700#endif
701        }
702        else
703        {
704          G4double p1a10=2.2e-28*a10;
705          G4double r4a16=6.e14/a16;
706          G4double s4a16=r4a16*r4a16;
707          // a24
708          // a36
709          // The main pre-exponent      (peh_sg)
710          lastPAR[ 9]=4.5*std::pow(a,1.15);                  // p1
711          lastPAR[10]=.06*std::pow(a,.6);                    // p2
712          lastPAR[11]=.6*a/(1.+2.e15/a16);                   // p3
713          lastPAR[12]=.17/(a+9.e5/a3+1.5e33/a32);            // p4
714          lastPAR[13]=(.001+7.e-11*a5)/(1.+4.4e-11*a5);      // p5
715          lastPAR[14]=(p1a10*p1a10+2.e-29)/(1.+2.e-22*a12);  // p6
716          // The main slope             (peh_sl)
717          lastPAR[15]=400./a12+2.e-22*a9;                    // p1
718          lastPAR[16]=1.e-32*a12/(1.+5.e22/a14);             // p2
719          lastPAR[17]=1000./a2+9.5*sa*ssa;                   // p3
720          lastPAR[18]=4.e-6*a*asa+1.e11/a16;                 // p4
721          lastPAR[19]=(120./a+.002*a2)/(1.+2.e14/a16);       // p5
722          lastPAR[20]=9.+100./a;                             // p6
723          // The main quadratic         (peh_sh)
724          lastPAR[21]=.002*a3+3.e7/a6;                       // p1
725          lastPAR[22]=7.e-15*a4*asa;                         // p2
726          lastPAR[23]=9000./a4;                              // p3
727          // The 1st max pre-exponent   (peh_qq)
728          lastPAR[24]=.0011*asa/(1.+3.e34/a32/a4);           // p1
729          lastPAR[25]=1.e-5*a2+2.e14/a16;                    // p2
730          lastPAR[26]=1.2e-11*a2/(1.+1.5e19/a12);            // p3
731          lastPAR[27]=.016*asa/(1.+5.e16/a16);               // p4
732          // The 1st max slope          (peh_qs)
733          lastPAR[28]=.002*a4/(1.+7.e7/std::pow(a-6.83,14)); // p1
734          lastPAR[29]=2.e6/a6+7.2/std::pow(a,.11);           // p2
735          lastPAR[30]=11.*a3/(1.+7.e23/a16/a8);              // p3
736          lastPAR[31]=100./asa;                              // p4
737          // The 2nd max pre-exponent   (peh_ss)
738          lastPAR[32]=(.1+4.4e-5*a2)/(1.+5.e5/a4);           // p1
739          lastPAR[33]=3.5e-4*a2/(1.+1.e8/a8);                // p2
740          lastPAR[34]=1.3+3.e5/a4;                           // p3
741          lastPAR[35]=500./(a2+50.)+3;                       // p4
742          lastPAR[36]=1.e-9/a+s4a16*s4a16;                   // p5
743          // The 2nd max slope          (peh_sb)
744          lastPAR[37]=.4*asa+3.e-9*a6;                       // p1
745          lastPAR[38]=.0005*a5;                              // p2
746          lastPAR[39]=.002*a5;                               // p3
747          lastPAR[40]=10.;                                   // p4
748          // The effective pre-exponent (peh_us)
749          lastPAR[41]=.05+.005*a;                            // p1
750          lastPAR[42]=7.e-8/sa;                              // p2
751          lastPAR[43]=.8*sa;                                 // p3
752          lastPAR[44]=.02*sa;                                // p4
753          lastPAR[45]=1.e8/a3;                               // p5
754          lastPAR[46]=3.e32/(a32+1.e32);                     // p6
755          // The effective slope        (peh_ub)
756          lastPAR[47]=24.;                                   // p1
757          lastPAR[48]=20./sa;                                // p2
758          lastPAR[49]=7.e3*a/(sa+1.);                        // p3
759          lastPAR[50]=900.*sa/(1.+500./a3);                  // p4
760#ifdef pdebug
761         G4cout<<"G4QElCS::CalcCS:ha "<<lastPAR[41]<<", "<<lastPAR[42]<<", "<<lastPAR[43]
762               <<", "<<lastPAR[44]<<", "<<lastPAR[45]<<", "<<lastPAR[46]<<G4endl;
763#endif
764        }
765        // Parameter for lowEnergyNeutrons
766        lastPAR[51]=1.e15+2.e27/a4/(1.+2.e-18*a16);
767      }
768      lastPAR[nLast]=pwd;
769      // and initialize the zero element of the table
770      G4double lp=lPMin;                                      // ln(momentum)
771      G4bool memCS=onlyCS;                                    // ??
772      onlyCS=false;
773      lastCST[0]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables
774      onlyCS=memCS;
775      lastSST[0]=theSS;
776      lastS1T[0]=theS1;
777      lastB1T[0]=theB1;
778      lastS2T[0]=theS2;
779      lastB2T[0]=theB2;
780      lastS3T[0]=theS3;
781      lastB3T[0]=theB3;
782      lastS4T[0]=theS4;
783      lastB4T[0]=theB4;
784#ifdef pdebug
785      G4cout<<"G4QElasticCrossSection::GetPTables:ip=0(init), lp="<<lp<<",S1="<<theS1
786            <<",B1="<<theB1<<",S2="<<theS2<<",B2="<<theB3<<",S3="<<theS3
787            <<",B3="<<theB3<<",S4="<<theS4<<",B4="<<theB4<<G4endl;
788#endif
789    }
790    if(LP>ILP)
791    {
792      G4int ini = static_cast<int>((ILP-lPMin+.000001)/dlnP)+1; // already inited till this
793      if(ini<0) ini=0;
794      if(ini<nPoints)
795      {
796        G4int fin = static_cast<int>((LP-lPMin)/dlnP)+1; // final bin of initialization
797        if(fin>=nPoints) fin=nLast;               // Limit of the tabular initialization
798        if(fin>=ini)
799        {
800          G4double lp=0.;
801          for(G4int ip=ini; ip<=fin; ip++)        // Calculate tabular CS,S1,B1,S2,B2,S3,B3
802          {
803            lp=lPMin+ip*dlnP;                     // ln(momentum)
804            G4bool memCS=onlyCS;
805            onlyCS=false;
806            lastCST[ip]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables (ret CS)
807            onlyCS=memCS;
808            lastSST[ip]=theSS;
809            lastS1T[ip]=theS1;
810            lastB1T[ip]=theB1;
811            lastS2T[ip]=theS2;
812            lastB2T[ip]=theB2;
813            lastS3T[ip]=theS3;
814            lastB3T[ip]=theB3;
815            lastS4T[ip]=theS4;
816            lastB4T[ip]=theB4;
817#ifdef pdebug
818            G4cout<<"G4QElasticCrossSection::GetPTables:ip="<<ip<<",lp="<<lp<<",S1="<<theS1
819                  <<",B1="<<theB1<<",S2="<<theS2<<",B2="<<theB2<<",S3="
820                  <<theS3<<",B3="<<theB3<<",S4="<<theS4<<",B4="<<theB4<<G4endl;
821#endif
822          }
823          return lp;
824        }
825        else G4cout<<"*Warning*G4QElasticCrossSection::GetPTables: PDG="<<PDG<<", Z="<<tgZ
826                   <<", N="<<tgN<<", i="<<ini<<" > fin="<<fin<<", LP="<<LP<<" > ILP="<<ILP
827                   <<" nothing is done!"<<G4endl;
828      }
829      else G4cout<<"*Warning*G4QElasticCrossSection::GetPTables: PDG="<<PDG<<", Z="<<tgZ
830                 <<", N="<<tgN<<", i="<<ini<<">= max="<<nPoints<<", LP="<<LP<<" > ILP="
831                 <<ILP<<", lPMax="<<lPMax<<" nothing is done!"<<G4endl;
832    }
833#ifdef pdebug
834    else G4cout<<"*Warning*G4QElasticCrossSection::GetPTables: PDG="<<PDG<<", Z="<<tgZ
835               <<", N="<<tgN<<", LP="<<LP<<" <= ILP="<<ILP<<" nothing is done!"<<G4endl;
836#endif
837  }
838  else
839  {
840    G4cout<<"*Error*G4QElasticCrossSection::GetPTables: PDG="<<PDG<<", Z="<<tgZ<<", N="
841          <<tgN<<", while it is defined only for PDG=2212(p)/2112(n)"<<G4endl;
842    throw G4QException("G4QElasticCrossSection::GetPTables: only nA & pA are implemented");
843  }
844  return ILP;
845}
846
847// Returns Q2=-t in independent units (MeV^2) (all internal calculations are in GeV)
848G4double G4QElasticCrossSection::GetExchangeT(G4int tgZ, G4int tgN, G4int PDG)
849{
850  static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
851  static const G4double third=1./3.;
852  static const G4double fifth=1./5.;
853  static const G4double sevth=1./7.;
854  if(tgZ==0 && tgN==1)                 // Temporary change for Quasi-Elastic
855  {
856    tgZ=1;
857    tgN=1;
858    if     (PDG==2212) PDG=2112;
859    else if(PDG==2112) PDG=2212;
860  }
861#ifdef tdebug
862  G4cout<<"G4QElasticCS::GetExcT:F="<<onlyCS<<",Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<G4endl;
863#endif
864  if(onlyCS) G4cout<<"**Warning*G4QElasticCrossSection::GetExchangeT: onlyCS=true"<<G4endl;
865  if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
866  G4double q2=0.;
867  if(PDG==2112 && tgZ==1 && tgN==0)                // ===> n+p=n+p
868  {
869#ifdef tdebug
870    G4cout<<"G4QElasticCS::GetExchangeT: TM="<<lastTM<<",S1="<<theS1<<",B1="<<theB1<<",S2="
871          <<theS2<<",B2="<<theB2<<",GeV2="<<GeVSQ<<G4endl;
872#endif
873    G4double E1=lastTM*theB1;
874    G4double R1=(1.-std::exp(-E1));
875#ifdef tdebug
876    G4double ts1=-std::log(1.-R1)/theB1;
877    G4double ds1=std::fabs(ts1-lastTM)/lastTM;
878    if(ds1>.0001)
879      G4cout<<"*Warning*G4QElCS::GetExT:1n "<<ts1<<"#"<<lastTM<<",d="<<ds1
880            <<",R1="<<R1<<",E1="<<E1<<G4endl;
881#endif
882    G4double E2=lastTM*theB2;
883    G4double R2=(1.-std::exp(-E2));
884#ifdef tdebug
885    G4double ts2=-std::log(1.-R2)/theB2;
886    G4double ds2=std::fabs(ts2-lastTM)/lastTM;
887    if(ds2>.0001)
888      G4cout<<"*Warning*G4QElCS::GetExT:2n "<<ts2<<"#"<<lastTM<<",d="<<ds2
889            <<",R2="<<R2<<",E2="<<E2<<G4endl;
890#endif
891    //G4double E3=lastTM*theB3;
892    //G4double R3=(1.-std::exp(-E3));
893#ifdef tdebug
894    //G4double ts3=-std::log(1.-R3)/theB3;
895    //G4double ds3=std::fabs(ts3-lastTM)/lastTM;
896    //if(ds3>.01)G4cout<<"*Warn*G4QElCS::GetExT:3n "<<ts3<<"#"<<lastTM<<",d="<<ds3<<G4endl;
897#endif
898    G4double I1=R1*theS1;
899    G4double I2=R2*theS2/theB2;
900    //G4double I3=R3*theS3/theB3;
901    G4double I12=I1+I2;
902    //G4double rand=(I12+I3)*G4UniformRand();
903    G4double rand=I12*G4UniformRand();
904    if     (rand<I1 )
905    {
906      G4double ran=R1*G4UniformRand();
907      if(ran>1.) ran=1.;
908      q2=-std::log(1.-ran)/theB1;       // t-chan
909    }
910    else
911    {
912      G4double ran=R2*G4UniformRand();
913      if(ran>1.) ran=1.;
914      q2=lastTM+std::log(1.-ran)/theB2; // u-chan (ChEx)
915    }
916  }
917  else if(PDG==2212 && tgZ==1 && tgN==0)                // ===> p+p=p+p
918  {
919#ifdef tdebug
920    G4cout<<"G4QElasticCS::GetExchangeT: TM="<<lastTM<<",S1="<<theS1<<",B1="<<theB1<<",S2="
921          <<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",GeV2="<<GeVSQ<<G4endl;
922#endif
923    G4double E1=lastTM*theB1;
924    G4double R1=(1.-std::exp(-E1));
925#ifdef tdebug
926    G4double ts1=-std::log(1.-R1)/theB1;
927    G4double ds1=std::fabs(ts1-lastTM)/lastTM;
928    if(ds1>.0001)
929      G4cout<<"*Warning*G4QElCS::GetExT:1p "<<ts1<<"#"<<lastTM<<",d="<<ds1
930            <<",R1="<<R1<<",E1="<<E1<<G4endl;
931#endif
932    G4double E2=lastTM*theB2;
933    G4double R2=(1.-std::exp(-E2*E2*E2));
934#ifdef tdebug
935    G4double ts2=std::pow(-std::log(1.-R2),.333333333)/theB2;
936    G4double ds2=std::fabs(ts2-lastTM)/lastTM;
937    if(ds2>.0001)
938      G4cout<<"*Warning*G4QElCS::GetExT:2p "<<ts2<<"#"<<lastTM<<",d="<<ds2
939            <<",R2="<<R2<<",E2="<<E2<<G4endl;
940#endif
941    G4double E3=lastTM*theB3;
942    G4double R3=(1.-std::exp(-E3));
943#ifdef tdebug
944    G4double ts3=-std::log(1.-R3)/theB3;
945    G4double ds3=std::fabs(ts3-lastTM)/lastTM;
946    if(ds3>.0001)
947      G4cout<<"*Warning*G4QElCS::GetExT:3p "<<ts3<<"#"<<lastTM<<",d="<<ds3
948            <<",R3="<<R1<<",E3="<<E3<<G4endl;
949#endif
950    G4double I1=R1*theS1/theB1;
951    G4double I2=R2*theS2;
952    G4double I3=R3*theS3;
953    G4double I12=I1+I2;
954    G4double rand=(I12+I3)*G4UniformRand();
955    if     (rand<I1 )
956    {
957      G4double ran=R1*G4UniformRand();
958      if(ran>1.) ran=1.;
959      q2=-std::log(1.-ran)/theB1;
960    }
961    else if(rand<I12)
962    {
963      G4double ran=R2*G4UniformRand();
964      if(ran>1.) ran=1.;
965      q2=-std::log(1.-ran);
966      if(q2<0.) q2=0.;
967      q2=std::pow(q2,third)/theB2;
968    }
969    else
970    {
971      G4double ran=R3*G4UniformRand();
972      if(ran>1.) ran=1.;
973      q2=-std::log(1.-ran)/theB3;
974    }
975  }
976  else if(PDG==2212 || PDG==2112)
977  {
978    G4double a=tgZ+tgN;
979#ifdef tdebug
980    G4cout<<"G4QElCS::GetExT: a="<<a<<",t="<<lastTM<<",S1="<<theS1<<",B1="<<theB1<<",SS="
981          <<theSS<<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",S4="
982          <<theS4<<",B4="<<theB4<<G4endl;
983#endif
984    G4double E1=lastTM*(theB1+lastTM*theSS);
985    G4double R1=(1.-std::exp(-E1));
986    G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
987#ifdef tdebug
988    G4double ts1=-std::log(1.-R1)/theB1;
989    if(std::fabs(tss)>1.e-7) ts1=(std::sqrt(theB1*(theB1+(tss+tss)*ts1))-theB1)/tss;
990    G4double ds1=(ts1-lastTM)/lastTM;
991    if(ds1>.0001)
992      G4cout<<"*Warning*G4QElCS::GetExT:1a "<<ts1<<"#"<<lastTM<<",d="<<ds1
993            <<",R1="<<R1<<",E1="<<E1<<G4endl;
994#endif
995    G4double tm2=lastTM*lastTM;
996    G4double E2=lastTM*tm2*theB2;                   // power 3 for lowA, 5 for HighA (1st)
997    if(a>6.5)E2*=tm2;                               // for heavy nuclei
998    G4double R2=(1.-std::exp(-E2));
999#ifdef tdebug
1000    G4double ts2=-std::log(1.-R2)/theB2;
1001    if(a<6.5)ts2=std::pow(ts2,third);
1002    else     ts2=std::pow(ts2,fifth);
1003    G4double ds2=std::fabs(ts2-lastTM)/lastTM;
1004    if(ds2>.0001)
1005      G4cout<<"*Warning*G4QElCS::GetExT:2a "<<ts2<<"#"<<lastTM<<",d="<<ds2
1006            <<",R2="<<R2<<",E2="<<E2<<G4endl;
1007#endif
1008    G4double E3=lastTM*theB3;
1009    if(a>6.5)E3*=tm2*tm2*tm2;                       // power 1 for lowA, 7 (2nd) for HighA
1010    G4double R3=(1.-std::exp(-E3));
1011#ifdef tdebug
1012    G4double ts3=-std::log(1.-R3)/theB3;
1013    if(a>6.5)ts3=std::pow(ts3,sevth);
1014    G4double ds3=std::fabs(ts3-lastTM)/lastTM;
1015    if(ds3>.0001)
1016      G4cout<<"*Warning*G4QElCS::GetExT:3a "<<ts3<<"#"<<lastTM<<",d="<<ds3
1017            <<",R3="<<R3<<",E3="<<E3<<G4endl;
1018#endif
1019    G4double E4=lastTM*theB4;
1020    G4double R4=(1.-std::exp(-E4));
1021#ifdef tdebug
1022    G4double ts4=-std::log(1.-R4)/theB4;
1023    G4double ds4=std::fabs(ts4-lastTM)/lastTM;
1024    if(ds4>.0001)
1025      G4cout<<"*Warning*G4QElCS::GetExT:4a "<<ts4<<"#"<<lastTM<<",d="<<ds4
1026            <<",R4="<<R4<<",E4="<<E4<<G4endl;
1027#endif
1028    G4double I1=R1*theS1;
1029    G4double I2=R2*theS2;
1030    G4double I3=R3*theS3;
1031    G4double I4=R4*theS4;
1032    G4double I12=I1+I2;
1033    G4double I13=I12+I3;
1034    G4double rand=(I13+I4)*G4UniformRand();
1035#ifdef tdebug
1036    G4cout<<"G4QElCS::GtExT:1="<<I1<<",2="<<I2<<",3="<<I3<<",4="<<I4<<",r="<<rand<<G4endl;
1037#endif
1038    if(rand<I1)
1039    {
1040      G4double ran=R1*G4UniformRand();
1041      if(ran>1.) ran=1.;
1042      q2=-std::log(1.-ran)/theB1;
1043      if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
1044#ifdef tdebug
1045      G4cout<<"G4QElCS::GetExT:Q2="<<q2<<",ss="<<tss/2<<",b1="<<theB1<<",t1="<<ts1<<G4endl;
1046#endif
1047    }
1048    else if(rand<I12)
1049    {
1050      G4double ran=R2*G4UniformRand();
1051      if(ran>1.) ran=1.;
1052      q2=-std::log(1.-ran)/theB2;
1053      if(q2<0.) q2=0.;
1054      if(a<6.5) q2=std::pow(q2,third);
1055      else      q2=std::pow(q2,fifth);
1056#ifdef tdebug
1057      G4cout<<"G4QElCS::GetExT: Q2="<<q2<<", r2="<<R2<<", b2="<<theB2<<",t2="<<ts2<<G4endl;
1058#endif
1059    }
1060    else if(rand<I13)
1061    {
1062      G4double ran=R3*G4UniformRand();
1063      if(ran>1.) ran=1.;
1064      q2=-std::log(1.-ran)/theB3;
1065      if(q2<0.) q2=0.;
1066      if(a>6.5) q2=std::pow(q2,sevth);
1067#ifdef tdebug
1068      G4cout<<"G4QElCS::GetExT:Q2="<<q2<<", r3="<<R2<<", b3="<<theB2<<",t3="<<ts2<<G4endl;
1069#endif
1070    }
1071    else
1072    {
1073      G4double ran=R4*G4UniformRand();
1074      if(ran>1.) ran=1.;
1075      q2=-std::log(1.-ran)/theB4;
1076      if(a<6.5) q2=lastTM-q2;                    // u reduced for lightA (starts from 0)
1077#ifdef tdebug
1078      G4cout<<"G4QElCS::GetExT:Q2="<<q2<<",m="<<lastTM<<",b4="<<theB3<<",t4="<<ts3<<G4endl;
1079#endif
1080    }
1081  }
1082  else
1083  {
1084    G4cout<<"*Error*G4QElasticCrossSection::GetExchangeT: PDG="<<PDG<<", Z="<<tgZ<<", N="
1085          <<tgN<<", while it is defined only for PDG=2212 & PDG==2112"<<G4endl;
1086    throw G4QException("G4QElasticCrossSection::GetExchangeT: nA and pA are implemented");
1087  }
1088  if(q2<0.) q2=0.;
1089  if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QElasticCrossSect::GetExchangeT: -t="<<q2<<G4endl;
1090  if(q2>lastTM)
1091  {
1092#ifdef tdebug
1093    G4cout<<"*Warning*G4QElasticCrossSect::GetExT:-t="<<q2<<">"<<lastTM<<G4endl;
1094#endif
1095    q2=lastTM;
1096  }
1097  return q2*GeVSQ;
1098}
1099
1100// Returns B in independent units (MeV^-2) (all internal calculations are in GeV) see ExT
1101G4double G4QElasticCrossSection::GetSlope(G4int tgZ, G4int tgN, G4int PDG)
1102{
1103  static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
1104#ifdef tdebug
1105  G4cout<<"G4QElasticCS::GetSlope:"<<onlyCS<<", Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<G4endl;
1106#endif
1107  if(onlyCS) G4cout<<"*Warning*G4QElasticCrossSection::GetSlope: onlyCS=true"<<G4endl;
1108  if(lastLP<-4.3) return 0.;          // S-wave for p<14 MeV/c (kinE<.1MeV)
1109  if(PDG!=2212 && PDG!=2112)
1110  {
1111    G4cout<<"*Error*G4QElasticCrossSection::GetSlope: PDG="<<PDG<<", Z="<<tgZ<<", N="
1112          <<tgN<<", while it is defined only for PDG=2212 & PDG=2112"<<G4endl;
1113    throw G4QException("G4QElasticCrossSection::GetSlope: nA and pA are implemented");
1114  }
1115  if(theB1<0.) theB1=0.;
1116  if(!(theB1>=-1.||theB1<=1.))G4cout<<"*NAN*G4QElasticCrossSect::Getslope:"<<theB1<<G4endl;
1117  return theB1/GeVSQ;
1118}
1119
1120// Returns half max(Q2=-t) in independent units (MeV^2)
1121G4double G4QElasticCrossSection::GetHMaxT()
1122{
1123  static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
1124  return lastTM*HGeVSQ;
1125}
1126
1127// lastLP is used, so calculating tables, one need to remember and then recover lastLP
1128G4double G4QElasticCrossSection::GetTabValues(G4double lp, G4int PDG, G4int tgZ,G4int tgN)
1129{
1130  if(tgZ==0 && tgN==1)                 // Temporary change for Quasi-Elastic
1131  {
1132    tgZ=1;
1133    tgN=1;
1134    if     (PDG==2212) PDG=2112;
1135    else if(PDG==2112) PDG=2212;
1136  }
1137  if(tgZ<0 || tgZ>92)
1138  {
1139    G4cout<<"*Warning*G4QElasticCS::GetTabValue: (1-92) No isotopes for Z="<<tgZ<<G4endl;
1140    return 0.;
1141  }
1142  G4int iZ=tgZ-1; // Z index
1143  if(iZ<0)
1144  {
1145    iZ=0;         // conversion of the neutron target to the proton target
1146    tgZ=1;
1147    tgN=0;
1148  }
1149  //if(nN[iZ][0] < 0)
1150  //{
1151#ifdef isodebug
1152  //  G4cout<<"*Warning*G4QElasticCS::GetTabValue: No isotopes for Z="<<tgZ<<G4endl;
1153#endif
1154  //  return 0.;
1155  //}
1156#ifdef pdebug
1157  G4cout<<"G4QElasticCS::GetTabVal: lp="<<lp<<",Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<G4endl;
1158#endif
1159  G4double p=std::exp(lp);              // momentum
1160  G4double sp=std::sqrt(p);             // sqrt(p)
1161  G4double p2=p*p;           
1162  G4double p3=p2*p;
1163  G4double p4=p3*p;
1164  if ( (PDG == 2112 && tgZ == 1 && tgN == 0) || 
1165       (PDG == 2212 && tgZ == 0 && tgN == 1) ) {       // np/pn
1166
1167    G4double ssp=std::sqrt(sp);           // sqrt(sqrt(p))=p^.25
1168    G4double p2s=p2*sp;
1169    G4double dl1=lp-lastPAR[3];
1170    theSS=lastPAR[21];
1171    theS1=(lastPAR[9]+lastPAR[10]*dl1*dl1+lastPAR[11]/p)/(1.+lastPAR[12]/p4)
1172          +lastPAR[13]/(p4+lastPAR[14]);
1173    theB1=(lastPAR[17]+lastPAR[18]/(p4*p4+lastPAR[19]*p3))/(1.+lastPAR[20]/p4);
1174    theS2=(lastPAR[15]+lastPAR[16]/p4/p)/p3;
1175    theB2=lastPAR[22]/(p*sp+lastPAR[23]); 
1176    theS3=0.;
1177    theB3=0.; 
1178    theS4=0.;
1179    theB4=0.; 
1180#ifdef tdebug
1181    G4cout<<"G4QElasticCS::GetTableValues:(np) TM="<<lastTM<<",S1="<<theS1<<",B1="<<theB1
1182          <<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS1<<",B3="<<theB1<<G4endl;
1183#endif
1184    // Returns the total elastic pp cross-section (to avoid spoiling lastSIG)
1185    return lastPAR[0]/(p2s+lastPAR[1]*p+lastPAR[2]/ssp)+lastPAR[4]/p
1186           +(lastPAR[5]+lastPAR[6]*dl1*dl1+lastPAR[7]/p)/(1.+lastPAR[8]/p4);
1187
1188  } else if ( (PDG == 2212 && tgZ == 1 && tgN == 0) || 
1189              (PDG == 2112 && tgZ == 0 && tgN == 1) ) { // pp/nn
1190
1191    G4double p2s=p2*sp;
1192    G4double dl1=lp-lastPAR[3];
1193    G4double dl2=lp-lastPAR[8];
1194    theSS=lastPAR[31];
1195    theS1=(lastPAR[9]+lastPAR[10]*dl2*dl2)/(1.+lastPAR[11]/p4/p)+
1196          (lastPAR[12]/p2+lastPAR[13]*p)/(p4+lastPAR[14]*sp);
1197    theB1=lastPAR[15]*std::pow(p,lastPAR[16])/(1.+lastPAR[17]/p3);
1198    theS2=lastPAR[18]+lastPAR[19]/(p4+lastPAR[20]*p);
1199    theB2=lastPAR[21]+lastPAR[22]/(p4+lastPAR[23]/sp); 
1200    theS3=lastPAR[24]+lastPAR[25]/(p4*p4+lastPAR[26]*p2+lastPAR[27]);
1201    theB3=lastPAR[28]+lastPAR[29]/(p4+lastPAR[30]); 
1202    theS4=0.;
1203    theB4=0.; 
1204#ifdef tdebug
1205    G4cout<<"G4QElasticCS::GetTableValues:(pp) TM="<<lastTM<<",S1="<<theS1<<",B1="<<theB1
1206          <<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS1<<",B3="<<theB1<<G4endl;
1207#endif
1208    // Returns the total elastic pp cross-section (to avoid spoiling lastSIG)
1209    return lastPAR[0]/p2s/(1.+lastPAR[7]/p2s)+(lastPAR[1]+lastPAR[2]*dl1*dl1+lastPAR[4]/p)
1210                                                   /(1.+lastPAR[5]*lp)/(1.+lastPAR[6]/p4);
1211  }
1212  else if(PDG==2212 || PDG==2112) // n/p + A
1213  {
1214    G4double p5=p4*p;
1215    G4double p6=p5*p;
1216    G4double p8=p6*p2;
1217    G4double p10=p8*p2;
1218    G4double p12=p10*p2;
1219    G4double p16=p8*p8;
1220    //G4double p24=p16*p8;
1221    G4double dl=lp-5.;
1222    G4double a=tgZ+tgN;
1223    G4double pah=std::pow(p,a/2);
1224    G4double pa=pah*pah;
1225    G4double pa2=pa*pa;
1226    if(a<6.5)
1227    {
1228      theS1=lastPAR[9]/(1.+lastPAR[10]*p4*pa)+lastPAR[11]/(p4+lastPAR[12]*p4/pa2)+
1229            (lastPAR[13]*dl*dl+lastPAR[14])/(1.+lastPAR[15]/p2);
1230      theB1=(lastPAR[16]+lastPAR[17]*p2)/(p4+lastPAR[18]/pah)+lastPAR[19];
1231      theSS=lastPAR[20]/(1.+lastPAR[21]/p2)+lastPAR[22]/(p6/pa+lastPAR[23]/p16);
1232      theS2=lastPAR[24]/(pa/p2+lastPAR[25]/p4)+lastPAR[26];
1233      theB2=lastPAR[27]*std::pow(p,lastPAR[28])+lastPAR[29]/(p8+lastPAR[30]/p16);
1234      theS3=lastPAR[31]/(pa*p+lastPAR[32]/pa)+lastPAR[33];
1235      theB3=lastPAR[34]/(p3+lastPAR[35]/p6)+lastPAR[36]/(1.+lastPAR[37]/p2);
1236      theS4=p2*(pah*lastPAR[38]*std::exp(-pah*lastPAR[39])+
1237                lastPAR[40]/(1.+lastPAR[41]*std::pow(p,lastPAR[42])));
1238      theB4=lastPAR[43]*pa/p2/(1.+pa*lastPAR[44]);
1239#ifdef tdebug
1240      G4cout<<"G4QElCS::GetTabV: lA, p="<<p<<",S1="<<theS1<<",B1="<<theB1<<",SS="<<theSS
1241            <<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",S4="<<theS4
1242            <<",B4="<<theB4<<G4endl;
1243#endif
1244    }
1245    else
1246    {
1247      theS1=lastPAR[9]/(1.+lastPAR[10]/p4)+lastPAR[11]/(p4+lastPAR[12]/p2)+
1248            lastPAR[13]/(p5+lastPAR[14]/p16);
1249      theB1=(lastPAR[15]/p8+lastPAR[19])/(p+lastPAR[16]/std::pow(p,lastPAR[20]))+
1250            lastPAR[17]/(1.+lastPAR[18]/p4);
1251      theSS=lastPAR[21]/(p4/std::pow(p,lastPAR[23])+lastPAR[22]/p4);
1252      theS2=lastPAR[24]/p4/(std::pow(p,lastPAR[25])+lastPAR[26]/p12)+lastPAR[27];
1253      theB2=lastPAR[28]/std::pow(p,lastPAR[29])+lastPAR[30]/std::pow(p,lastPAR[31]);
1254      theS3=lastPAR[32]/std::pow(p,lastPAR[35])/(1.+lastPAR[36]/p12)+
1255            lastPAR[33]/(1.+lastPAR[34]/p6);
1256      theB3=lastPAR[37]/p8+lastPAR[38]/p2+lastPAR[39]/(1.+lastPAR[40]/p8);
1257      theS4=(lastPAR[41]/p4+lastPAR[46]/p)/(1.+lastPAR[42]/p10)+
1258            (lastPAR[43]+lastPAR[44]*dl*dl)/(1.+lastPAR[45]/p12);
1259      theB4=lastPAR[47]/(1.+lastPAR[48]/p)+lastPAR[49]*p4/(1.+lastPAR[50]*p5);
1260#ifdef tdebug
1261      G4cout<<"G4QElCS::GetTabV: hA, p="<<p<<",S1="<<theS1<<",B1="<<theB1<<",SS="<<theSS
1262            <<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",S4="<<theS4
1263            <<",B4="<<theB4<<G4endl;
1264#endif
1265    }
1266    // Returns the total elastic (n/p)A cross-section (to avoid spoiling lastSIG)
1267#ifdef tdebug
1268    G4cout<<"G4QElCS::GetTabV: PDG="<<PDG<<",P="<<p<<",N="<<tgN<<",Z="<<tgZ<<G4endl;
1269#endif
1270    //if(PDG!=2112 || !tgN || p>.1)      // No Low Energy neutron addition
1271    if(PDG==2212)                        // Unique formula for protons
1272    {
1273      G4double rp16=lastPAR[6]/p16;
1274      return (lastPAR[0]*dl*dl+lastPAR[1])/(1.+lastPAR[2]/p)+lastPAR[3]/(p4+lastPAR[4]/p2)
1275             + lastPAR[5]/(p5+rp16*rp16)+lastPAR[7]/(p4+lastPAR[8]/p4);
1276    }
1277    else if(a<6.5)                       // neutrons on light nuclei
1278     return (lastPAR[0]*dl*dl+lastPAR[1])/(1.+lastPAR[2]/p)+lastPAR[3]/(p4+lastPAR[4]/p2)+
1279            lastPAR[5]/(p5+lastPAR[6]/p8)+lastPAR[7]/(p4+lastPAR[8]);
1280    else                                 // neutrons on heavy nuclei
1281     return (lastPAR[0]*dl*dl+lastPAR[1])/(1.+lastPAR[2]/p)+lastPAR[3]/(p4+lastPAR[4]/p2)+
1282            lastPAR[5]/(p8+lastPAR[6]/p8)+lastPAR[7]/(p2+lastPAR[8]);
1283    //else
1284    //{
1285    //  G4int nI=nn[iZ];
1286    //  for(G4int i=0; i<nI; i++) if(tgN==nN[iZ][i])
1287    //  {
1288    //    G4double kE=p2/dNM;
1289#ifdef tdebug
1290    //  G4cout<<"G4QECS::GTV:P="<<p<<",E="<<kE<<",X="<<nX[iZ][i]<<",T="<<nT[iZ][i]<<G4endl;
1291#endif
1292    //    if(kE<nT[iZ][i]) return 0.;    // 0 below the threshold (in MeV)
1293    //   if(p<2.) return nX[iZ][i];     // At very low momentum(<2MeV/c) -> only LECS
1294    //    return lastPAR[3]/(p4+lastPAR[4]/p2)+lastPAR[5]/(p5+rp16*rp16)+
1295    //           lastPAR[7]/(p4+lastPAR[8]/p4)+nX[iZ][i]/(1.+lastPAR[51]*p16);
1296    //  }
1297    //G4cerr<<"*G4QElCS::GTV:PDG="<<PDG<<",Z="<<tgZ<<",N="<<tgN<<": Sig not found"<<G4endl;
1298    //}
1299  }
1300  else
1301  {
1302    G4cout<<"*Error*G4QElasticCrossSection::GetTabValues: PDG="<<PDG<<", Z="<<tgZ<<", N="
1303          <<tgN<<", while it is defined only for PDG=2212/2112, Z>0, N>0"<<G4endl;
1304    throw G4QException("G4QElasticCrossSection::GetTabValues: only NA is implemented");
1305  }
1306  return 0.;
1307} // End of GetTableValues
1308
1309// Returns max -t=Q2 (GeV^2) for the momentum pP(GeV) and the target nucleus (tgN,tgZ)
1310G4double G4QElasticCrossSection::GetQ2max(G4int PDG, G4int tgZ, G4int tgN, G4double pP)
1311{
1312  static const G4double mNeut= G4QPDGCode(2112).GetMass()*.001; // MeV to GeV
1313  static const G4double mProt= G4QPDGCode(2212).GetMass()*.001; // MeV to GeV
1314  //static const G4double mLamb= G4QPDGCode(3122).GetMass()*.001; // MeV to GeV
1315  //static const G4double mHe3 = G4QPDGCode(2112).GetNuclMass(2,1,0)*.001; // MeV to GeV
1316  //static const G4double mAlph = G4QPDGCode(2112).GetNuclMass(2,2,0)*.001; // MeV to GeV
1317  //static const G4double mDeut = G4QPDGCode(2112).GetNuclMass(1,1,0)*.001; // MeV to GeV
1318  static const G4double mProt2= mProt*mProt;
1319  static const G4double mNeut2= mNeut*mNeut;
1320  //static const G4double mDeut2= mDeut*mDeut;
1321  if(tgZ==0 && tgN==1)                 // Temporary change for Quasi-Elastic
1322  {
1323    tgZ=1;
1324    tgN=1;
1325    if     (PDG==2212) PDG=2112;
1326    else if(PDG==2112) PDG=2212;
1327  }
1328  G4double pP2=pP*pP;                                 // squared momentum of the projectile
1329  if(PDG==2212 && tgZ==1 && tgN==0)
1330  {
1331    G4double tMid=std::sqrt(pP2+mProt2)*mProt-mProt2; // CMS 90deg value of -t=Q2 (GeV^2)
1332    return tMid+tMid;
1333  }
1334  else if(PDG==2112 && tgZ==0 && tgN==1)
1335  {
1336    G4double tMid=std::sqrt(pP2+mNeut2)*mNeut-mNeut2;  // CMS 90deg value of -t=Q2 (GeV^2)
1337    return tMid+tMid;
1338  }
1339  else if(PDG==2112 && (tgZ || tgN))                   // ---> nA
1340  {
1341    G4double mt=mProt;                                 // Target mass in GeV
1342    if(tgN||tgZ>1) mt=G4QPDGCode(90000000+tgZ*1000+tgN).GetMass()*.001; // Target mass GeV
1343    G4double dmt=mt+mt;
1344    G4double s=dmt*std::sqrt(pP2+mNeut2)+mNeut2+mt*mt; // Mondelstam s (in GeV^2)
1345    return dmt*dmt*pP2/s;
1346  }
1347  else if(PDG==2212 && (tgZ || tgN))                   // ---> pA
1348  {
1349    G4double mt=G4QPDGCode(90000000+tgZ*1000+tgN).GetMass()*.001; // Target mass in GeV
1350    G4double dmt=mt+mt;
1351    G4double s=dmt*std::sqrt(pP2+mProt2)+mProt2+mt*mt; // Mondelstam s
1352    return dmt*dmt*pP2/s;
1353  }
1354  else
1355  {
1356    G4cout<<"*Error*G4QElasticCrossSection::GetQ2max: PDG="<<PDG<<", Z="<<tgZ<<", N="
1357          <<tgN<<", while it is defined only for p,n projectiles & Z_target>0"<<G4endl;
1358    throw G4QException("G4QElasticCrossSection::GetQ2max: only nA & pA are implemented");
1359  }
1360}
Note: See TracBrowser for help on using the repository browser.