G4ChipsAntiBaryonElasticXS Class Reference

#include <G4ChipsAntiBaryonElasticXS.hh>

Inheritance diagram for G4ChipsAntiBaryonElasticXS:

G4VCrossSectionDataSet

Public Member Functions

 G4ChipsAntiBaryonElasticXS ()
 ~G4ChipsAntiBaryonElasticXS ()
virtual G4bool IsIsoApplicable (const G4DynamicParticle *Pt, G4int Z, G4int A, const G4Element *elm, const G4Material *mat)
virtual G4double GetIsoCrossSection (const G4DynamicParticle *, G4int tgZ, G4int A, const G4Isotope *iso=0, const G4Element *elm=0, const G4Material *mat=0)
virtual G4double GetChipsCrossSection (G4double momentum, G4int Z, G4int N, G4int pdg)
G4double GetExchangeT (G4int tZ, G4int tN, G4int pPDG)

Static Public Member Functions

static const char * Default_Name ()

Detailed Description

Definition at line 47 of file G4ChipsAntiBaryonElasticXS.hh.


Constructor & Destructor Documentation

G4ChipsAntiBaryonElasticXS::G4ChipsAntiBaryonElasticXS (  ) 

Definition at line 55 of file G4ChipsAntiBaryonElasticXS.cc.

00055                                                       :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
00056 {
00057   lPMin=-8.;  //Min tabulatedLogarithmMomentum(D)
00058   lPMax= 8.;  //Max tabulatedLogarithmMomentum(D)
00059   dlnP=(lPMax-lPMin)/nLast;// LogStep inTable (D)
00060   onlyCS=true;//Flag toCalculOnlyCS(not Si/Bi)(L)
00061   lastSIG=0.; //Last calculated cross section (L)
00062   lastLP=-10.;//LastLog(mom_of IncidentHadron)(L)
00063   lastTM=0.;  //Last t_maximum                (L)
00064   theSS=0.;   //TheLastSqSlope of 1st difr.Max(L)
00065   theS1=0.;   //TheLastMantissa of 1st difrMax(L)
00066   theB1=0.;   //TheLastSlope of 1st difructMax(L)
00067   theS2=0.;   //TheLastMantissa of 2nd difrMax(L)
00068   theB2=0.;   //TheLastSlope of 2nd difructMax(L)
00069   theS3=0.;   //TheLastMantissa of 3d difr.Max(L)
00070   theB3=0.;   //TheLastSlope of 3d difruct.Max(L)
00071   theS4=0.;   //TheLastMantissa of 4th difrMax(L)
00072   theB4=0.;   //TheLastSlope of 4th difructMax(L)
00073   lastTZ=0;   // Last atomic number of the target
00074   lastTN=0;   // Last # of neutrons in the target
00075   lastPIN=0.; // Last initialized max momentum
00076   lastCST=0;  // Elastic cross-section table
00077   lastPAR=0;  // ParametersForFunctionCalculation
00078   lastSST=0;  // E-dep ofSqardSlope of 1st difMax
00079   lastS1T=0;  // E-dep of mantissa of 1st dif.Max
00080   lastB1T=0;  // E-dep of the slope of 1st difMax
00081   lastS2T=0;  // E-dep of mantissa of 2nd difrMax
00082   lastB2T=0;  // E-dep of the slope of 2nd difMax
00083   lastS3T=0;  // E-dep of mantissa of 3d difr.Max
00084   lastB3T=0;  // E-dep of the slope of 3d difrMax
00085   lastS4T=0;  // E-dep of mantissa of 4th difrMax
00086   lastB4T=0;  // E-dep of the slope of 4th difMax
00087   lastN=0;    // The last N of calculated nucleus
00088   lastZ=0;    // The last Z of calculated nucleus
00089   lastP=0.;   // LastUsed inCrossSection Momentum
00090   lastTH=0.;  // Last threshold momentum
00091   lastCS=0.;  // Last value of the Cross Section
00092   lastI=0;    // The last position in the DAMDB
00093 }

G4ChipsAntiBaryonElasticXS::~G4ChipsAntiBaryonElasticXS (  ) 

Definition at line 95 of file G4ChipsAntiBaryonElasticXS.cc.

