source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/src/G4QElasticCrossSection.cc @ 962

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

update processes

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