G4ChipsProtonElasticXS.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: G4ChipsProtonElasticXS 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 elastic process. 
00036 // Class extracted from CHIPS and integrated in Geant4 by W.Pokorski
00037 // -------------------------------------------------------------------------------
00038 
00039 
00040 #include "G4ChipsProtonElasticXS.hh"
00041 #include "G4SystemOfUnits.hh"
00042 #include "G4DynamicParticle.hh"
00043 #include "G4ParticleDefinition.hh"
00044 #include "G4Proton.hh"
00045 #include "G4Nucleus.hh"
00046 #include "G4ParticleTable.hh"
00047 #include "G4NucleiProperties.hh"
00048 
00049 // factory
00050 #include "G4CrossSectionFactory.hh"
00051 //
00052 G4_DECLARE_XS_FACTORY(G4ChipsProtonElasticXS);
00053 
00054 
00055 G4ChipsProtonElasticXS::G4ChipsProtonElasticXS():G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
00056 {
00057   // Initialization of the parameters
00058   lPMin=-8.;  // Min tabulated logarithmicMomentum(D)
00059   lPMax= 8.;  // Max tabulated logarithmicMomentum(D)
00060   dlnP=(lPMax-lPMin)/nLast;// LogStep in the table(D)
00061   onlyCS=false;// Flag toCalculateOnlyCS(not Si/Bi)(L)
00062   lastSIG=0.; // Last calculated cross section    (L)
00063   lastLP=-10.;// Last log(mom_ofTheIncidentHadron)(L)
00064   lastTM=0.;  // Last t_maximum                   (L)
00065   theSS=0.;   // The Last sq.slope of 1st difr.Max(L)
00066   theS1=0.;   // The Last mantissa of 1st difr.Max(L)
00067   theB1=0.;   // The Last slope of 1st difruct.Max(L)
00068   theS2=0.;   // The Last mantissa of 2nd difr.Max(L)
00069   theB2=0.;   // The Last slope of 2nd difruct.Max(L)
00070   theS3=0.;   // The Last mantissa of 3d difr. Max(L)
00071   theB3=0.;   // The Last slope of 3d difruct. Max(L)
00072   theS4=0.;   // The Last mantissa of 4th difr.Max(L)
00073   theB4=0.;   // The Last slope of 4th difruct.Max(L)
00074   lastTZ=0;   // Last atomic number of the target
00075   lastTN=0;   // Last # of neutrons in the target
00076   lastPIN=0.; // Last initialized max momentum
00077   lastCST=0;  // Elastic cross-section table
00078   lastPAR=0;  // Parameters for FunctionalCalculation
00079   lastSST=0;  // E-dep of sq.slope of the 1st dif.Max
00080   lastS1T=0;  // E-dep of mantissa of the 1st dif.Max
00081   lastB1T=0;  // E-dep of the slope of the 1st difMax
00082   lastS2T=0;  // E-dep of mantissa of the 2nd difrMax
00083   lastB2T=0;  // E-dep of the slope of the 2nd difMax
00084   lastS3T=0;  // E-dep of mantissa of the 3d difr.Max
00085   lastB3T=0;  // E-dep of the slope of the 3d difrMax
00086   lastS4T=0;  // E-dep of mantissa of the 4th difrMax
00087   lastB4T=0;  // E-dep of the slope of the 4th difMax
00088   lastN=0;    // The last N of calculated nucleus
00089   lastZ=0;    // The last Z of calculated nucleus
00090   lastP=0.;   // Last used in cross section Momentum
00091   lastTH=0.;  // Last threshold momentum
00092   lastCS=0.;  // Last value of the Cross Section
00093   lastI=0;    // The last position in the DAMDB
00094 }
00095 
00096 
00097 G4ChipsProtonElasticXS::~G4ChipsProtonElasticXS()
00098 {
00099   std::vector<G4double*>::iterator pos;
00100   for (pos=CST.begin(); pos<CST.end(); pos++)
00101   { delete [] *pos; }
00102   CST.clear();
00103   for (pos=PAR.begin(); pos<PAR.end(); pos++)
00104   { delete [] *pos; }
00105   PAR.clear();
00106   for (pos=SST.begin(); pos<SST.end(); pos++)
00107   { delete [] *pos; }
00108   SST.clear();
00109   for (pos=S1T.begin(); pos<S1T.end(); pos++)
00110   { delete [] *pos; }
00111   S1T.clear();
00112   for (pos=B1T.begin(); pos<B1T.end(); pos++)
00113   { delete [] *pos; }
00114   B1T.clear();
00115   for (pos=S2T.begin(); pos<S2T.end(); pos++)
00116   { delete [] *pos; }
00117   S2T.clear();
00118   for (pos=B2T.begin(); pos<B2T.end(); pos++)
00119   { delete [] *pos; }
00120   B2T.clear();
00121   for (pos=S3T.begin(); pos<S3T.end(); pos++)
00122   { delete [] *pos; }
00123   S3T.clear();
00124   for (pos=B3T.begin(); pos<B3T.end(); pos++)
00125   { delete [] *pos; }
00126   B3T.clear();
00127   for (pos=S4T.begin(); pos<S4T.end(); pos++)
00128   { delete [] *pos; }
00129   S4T.clear();
00130   for (pos=B4T.begin(); pos<B4T.end(); pos++)
00131   { delete [] *pos; }
00132   B4T.clear();
00133  
00134 }
00135 
00136 G4bool G4ChipsProtonElasticXS::IsIsoApplicable(const G4DynamicParticle* Pt, G4int, G4int,    
00137                                  const G4Element*,
00138                                  const G4Material*)
00139 {
00140   G4ParticleDefinition* particle = Pt->GetDefinition();
00141   if (particle == G4Proton::Proton()      ) return true;
00142   return false;
00143 }
00144 
00145 
00146 G4double G4ChipsProtonElasticXS::GetIsoCrossSection(const G4DynamicParticle* Pt, G4int tgZ, G4int A,  
00147                                                         const G4Isotope*,
00148                                                         const G4Element*,
00149                                                         const G4Material*)
00150 {
00151   G4double pMom=Pt->GetTotalMomentum();
00152   G4int tgN = A - tgZ;
00153 
00154   return GetChipsCrossSection(pMom, tgZ, tgN, 2212);
00155 }
00156 
00157 
00158 // The main member function giving the collision cross section (P is in IU, CS is in mb)
00159 // Make pMom in independent units ! (Now it is MeV)
00160 G4double G4ChipsProtonElasticXS::GetChipsCrossSection(G4double pMom, G4int tgZ, G4int tgN, G4int)
00161 {
00162   static std::vector <G4int>    colN;  // Vector of N for calculated nuclei (isotops)
00163   static std::vector <G4int>    colZ;  // Vector of Z for calculated nuclei (isotops)
00164   static std::vector <G4double> colP;  // Vector of last momenta for the reaction
00165   static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
00166   static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
00167   // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
00168 
00169   G4double pEn=pMom;
00170   onlyCS=false;
00171 
00172   G4bool in=false;                   // By default the isotope must be found in the AMDB
00173   lastP   = 0.;                      // New momentum history (nothing to compare with)
00174   lastN   = tgN;                     // The last N of the calculated nucleus
00175   lastZ   = tgZ;                     // The last Z of the calculated nucleus
00176   lastI   = colN.size();             // Size of the Associative Memory DB in the heap
00177   if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
00178   {                                  // The nucleus with projPDG is found in AMDB
00179     if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
00180     {
00181       lastI=i;
00182       lastTH =colTH[i];              // Last THreshold (A-dependent)
00183       if(pEn<=lastTH)
00184       {
00185         return 0.;                   // Energy is below the Threshold value
00186       }
00187       lastP  =colP [i];              // Last Momentum  (A-dependent)
00188       lastCS =colCS[i];              // Last CrossSect (A-dependent)
00189       if(lastP == pMom)              // Do not recalculate
00190       {
00191         CalculateCrossSection(onlyCS,-1,i,2212,lastZ,lastN,pMom); // Update param's only
00192         return lastCS*millibarn;     // Use theLastCS
00193       }
00194       in = true;                       // This is the case when the isotop is found in DB
00195       // Momentum pMom is in IU ! @@ Units
00196       lastCS=CalculateCrossSection(onlyCS,-1,i,2212,lastZ,lastN,pMom); // read & update
00197       if(lastCS<=0. && pEn>lastTH)    // Correct the threshold
00198       {
00199         lastTH=pEn;
00200       }
00201       break;                           // Go out of the LOOP with found lastI
00202     }
00203   } // End of attampt to find the nucleus in DB
00204   if(!in)                            // This nucleus has not been calculated previously
00205   {
00207     lastCS=CalculateCrossSection(onlyCS,0,lastI,2212,lastZ,lastN,pMom);//calculate&create
00208     if(lastCS<=0.)
00209     {
00210       lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
00211       if(pEn>lastTH)
00212       {
00213         lastTH=pEn;
00214       }
00215     }
00216     colN.push_back(tgN);
00217     colZ.push_back(tgZ);
00218     colP.push_back(pMom);
00219     colTH.push_back(lastTH);
00220     colCS.push_back(lastCS);
00221     return lastCS*millibarn;
00222   } // End of creation of the new set of parameters
00223   else
00224   {
00225     colP[lastI]=pMom;
00226     colCS[lastI]=lastCS;
00227   }
00228   return lastCS*millibarn;
00229 }
00230 
00231 // Calculation of total elastic cross section (p in IU, CS in mb) @@ Units (?)
00232 // F=0 - create AMDB, F=-1 - read&update AMDB, F=1 - update AMDB (sinchro with higher AMDB)
00233 G4double G4ChipsProtonElasticXS::CalculateCrossSection(G4bool CS, G4int F, G4int I,
00234                                              G4int PDG, G4int tgZ, G4int tgN, G4double pIU)
00235 {
00236   // *** Begin of Associative Memory DB for acceleration of the cross section calculations
00237   static std::vector <G4double>  PIN;   // Vector of max initialized log(P) in the table
00238   // *** End of Static Definitions (Associative Memory Data Base) ***
00239   G4double pMom=pIU/GeV;                // All calculations are in GeV
00240   onlyCS=CS;                            // Flag to calculate only CS (not Si/Bi)
00241   lastLP=std::log(pMom);                // Make a logarithm of the momentum for calculation
00242   if(F)                                 // This isotope was found in AMDB =>RETRIEVE/UPDATE
00243   {
00244     if(F<0)                             // the AMDB must be loded
00245     {
00246       lastPIN = PIN[I];                 // Max log(P) initialised for this table set
00247       lastPAR = PAR[I];                 // Pointer to the parameter set
00248       lastCST = CST[I];                 // Pointer to the total sross-section table
00249       lastSST = SST[I];                 // Pointer to the first squared slope
00250       lastS1T = S1T[I];                 // Pointer to the first mantissa
00251       lastB1T = B1T[I];                 // Pointer to the first slope
00252       lastS2T = S2T[I];                 // Pointer to the second mantissa
00253       lastB2T = B2T[I];                 // Pointer to the second slope
00254       lastS3T = S3T[I];                 // Pointer to the third mantissa
00255       lastB3T = B3T[I];                 // Pointer to the rhird slope
00256       lastS4T = S4T[I];                 // Pointer to the 4-th mantissa
00257       lastB4T = B4T[I];                 // Pointer to the 4-th slope
00258     }
00259     if(lastLP>lastPIN && lastLP<lPMax)
00260     {
00261       lastPIN=GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);// Can update upper logP-Limit in tabs
00262       PIN[I]=lastPIN;                   // Remember the new P-Limit of the tables
00263     }
00264   }
00265   else                                  // This isotope wasn't initialized => CREATE
00266   {
00267     lastPAR = new G4double[nPoints];    // Allocate memory for parameters of CS function
00268     lastPAR[nLast]=0;                   // Initialization for VALGRIND
00269     lastCST = new G4double[nPoints];    // Allocate memory for Tabulated CS function    
00270     lastSST = new G4double[nPoints];    // Allocate memory for Tabulated first sqaredSlope 
00271     lastS1T = new G4double[nPoints];    // Allocate memory for Tabulated first mantissa 
00272     lastB1T = new G4double[nPoints];    // Allocate memory for Tabulated first slope    
00273     lastS2T = new G4double[nPoints];    // Allocate memory for Tabulated second mantissa
00274     lastB2T = new G4double[nPoints];    // Allocate memory for Tabulated second slope   
00275     lastS3T = new G4double[nPoints];    // Allocate memory for Tabulated third mantissa 
00276     lastB3T = new G4double[nPoints];    // Allocate memory for Tabulated third slope    
00277     lastS4T = new G4double[nPoints];    // Allocate memory for Tabulated 4-th mantissa 
00278     lastB4T = new G4double[nPoints];    // Allocate memory for Tabulated 4-th slope    
00279     lastPIN = GetPTables(lastLP,lPMin,PDG,tgZ,tgN); // Returns the new P-limit for tables
00280     PIN.push_back(lastPIN);             // Fill parameters of CS function to AMDB
00281     PAR.push_back(lastPAR);             // Fill parameters of CS function to AMDB
00282     CST.push_back(lastCST);             // Fill Tabulated CS function to AMDB    
00283     SST.push_back(lastSST);             // Fill Tabulated first sq.slope to AMDB 
00284     S1T.push_back(lastS1T);             // Fill Tabulated first mantissa to AMDB 
00285     B1T.push_back(lastB1T);             // Fill Tabulated first slope to AMDB    
00286     S2T.push_back(lastS2T);             // Fill Tabulated second mantissa to AMDB 
00287     B2T.push_back(lastB2T);             // Fill Tabulated second slope to AMDB    
00288     S3T.push_back(lastS3T);             // Fill Tabulated third mantissa to AMDB 
00289     B3T.push_back(lastB3T);             // Fill Tabulated third slope to AMDB    
00290     S4T.push_back(lastS4T);             // Fill Tabulated 4-th mantissa to AMDB 
00291     B4T.push_back(lastB4T);             // Fill Tabulated 4-th slope to AMDB    
00292   } // End of creation/update of the new set of parameters and tables
00293   // =--------= NOW Update (if necessary) and Calculate the Cross Section =------------=
00294   if(lastLP>lastPIN && lastLP<lPMax)
00295   {
00296     lastPIN = GetPTables(lastLP,lastPIN,PDG,tgZ,tgN);
00297   }
00298   if(!onlyCS) lastTM=GetQ2max(PDG, tgZ, tgN, pMom); // Calculate (-t)_max=Q2_max (GeV2)
00299   if(lastLP>lPMin && lastLP<=lastPIN)   // Linear fit is made using precalculated tables
00300   {
00301     if(lastLP==lastPIN)
00302     {
00303       G4double shift=(lastLP-lPMin)/dlnP+.000001; // Log distance from lPMin
00304       G4int    blast=static_cast<int>(shift); // this is a bin number of the lower edge (0)
00305       if(blast<0 || blast>=nLast) G4cout<<"G4QEleastCS::CCS:b="<<blast<<","<<nLast<<G4endl;
00306       lastSIG = lastCST[blast];
00307       if(!onlyCS)                       // Skip the differential cross-section parameters
00308       {
00309         theSS  = lastSST[blast];
00310         theS1  = lastS1T[blast];
00311         theB1  = lastB1T[blast];
00312         theS2  = lastS2T[blast];
00313         theB2  = lastB2T[blast];
00314         theS3  = lastS3T[blast];
00315         theB3  = lastB3T[blast];
00316         theS4  = lastS4T[blast];
00317         theB4  = lastB4T[blast];
00318       }
00319     }
00320     else
00321     {
00322       G4double shift=(lastLP-lPMin)/dlnP;        // a shift from the beginning of the table
00323       G4int    blast=static_cast<int>(shift);    // the lower bin number
00324       if(blast<0)   blast=0;
00325       if(blast>=nLast) blast=nLast-1;            // low edge of the last bin
00326       shift-=blast;                              // step inside the unit bin
00327       G4int lastL=blast+1;                       // the upper bin number
00328       G4double SIGL=lastCST[blast];              // the basic value of the cross-section
00329       lastSIG= SIGL+shift*(lastCST[lastL]-SIGL); // calculated total elastic cross-section
00330       if(!onlyCS)                       // Skip the differential cross-section parameters
00331       {
00332         G4double SSTL=lastSST[blast];           // the low bin of the first squared slope
00333         theSS=SSTL+shift*(lastSST[lastL]-SSTL); // the basic value of the first sq.slope
00334         G4double S1TL=lastS1T[blast];           // the low bin of the first mantissa
00335         theS1=S1TL+shift*(lastS1T[lastL]-S1TL); // the basic value of the first mantissa
00336         G4double B1TL=lastB1T[blast];           // the low bin of the first slope
00337         theB1=B1TL+shift*(lastB1T[lastL]-B1TL); // the basic value of the first slope
00338         G4double S2TL=lastS2T[blast];           // the low bin of the second mantissa
00339         theS2=S2TL+shift*(lastS2T[lastL]-S2TL); // the basic value of the second mantissa
00340         G4double B2TL=lastB2T[blast];           // the low bin of the second slope
00341         theB2=B2TL+shift*(lastB2T[lastL]-B2TL); // the basic value of the second slope
00342         G4double S3TL=lastS3T[blast];           // the low bin of the third mantissa
00343         theS3=S3TL+shift*(lastS3T[lastL]-S3TL); // the basic value of the third mantissa
00344         G4double B3TL=lastB3T[blast];           // the low bin of the third slope
00345         theB3=B3TL+shift*(lastB3T[lastL]-B3TL); // the basic value of the third slope
00346         G4double S4TL=lastS4T[blast];           // the low bin of the 4-th mantissa
00347         theS4=S4TL+shift*(lastS4T[lastL]-S4TL); // the basic value of the 4-th mantissa
00348         G4double B4TL=lastB4T[blast];           // the low bin of the 4-th slope
00349         theB4=B4TL+shift*(lastB4T[lastL]-B4TL); // the basic value of the 4-th slope
00350       }
00351     }
00352   }
00353   else lastSIG=GetTabValues(lastLP, PDG, tgZ, tgN); // Direct calculation beyond the table
00354   if(lastSIG<0.) lastSIG = 0.;                   // @@ a Warning print can be added
00355   return lastSIG;
00356 }
00357 
00358 // It has parameter sets for all tZ/tN/PDG, using them the tables can be created/updated
00359 G4double G4ChipsProtonElasticXS::GetPTables(G4double LP, G4double ILP, G4int PDG,
00360                                                   G4int tgZ, G4int tgN)
00361 {
00362   // @@ At present all nA==pA ---------> Each neucleus can have not more than 51 parameters
00363   static const G4double pwd=2727;
00364   const G4int n_npel=24;                // #of parameters for np-elastic (<nPoints=128)
00365   const G4int n_ppel=32;                // #of parameters for pp-elastic (<nPoints=128)
00366   //                      -0- -1-  -2- -3- -4-  -5- -6- -7- -8- -9--10--11--12--13- -14-
00367   G4double np_el[n_npel]={12.,.05,.0001,5.,.35,6.75,.14,19.,.6,6.75,.14,13.,.14,.6,.00013,
00368                           75.,.001,7.2,4.32,.012,2.5,0.0,12.,.34};
00369   //                      -15--16--17- -18- -19--20--21--22--23-
00370   //                       -0-   -1-  -2- -3- -4- -5-  -6-  -7-  -8--9--10--11--12--13-
00371   G4double pp_el[n_ppel]={2.865,18.9,.6461,3.,9.,.425,.4276,.0022,5.,74.,3.,3.4,.2,.17,
00372                           .001,8.,.055,3.64,5.e-5,4000.,1500.,.46,1.2e6,3.5e6,5.e-5,1.e10,
00373                           8.5e8,1.e10,1.1,3.4e6,6.8e6,0.};
00374   //                      -14--15- -16- -17- -18-  -19- -20- -21- -22-  -23-   -24-  -25-
00375   //                       -26- -27-  -28- -29- -30- -31-
00376   if(PDG==2212)
00377   {
00378     // -- Total pp elastic cross section cs & s1/b1 (main), s2/b2 (tail1), s3/b3 (tail2) --
00379     //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|
00380     //CS=2.865/p2s/(1+.0022/p2s)+(18.9+.6461*dl1*dl1+9./p)/(1.+.425*lp)/(1.+.4276/p4);
00381     //   par(0)       par(7)     par(1) par(2)      par(4)      par(5)         par(6)
00382     //dl2=lp-5., s1=(74.+3.*dl2*dl2)/(1+3.4/p4/p)+(.2/p2+17.*p)/(p4+.001*sp),
00383     //     par(8) par(9) par(10)        par(11)   par(12)par(13)    par(14)
00384     // b1=8.*p**.055/(1.+3.64/p3); s2=5.e-5+4000./(p4+1500.*p); b2=.46+1.2e6/(p4+3.5e6/sp);
00385     // par(15) par(16)  par(17)     par(18) par(19)  par(20)   par(21) par(22)  par(23)
00386     // s3=5.e-5+1.e10/(p4*p4+8.5e8*p2+1.e10); b3=1.1+3.4e6/(p4+6.8e6); ss=0.
00387     //  par(24) par(25)     par(26)  par(27) par(28) par(29)  par(30)   par(31)
00388     //
00389     if(lastPAR[nLast]!=pwd) // A unique flag to avoid the repeatable definition
00390     {
00391       if ( tgZ == 0 && tgN == 1 )
00392       {
00393         for (G4int ip=0; ip<n_npel; ip++) lastPAR[ip]=np_el[ip]; // pn
00394 
00395       }
00396       else if ( tgZ == 1 && tgN == 0 )
00397       {
00398         for (G4int ip=0; ip<n_ppel; ip++) lastPAR[ip]=pp_el[ip]; // pp
00399       }
00400       else
00401       {
00402         G4double a=tgZ+tgN;
00403         G4double sa=std::sqrt(a);
00404         G4double ssa=std::sqrt(sa);
00405         G4double asa=a*sa;
00406         G4double a2=a*a;
00407         G4double a3=a2*a;
00408         G4double a4=a3*a;
00409         G4double a5=a4*a;
00410         G4double a6=a4*a2;
00411         G4double a7=a6*a;
00412         G4double a8=a7*a;
00413         G4double a9=a8*a;
00414         G4double a10=a5*a5;
00415         G4double a12=a6*a6;
00416         G4double a14=a7*a7;
00417         G4double a16=a8*a8;
00418         G4double a17=a16*a;
00419         G4double a20=a16*a4;
00420         G4double a32=a16*a16;
00421         // Reaction cross-section parameters (pel=peh_fit.f)
00422         lastPAR[0]=5./(1.+22./asa);                                          // p1
00423         lastPAR[1]=4.8*std::pow(a,1.14)/(1.+3.6/a3);                         // p2
00424         lastPAR[2]=1./(1.+4.E-3*a4)+2.E-6*a3/(1.+1.3E-6*a3);                 // p3
00425         lastPAR[3]=1.3*a;                                                    // p4
00426         lastPAR[4]=3.E-8*a3/(1.+4.E-7*a4);                                   // p5
00427         lastPAR[5]=.07*asa/(1.+.009*a2);                                     // p6
00428         lastPAR[6]=(3.+3.E-16*a20)/(1.+a20*(2.E-16/a+3.E-19*a));             // p7 (11)
00429         lastPAR[7]=(5.E-9*a4*sa+.27/a)/(1.+5.E16/a20)/(1.+6.E-9*a4)+.015/a2; // p8
00430         lastPAR[8]=(.001*a+.07/a)/(1.+5.E13/a16+5.E-7*a3)+.0003/sa;          // p9 (10)
00431         // @@ the differential cross-section is parameterized separately for A>6 & A<7
00432         if(a<6.5)
00433         {
00434           G4double a28=a16*a12;
00435           // The main pre-exponent      (pel_sg)
00436           lastPAR[ 9]=4000*a;                                // p1
00437           lastPAR[10]=1.2e7*a8+380*a17;                      // p2
00438           lastPAR[11]=.7/(1.+4.e-12*a16);                    // p3
00439           lastPAR[12]=2.5/a8/(a4+1.e-16*a32);                // p4
00440           lastPAR[13]=.28*a;                                 // p5
00441           lastPAR[14]=1.2*a2+2.3;                            // p6
00442           lastPAR[15]=3.8/a;                                 // p7
00443           // The main slope             (pel_sl)
00444           lastPAR[16]=.01/(1.+.0024*a5);                     // p1
00445           lastPAR[17]=.2*a;                                  // p2
00446           lastPAR[18]=9.e-7/(1.+.035*a5);                    // p3
00447           lastPAR[19]=(42.+2.7e-11*a16)/(1.+.14*a);          // p4
00448           // The main quadratic         (pel_sh)
00449           lastPAR[20]=2.25*a3;                               // p1
00450           lastPAR[21]=18.;                                   // p2
00451           lastPAR[22]=2.4e-3*a8/(1.+2.6e-4*a7);              // p3
00452           lastPAR[23]=3.5e-36*a32*a8/(1.+5.e-15*a32/a);      // p4
00453           // The 1st max pre-exponent   (pel_qq)
00454           lastPAR[24]=1.e5/(a8+2.5e12/a16);                  // p1
00455           lastPAR[25]=8.e7/(a12+1.e-27*a28*a28);             // p2 
00456           lastPAR[26]=.0006*a3;                              // p3
00457           // The 1st max slope          (pel_qs)
00458           lastPAR[27]=10.+4.e-8*a12*a;                       // p1
00459           lastPAR[28]=.114;                                  // p2
00460           lastPAR[29]=.003;                                  // p3
00461           lastPAR[30]=2.e-23;                                // p4
00462           // The effective pre-exponent (pel_ss)
00463           lastPAR[31]=1./(1.+.0001*a8);                      // p1
00464           lastPAR[32]=1.5e-4/(1.+5.e-6*a12);                 // p2
00465           lastPAR[33]=.03;                                   // p3
00466           // The effective slope        (pel_sb)
00467           lastPAR[34]=a/2;                                   // p1
00468           lastPAR[35]=2.e-7*a4;                              // p2
00469           lastPAR[36]=4.;                                    // p3
00470           lastPAR[37]=64./a3;                                // p4
00471           // The gloria pre-exponent    (pel_us)
00472           lastPAR[38]=1.e8*std::exp(.32*asa);                // p1
00473           lastPAR[39]=20.*std::exp(.45*asa);                 // p2
00474           lastPAR[40]=7.e3+2.4e6/a5;                         // p3
00475           lastPAR[41]=2.5e5*std::exp(.085*a3);               // p4
00476           lastPAR[42]=2.5*a;                                 // p5
00477           // The gloria slope           (pel_ub)
00478           lastPAR[43]=920.+.03*a8*a3;                        // p1
00479           lastPAR[44]=93.+.0023*a12;                         // p2
00480         }
00481         else
00482         {
00483           G4double p1a10=2.2e-28*a10;
00484           G4double r4a16=6.e14/a16;
00485           G4double s4a16=r4a16*r4a16;
00486           // a24
00487           // a36
00488           // The main pre-exponent      (peh_sg)
00489           lastPAR[ 9]=4.5*std::pow(a,1.15);                  // p1
00490           lastPAR[10]=.06*std::pow(a,.6);                    // p2
00491           lastPAR[11]=.6*a/(1.+2.e15/a16);                   // p3
00492           lastPAR[12]=.17/(a+9.e5/a3+1.5e33/a32);            // p4
00493           lastPAR[13]=(.001+7.e-11*a5)/(1.+4.4e-11*a5);      // p5
00494           lastPAR[14]=(p1a10*p1a10+2.e-29)/(1.+2.e-22*a12);  // p6
00495           // The main slope             (peh_sl)
00496           lastPAR[15]=400./a12+2.e-22*a9;                    // p1
00497           lastPAR[16]=1.e-32*a12/(1.+5.e22/a14);             // p2
00498           lastPAR[17]=1000./a2+9.5*sa*ssa;                   // p3
00499           lastPAR[18]=4.e-6*a*asa+1.e11/a16;                 // p4
00500           lastPAR[19]=(120./a+.002*a2)/(1.+2.e14/a16);       // p5
00501           lastPAR[20]=9.+100./a;                             // p6
00502           // The main quadratic         (peh_sh)
00503           lastPAR[21]=.002*a3+3.e7/a6;                       // p1
00504           lastPAR[22]=7.e-15*a4*asa;                         // p2
00505           lastPAR[23]=9000./a4;                              // p3
00506           // The 1st max pre-exponent   (peh_qq)
00507           lastPAR[24]=.0011*asa/(1.+3.e34/a32/a4);           // p1
00508           lastPAR[25]=1.e-5*a2+2.e14/a16;                    // p2
00509           lastPAR[26]=1.2e-11*a2/(1.+1.5e19/a12);            // p3
00510           lastPAR[27]=.016*asa/(1.+5.e16/a16);               // p4
00511           // The 1st max slope          (peh_qs)
00512           lastPAR[28]=.002*a4/(1.+7.e7/std::pow(a-6.83,14)); // p1
00513           lastPAR[29]=2.e6/a6+7.2/std::pow(a,.11);           // p2
00514           lastPAR[30]=11.*a3/(1.+7.e23/a16/a8);              // p3
00515           lastPAR[31]=100./asa;                              // p4
00516           // The 2nd max pre-exponent   (peh_ss)
00517           lastPAR[32]=(.1+4.4e-5*a2)/(1.+5.e5/a4);           // p1
00518           lastPAR[33]=3.5e-4*a2/(1.+1.e8/a8);                // p2
00519           lastPAR[34]=1.3+3.e5/a4;                           // p3
00520           lastPAR[35]=500./(a2+50.)+3;                       // p4
00521           lastPAR[36]=1.e-9/a+s4a16*s4a16;                   // p5
00522           // The 2nd max slope          (peh_sb)
00523           lastPAR[37]=.4*asa+3.e-9*a6;                       // p1
00524           lastPAR[38]=.0005*a5;                              // p2
00525           lastPAR[39]=.002*a5;                               // p3
00526           lastPAR[40]=10.;                                   // p4
00527           // The effective pre-exponent (peh_us)
00528           lastPAR[41]=.05+.005*a;                            // p1
00529           lastPAR[42]=7.e-8/sa;                              // p2
00530           lastPAR[43]=.8*sa;                                 // p3
00531           lastPAR[44]=.02*sa;                                // p4
00532           lastPAR[45]=1.e8/a3;                               // p5
00533           lastPAR[46]=3.e32/(a32+1.e32);                     // p6
00534           // The effective slope        (peh_ub)
00535           lastPAR[47]=24.;                                   // p1
00536           lastPAR[48]=20./sa;                                // p2
00537           lastPAR[49]=7.e3*a/(sa+1.);                        // p3
00538           lastPAR[50]=900.*sa/(1.+500./a3);                  // p4
00539         }
00540         // Parameter for lowEnergyNeutrons
00541         lastPAR[51]=1.e15+2.e27/a4/(1.+2.e-18*a16);
00542       }
00543       lastPAR[nLast]=pwd;
00544       // and initialize the zero element of the table
00545       G4double lp=lPMin;                                      // ln(momentum)
00546       G4bool memCS=onlyCS;                                    // ??
00547       onlyCS=false;
00548       lastCST[0]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables
00549       onlyCS=memCS;
00550       lastSST[0]=theSS;
00551       lastS1T[0]=theS1;
00552       lastB1T[0]=theB1;
00553       lastS2T[0]=theS2;
00554       lastB2T[0]=theB2;
00555       lastS3T[0]=theS3;
00556       lastB3T[0]=theB3;
00557       lastS4T[0]=theS4;
00558       lastB4T[0]=theB4;
00559     }
00560     if(LP>ILP)
00561     {
00562       G4int ini = static_cast<int>((ILP-lPMin+.000001)/dlnP)+1; // already inited till this
00563       if(ini<0) ini=0;
00564       if(ini<nPoints)
00565       {
00566         G4int fin = static_cast<int>((LP-lPMin)/dlnP)+1; // final bin of initialization
00567         if(fin>=nPoints) fin=nLast;               // Limit of the tabular initialization
00568         if(fin>=ini)
00569         {
00570           G4double lp=0.;
00571           for(G4int ip=ini; ip<=fin; ip++)        // Calculate tabular CS,S1,B1,S2,B2,S3,B3
00572           {
00573             lp=lPMin+ip*dlnP;                     // ln(momentum)
00574             G4bool memCS=onlyCS;
00575             onlyCS=false;
00576             lastCST[ip]=GetTabValues(lp, PDG, tgZ, tgN); // Calculate AMDB tables (ret CS)
00577             onlyCS=memCS;
00578             lastSST[ip]=theSS;
00579             lastS1T[ip]=theS1;
00580             lastB1T[ip]=theB1;
00581             lastS2T[ip]=theS2;
00582             lastB2T[ip]=theB2;
00583             lastS3T[ip]=theS3;
00584             lastB3T[ip]=theB3;
00585             lastS4T[ip]=theS4;
00586             lastB4T[ip]=theB4;
00587           }
00588           return lp;
00589         }
00590         else G4cout<<"*Warning*G4ChipsProtonElasticXS::GetPTables: PDG="<<PDG<<", Z="
00591                    <<tgZ<<", N="<<tgN<<", i="<<ini<<" > fin="<<fin<<", LP="<<LP<<" > ILP="
00592                    <<ILP<<" nothing is done!"<<G4endl;
00593       }
00594       else G4cout<<"*Warning*G4ChipsProtonElasticXS::GetPTables: PDG="<<PDG<<", Z="
00595                  <<tgZ<<", N="<<tgN<<", i="<<ini<<">= max="<<nPoints<<", LP="<<LP
00596                  <<" > ILP="<<ILP<<", lPMax="<<lPMax<<" nothing is done!"<<G4endl;
00597     }
00598   }
00599   else
00600   {
00601     // G4cout<<"*Error*G4ChipsProtonElasticXS::GetPTables: PDG="<<PDG<<", Z="<<tgZ
00602     //       <<", N="<<tgN<<", while it is defined only for PDG=2212"<<G4endl;
00603     // throw G4QException("G4ChipsProtonElasticXS::GetPTables: only pA're implemented");
00604     G4ExceptionDescription ed;
00605     ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
00606        << ", while it is defined only for PDG=2212 (p)" << G4endl;
00607     G4Exception("G4ChipsProtonElasticXS::GetPTables()", "HAD_CHPS_0000",
00608                 FatalException, ed);
00609   }
00610   return ILP;
00611 }
00612 
00613 // Returns Q2=-t in independent units (MeV^2) (all internal calculations are in GeV)
00614 G4double G4ChipsProtonElasticXS::GetExchangeT(G4int tgZ, G4int tgN, G4int PDG)
00615 {
00616   static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
00617   static const G4double third=1./3.;
00618   static const G4double fifth=1./5.;
00619   static const G4double sevth=1./7.;
00620   if(PDG!=2212) G4cout<<"**Warning*G4ChipsProtonElasticXS::GetExT:PDG="<<PDG<<G4endl;
00621   if(onlyCS) G4cout<<"**Warning*G4ChipsProtonElasticXS::GetExchanT:onlyCS=1"<<G4endl;
00622   if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
00623   G4double q2=0.;
00624   if(tgZ==1 && tgN==0)                // ===> p+p=p+p
00625   {
00626     G4double E1=lastTM*theB1;
00627     G4double R1=(1.-std::exp(-E1));
00628     G4double E2=lastTM*theB2;
00629     G4double R2=(1.-std::exp(-E2*E2*E2));
00630     G4double E3=lastTM*theB3;
00631     G4double R3=(1.-std::exp(-E3));
00632     G4double I1=R1*theS1/theB1;
00633     G4double I2=R2*theS2;
00634     G4double I3=R3*theS3;
00635     G4double I12=I1+I2;
00636     G4double rand=(I12+I3)*G4UniformRand();
00637     if     (rand<I1 )
00638     {
00639       G4double ran=R1*G4UniformRand();
00640       if(ran>1.) ran=1.;
00641       q2=-std::log(1.-ran)/theB1;
00642     }
00643     else if(rand<I12)
00644     {
00645       G4double ran=R2*G4UniformRand();
00646       if(ran>1.) ran=1.;
00647       q2=-std::log(1.-ran);
00648       if(q2<0.) q2=0.;
00649       q2=std::pow(q2,third)/theB2;
00650     }
00651     else
00652     {
00653       G4double ran=R3*G4UniformRand();
00654       if(ran>1.) ran=1.;
00655       q2=-std::log(1.-ran)/theB3;
00656     }
00657   }
00658   else
00659   {
00660     G4double a=tgZ+tgN;
00661     G4double E1=lastTM*(theB1+lastTM*theSS);
00662     G4double R1=(1.-std::exp(-E1));
00663     G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
00664     G4double tm2=lastTM*lastTM;
00665     G4double E2=lastTM*tm2*theB2;                   // power 3 for lowA, 5 for HighA (1st)
00666     if(a>6.5)E2*=tm2;                               // for heavy nuclei
00667     G4double R2=(1.-std::exp(-E2));
00668     G4double E3=lastTM*theB3;
00669     if(a>6.5)E3*=tm2*tm2*tm2;                       // power 1 for lowA, 7 (2nd) for HighA
00670     G4double R3=(1.-std::exp(-E3));
00671     G4double E4=lastTM*theB4;
00672     G4double R4=(1.-std::exp(-E4));
00673     G4double I1=R1*theS1;
00674     G4double I2=R2*theS2;
00675     G4double I3=R3*theS3;
00676     G4double I4=R4*theS4;
00677     G4double I12=I1+I2;
00678     G4double I13=I12+I3;
00679     G4double rand=(I13+I4)*G4UniformRand();
00680     if(rand<I1)
00681     {
00682       G4double ran=R1*G4UniformRand();
00683       if(ran>1.) ran=1.;
00684       q2=-std::log(1.-ran)/theB1;
00685       if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
00686     }
00687     else if(rand<I12)
00688     {
00689       G4double ran=R2*G4UniformRand();
00690       if(ran>1.) ran=1.;
00691       q2=-std::log(1.-ran)/theB2;
00692       if(q2<0.) q2=0.;
00693       if(a<6.5) q2=std::pow(q2,third);
00694       else      q2=std::pow(q2,fifth);
00695     }
00696     else if(rand<I13)
00697     {
00698       G4double ran=R3*G4UniformRand();
00699       if(ran>1.) ran=1.;
00700       q2=-std::log(1.-ran)/theB3;
00701       if(q2<0.) q2=0.;
00702       if(a>6.5) q2=std::pow(q2,sevth);
00703     }
00704     else
00705     {
00706       G4double ran=R4*G4UniformRand();
00707       if(ran>1.) ran=1.;
00708       q2=-std::log(1.-ran)/theB4;
00709       if(a<6.5) q2=lastTM-q2;                    // u reduced for lightA (starts from 0)
00710     }
00711   }
00712   if(q2<0.) q2=0.;
00713   if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QElasticCrossSect::GetExchangeT: -t="<<q2<<G4endl;
00714   if(q2>lastTM)
00715   {
00716     q2=lastTM;
00717   }
00718   return q2*GeVSQ;
00719 }
00720 
00721 // Returns B in independent units (MeV^-2) (all internal calculations are in GeV) see ExT
00722 G4double G4ChipsProtonElasticXS::GetSlope(G4int tgZ, G4int tgN, G4int PDG)
00723 {
00724   static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
00725   if(onlyCS) G4cout<<"*Warning*G4ChipsProtonElasticXS::GetSlope:onlyCS=true"<<G4endl;
00726   if(lastLP<-4.3) return 0.;          // S-wave for p<14 MeV/c (kinE<.1MeV)
00727   if(PDG!=2212)
00728   {
00729     // G4cout<<"*Error*G4ChipsProtonElasticXS::GetSlope: PDG="<<PDG<<", Z="<<tgZ<<", N="
00730     //       <<tgN<<", while it is defined only for PDG=2212"<<G4endl;
00731     // throw G4QException("G4ChipsProtonElasticXS::GetSlope: pA are implemented");
00732     G4ExceptionDescription ed;
00733     ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
00734        << ", while it is defined only for PDG=2212 (p)" << G4endl;
00735     G4Exception("G4ChipsProtonElasticXS::GetSlope()", "HAD_CHPS_0000",
00736                 FatalException, ed);
00737   }
00738   if(theB1<0.) theB1=0.;
00739   if(!(theB1>=-1.||theB1<=1.))G4cout<<"*NAN*G4QElasticCrossSect::Getslope:"<<theB1<<G4endl;
00740   return theB1/GeVSQ;
00741 }
00742 
00743 // Returns half max(Q2=-t) in independent units (MeV^2)
00744 G4double G4ChipsProtonElasticXS::GetHMaxT()
00745 {
00746   static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
00747   return lastTM*HGeVSQ;
00748 }
00749 
00750 // lastLP is used, so calculating tables, one need to remember and then recover lastLP
00751 G4double G4ChipsProtonElasticXS::GetTabValues(G4double lp, G4int PDG, G4int tgZ,
00752                                                     G4int tgN)
00753 {
00754   if(PDG!=2212) G4cout<<"*Warning*G4ChipsProtonElasticXS::GetTabV:PDG="<<PDG<<G4endl;
00755   if(tgZ<0 || tgZ>92)
00756   {
00757     G4cout<<"*Warning*G4QProtonElCS::GetTabValue: (1-92) No isotopes for Z="<<tgZ<<G4endl;
00758     return 0.;
00759   }
00760   G4int iZ=tgZ-1; // Z index
00761   if(iZ<0)
00762   {
00763     iZ=0;         // conversion of the neutron target to the proton target
00764     tgZ=1;
00765     tgN=0;
00766   }
00767   G4double p=std::exp(lp);              // momentum
00768   G4double sp=std::sqrt(p);             // sqrt(p)
00769   G4double p2=p*p;            
00770   G4double p3=p2*p;
00771   G4double p4=p3*p;
00772   if ( tgZ == 1 && tgN == 0 ) // pp/nn
00773   {
00774     G4double p2s=p2*sp;
00775     G4double dl2=lp-lastPAR[8];
00776     theSS=lastPAR[31];
00777     theS1=(lastPAR[9]+lastPAR[10]*dl2*dl2)/(1.+lastPAR[11]/p4/p)+
00778           (lastPAR[12]/p2+lastPAR[13]*p)/(p4+lastPAR[14]*sp);
00779     theB1=lastPAR[15]*std::pow(p,lastPAR[16])/(1.+lastPAR[17]/p3);
00780     theS2=lastPAR[18]+lastPAR[19]/(p4+lastPAR[20]*p);
00781     theB2=lastPAR[21]+lastPAR[22]/(p4+lastPAR[23]/sp); 
00782     theS3=lastPAR[24]+lastPAR[25]/(p4*p4+lastPAR[26]*p2+lastPAR[27]);
00783     theB3=lastPAR[28]+lastPAR[29]/(p4+lastPAR[30]); 
00784     theS4=0.;
00785     theB4=0.; 
00786     // Returns the total elastic pp cross-section (to avoid spoiling lastSIG)
00787     G4double dl1=lp-lastPAR[3];
00788     return lastPAR[0]/p2s/(1.+lastPAR[7]/p2s)+(lastPAR[1]+lastPAR[2]*dl1*dl1+lastPAR[4]/p)
00789                                                    /(1.+lastPAR[5]*lp)/(1.+lastPAR[6]/p4);
00790   }
00791   else
00792   {
00793     G4double p5=p4*p;
00794     G4double p6=p5*p;
00795     G4double p8=p6*p2;
00796     G4double p10=p8*p2;
00797     G4double p12=p10*p2;
00798     G4double p16=p8*p8;
00799     //G4double p24=p16*p8;
00800     G4double dl=lp-5.;
00801     G4double a=tgZ+tgN;
00802     G4double pah=std::pow(p,a/2);
00803     G4double pa=pah*pah;
00804     G4double pa2=pa*pa;
00805     if(a<6.5)
00806     {
00807       theS1=lastPAR[9]/(1.+lastPAR[10]*p4*pa)+lastPAR[11]/(p4+lastPAR[12]*p4/pa2)+
00808             (lastPAR[13]*dl*dl+lastPAR[14])/(1.+lastPAR[15]/p2);
00809       theB1=(lastPAR[16]+lastPAR[17]*p2)/(p4+lastPAR[18]/pah)+lastPAR[19];
00810       theSS=lastPAR[20]/(1.+lastPAR[21]/p2)+lastPAR[22]/(p6/pa+lastPAR[23]/p16);
00811       theS2=lastPAR[24]/(pa/p2+lastPAR[25]/p4)+lastPAR[26];
00812       theB2=lastPAR[27]*std::pow(p,lastPAR[28])+lastPAR[29]/(p8+lastPAR[30]/p16);
00813       theS3=lastPAR[31]/(pa*p+lastPAR[32]/pa)+lastPAR[33];
00814       theB3=lastPAR[34]/(p3+lastPAR[35]/p6)+lastPAR[36]/(1.+lastPAR[37]/p2);
00815       theS4=p2*(pah*lastPAR[38]*std::exp(-pah*lastPAR[39])+
00816                 lastPAR[40]/(1.+lastPAR[41]*std::pow(p,lastPAR[42])));
00817       theB4=lastPAR[43]*pa/p2/(1.+pa*lastPAR[44]);
00818     }
00819     else
00820     {
00821       theS1=lastPAR[9]/(1.+lastPAR[10]/p4)+lastPAR[11]/(p4+lastPAR[12]/p2)+
00822             lastPAR[13]/(p5+lastPAR[14]/p16);
00823       theB1=(lastPAR[15]/p8+lastPAR[19])/(p+lastPAR[16]/std::pow(p,lastPAR[20]))+
00824             lastPAR[17]/(1.+lastPAR[18]/p4);
00825       theSS=lastPAR[21]/(p4/std::pow(p,lastPAR[23])+lastPAR[22]/p4);
00826       theS2=lastPAR[24]/p4/(std::pow(p,lastPAR[25])+lastPAR[26]/p12)+lastPAR[27];
00827       theB2=lastPAR[28]/std::pow(p,lastPAR[29])+lastPAR[30]/std::pow(p,lastPAR[31]);
00828       theS3=lastPAR[32]/std::pow(p,lastPAR[35])/(1.+lastPAR[36]/p12)+
00829             lastPAR[33]/(1.+lastPAR[34]/p6);
00830       theB3=lastPAR[37]/p8+lastPAR[38]/p2+lastPAR[39]/(1.+lastPAR[40]/p8);
00831       theS4=(lastPAR[41]/p4+lastPAR[46]/p)/(1.+lastPAR[42]/p10)+
00832             (lastPAR[43]+lastPAR[44]*dl*dl)/(1.+lastPAR[45]/p12);
00833       theB4=lastPAR[47]/(1.+lastPAR[48]/p)+lastPAR[49]*p4/(1.+lastPAR[50]*p5);
00834     }
00835     // Returns the total elastic (n/p)A cross-section (to avoid spoiling lastSIG)
00836     //         p1               p2              p3            p6
00837     return (lastPAR[0]*dl*dl+lastPAR[1])/(1.+lastPAR[2]/p+lastPAR[5]/p6)+
00838      lastPAR[3]/(p3+lastPAR[4]/p3)+lastPAR[7]/(p4+std::pow((lastPAR[8]/p),lastPAR[6]));
00839     //   p4            p5               p8                 p9             p7
00840   }
00841   return 0.;
00842 } // End of GetTableValues
00843 
00844 // Returns max -t=Q2 (GeV^2) for the momentum pP(GeV) and the target nucleus (tgN,tgZ)
00845 G4double G4ChipsProtonElasticXS::GetQ2max(G4int PDG, G4int tgZ, G4int tgN,
00846                                                 G4double pP)
00847 {
00848   static const G4double mProt= G4Proton::Proton()->GetPDGMass()*.001; // MeV to GeV
00849   static const G4double mProt2= mProt*mProt;
00850 
00851   G4double pP2=pP*pP;                                 // squared momentum of the projectile
00852   if(tgZ==1 && tgN==0)
00853   {
00854     G4double tMid=std::sqrt(pP2+mProt2)*mProt-mProt2; // CMS 90deg value of -t=Q2 (GeV^2)
00855     return tMid+tMid;
00856   }
00857   else if(tgZ || tgN)                                 // ---> pA
00858   {
00859     G4double mt=G4ParticleTable::GetParticleTable()->FindIon(tgZ,tgZ+tgN,0,tgZ)->GetPDGMass()*.001; // Target mass in GeV
00860     G4double dmt=mt+mt;
00861     G4double mds=dmt*std::sqrt(pP2+mProt2)+mProt2+mt*mt;// Mondelstam mds
00862     return dmt*dmt*pP2/mds;
00863   }
00864   else
00865   {
00866     // G4cout<<"*Error*G4ChipsProtonElasticXS::GetQ2max: PDG="<<PDG<<", Z="<<tgZ<<", N="
00867     //       <<tgN<<", while it is defined only for p projectiles & Z_target>0"<<G4endl;
00868     // throw G4QException("G4ChipsProtonElasticXS::GetQ2max: only pA are implemented");
00869     G4ExceptionDescription ed;
00870     ed << "PDG = " << PDG << ", Z = " << tgZ << ", N = " << tgN
00871        << ", while it is defined only for p projectiles & Z_target>0" << G4endl;
00872     G4Exception("G4ChipsProtonElasticXS::GetQ2max()", "HAD_CHPS_0000",
00873                 FatalException, ed);
00874     return 0;
00875   }
00876 }

Generated on Mon May 27 17:47:53 2013 for Geant4 by  doxygen 1.4.7