00096 {
00097   std::vector<G4double*>::iterator pos;
00098   for (pos=CST.begin(); pos<CST.end(); pos++)
00099   { delete [] *pos; }
00100   CST.clear();
00101   for (pos=PAR.begin(); pos<PAR.end(); pos++)
00102   { delete [] *pos; }
00103   PAR.clear();
00104   for (pos=SST.begin(); pos<SST.end(); pos++)
00105   { delete [] *pos; }
00106   SST.clear();
00107   for (pos=S1T.begin(); pos<S1T.end(); pos++)
00108   { delete [] *pos; }
00109   S1T.clear();
00110   for (pos=B1T.begin(); pos<B1T.end(); pos++)
00111   { delete [] *pos; }
00112   B1T.clear();
00113   for (pos=S2T.begin(); pos<S2T.end(); pos++)
00114   { delete [] *pos; }
00115   S2T.clear();
00116   for (pos=B2T.begin(); pos<B2T.end(); pos++)
00117   { delete [] *pos; }
00118   B2T.clear();
00119   for (pos=S3T.begin(); pos<S3T.end(); pos++)
00120   { delete [] *pos; }
00121   S3T.clear();
00122   for (pos=B3T.begin(); pos<B3T.end(); pos++)
00123   { delete [] *pos; }
00124   B3T.clear();
00125   for (pos=S4T.begin(); pos<S4T.end(); pos++)
00126   { delete [] *pos; }
00127   S4T.clear();
00128   for (pos=B4T.begin(); pos<B4T.end(); pos++)
00129   { delete [] *pos; }
00130   B4T.clear();
00131 }


Member Function Documentation

static const char* G4ChipsAntiBaryonElasticXS::Default_Name (  )  [inline, static]

Definition at line 56 of file G4ChipsAntiBaryonElasticXS.hh.

Referenced by G4ChipsComponentXS::G4ChipsComponentXS(), and G4ChipsElasticModel::G4ChipsElasticModel().

00056 {return "ChipsAntiBaryonElasticXS";}

G4double G4ChipsAntiBaryonElasticXS::GetChipsCrossSection ( G4double  momentum,
G4int  Z,
G4int  N,
G4int  pdg 
) [virtual]

Definition at line 193 of file G4ChipsAntiBaryonElasticXS.cc.

Referenced by GetIsoCrossSection(), and G4ChipsElasticModel::SampleInvariantT().

00194 {
00195   static std::vector <G4int>    colN;  // Vector of N for calculated nuclei (isotops)
00196   static std::vector <G4int>    colZ;  // Vector of Z for calculated nuclei (isotops)
00197   static std::vector <G4double> colP;  // Vector of last momenta for the reaction
00198   static std::vector <G4double> colTH; // Vector of energy thresholds for the reaction
00199   static std::vector <G4double> colCS; // Vector of last cross sections for the reaction
00200   // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
00201 
00202   G4bool fCS = false;
00203 
00204   G4double pEn=pMom;
00205   onlyCS=fCS;
00206 
00207   G4bool in=false;                   // By default the isotope must be found in the AMDB
00208   lastP   = 0.;                      // New momentum history (nothing to compare with)
00209   lastN   = tgN;                     // The last N of the calculated nucleus
00210   lastZ   = tgZ;                     // The last Z of the calculated nucleus
00211   lastI   = colN.size();             // Size of the Associative Memory DB in the heap
00212   if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
00213   {                                  // The nucleus with projPDG is found in AMDB
00214     if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
00215     {
00216       lastI=i;
00217       lastTH =colTH[i];              // Last THreshold (A-dependent)
00218       if(pEn<=lastTH)
00219       {
00220         return 0.;                   // Energy is below the Threshold value
00221       }
00222       lastP  =colP [i];              // Last Momentum  (A-dependent)
00223       lastCS =colCS[i];              // Last CrossSect (A-dependent)
00224       //  if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
00225       if(lastP == pMom)              // Do not recalculate
00226       {
00227         CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // Update param's only
00228         return lastCS*millibarn;     // Use theLastCS
00229       }
00230       in = true;                       // This is the case when the isotop is found in DB
00231       // Momentum pMom is in IU ! @@ Units
00232       lastCS=CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // read & update
00233       if(lastCS<=0. && pEn>lastTH)    // Correct the threshold
00234       {
00235         lastTH=pEn;
00236       }
00237       break;                           // Go out of the LOOP with found lastI
00238     }
00239   } // End of attampt to find the nucleus in DB
00240   if(!in)                            // This nucleus has not been calculated previously
00241   {
00243     lastCS=CalculateCrossSection(fCS,0,lastI,pPDG,lastZ,lastN,pMom);//calculate&create
00244     if(lastCS<=0.)
00245     {
00246       lastTH = 0; // ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
00247       if(pEn>lastTH)
00248       {
00249         lastTH=pEn;
00250       }
00251     }
00252     colN.push_back(tgN);
00253     colZ.push_back(tgZ);
00254     colP.push_back(pMom);
00255     colTH.push_back(lastTH);
00256     colCS.push_back(lastCS);
00257     return lastCS*millibarn;
00258   } // End of creation of the new set of parameters
00259   else
00260   {
00261     colP[lastI]=pMom;
00262     colCS[lastI]=lastCS;
00263   }
00264   return lastCS*millibarn;
00265 }

