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

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