G4QProtonElasticCrossSection.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 //
00030 // G4 Physics class: G4QProtonElasticCrossSection for pA elastic cross sections
00031 // Created: M.V. Kossov, CERN/ITEP(Moscow), 10-OCT-01
00032 // The last update: M.V. Kossov, CERN/ITEP (Moscow) 12-Jan-10 (from G4QElCrSect)
00033 //
00034 // -------------------------------------------------------------------------------
00035 // Short description: Interaction cross-sections for the G4QElastic process
00036 // -------------------------------------------------------------------------------
00037 
00038 //#define debug
00039 //#define isodebug
00040 //#define pdebug
00041 //#define ppdebug
00042 //#define tdebug
00043 //#define sdebug
00044 
00045 #include "G4QProtonElasticCrossSection.hh"
00046 #include "G4SystemOfUnits.hh"
00047 
00048 // Initialization of the static parameters
00049 const G4int G4QProtonElasticCrossSection::nPoints=128;//#ofPt in the AMDB table(>anyPar)(D)
00050 const G4int G4QProtonElasticCrossSection::nLast=nPoints-1;// theLastElement in the table(D)
00051 G4double  G4QProtonElasticCrossSection::lPMin=-8.;  // Min tabulated logarithmicMomentum(D)
00052 G4double  G4QProtonElasticCrossSection::lPMax= 8.;  // Max tabulated logarithmicMomentum(D)
00053 G4double  G4QProtonElasticCrossSection::dlnP=(lPMax-lPMin)/nLast;// LogStep in the table(D)
00054 G4bool    G4QProtonElasticCrossSection::onlyCS=true;// Flag toCalculateOnlyCS(not Si/Bi)(L)
00055 G4double  G4QProtonElasticCrossSection::lastSIG=0.; // Last calculated cross section    (L)
00056 G4double  G4QProtonElasticCrossSection::lastLP=-10.;// Last log(mom_ofTheIncidentHadron)(L)
00057 G4double  G4QProtonElasticCrossSection::lastTM=0.;  // Last t_maximum                   (L)
00058 G4double  G4QProtonElasticCrossSection::theSS=0.;   // The Last sq.slope of 1st difr.Max(L)
00059 G4double  G4QProtonElasticCrossSection::theS1=0.;   // The Last mantissa of 1st difr.Max(L)
00060 G4double  G4QProtonElasticCrossSection::theB1=0.;   // The Last slope of 1st difruct.Max(L)
00061 G4double  G4QProtonElasticCrossSection::theS2=0.;   // The Last mantissa of 2nd difr.Max(L)
00062 G4double  G4QProtonElasticCrossSection::theB2=0.;   // The Last slope of 2nd difruct.Max(L)
00063 G4double  G4QProtonElasticCrossSection::theS3=0.;   // The Last mantissa of 3d difr. Max(L)
00064 G4double  G4QProtonElasticCrossSection::theB3=0.;   // The Last slope of 3d difruct. Max(L)
00065 G4double  G4QProtonElasticCrossSection::theS4=0.;   // The Last mantissa of 4th difr.Max(L)
00066 G4double  G4QProtonElasticCrossSection::theB4=0.;   // The Last slope of 4th difruct.Max(L)
00067 G4int     G4QProtonElasticCrossSection::lastTZ=0;   // Last atomic number of the target
00068 G4int     G4QProtonElasticCrossSection::lastTN=0;   // Last # of neutrons in the target
00069 G4double  G4QProtonElasticCrossSection::lastPIN=0.; // Last initialized max momentum
00070 G4double* G4QProtonElasticCrossSection::lastCST=0;  // Elastic cross-section table
00071 G4double* G4QProtonElasticCrossSection::lastPAR=0;  // Parameters for FunctionalCalculation
00072 G4double* G4QProtonElasticCrossSection::lastSST=0;  // E-dep of sq.slope of the 1st dif.Max
00073 G4double* G4QProtonElasticCrossSection::lastS1T=0;  // E-dep of mantissa of the 1st dif.Max
00074 G4double* G4QProtonElasticCrossSection::lastB1T=0;  // E-dep of the slope of the 1st difMax
00075 G4double* G4QProtonElasticCrossSection::lastS2T=0;  // E-dep of mantissa of the 2nd difrMax
00076 G4double* G4QProtonElasticCrossSection::lastB2T=0;  // E-dep of the slope of the 2nd difMax
00077 G4double* G4QProtonElasticCrossSection::lastS3T=0;  // E-dep of mantissa of the 3d difr.Max
00078 G4double* G4QProtonElasticCrossSection::lastB3T=0;  // E-dep of the slope of the 3d difrMax
00079 G4double* G4QProtonElasticCrossSection::lastS4T=0;  // E-dep of mantissa of the 4th difrMax
00080 G4double* G4QProtonElasticCrossSection::lastB4T=0;  // E-dep of the slope of the 4th difMax
00081 G4int     G4QProtonElasticCrossSection::lastN=0;    // The last N of calculated nucleus
00082 G4int     G4QProtonElasticCrossSection::lastZ=0;    // The last Z of calculated nucleus
00083 G4double  G4QProtonElasticCrossSection::lastP=0.;   // Last used in cross section Momentum
00084 G4double  G4QProtonElasticCrossSection::lastTH=0.;  // Last threshold momentum
00085 G4double  G4QProtonElasticCrossSection::lastCS=0.;  // Last value of the Cross Section
00086 G4int     G4QProtonElasticCrossSection::lastI=0;    // The last position in the DAMDB
00087 
00088 std::vector<G4double*> G4QProtonElasticCrossSection::PAR; // Vector of pars for functCalcul
00089 std::vector<G4double*> G4QProtonElasticCrossSection::CST; // Vector of cross-section table
00090 std::vector<G4double*> G4QProtonElasticCrossSection::SST; // Vector of the 1st SquaredSlope
00091 std::vector<G4double*> G4QProtonElasticCrossSection::S1T; // Vector of the 1st mantissa
00092 std::vector<G4double*> G4QProtonElasticCrossSection::B1T; // Vector of the 1st slope
00093 std::vector<G4double*> G4QProtonElasticCrossSection::S2T; // Vector of the 2nd mantissa
00094 std::vector<G4double*> G4QProtonElasticCrossSection::B2T; // Vector of the 2nd slope
00095 std::vector<G4double*> G4QProtonElasticCrossSection::S3T; // Vector of the 3d mantissa
00096 std::vector<G4double*> G4QProtonElasticCrossSection::B3T; // Vector of the 3d slope
00097 std::vector<G4double*> G4QProtonElasticCrossSection::S4T; // Vector of the 4th mantissa (g)
00098 std::vector<G4double*> G4QProtonElasticCrossSection::B4T; // Vector of the 4th slope (glor)
00099 
00100 G4QProtonElasticCrossSection::G4QProtonElasticCrossSection()
00101 {
00102 }
00103 
00104 G4QProtonElasticCrossSection::~G4QProtonElasticCrossSection()
00105 {
00106   std::vector<G4double*>::iterator pos;
00107   for (pos=CST.begin(); pos<CST.end(); pos++)
00108   { delete [] *pos; }
00109   CST.clear();
00110   for (pos=PAR.begin(); pos<PAR.end(); pos++)
00111   { delete [] *pos; }
00112   PAR.clear();
00113   for (pos=SST.begin(); pos<SST.end(); pos++)
00114   { delete [] *pos; }
00115   SST.clear();
00116   for (pos=S1T.begin(); pos<S1T.end(); pos++)
00117   { delete [] *pos; }
00118   S1T.clear();
00119   for (pos=B1T.begin(); pos<B1T.end(); pos++)
00120   { delete [] *pos; }
00121   B1T.clear();
00122   for (pos=S2T.begin(); pos<S2T.end(); pos++)
00123   { delete [] *pos; }
00124   S2T.clear();
00125   for (pos=B2T.begin(); pos<B2T.end(); pos++)
00126   { delete [] *pos; }
00127   B2T.clear();
00128   for (pos=S3T.begin(); pos<S3T.end(); pos++)
00129   { delete [] *pos; }
00130   S3T.clear();
00131   for (pos=B3T.begin(); pos<B3T.end(); pos++)
00132   { delete [] *pos; }
00133   B3T.clear();
00134   for (pos=S4T.begin(); pos<S4T.end(); pos++)
00135   { delete [] *pos; }
00136   S4T.clear();
00137   for (pos=B4T.begin(); pos<B4T.end(); pos++)
00138   { delete [] *pos; }
00139   B4T.clear();
00140 }
00141 
00142 // Returns Pointer to the G4VQCrossSection class
00143 G4VQCrossSection* G4QProtonElasticCrossSection::GetPointer()
00144 {
00145   static G4QProtonElasticCrossSection theCrossSection;//*StaticBody of theQEl CrossSection*
00146   return &theCrossSection;
00147 }
00148 
00149 // The main member function giving the collision cross section (P is in IU, CS is in mb)
00150 // Make pMom in independent units ! (Now it is MeV)
00151 G4double G4QProtonElasticCrossSection::GetCrossSection(G4bool fCS,G4double pMom, G4int tgZ,
00152                                                        G4int tgN, G4int pPDG)
00153 {
00154   static std::vector <G4int>    colN;  // Vector of N for calculated nuclei (isotops)
00155   static std::vector <G4int>    colZ;  // Vector of Z for calculated nuclei (isotops)
00156   static std::vector <G4double> colP;  // Vector of last momenta for the reaction
00157   static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
00158   static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
00159   // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
00160   G4double pEn=pMom;
00161   onlyCS=fCS;
00162 #ifdef pdebug
00163   G4cout<<"G4QPElCS::GetCS:>>> f="<<fCS<<", p="<<pMom<<", Z="<<tgZ<<"("<<lastZ<<") ,N="
00164         <<tgN<<"("<<lastN<<"), T="<<pEn<<"("<<lastTH<<")"<<",Sz="<<colN.size()<<G4endl;
00165   //CalculateCrossSection(fCS,-27,j,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00166 #endif
00167   if(pPDG!=2212)
00168   {
00169     G4cout<<"*Warning*G4QProtonElCS::GetCS: *** Found pPDG="<<pPDG<<" =--=> CS=0"<<G4endl;
00170     //CalculateCrossSection(fCS,-27,j,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00171     return 0.;                         // projectile PDG=0 is a mistake (?!) @@
00172   }
00173   G4bool in=false;                   // By default the isotope must be found in the AMDB
00174   lastP   = 0.;                      // New momentum history (nothing to compare with)
00175   lastN   = tgN;                     // The last N of the calculated nucleus
00176   lastZ   = tgZ;                     // The last Z of the calculated nucleus
00177   lastI   = colN.size();             // Size of the Associative Memory DB in the heap
00178   if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
00179   {                                  // The nucleus with projPDG is found in AMDB
00180     if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
00181     {
00182       lastI=i;
00183       lastTH =colTH[i];              // Last THreshold (A-dependent)
00184 #ifdef pdebug
00185       G4cout<<"G4QElCS::GetCS:*Found* P="<<pMom<<",Threshold="<<lastTH<<",i="<<i<<G4endl;
00186       //CalculateCrossSection(fCS,-27,i,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00187 #endif
00188       if(pEn<=lastTH)
00189       {
00190 #ifdef pdebug
00191         G4cout<<"G4QElCS::GetCS:Found T="<<pEn<<" < Threshold="<<lastTH<<",CS=0"<<G4endl;
00192         //CalculateCrossSection(fCS,-27,i,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00193 #endif
00194         return 0.;                   // Energy is below the Threshold value
00195       }
00196       lastP  =colP [i];              // Last Momentum  (A-dependent)
00197       lastCS =colCS[i];              // Last CrossSect (A-dependent)
00198       //  if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
00199       if(lastP == pMom)              // Do not recalculate
00200       {
00201 #ifdef pdebug
00202         G4cout<<"G4QElCS::GetCS:P="<<pMom<<",CS="<<lastCS*millibarn<<G4endl;
00203 #endif
00204         CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // Update param's only
00205         return lastCS*millibarn;     // Use theLastCS
00206       }
00207       in = true;                       // This is the case when the isotop is found in DB
00208       // Momentum pMom is in IU ! @@ Units
00209 #ifdef pdebug
00210       G4cout<<"G4QElCS::G:UpdateDB P="<<pMom<<",f="<<fCS<<",I="<<lastI<<",i="<<i<<G4endl;
00211 #endif
00212       lastCS=CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // read & update
00213 #ifdef pdebug
00214       G4cout<<"G4QElCS::GetCrosSec: *****> New (inDB) Calculated CS="<<lastCS<<G4endl;
00215       //CalculateCrossSection(fCS,-27,i,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00216 #endif
00217       if(lastCS<=0. && pEn>lastTH)    // Correct the threshold
00218       {
00219 #ifdef pdebug
00220         G4cout<<"G4QElCS::GetCS: New T="<<pEn<<"(CS=0) > Threshold="<<lastTH<<G4endl;
00221 #endif
00222         lastTH=pEn;
00223       }
00224       break;                           // Go out of the LOOP with found lastI
00225     }
00226 #ifdef pdebug
00227     G4cout<<"---G4QElCrossSec::GetCrosSec:pPDG="<<pPDG<<",i="<<i<<",N="<<colN[i]
00228           <<",Z["<<i<<"]="<<colZ[i]<<G4endl;
00229     //CalculateCrossSection(fCS,-27,i,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00230 #endif
00231   } // End of attampt to find the nucleus in DB
00232   if(!in)                            // This nucleus has not been calculated previously
00233   {
00234 #ifdef pdebug
00235     G4cout<<"G4QElCS::GetCrosSec:CalcNew P="<<pMom<<",f="<<fCS<<",lastI="<<lastI<<G4endl;
00236 #endif
00238     lastCS=CalculateCrossSection(fCS,0,lastI,pPDG,lastZ,lastN,pMom);//calculate&create
00239     if(lastCS<=0.)
00240     {
00241       lastTH = ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
00242 #ifdef pdebug
00243       G4cout<<"G4QElCrossSection::GetCrossSect: NewThresh="<<lastTH<<",T="<<pEn<<G4endl;
00244 #endif
00245       if(pEn>lastTH)
00246       {
00247 #ifdef pdebug
00248         G4cout<<"G4QElCS::GetCS: First T="<<pEn<<"(CS=0) > Threshold="<<lastTH<<G4endl;
00249 #endif
00250         lastTH=pEn;
00251       }
00252     }
00253 #ifdef pdebug
00254     G4cout<<"G4QElCS::GetCrosSec: New CS="<<lastCS<<",lZ="<<lastN<<",lN="<<lastZ<<G4endl;
00255     //CalculateCrossSection(fCS,-27,lastI,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00256 #endif
00257     colN.push_back(tgN);
00258     colZ.push_back(tgZ);
00259     colP.push_back(pMom);
00260     colTH.push_back(lastTH);
00261     colCS.push_back(lastCS);
00262 #ifdef pdebug
00263     G4cout<<"G4QElCS::GetCS:1st,P="<<pMom<<"(MeV),CS="<<lastCS*millibarn<<"(mb)"<<G4endl;
00264     //CalculateCrossSection(fCS,-27,lastI,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00265 #endif
00266     return lastCS*millibarn;
00267   } // End of creation of the new set of parameters
00268   else
00269   {
00270 #ifdef pdebug
00271     G4cout<<"G4QElCS::GetCS: Update lastI="<<lastI<<G4endl;
00272 #endif
00273     colP[lastI]=pMom;
00274     colCS[lastI]=lastCS;
00275   }
00276 #ifdef pdebug
00277   G4cout<<"G4QElCS::GetCrSec:End,P="<<pMom<<"(MeV),CS="<<lastCS*millibarn<<"(mb)"<<G4endl;
00278   //CalculateCrossSection(fCS,-27,lastI,pPDG,lastZ,lastN,pMom); // DUMMY TEST
00279   G4cout<<"G4QElCS::GetCrSec:***End***, onlyCS="<<onlyCS<<G4endl;
00280 #endif
00281   return lastCS*millibarn;
00282 }
00283 
00284 // Calculation of total elastic cross section (p in IU, CS in mb) @@ Units (?)
00285 // F=0 - create AMDB, F=-1 - read&update AMDB, F=1 - update AMDB (sinchro with higher AMDB)
00286 G4double G4QProtonElasticCrossSection::CalculateCrossSection(G4bool CS, G4int F, G4int I,
00287                                              G4int PDG, G4int tgZ, G4int tgN, G4double pIU)
00288 {
00289   // *** Begin of Associative Memory DB for acceleration of the cross section calculations
00290   static std::vector <G4double>  PIN;   // Vector of max initialized log(P) in the table
00291   // *** End of Static Definitions (Associative Memory Data Base) ***
00292   G4double pMom=pIU/GeV;                // All calculations are in GeV
00293   onlyCS=CS;                            // Flag to calculate only CS (not Si/Bi)
00294 #ifdef pdebug
00295   G4cout<<"G4QProtonElasticCrossS::CalcCS:->onlyCS="<<onlyCS<<",F="<<F<<",p="<<pIU<<G4endl;
00296 #endif
00297   lastLP=std::log(pMom);                // Make a logarithm of the momentum for calculation
00298   if(F)                                 // This isotope was found in AMDB =>RETRIEVE/UPDATE
00299   {
00300     if(F<0)                             // the AMDB must be loded
00301     {
00302       lastPIN = PIN[I];                 // Max log(P) initialised for this table set
00303       lastPAR = PAR[I];                 // Pointer to the parameter set
00304       lastCST = CST[I];                 // Pointer to the total sross-section table
00305       lastSST = SST[I];                 // Pointer to the first squared slope
00306       lastS1T = S1T[I];                 // Pointer to the first mantissa
00307       lastB1T = B1T[I];                 // Pointer to the first slope
00308       lastS2T = S2T[I];                 // Pointer to the second mantissa
00309       lastB2T = B2T[I];                 // Pointer to the second slope
00310       lastS3T = S3T[I];                 // Pointer to the third mantissa
00311       lastB3T = B3T[I];                 // Pointer to the rhird slope
00312       lastS4T = S4T[I];                 // Pointer to the 4-th mantissa
00313       lastB4T = B4T[I];                 // Pointer to the 4-th slope
00314 #ifdef pdebug
00315       G4cout<<"G4QElasticCS::CalcCS: DB is updated for I="<<I<<",*,PIN4="<<PIN[4]<<G4endl;
00316 #endif
00317     }
00318 #ifdef pdebug
00319     G4cout<<"G4QProtonElasticCrossS::CalcCS:*read*, LP="<<lastLP<<",PIN="<<lastPIN<<G4endl;
00320 #endif
00321     if(lastLP>lastPIN && lastLP<lPMax)
00322     {
00323       lastPIN=GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);// Can update upper logP-Limit in tabs
00324 #ifdef pdebug
00325       G4cout<<"G4QElCS::CalcCS:*updated(I)*,LP="<<lastLP<<"<IN["<<I<<"]="<<lastPIN<<G4endl;
00326 #endif
00327       PIN[I]=lastPIN;                   // Remember the new P-Limit of the tables
00328     }
00329   }
00330   else                                  // This isotope wasn't initialized => CREATE
00331   {
00332     lastPAR = new G4double[nPoints];    // Allocate memory for parameters of CS function
00333     lastPAR[nLast]=0;                   // Initialization for VALGRIND
00334     lastCST = new G4double[nPoints];    // Allocate memory for Tabulated CS function    
00335     lastSST = new G4double[nPoints];    // Allocate memory for Tabulated first sqaredSlope 
00336     lastS1T = new G4double[nPoints];    // Allocate memory for Tabulated first mantissa 
00337     lastB1T = new G4double[nPoints];    // Allocate memory for Tabulated first slope    
00338     lastS2T = new G4double[nPoints];    // Allocate memory for Tabulated second mantissa
00339     lastB2T = new G4double[nPoints];    // Allocate memory for Tabulated second slope   
00340     lastS3T = new G4double[nPoints];    // Allocate memory for Tabulated third mantissa 
00341     lastB3T = new G4double[nPoints];    // Allocate memory for Tabulated third slope    
00342     lastS4T = new G4double[nPoints];    // Allocate memory for Tabulated 4-th mantissa 
00343     lastB4T = new G4double[nPoints];    // Allocate memory for Tabulated 4-th slope    
00344 #ifdef pdebug
00345     G4cout<<"G4QProtonElasticCrossS::CalcCS:*ini*,lastLP="<<lastLP<<",min="<<lPMin<<G4endl;
00346 #endif
00347     lastPIN = GetPTables(lastLP,lPMin,PDG,tgZ,tgN); // Returns the new P-limit for tables
00348 #ifdef pdebug
00349     G4cout<<"G4QProtElCS::CCS:i,Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<",LP"<<lastPIN<<G4endl;
00350 #endif
00351     PIN.push_back(lastPIN);             // Fill parameters of CS function to AMDB
00352     PAR.push_back(lastPAR);             // Fill parameters of CS function to AMDB
00353     CST.push_back(lastCST);             // Fill Tabulated CS function to AMDB    
00354     SST.push_back(lastSST);             // Fill Tabulated first sq.slope to AMDB 
00355     S1T.push_back(lastS1T);             // Fill Tabulated first mantissa to AMDB 
00356     B1T.push_back(lastB1T);             // Fill Tabulated first slope to AMDB    
00357     S2T.push_back(lastS2T);             // Fill Tabulated second mantissa to AMDB 
00358     B2T.push_back(lastB2T);             // Fill Tabulated second slope to AMDB    
00359     S3T.push_back(lastS3T);             // Fill Tabulated third mantissa to AMDB 
00360     B3T.push_back(lastB3T);             // Fill Tabulated third slope to AMDB    
00361     S4T.push_back(lastS4T);             // Fill Tabulated 4-th mantissa to AMDB 
00362     B4T.push_back(lastB4T);             // Fill Tabulated 4-th slope to AMDB    
00363   } // End of creation/update of the new set of parameters and tables
00364   // =--------= NOW Update (if necessary) and Calculate the Cross Section =------------=
00365 #ifdef pdebug
00366   G4cout<<"G4QElCS::CalcCS:?update?,LP="<<lastLP<<",IN="<<lastPIN<<",ML="<<lPMax<<G4endl;
00367 #endif
00368   if(lastLP>lastPIN && lastLP<lPMax)
00369   {
00370     lastPIN = GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);
00371 #ifdef pdebug
00372     G4cout<<"G4QElCS::CalcCS: *updated(O)*, LP="<<lastLP<<" < IN="<<lastPIN<<G4endl;
00373 #endif
00374   }
00375 #ifdef pdebug
00376   G4cout<<"G4QElastCS::CalcCS: lastLP="<<lastLP<<",lPM="<<lPMin<<",lPIN="<<lastPIN<<G4endl;
00377 #endif
00378   if(!onlyCS) lastTM=GetQ2max(PDG, tgZ, tgN, pMom); // Calculate (-t)_max=Q2_max (GeV2)
00379 #ifdef pdebug
00380   G4cout<<"G4QElasticCrosSec::CalcCS:oCS="<<onlyCS<<",-t="<<lastTM<<", p="<<lastLP<<G4endl;
00381 #endif
00382   if(lastLP>lPMin && lastLP<=lastPIN)   // Linear fit is made using precalculated tables
00383   {
00384     if(lastLP==lastPIN)
00385     {
00386       G4double shift=(lastLP-lPMin)/dlnP+.000001; // Log distance from lPMin
00387       G4int    blast=static_cast<int>(shift); // this is a bin number of the lower edge (0)
00388       if(blast<0 || blast>=nLast) G4cout<<"G4QEleastCS::CCS:b="<<blast<<","<<nLast<<G4endl;
00389       lastSIG = lastCST[blast];
00390       if(!onlyCS)                       // Skip the differential cross-section parameters
00391       {
00392         theSS  = lastSST[blast];
00393         theS1  = lastS1T[blast];
00394         theB1  = lastB1T[blast];
00395         theS2  = lastS2T[blast];
00396         theB2  = lastB2T[blast];
00397         theS3  = lastS3T[blast];
00398         theB3  = lastB3T[blast];
00399         theS4  = lastS4T[blast];
00400         theB4  = lastB4T[blast];
00401       }
00402 #ifdef pdebug
00403       G4cout<<"G4QProtonElasticCrossS::CalculateCS:(E) S1="<<theS1<<", B1="<<theB1<<G4endl;
00404 #endif
00405     }
00406     else
00407     {
00408       G4double shift=(lastLP-lPMin)/dlnP;        // a shift from the beginning of the table
00409       G4int    blast=static_cast<int>(shift);    // the lower bin number
00410       if(blast<0)   blast=0;
00411       if(blast>=nLast) blast=nLast-1;            // low edge of the last bin
00412       shift-=blast;                              // step inside the unit bin
00413       G4int lastL=blast+1;                       // the upper bin number
00414       G4double SIGL=lastCST[blast];              // the basic value of the cross-section
00415       lastSIG= SIGL+shift*(lastCST[lastL]-SIGL); // calculated total elastic cross-section
00416 #ifdef pdebug
00417       G4cout<<"G4QElCS::CalcCrossSection: Sig="<<lastSIG<<", P="<<pMom<<", Z="<<tgZ<<", N="
00418             <<tgN<<", PDG="<<PDG<<", onlyCS="<<onlyCS<<G4endl;
00419 #endif
00420       if(!onlyCS)                       // Skip the differential cross-section parameters
00421       {
00422         G4double SSTL=lastSST[blast];           // the low bin of the first squared slope
00423         theSS=SSTL+shift*(lastSST[lastL]-SSTL); // the basic value of the first sq.slope
00424         G4double S1TL=lastS1T[blast];           // the low bin of the first mantissa
00425         theS1=S1TL+shift*(lastS1T[lastL]-S1TL); // the basic value of the first mantissa
00426         G4double B1TL=lastB1T[blast];           // the low bin of the first slope
00427 #ifdef pdebug
00428         G4cout<<"G4QElCS::CalcCrossSection:bl="<<blast<<",ls="<<lastL<<",SL="<<S1TL<<",SU="
00429               <<lastS1T[lastL]<<",BL="<<B1TL<<",BU="<<lastB1T[lastL]<<G4endl;
00430 #endif
00431         theB1=B1TL+shift*(lastB1T[lastL]-B1TL); // the basic value of the first slope
00432         G4double S2TL=lastS2T[blast];           // the low bin of the second mantissa
00433         theS2=S2TL+shift*(lastS2T[lastL]-S2TL); // the basic value of the second mantissa
00434         G4double B2TL=lastB2T[blast];           // the low bin of the second slope
00435         theB2=B2TL+shift*(lastB2T[lastL]-B2TL); // the basic value of the second slope
00436         G4double S3TL=lastS3T[blast];           // the low bin of the third mantissa
00437         theS3=S3TL+shift*(lastS3T[lastL]-S3TL); // the basic value of the third mantissa
00438 #ifdef pdebug
00439         G4cout<<"G4QElCS::CCS: s3l="<<S3TL<<",sh3="<<shift<<",s3h="<<lastS3T[lastL]<<",b="
00440               <<blast<<",l="<<lastL<<G4endl;
00441 #endif
00442         G4double B3TL=lastB3T[blast];           // the low bin of the third slope
00443         theB3=B3TL+shift*(lastB3T[lastL]-B3TL); // the basic value of the third slope
00444         G4double S4TL=lastS4T[blast];           // the low bin of the 4-th mantissa
00445         theS4=S4TL+shift*(lastS4T[lastL]-S4TL); // the basic value of the 4-th mantissa
00446 #ifdef pdebug
00447         G4cout<<"G4QElCS::CCS: s4l="<<S4TL<<",sh4="<<shift<<",s4h="<<lastS4T[lastL]<<",b="
00448               <<blast<<",l="<<lastL<<G4endl;
00449 #endif
00450         G4double B4TL=lastB4T[blast];           // the low bin of the 4-th slope
00451         theB4=B4TL+shift*(lastB4T[lastL]-B4TL); // the basic value of the 4-th slope
00452       }
00453 #ifdef pdebug
00454       G4cout<<"G4QProtonElasticCrossS::CalculateCS:(I) S1="<<theS1<<", B1="<<theB1<<G4endl;
00455 #endif
00456     }
00457   }
00458   else lastSIG=GetTabValues(lastLP, PDG, tgZ, tgN); // Direct calculation beyond the table
00459   if(lastSIG<0.) lastSIG = 0.;                   // @@ a Warning print can be added
00460 #ifdef pdebug
00461   G4cout<<"G4QProtonElasticCrossSection::CalculateCS: END, onlyCS="<<onlyCS<<G4endl;
00462 #endif
00463   return lastSIG;
00464 }
00465 
00466 // It has parameter sets for all tZ/tN/PDG, using them the tables can be created/updated
00467 G4double G4QProtonElasticCrossSection::GetPTables(G4double LP, G4double ILP, G4int PDG,
00468                                                   G4int tgZ, G4int tgN)
00469 {
00470   // @@ At present all nA==pA ---------> Each neucleus can have not more than 51 parameters
00471   static const G4double pwd=2727;
00472   const G4int n_npel=24;                // #of parameters for np-elastic (<nPoints=128)
00473   const G4int n_ppel=32;                // #of parameters for pp-elastic (<nPoints=128)
00474   //                      -0- -1-  -2- -3- -4-  -5- -6- -7- -8- -9--10--11--12--13- -14-
00475   G4double np_el[n_npel]={12.,.05,.0001,5.,.35,6.75,.14,19.,.6,6.75,.14,13.,.14,.6,.00013,
00476                           75.,.001,7.2,4.32,.012,2.5,0.0,12.,.34};
00477   //                      -15--16--17- -18- -19--20--21--22--23-
00478   //                       -0-   -1-  -2- -3- -4- -5-  -6-  -7-  -8--9--10--11--12--13-
00479   G4double pp_el[n_ppel]={2.865,18.9,.6461,3.,9.,.425,.4276,.0022,5.,74.,3.,3.4,.2,.17,
00480                           .001,8.,.055,3.64,5.e-5,4000.,1500.,.46,1.2e6,3.5e6,5.e-5,1.e10,
00481                           8.5e8,1.e10,1.1,3.4e6,6.8e6,0.};
00482   //                      -14--15- -16- -17- -18-  -19- -20- -21- -22-  -23-   -24-  -25-
00483   //                       -26- -27-  -28- -29- -30- -31-
00484   if(PDG==2212)
00485   {
00486     // -- Total pp elastic cross section cs & s1/b1 (main), s2/b2 (tail1), s3/b3 (tail2) --
00487     //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|
00488     //CS=2.865/p2s/(1+.0022/p2s)+(18.9+.6461*dl1*dl1+9./p)/(1.+.425*lp)/(1.+.4276/p4);
00489     //   par(0)       par(7)     par(1) par(2)      par(4)      par(5)         par(6)
00490     //dl2=lp-5., s1=(74.+3.*dl2*dl2)/(1+3.4/p4/p)+(.2/p2+17.*p)/(p4+.001*sp),
00491     //     par(8) par(9) par(10)        par(11)   par(12)par(13)    par(14)
00492     // b1=8.*p**.055/(1.+3.64/p3); s2=5.e-5+4000./(p4+1500.*p); b2=.46+1.2e6/(p4+3.5e6/sp);
00493     // par(15) par(16)  par(17)     par(18) par(19)  par(20)   par(21) par(22)  par(23)
00494     // s3=5.e-5+1.e10/(p4*p4+8.5e8*p2+1.e10); b3=1.1+3.4e6/(p4+6.8e6); ss=0.
00495     //  par(24) par(25)     par(26)  par(27) par(28) par(29)  par(30)   par(31)
00496     //
00497     if(lastPAR[nLast]!=pwd) // A unique flag to avoid the repeatable definition
00498     {
00499       if ( tgZ == 0 && tgN == 1 )
00500       {
00501         for (G4int ip=0; ip<n_npel; ip++) lastPAR[ip]=np_el[ip]; // pn
00502 
00503       }
00504       else if ( tgZ == 1 && tgN == 0 )
00505       {
00506         for (G4int ip=0; ip<n_ppel; ip++) lastPAR[ip]=pp_el[ip]; // pp
00507       }
00508       else
00509       {
00510         G4double a=tgZ+tgN;
00511         G4double sa=std::sqrt(a);
00512         G4double ssa=std::sqrt(sa);
00513         G4double asa=a*sa;
00514         G4double a2=a*a;
00515         G4double a3=a2*a;
00516         G4double a4=a3*a;
00517         G4double a5=a4*a;
00518         G4double a6=a4*a2;
00519         G4double a7=a6*a;
00520         G4double a8=a7*a;
00521         G4double a9=a8*a;
00522         G4double a10=a5*a5;
00523         G4double a12=a6*a6;
00524         G4double a14=a7*a7;
00525         G4double a16=a8*a8;
00526         G4double a17=a16*a;
00527         G4double a20=a16*a4;
00528         G4double a32=a16*a16;
00529         // Reaction cross-section parameters (pel=peh_fit.f)
00530         lastPAR[0]=5./(1.+22./asa);                                          // p1
00531         lastPAR[1]=4.8*std::pow(a,1.14)/(1.+3.6/a3);                         // p2
00532         lastPAR[2]=1./(1.+4.E-3*a4)+2.E-6*a3/(1.+1.3E-6*a3);                 // p3
00533         lastPAR[3]=1.3*a;                                                    // p4
00534         lastPAR[4]=3.E-8*a3/(1.+4.E-7*a4);                                   // p5
00535         lastPAR[5]=.07*asa/(1.+.009*a2);                                     // p6
00536         lastPAR[6]=(3.+3.E-16*a20)/(1.+a20*(2.E-16/a+3.E-19*a));             // p7 (11)
00537         lastPAR[7]=(5.E-9*a4*sa+.27/a)/(1.+5.E16/a20)/(1.+6.E-9*a4)+.015/a2; // p8
00538         lastPAR[8]=(.001*a+.07/a)/(1.+5.E13/a16+5.E-7*a3)+.0003/sa;          // p9 (10)
00539         // @@ the differential cross-section is parameterized separately for A>6 & A<7
00540         if(a<6.5)
00541         {
00542           G4double a28=a16*a12;
00543           // The main pre-exponent      (pel_sg)
00544           lastPAR[ 9]=4000*a;                                // p1
00545           lastPAR[10]=1.2e7*a8+380*a17;                      // p2
00546           lastPAR[11]=.7/(1.+4.e-12*a16);                    // p3
00547           lastPAR[12]=2.5/a8/(a4+1.e-16*a32);                // p4
00548           lastPAR[13]=.28*a;                                 // p5
00549           lastPAR[14]=1.2*a2+2.3;                            // p6
00550           lastPAR[15]=3.8/a;                                 // p7
00551           // The main slope             (pel_sl)
00552           lastPAR[16]=.01/(1.+.0024*a5);                     // p1
00553           lastPAR[17]=.2*a;                                  // p2
00554           lastPAR[18]=9.e-7/(1.+.035*a5);                    // p3
00555           lastPAR[19]=(42.+2.7e-11*a16)/(1.+.14*a);          // p4
00556           // The main quadratic         (pel_sh)
00557           lastPAR[20]=2.25*a3;                               // p1
00558           lastPAR[21]=18.;                                   // p2
00559           lastPAR[22]=2.4e-3*a8/(1.+2.6e-4*a7);              // p3
00560           lastPAR[23]=3.5e-36*a32*a8/(1.+5.e-15*a32/a);      // p4
00561           // The 1st max pre-exponent   (pel_qq)
00562           lastPAR[24]=1.e5/(a8+2.5e12/a16);                  // p1
00563           lastPAR[25]=8.e7/(a12+1.e-27*a28*a28);             // p2 
00564           lastPAR[26]=.0006*a3;                              // p3
00565           // The 1st max slope          (pel_qs)
00566           lastPAR[27]=10.+4.e-8*a12*a;                       // p1
00567           lastPAR[28]=.114;                                  // p2
00568           lastPAR[29]=.003;                                  // p3
00569           lastPAR[30]=2.e-23;                                // p4
00570           // The effective pre-exponent (pel_ss)
00571           lastPAR[31]=1./(1.+.0001*a8);                      // p1
00572           lastPAR[32]=1.5e-4/(1.+5.e-6*a12);                 // p2
00573           lastPAR[33]=.03;                                   // p3
00574           // The effective slope        (pel_sb)
00575           lastPAR[34]=a/2;                                   // p1
00576           lastPAR[35]=2.e-7*a4;                              // p2
00577           lastPAR[36]=4.;                                    // p3
00578           lastPAR[37]=64./a3;                                // p4
00579           // The gloria pre-exponent    (pel_us)
00580           lastPAR[38]=1.e8*std::exp(.32*asa);                // p1
00581           lastPAR[39]=20.*std::exp(.45*asa);                 // p2
00582           lastPAR[40]=7.e3+2.4e6/a5;                         // p3
00583           lastPAR[41]=2.5e5*std::exp(.085*a3);               // p4
00584           lastPAR[42]=2.5*a;                                 // p5
00585           // The gloria slope           (pel_ub)
00586           lastPAR[43]=920.+.03*a8*a3;                        // p1
00587           lastPAR[44]=93.+.0023*a12;                         // p2
00588 #ifdef pdebug
00589          G4cout<<"G4QElCS::CalcCS:la "<<lastPAR[38]<<", "<<lastPAR[39]<<", "<<lastPAR[40]
00590                <<", "<<lastPAR[42]<<", "<<lastPAR[43]<<", "<<lastPAR[44]<<G4endl;
00591 #endif
00592         }
00593         else
00594         {
00595           G4double p1a10=2.2e-28*a10;
00596           G4double r4a16=6.e14/a16;
00597           G4double s4a16=r4a16*r4a16;
00598           // a24
00599           // a36
00600           // The main pre-exponent      (peh_sg)
00601           lastPAR[ 9]=4.5*std::pow(a,1.15);                  // p1
00602           lastPAR[10]=.06*std::pow(a,.6);                    // p2
00603           lastPAR[11]=.6*a/(1.+2.e15/a16);                   // p3
00604           lastPAR[12]=.17/(a+9.e5/a3+1.5e33/a32);            // p4
00605           lastPAR[13]=(.001+7.e-11*a5)/(1.+4.4e-11*a5);      // p5
00606           lastPAR[14]=(p1a10*p1a10+2.e-29)/(1.+2.e-22*a12);  // p6
00607           // The main slope             (peh_sl)
00608           lastPAR[15]=400./a12+2.e-22*a9;                    // p1
00609           lastPAR[16]=1.e-32*a12/(1.+5.e22/a14);             // p2
00610           lastPAR[17]=1000./a2+9.5*sa*ssa;                   // p3
00611           lastPAR[18]=4.e-6*a*asa+1.e11/a16;                 // p4
00612           lastPAR[19]=(120./a+.002*a2)/(1.+2.e14/a16);       // p5
00613           lastPAR[20]=9.+100./a;                             // p6
00614           // The main quadratic         (peh_sh)
00615           lastPAR[21]=.002*a3+3.e7/a6;                       // p1
00616           lastPAR[22]=7.e-15*a4*asa;                         // p2
00617           lastPAR[23]=9000./a4;                              // p3
00618           // The 1st max pre-exponent   (peh_qq)
00619           lastPAR[24]=.0011*asa/(1.+3.e34/a32/a4);           // p1
00620           lastPAR[25]=1.e-5*a2+2.e14/a16;                    // p2
00621           lastPAR[26]=1.2e-11*a2/(1.+1.5e19/a12);            // p3
00622           lastPAR[27]=.016*asa/(1.+5.e16/a16);               // p4
00623           // The 1st max slope          (peh_qs)
00624           lastPAR[28]=.002*a4/(1.+7.e7/std::pow(a-6.83,14)); // p1
00625           lastPAR[29]=2.e6/a6+7.2/std::pow(a,.11);           // p2
00626           lastPAR[30]=11.*a3/(1.+7.e23/a16/a8);              // p3
00627           lastPAR[31]=100./asa;                              // p4
00628           // The 2nd max pre-exponent   (peh_ss)
00629           lastPAR[32]=(.1+4.4e-5*a2)/(1.+5.e5/a4);           // p1
00630           lastPAR[33]=3.5e-4*a2/(1.+1.e8/a8);                // p2
00631           lastPAR[34]=1.3+3.e5/a4;                           // p3
00632           lastPAR[35]=500./(a2+50.)+3;                       // p4
00633           lastPAR[36]=1.e-9/a+s4a16*s4a16;                   // p5
00634           // The 2nd max slope          (peh_sb)
00635           lastPAR[37]=.4*asa+3.e-9*a6;                       // p1
00636           lastPAR[38]=.0005*a5;                              // p2
00637           lastPAR[39]=.002*a5;                               // p3
00638           lastPAR[40]=10.;                                   // p4
00639           // The effective pre-exponent (peh_us)
00640           lastPAR[41]=.05+.005*a;                            // p1
00641           lastPAR[42]=7.e-8/sa;                              // p2
00642           lastPAR[43]=.8*sa;                                 // p3
00643           lastPAR[44]=.02*sa;                                // p4
00644           lastPAR[45]=1.e8/a3;                               // p5
00645           lastPAR[46]=3.e32/(a32+1.e32);                     // p6
00646           // The effective slope        (peh_ub)
00647           lastPAR[47]=24.;                                   // p1
00648           lastPAR[48]=20./sa;                                // p2
00649           lastPAR[49]=7.e3*a/(sa+1.);                        // p3
00650           lastPAR[50]=900.*sa/(1.+500./a3);                  // p4
00651 #ifdef pdebug
00652          G4cout<<"G4QElCS::CalcCS:ha "<<lastPAR[41]<<", "<<lastPAR[42]<<", "<<lastPAR[43]
00653                <<", "<<lastPAR[44]<<", "<<lastPAR[45]<<", "<<lastPAR[46]<<G4endl;
00654 #endif
00655         }
00656         // Parameter for lowEnergyNeutrons
00657         lastPAR[51]=1.e15+2.e27/a4/(1.+2.e-18*a16);
00658       }
00659       lastPAR[nLast]=pwd;
00660       // and initialize the zero element of the table
00661       G4double lp=lPMin;                                      // ln(momentum)
00662       G4bool memCS=onlyCS;                                    // ??
00663       onlyCS=false;
00664       lastCST[0]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables
00665       onlyCS=memCS;
00666       lastSST[0]=theSS;
00667       lastS1T[0]=theS1;
00668       lastB1T[0]=theB1;
00669       lastS2T[0]=theS2;
00670       lastB2T[0]=theB2;
00671       lastS3T[0]=theS3;
00672       lastB3T[0]=theB3;
00673       lastS4T[0]=theS4;
00674       lastB4T[0]=theB4;
00675 #ifdef pdebug
00676       G4cout<<"G4QProtonElasticCrossSection::GetPTables:ip=0(init), lp="<<lp<<",S1="<<theS1
00677             <<",B1="<<theB1<<",S2="<<theS2<<",B2="<<theB3<<",S3="<<theS3
00678             <<",B3="<<theB3<<",S4="<<theS4<<",B4="<<theB4<<G4endl;
00679 #endif
00680     }
00681     if(LP>ILP)
00682     {
00683       G4int ini = static_cast<int>((ILP-lPMin+.000001)/dlnP)+1; // already inited till this
00684       if(ini<0) ini=0;
00685       if(ini<nPoints)
00686       {
00687         G4int fin = static_cast<int>((LP-lPMin)/dlnP)+1; // final bin of initialization
00688         if(fin>=nPoints) fin=nLast;               // Limit of the tabular initialization
00689         if(fin>=ini)
00690         {
00691           G4double lp=0.;
00692           for(G4int ip=ini; ip<=fin; ip++)        // Calculate tabular CS,S1,B1,S2,B2,S3,B3
00693           {
00694             lp=lPMin+ip*dlnP;                     // ln(momentum)
00695             G4bool memCS=onlyCS;
00696             onlyCS=false;
00697             lastCST[ip]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables (ret CS)
00698             onlyCS=memCS;
00699             lastSST[ip]=theSS;
00700             lastS1T[ip]=theS1;
00701             lastB1T[ip]=theB1;
00702             lastS2T[ip]=theS2;
00703             lastB2T[ip]=theB2;
00704             lastS3T[ip]=theS3;
00705             lastB3T[ip]=theB3;
00706             lastS4T[ip]=theS4;
00707             lastB4T[ip]=theB4;
00708 #ifdef pdebug
00709             G4cout<<"G4QProtonElasticCrossSection::GetPTables:ip="<<ip<<",lp="<<lp<<",S1="
00710                   <<theS1<<",B1="<<theB1<<",S2="<<theS2<<",B2="<<theB2<<",S3="
00711                   <<theS3<<",B3="<<theB3<<",S4="<<theS4<<",B4="<<theB4<<G4endl;
00712 #endif
00713           }
00714           return lp;
00715         }
00716         else G4cout<<"*Warning*G4QProtonElasticCrossSection::GetPTables: PDG="<<PDG<<", Z="
00717                    <<tgZ<<", N="<<tgN<<", i="<<ini<<" > fin="<<fin<<", LP="<<LP<<" > ILP="
00718                    <<ILP<<" nothing is done!"<<G4endl;
00719       }
00720       else G4cout<<"*Warning*G4QProtonElasticCrossSection::GetPTables: PDG="<<PDG<<", Z="
00721                  <<tgZ<<", N="<<tgN<<", i="<<ini<<">= max="<<nPoints<<", LP="<<LP
00722                  <<" > ILP="<<ILP<<", lPMax="<<lPMax<<" nothing is done!"<<G4endl;
00723     }
00724 #ifdef pdebug
00725     else G4cout<<"*Warning*G4QProtonElasticCrossSection::GetPTables:PDG="<<PDG<<", Z="<<tgZ
00726                <<", N="<<tgN<<", LP="<<LP<<" <= ILP="<<ILP<<" nothing is done!"<<G4endl;
00727 #endif
00728   }
00729   else
00730   {
00731     // G4cout<<"*Error*G4QProtonElasticCrossSection::GetPTables: PDG="<<PDG<<", Z="<<tgZ
00732     //       <<", N="<<tgN<<", while it is defined only for PDG=2212"<<G4endl;
00733     // throw G4QException("G4QProtonElasticCrossSection::GetPTables: only pA're implemented");
00734     G4ExceptionDescription ed;
00735     ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
00736        << ", while it is defined only for PDG=2212 (p)" << G4endl;
00737     G4Exception("G4QProtonElasticCrossSection::GetPTables()", "HAD_CHPS_0000",
00738                 FatalException, ed);
00739   }
00740   return ILP;
00741 }
00742 
00743 // Returns Q2=-t in independent units (MeV^2) (all internal calculations are in GeV)
00744 G4double G4QProtonElasticCrossSection::GetExchangeT(G4int tgZ, G4int tgN, G4int PDG)
00745 {
00746   static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
00747   static const G4double third=1./3.;
00748   static const G4double fifth=1./5.;
00749   static const G4double sevth=1./7.;
00750 #ifdef tdebug
00751   G4cout<<"G4QProtElCS::GetExcT: F="<<onlyCS<<",Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<G4endl;
00752 #endif
00753   if(PDG!=2212) G4cout<<"**Warning*G4QProtonElasticCrossSection::GetExT:PDG="<<PDG<<G4endl;
00754   if(onlyCS) G4cout<<"**Warning*G4QProtonElasticCrossSection::GetExchanT:onlyCS=1"<<G4endl;
00755   if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
00756   G4double q2=0.;
00757   if(tgZ==1 && tgN==0)                // ===> p+p=p+p
00758   {
00759 #ifdef tdebug
00760     G4cout<<"G4QElasticCS::GetExchangeT: TM="<<lastTM<<",S1="<<theS1<<",B1="<<theB1<<",S2="
00761           <<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",GeV2="<<GeVSQ<<G4endl;
00762 #endif
00763     G4double E1=lastTM*theB1;
00764     G4double R1=(1.-std::exp(-E1));
00765 #ifdef tdebug
00766     G4double ts1=-std::log(1.-R1)/theB1;
00767     G4double ds1=std::fabs(ts1-lastTM)/lastTM;
00768     if(ds1>.0001)
00769       G4cout<<"*Warning*G4QElCS::GetExT:1p "<<ts1<<"#"<<lastTM<<",d="<<ds1
00770             <<",R1="<<R1<<",E1="<<E1<<G4endl;
00771 #endif
00772     G4double E2=lastTM*theB2;
00773     G4double R2=(1.-std::exp(-E2*E2*E2));
00774 #ifdef tdebug
00775     G4double ts2=std::pow(-std::log(1.-R2),.333333333)/theB2;
00776     G4double ds2=std::fabs(ts2-lastTM)/lastTM;
00777     if(ds2>.0001)
00778       G4cout<<"*Warning*G4QElCS::GetExT:2p "<<ts2<<"#"<<lastTM<<",d="<<ds2
00779             <<",R2="<<R2<<",E2="<<E2<<G4endl;
00780 #endif
00781     G4double E3=lastTM*theB3;
00782     G4double R3=(1.-std::exp(-E3));
00783 #ifdef tdebug
00784     G4double ts3=-std::log(1.-R3)/theB3;
00785     G4double ds3=std::fabs(ts3-lastTM)/lastTM;
00786     if(ds3>.0001)
00787       G4cout<<"*Warning*G4QElCS::GetExT:3p "<<ts3<<"#"<<lastTM<<",d="<<ds3
00788             <<",R3="<<R1<<",E3="<<E3<<G4endl;
00789 #endif
00790     G4double I1=R1*theS1/theB1;
00791     G4double I2=R2*theS2;
00792     G4double I3=R3*theS3;
00793     G4double I12=I1+I2;
00794     G4double rand=(I12+I3)*G4UniformRand();
00795     if     (rand<I1 )
00796     {
00797       G4double ran=R1*G4UniformRand();
00798       if(ran>1.) ran=1.;
00799       q2=-std::log(1.-ran)/theB1;
00800     }
00801     else if(rand<I12)
00802     {
00803       G4double ran=R2*G4UniformRand();
00804       if(ran>1.) ran=1.;
00805       q2=-std::log(1.-ran);
00806       if(q2<0.) q2=0.;
00807       q2=std::pow(q2,third)/theB2;
00808     }
00809     else
00810     {
00811       G4double ran=R3*G4UniformRand();
00812       if(ran>1.) ran=1.;
00813       q2=-std::log(1.-ran)/theB3;
00814     }
00815   }
00816   else
00817   {
00818     G4double a=tgZ+tgN;
00819 #ifdef tdebug
00820     G4cout<<"G4QElCS::GetExT: a="<<a<<",t="<<lastTM<<",S1="<<theS1<<",B1="<<theB1<<",SS="
00821           <<theSS<<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",S4="
00822           <<theS4<<",B4="<<theB4<<G4endl;
00823 #endif
00824     G4double E1=lastTM*(theB1+lastTM*theSS);
00825     G4double R1=(1.-std::exp(-E1));
00826     G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
00827 #ifdef tdebug
00828     G4double ts1=-std::log(1.-R1)/theB1;
00829     if(std::fabs(tss)>1.e-7) ts1=(std::sqrt(theB1*(theB1+(tss+tss)*ts1))-theB1)/tss;
00830     G4double ds1=(ts1-lastTM)/lastTM;
00831     if(ds1>.0001)
00832       G4cout<<"*Warning*G4QElCS::GetExT:1a "<<ts1<<"#"<<lastTM<<",d="<<ds1
00833             <<",R1="<<R1<<",E1="<<E1<<G4endl;
00834 #endif
00835     G4double tm2=lastTM*lastTM;
00836     G4double E2=lastTM*tm2*theB2;                   // power 3 for lowA, 5 for HighA (1st)
00837     if(a>6.5)E2*=tm2;                               // for heavy nuclei
00838     G4double R2=(1.-std::exp(-E2));
00839 #ifdef tdebug
00840     G4double ts2=-std::log(1.-R2)/theB2;
00841     if(a<6.5)ts2=std::pow(ts2,third);
00842     else     ts2=std::pow(ts2,fifth);
00843     G4double ds2=std::fabs(ts2-lastTM)/lastTM;
00844     if(ds2>.0001)
00845       G4cout<<"*Warning*G4QElCS::GetExT:2a "<<ts2<<"#"<<lastTM<<",d="<<ds2
00846             <<",R2="<<R2<<",E2="<<E2<<G4endl;
00847 #endif
00848     G4double E3=lastTM*theB3;
00849     if(a>6.5)E3*=tm2*tm2*tm2;                       // power 1 for lowA, 7 (2nd) for HighA
00850     G4double R3=(1.-std::exp(-E3));
00851 #ifdef tdebug
00852     G4double ts3=-std::log(1.-R3)/theB3;
00853     if(a>6.5)ts3=std::pow(ts3,sevth);
00854     G4double ds3=std::fabs(ts3-lastTM)/lastTM;
00855     if(ds3>.0001)
00856       G4cout<<"*Warning*G4QElCS::GetExT:3a "<<ts3<<"#"<<lastTM<<",d="<<ds3
00857             <<",R3="<<R3<<",E3="<<E3<<G4endl;
00858 #endif
00859     G4double E4=lastTM*theB4;
00860     G4double R4=(1.-std::exp(-E4));
00861 #ifdef tdebug
00862     G4double ts4=-std::log(1.-R4)/theB4;
00863     G4double ds4=std::fabs(ts4-lastTM)/lastTM;
00864     if(ds4>.0001)
00865       G4cout<<"*Warning*G4QElCS::GetExT:4a "<<ts4<<"#"<<lastTM<<",d="<<ds4
00866             <<",R4="<<R4<<",E4="<<E4<<G4endl;
00867 #endif
00868     G4double I1=R1*theS1;
00869     G4double I2=R2*theS2;
00870     G4double I3=R3*theS3;
00871     G4double I4=R4*theS4;
00872     G4double I12=I1+I2;
00873     G4double I13=I12+I3;
00874     G4double rand=(I13+I4)*G4UniformRand();
00875 #ifdef tdebug
00876     G4cout<<"G4QElCS::GtExT:1="<<I1<<",2="<<I2<<",3="<<I3<<",4="<<I4<<",r="<<rand<<G4endl;
00877 #endif
00878     if(rand<I1)
00879     {
00880       G4double ran=R1*G4UniformRand();
00881       if(ran>1.) ran=1.;
00882       q2=-std::log(1.-ran)/theB1;
00883       if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
00884 #ifdef tdebug
00885       G4cout<<"G4QElCS::GetExT:Q2="<<q2<<",ss="<<tss/2<<",b1="<<theB1<<",t1="<<ts1<<G4endl;
00886 #endif
00887     }
00888     else if(rand<I12)
00889     {
00890       G4double ran=R2*G4UniformRand();
00891       if(ran>1.) ran=1.;
00892       q2=-std::log(1.-ran)/theB2;
00893       if(q2<0.) q2=0.;
00894       if(a<6.5) q2=std::pow(q2,third);
00895       else      q2=std::pow(q2,fifth);
00896 #ifdef tdebug
00897       G4cout<<"G4QElCS::GetExT: Q2="<<q2<<", r2="<<R2<<", b2="<<theB2<<",t2="<<ts2<<G4endl;
00898 #endif
00899     }
00900     else if(rand<I13)
00901     {
00902       G4double ran=R3*G4UniformRand();
00903       if(ran>1.) ran=1.;
00904       q2=-std::log(1.-ran)/theB3;
00905       if(q2<0.) q2=0.;
00906       if(a>6.5) q2=std::pow(q2,sevth);
00907 #ifdef tdebug
00908       G4cout<<"G4QElCS::GetExT:Q2="<<q2<<", r3="<<R2<<", b3="<<theB2<<",t3="<<ts2<<G4endl;
00909 #endif
00910     }
00911     else
00912     {
00913       G4double ran=R4*G4UniformRand();
00914       if(ran>1.) ran=1.;
00915       q2=-std::log(1.-ran)/theB4;
00916       if(a<6.5) q2=lastTM-q2;                    // u reduced for lightA (starts from 0)
00917 #ifdef tdebug
00918       G4cout<<"G4QElCS::GetExT:Q2="<<q2<<",m="<<lastTM<<",b4="<<theB3<<",t4="<<ts3<<G4endl;
00919 #endif
00920     }
00921   }
00922   if(q2<0.) q2=0.;
00923   if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QElasticCrossSect::GetExchangeT: -t="<<q2<<G4endl;
00924   if(q2>lastTM)
00925   {
00926 #ifdef tdebug
00927     G4cout<<"*Warning*G4QElasticCrossSect::GetExT:-t="<<q2<<">"<<lastTM<<G4endl;
00928 #endif
00929     q2=lastTM;
00930   }
00931   return q2*GeVSQ;
00932 }
00933 
00934 // Returns B in independent units (MeV^-2) (all internal calculations are in GeV) see ExT
00935 G4double G4QProtonElasticCrossSection::GetSlope(G4int tgZ, G4int tgN, G4int PDG)
00936 {
00937   static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
00938 #ifdef tdebug
00939   G4cout<<"G4QElasticCS::GetSlope:"<<onlyCS<<", Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<G4endl;
00940 #endif
00941   if(onlyCS) G4cout<<"*Warning*G4QProtonElasticCrossSection::GetSlope:onlyCS=true"<<G4endl;
00942   if(lastLP<-4.3) return 0.;          // S-wave for p<14 MeV/c (kinE<.1MeV)
00943   if(PDG!=2212)
00944   {
00945     // G4cout<<"*Error*G4QProtonElasticCrossSection::GetSlope: PDG="<<PDG<<", Z="<<tgZ<<", N="
00946     //       <<tgN<<", while it is defined only for PDG=2212"<<G4endl;
00947     // throw G4QException("G4QProtonElasticCrossSection::GetSlope: pA are implemented");
00948     G4ExceptionDescription ed;
00949     ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
00950        << ", while it is defined only for PDG=2212 (p)" << G4endl;
00951     G4Exception("G4QProtonElasticCrossSection::GetSlope()", "HAD_CHPS_0000",
00952                 FatalException, ed);
00953   }
00954   if(theB1<0.) theB1=0.;
00955   if(!(theB1>=-1.||theB1<=1.))G4cout<<"*NAN*G4QElasticCrossSect::Getslope:"<<theB1<<G4endl;
00956   return theB1/GeVSQ;
00957 }
00958 
00959 // Returns half max(Q2=-t) in independent units (MeV^2)
00960 G4double G4QProtonElasticCrossSection::GetHMaxT()
00961 {
00962   static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
00963   return lastTM*HGeVSQ;
00964 }
00965 
00966 // lastLP is used, so calculating tables, one need to remember and then recover lastLP
00967 G4double G4QProtonElasticCrossSection::GetTabValues(G4double lp, G4int PDG, G4int tgZ,
00968                                                     G4int tgN)
00969 {
00970   if(PDG!=2212) G4cout<<"*Warning*G4QProtonElasticCrossSection::GetTabV:PDG="<<PDG<<G4endl;
00971   if(tgZ<0 || tgZ>92)
00972   {
00973     G4cout<<"*Warning*G4QProtonElCS::GetTabValue: (1-92) No isotopes for Z="<<tgZ<<G4endl;
00974     return 0.;
00975   }
00976   G4int iZ=tgZ-1; // Z index
00977   if(iZ<0)
00978   {
00979     iZ=0;         // conversion of the neutron target to the proton target
00980     tgZ=1;
00981     tgN=0;
00982   }
00983   //if(nN[iZ][0] < 0)
00984   //{
00985 #ifdef isodebug
00986   //  G4cout<<"*Warning*G4QElasticCS::GetTabValue: No isotopes for Z="<<tgZ<<G4endl;
00987 #endif
00988   //  return 0.;
00989   //}
00990 #ifdef pdebug
00991   G4cout<<"G4QElasticCS::GetTabVal: lp="<<lp<<",Z="<<tgZ<<",N="<<tgN<<",PDG="<<PDG<<G4endl;
00992 #endif
00993   G4double p=std::exp(lp);              // momentum
00994   G4double sp=std::sqrt(p);             // sqrt(p)
00995   G4double p2=p*p;            
00996   G4double p3=p2*p;
00997   G4double p4=p3*p;
00998   if ( tgZ == 1 && tgN == 0 ) // pp/nn
00999   {
01000     G4double p2s=p2*sp;
01001     G4double dl2=lp-lastPAR[8];
01002     theSS=lastPAR[31];
01003     theS1=(lastPAR[9]+lastPAR[10]*dl2*dl2)/(1.+lastPAR[11]/p4/p)+
01004           (lastPAR[12]/p2+lastPAR[13]*p)/(p4+lastPAR[14]*sp);
01005     theB1=lastPAR[15]*std::pow(p,lastPAR[16])/(1.+lastPAR[17]/p3);
01006     theS2=lastPAR[18]+lastPAR[19]/(p4+lastPAR[20]*p);
01007     theB2=lastPAR[21]+lastPAR[22]/(p4+lastPAR[23]/sp); 
01008     theS3=lastPAR[24]+lastPAR[25]/(p4*p4+lastPAR[26]*p2+lastPAR[27]);
01009     theB3=lastPAR[28]+lastPAR[29]/(p4+lastPAR[30]); 
01010     theS4=0.;
01011     theB4=0.; 
01012 #ifdef tdebug
01013     G4cout<<"G4QElasticCS::GetTableValues:(pp) TM="<<lastTM<<",S1="<<theS1<<",B1="<<theB1
01014           <<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS1<<",B3="<<theB1<<G4endl;
01015 #endif
01016     // Returns the total elastic pp cross-section (to avoid spoiling lastSIG)
01017     G4double dl1=lp-lastPAR[3];
01018     return lastPAR[0]/p2s/(1.+lastPAR[7]/p2s)+(lastPAR[1]+lastPAR[2]*dl1*dl1+lastPAR[4]/p)
01019                                                    /(1.+lastPAR[5]*lp)/(1.+lastPAR[6]/p4);
01020   }
01021   else
01022   {
01023     G4double p5=p4*p;
01024     G4double p6=p5*p;
01025     G4double p8=p6*p2;
01026     G4double p10=p8*p2;
01027     G4double p12=p10*p2;
01028     G4double p16=p8*p8;
01029     //G4double p24=p16*p8;
01030     G4double dl=lp-5.;
01031     G4double a=tgZ+tgN;
01032     G4double pah=std::pow(p,a/2);
01033     G4double pa=pah*pah;
01034     G4double pa2=pa*pa;
01035     if(a<6.5)
01036     {
01037       theS1=lastPAR[9]/(1.+lastPAR[10]*p4*pa)+lastPAR[11]/(p4+lastPAR[12]*p4/pa2)+
01038             (lastPAR[13]*dl*dl+lastPAR[14])/(1.+lastPAR[15]/p2);
01039       theB1=(lastPAR[16]+lastPAR[17]*p2)/(p4+lastPAR[18]/pah)+lastPAR[19];
01040       theSS=lastPAR[20]/(1.+lastPAR[21]/p2)+lastPAR[22]/(p6/pa+lastPAR[23]/p16);
01041       theS2=lastPAR[24]/(pa/p2+lastPAR[25]/p4)+lastPAR[26];
01042       theB2=lastPAR[27]*std::pow(p,lastPAR[28])+lastPAR[29]/(p8+lastPAR[30]/p16);
01043       theS3=lastPAR[31]/(pa*p+lastPAR[32]/pa)+lastPAR[33];
01044       theB3=lastPAR[34]/(p3+lastPAR[35]/p6)+lastPAR[36]/(1.+lastPAR[37]/p2);
01045       theS4=p2*(pah*lastPAR[38]*std::exp(-pah*lastPAR[39])+
01046                 lastPAR[40]/(1.+lastPAR[41]*std::pow(p,lastPAR[42])));
01047       theB4=lastPAR[43]*pa/p2/(1.+pa*lastPAR[44]);
01048 #ifdef tdebug
01049       G4cout<<"G4QElCS::GetTabV: lA, p="<<p<<",S1="<<theS1<<",B1="<<theB1<<",SS="<<theSS
01050             <<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",S4="<<theS4
01051             <<",B4="<<theB4<<G4endl;
01052 #endif
01053     }
01054     else
01055     {
01056       theS1=lastPAR[9]/(1.+lastPAR[10]/p4)+lastPAR[11]/(p4+lastPAR[12]/p2)+
01057             lastPAR[13]/(p5+lastPAR[14]/p16);
01058       theB1=(lastPAR[15]/p8+lastPAR[19])/(p+lastPAR[16]/std::pow(p,lastPAR[20]))+
01059             lastPAR[17]/(1.+lastPAR[18]/p4);
01060       theSS=lastPAR[21]/(p4/std::pow(p,lastPAR[23])+lastPAR[22]/p4);
01061       theS2=lastPAR[24]/p4/(std::pow(p,lastPAR[25])+lastPAR[26]/p12)+lastPAR[27];
01062       theB2=lastPAR[28]/std::pow(p,lastPAR[29])+lastPAR[30]/std::pow(p,lastPAR[31]);
01063       theS3=lastPAR[32]/std::pow(p,lastPAR[35])/(1.+lastPAR[36]/p12)+
01064             lastPAR[33]/(1.+lastPAR[34]/p6);
01065       theB3=lastPAR[37]/p8+lastPAR[38]/p2+lastPAR[39]/(1.+lastPAR[40]/p8);
01066       theS4=(lastPAR[41]/p4+lastPAR[46]/p)/(1.+lastPAR[42]/p10)+
01067             (lastPAR[43]+lastPAR[44]*dl*dl)/(1.+lastPAR[45]/p12);
01068       theB4=lastPAR[47]/(1.+lastPAR[48]/p)+lastPAR[49]*p4/(1.+lastPAR[50]*p5);
01069 #ifdef tdebug
01070       G4cout<<"G4QElCS::GetTabV: hA, p="<<p<<",S1="<<theS1<<",B1="<<theB1<<",SS="<<theSS
01071             <<",S2="<<theS2<<",B2="<<theB2<<",S3="<<theS3<<",B3="<<theB3<<",S4="<<theS4
01072             <<",B4="<<theB4<<G4endl;
01073 #endif
01074     }
01075     // Returns the total elastic (n/p)A cross-section (to avoid spoiling lastSIG)
01076 #ifdef tdebug
01077     G4cout<<"G4QElCS::GetTabV: PDG="<<PDG<<",P="<<p<<",N="<<tgN<<",Z="<<tgZ<<G4endl;
01078 #endif
01079     //         p1               p2              p3            p6
01080     return (lastPAR[0]*dl*dl+lastPAR[1])/(1.+lastPAR[2]/p+lastPAR[5]/p6)+
01081      lastPAR[3]/(p3+lastPAR[4]/p3)+lastPAR[7]/(p4+std::pow((lastPAR[8]/p),lastPAR[6]));
01082     //   p4            p5               p8                 p9             p7
01083   }
01084   return 0.;
01085 } // End of GetTableValues
01086 
01087 // Returns max -t=Q2 (GeV^2) for the momentum pP(GeV) and the target nucleus (tgN,tgZ)
01088 G4double G4QProtonElasticCrossSection::GetQ2max(G4int PDG, G4int tgZ, G4int tgN,
01089                                                 G4double pP)
01090 {
01091   //static const G4double mNeut= G4QPDGCode(2112).GetMass()*.001; // MeV to GeV
01092   static const G4double mProt= G4QPDGCode(2212).GetMass()*.001; // MeV to GeV
01093   //static const G4double mLamb= G4QPDGCode(3122).GetMass()*.001; // MeV to GeV
01094   //static const G4double mHe3 = G4QPDGCode(2112).GetNuclMass(2,1,0)*.001; // MeV to GeV
01095   //static const G4double mAlph = G4QPDGCode(2112).GetNuclMass(2,2,0)*.001; // MeV to GeV
01096   //static const G4double mDeut = G4QPDGCode(2112).GetNuclMass(1,1,0)*.001; // MeV to GeV
01097   static const G4double mProt2= mProt*mProt;
01098   //static const G4double mNeut2= mNeut*mNeut;
01099   //static const G4double mDeut2= mDeut*mDeut;
01100   G4double pP2=pP*pP;                                 // squared momentum of the projectile
01101   if(tgZ==1 && tgN==0)
01102   {
01103     G4double tMid=std::sqrt(pP2+mProt2)*mProt-mProt2; // CMS 90deg value of -t=Q2 (GeV^2)
01104     return tMid+tMid;
01105   }
01106   else if(tgZ || tgN)                                 // ---> pA
01107   {
01108     G4double mt=G4QPDGCode(90000000+tgZ*1000+tgN).GetMass()*.001; // Target mass in GeV
01109     G4double dmt=mt+mt;
01110     G4double s_value=dmt*std::sqrt(pP2+mProt2)+mProt2+mt*mt;// Mondelstam s
01111     return dmt*dmt*pP2/s_value;
01112   }
01113   else
01114   {
01115     // G4cout<<"*Error*G4QProtonElasticCrossSection::GetQ2max: PDG="<<PDG<<", Z="<<tgZ<<", N="
01116     //       <<tgN<<", while it is defined only for p projectiles & Z_target>0"<<G4endl;
01117     // throw G4QException("G4QProtonElasticCrossSection::GetQ2max: only pA are implemented");
01118     G4ExceptionDescription ed;
01119     ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
01120        << ", while it is defined only for p projectiles & Z_target>0" << G4endl;
01121     G4Exception("G4QProtonElasticCrossSection::GetQ2max()", "HAD_CHPS_0000",
01122                 FatalException, ed);
01123     return 0;
01124   }
01125 }
01126 

Generated on Mon May 27 17:49:39 2013 for Geant4 by  doxygen 1.4.7