G4double G4ChipsAntiBaryonElasticXS::GetExchangeT ( G4int  tZ,
G4int  tN,
G4int  pPDG 
)

Definition at line 640 of file G4ChipsAntiBaryonElasticXS.cc.

References G4cout, and G4UniformRand.

Referenced by G4ChipsElasticModel::SampleInvariantT().

00641 {
00642   static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
00643   static const G4double third=1./3.;
00644   static const G4double fifth=1./5.;
00645   static const G4double sevth=1./7.;
00646 
00647   if(PDG<-3334 || PDG>-1111)G4cout<<"*Warning*G4QAntiBaryonElCS::GetExT:PDG="<<PDG<<G4endl;
00648   if(onlyCS)G4cout<<"WarningG4ChipsAntiBaryonElasticXS::GetExchanT:onlyCS=1"<<G4endl;
00649   if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
00650   G4double q2=0.;
00651   if(tgZ==1 && tgN==0)                // ===> p+p=p+p
00652   {
00653     G4double E1=lastTM*theB1;
00654     G4double R1=(1.-std::exp(-E1));
00655     G4double E2=lastTM*theB2;
00656     G4double R2=(1.-std::exp(-E2*E2*E2));
00657     G4double E3=lastTM*theB3;
00658     G4double R3=(1.-std::exp(-E3));
00659     G4double I1=R1*theS1/theB1;
00660     G4double I2=R2*theS2;
00661     G4double I3=R3*theS3;
00662     G4double I12=I1+I2;
00663     G4double rand=(I12+I3)*G4UniformRand();
00664     if     (rand<I1 )
00665     {
00666       G4double ran=R1*G4UniformRand();
00667       if(ran>1.) ran=1.;
00668       q2=-std::log(1.-ran)/theB1;
00669     }
00670     else if(rand<I12)
00671     {
00672       G4double ran=R2*G4UniformRand();
00673       if(ran>1.) ran=1.;
00674       q2=-std::log(1.-ran);
00675       if(q2<0.) q2=0.;
00676       q2=std::pow(q2,third)/theB2;
00677     }
00678     else
00679     {
00680       G4double ran=R3*G4UniformRand();
00681       if(ran>1.) ran=1.;
00682       q2=-std::log(1.-ran)/theB3;
00683     }
00684   }
00685   else
00686   {
00687     G4double a=tgZ+tgN;
00688     G4double E1=lastTM*(theB1+lastTM*theSS);
00689     G4double R1=(1.-std::exp(-E1));
00690     G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
00691     G4double tm2=lastTM*lastTM;
00692     G4double E2=lastTM*tm2*theB2;                   // power 3 for lowA, 5 for HighA (1st)
00693     if(a>6.5)E2*=tm2;                               // for heavy nuclei
00694     G4double R2=(1.-std::exp(-E2));
00695     G4double E3=lastTM*theB3;
00696     if(a>6.5)E3*=tm2*tm2*tm2;                       // power 1 for lowA, 7 (2nd) for HighA
00697     G4double R3=(1.-std::exp(-E3));
00698     G4double E4=lastTM*theB4;
00699     G4double R4=(1.-std::exp(-E4));
00700     G4double I1=R1*theS1;
00701     G4double I2=R2*theS2;
00702     G4double I3=R3*theS3;
00703     G4double I4=R4*theS4;
00704     G4double I12=I1+I2;
00705     G4double I13=I12+I3;
00706     G4double rand=(I13+I4)*G4UniformRand();
00707     if(rand<I1)
00708     {
00709       G4double ran=R1*G4UniformRand();
00710       if(ran>1.) ran=1.;
00711       q2=-std::log(1.-ran)/theB1;
00712       if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
00713     }
00714     else if(rand<I12)
00715     {
00716       G4double ran=R2*G4UniformRand();
00717       if(ran>1.) ran=1.;
00718       q2=-std::log(1.-ran)/theB2;
00719       if(q2<0.) q2=0.;
00720       if(a<6.5) q2=std::pow(q2,third);
00721       else      q2=std::pow(q2,fifth);
00722     }
00723     else if(rand<I13)
00724     {
00725       G4double ran=R3*G4UniformRand();
00726       if(ran>1.) ran=1.;
00727       q2=-std::log(1.-ran)/theB3;
00728       if(q2<0.) q2=0.;
00729       if(a>6.5) q2=std::pow(q2,sevth);
00730     }
00731     else
00732     {
00733       G4double ran=R4*G4UniformRand();
00734       if(ran>1.) ran=1.;
00735       q2=-std::log(1.-ran)/theB4;
00736       if(a<6.5) q2=lastTM-q2;                    // u reduced for lightA (starts from 0)
00737     }
00738   }
00739   if(q2<0.) q2=0.;
00740   if(!(q2>=-1.||q2<=1.))G4cout<<"*NAN*G4QaBElasticCrossSect::GetExchangeT:-t="<<q2<<G4endl;
00741   if(q2>lastTM)
00742   {
00743     q2=lastTM;
00744   }
00745   return q2*GeVSQ;
00746 }

