source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/cross_sections/src/G4QProtonElasticCrossSection.cc @ 1316

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

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