G4double G4ChipsAntiBaryonElasticXS::GetIsoCrossSection ( const G4DynamicParticle ,
G4int  tgZ,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
) [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 181 of file G4ChipsAntiBaryonElasticXS.cc.

References GetChipsCrossSection(), G4DynamicParticle::GetDefinition(), G4ParticleDefinition::GetPDGEncoding(), and G4DynamicParticle::GetTotalMomentum().

00185 {
00186   G4double pMom=Pt->GetTotalMomentum();
00187   G4int tgN = A - tgZ;
00188   G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
00189   
00190   return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
00191 }

G4bool G4ChipsAntiBaryonElasticXS::IsIsoApplicable ( const G4DynamicParticle Pt,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
) [virtual]

Reimplemented from G4VCrossSectionDataSet.

Definition at line 134 of file G4ChipsAntiBaryonElasticXS.cc.

References G4AntiLambda::AntiLambda(), G4AntiNeutron::AntiNeutron(), G4AntiOmegaMinus::AntiOmegaMinus(), G4AntiProton::AntiProton(), G4AntiSigmaMinus::AntiSigmaMinus(), G4AntiSigmaPlus::AntiSigmaPlus(), G4AntiSigmaZero::AntiSigmaZero(), G4AntiXiMinus::AntiXiMinus(), G4AntiXiZero::AntiXiZero(), and G4DynamicParticle::GetDefinition().

00137 {
00138   G4ParticleDefinition* particle = Pt->GetDefinition();
00139 
00140   if(particle == G4AntiNeutron::AntiNeutron())
00141   {
00142     return true;
00143   }
00144   else if(particle == G4AntiProton::AntiProton())
00145   {
00146     return true;
00147   }
00148   else if(particle == G4AntiLambda::AntiLambda())
00149   {
00150     return true;
00151   }
00152   else if(particle == G4AntiSigmaPlus::AntiSigmaPlus())
00153   {
00154     return true;
00155   }
00156   else if(particle == G4AntiSigmaMinus::AntiSigmaMinus())
00157   {
00158     return true;
00159   }
00160   else if(particle == G4AntiSigmaZero::AntiSigmaZero())
00161   {
00162     return true;
00163   }
00164   else if(particle == G4AntiXiMinus::AntiXiMinus())
00165   {
00166     return true;
00167   }
00168   else if(particle == G4AntiXiZero::AntiXiZero())
00169   {
00170     return true;
00171   }
00172   else if(particle == G4AntiOmegaMinus::AntiOmegaMinus())
00173   {
00174     return true;
00175   }
00176   return false;
00177 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:37 2013 for Geant4 by  doxygen 1.4.7