Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Static Public Member Functions
G4ChipsProtonElasticXS Class Reference

#include <G4ChipsProtonElasticXS.hh>

Inheritance diagram for G4ChipsProtonElasticXS:
G4VCrossSectionDataSet

Public Member Functions

 G4ChipsProtonElasticXS ()
 
 ~G4ChipsProtonElasticXS ()
 
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)
 
G4double GetHMaxT ()
 
- Public Member Functions inherited from G4VCrossSectionDataSet
 G4VCrossSectionDataSet (const G4String &nam="")
 
virtual ~G4VCrossSectionDataSet ()
 
virtual G4bool IsElementApplicable (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
G4double GetCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
G4double ComputeCrossSection (const G4DynamicParticle *, const G4Element *, const G4Material *mat=0)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
 
virtual G4IsotopeSelectIsotope (const G4Element *, G4double kinEnergy)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &)
 
virtual void CrossSectionDescription (std::ostream &) const
 
virtual G4int GetVerboseLevel () const
 
virtual void SetVerboseLevel (G4int value)
 
G4double GetMinKinEnergy () const
 
void SetMinKinEnergy (G4double value)
 
G4double GetMaxKinEnergy () const
 
void SetMaxKinEnergy (G4double value)
 
const G4StringGetName () const
 

Static Public Member Functions

static const char * Default_Name ()
 

Additional Inherited Members

- Protected Member Functions inherited from G4VCrossSectionDataSet
void SetName (const G4String &)
 
- Protected Attributes inherited from G4VCrossSectionDataSet
G4int verboseLevel
 

Detailed Description

Definition at line 45 of file G4ChipsProtonElasticXS.hh.

Constructor & Destructor Documentation

G4ChipsProtonElasticXS::G4ChipsProtonElasticXS ( )

Definition at line 56 of file G4ChipsProtonElasticXS.cc.

56  :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
57 {
58  // Initialization of the parameters
59  lPMin=-8.; // Min tabulated logarithmicMomentum(D)
60  lPMax= 8.; // Max tabulated logarithmicMomentum(D)
61  dlnP=(lPMax-lPMin)/nLast;// LogStep in the table(D)
62  onlyCS=false;// Flag toCalculateOnlyCS(not Si/Bi)(L)
63  lastSIG=0.; // Last calculated cross section (L)
64  lastLP=-10.;// Last log(mom_ofTheIncidentHadron)(L)
65  lastTM=0.; // Last t_maximum (L)
66  theSS=0.; // The Last sq.slope of 1st difr.Max(L)
67  theS1=0.; // The Last mantissa of 1st difr.Max(L)
68  theB1=0.; // The Last slope of 1st difruct.Max(L)
69  theS2=0.; // The Last mantissa of 2nd difr.Max(L)
70  theB2=0.; // The Last slope of 2nd difruct.Max(L)
71  theS3=0.; // The Last mantissa of 3d difr. Max(L)
72  theB3=0.; // The Last slope of 3d difruct. Max(L)
73  theS4=0.; // The Last mantissa of 4th difr.Max(L)
74  theB4=0.; // The Last slope of 4th difruct.Max(L)
75  lastTZ=0; // Last atomic number of the target
76  lastTN=0; // Last # of neutrons in the target
77  lastPIN=0.; // Last initialized max momentum
78  lastCST=0; // Elastic cross-section table
79  lastPAR=0; // Parameters for FunctionalCalculation
80  lastSST=0; // E-dep of sq.slope of the 1st dif.Max
81  lastS1T=0; // E-dep of mantissa of the 1st dif.Max
82  lastB1T=0; // E-dep of the slope of the 1st difMax
83  lastS2T=0; // E-dep of mantissa of the 2nd difrMax
84  lastB2T=0; // E-dep of the slope of the 2nd difMax
85  lastS3T=0; // E-dep of mantissa of the 3d difr.Max
86  lastB3T=0; // E-dep of the slope of the 3d difrMax
87  lastS4T=0; // E-dep of mantissa of the 4th difrMax
88  lastB4T=0; // E-dep of the slope of the 4th difMax
89  lastN=0; // The last N of calculated nucleus
90  lastZ=0; // The last Z of calculated nucleus
91  lastP=0.; // Last used in cross section Momentum
92  lastTH=0.; // Last threshold momentum
93  lastCS=0.; // Last value of the Cross Section
94  lastI=0; // The last position in the DAMDB
95 }
G4VCrossSectionDataSet(const G4String &nam="")
static const char * Default_Name()
G4ChipsProtonElasticXS::~G4ChipsProtonElasticXS ( )

Definition at line 98 of file G4ChipsProtonElasticXS.cc.

99 {
100  std::vector<G4double*>::iterator pos;
101  for (pos=CST.begin(); pos<CST.end(); pos++)
102  { delete [] *pos; }
103  CST.clear();
104  for (pos=PAR.begin(); pos<PAR.end(); pos++)
105  { delete [] *pos; }
106  PAR.clear();
107  for (pos=SST.begin(); pos<SST.end(); pos++)
108  { delete [] *pos; }
109  SST.clear();
110  for (pos=S1T.begin(); pos<S1T.end(); pos++)
111  { delete [] *pos; }
112  S1T.clear();
113  for (pos=B1T.begin(); pos<B1T.end(); pos++)
114  { delete [] *pos; }
115  B1T.clear();
116  for (pos=S2T.begin(); pos<S2T.end(); pos++)
117  { delete [] *pos; }
118  S2T.clear();
119  for (pos=B2T.begin(); pos<B2T.end(); pos++)
120  { delete [] *pos; }
121  B2T.clear();
122  for (pos=S3T.begin(); pos<S3T.end(); pos++)
123  { delete [] *pos; }
124  S3T.clear();
125  for (pos=B3T.begin(); pos<B3T.end(); pos++)
126  { delete [] *pos; }
127  B3T.clear();
128  for (pos=S4T.begin(); pos<S4T.end(); pos++)
129  { delete [] *pos; }
130  S4T.clear();
131  for (pos=B4T.begin(); pos<B4T.end(); pos++)
132  { delete [] *pos; }
133  B4T.clear();
134 
135 }

Member Function Documentation

static const char* G4ChipsProtonElasticXS::Default_Name ( )
inlinestatic
G4double G4ChipsProtonElasticXS::GetChipsCrossSection ( G4double  momentum,
G4int  Z,
G4int  N,
G4int  pdg 
)
virtual

!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)

Definition at line 161 of file G4ChipsProtonElasticXS.cc.

References G4ThreadLocal, and python.hepunit::millibarn.

Referenced by G4QuasiElRatios::ChExer(), GetIsoCrossSection(), G4ChipsElasticModel::SampleInvariantT(), and G4QuasiElRatios::Scatter().

162 {
163  static G4ThreadLocal std::vector <G4int> *colN_G4MT_TLS_ = 0 ; if (!colN_G4MT_TLS_) colN_G4MT_TLS_ = new std::vector <G4int> ; std::vector <G4int> &colN = *colN_G4MT_TLS_; // Vector of N for calculated nuclei (isotops)
164  static G4ThreadLocal std::vector <G4int> *colZ_G4MT_TLS_ = 0 ; if (!colZ_G4MT_TLS_) colZ_G4MT_TLS_ = new std::vector <G4int> ; std::vector <G4int> &colZ = *colZ_G4MT_TLS_; // Vector of Z for calculated nuclei (isotops)
165  static G4ThreadLocal std::vector <G4double> *colP_G4MT_TLS_ = 0 ; if (!colP_G4MT_TLS_) colP_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colP = *colP_G4MT_TLS_; // Vector of last momenta for the reaction
166  static G4ThreadLocal std::vector <G4double> *colTH_G4MT_TLS_ = 0 ; if (!colTH_G4MT_TLS_) colTH_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colTH = *colTH_G4MT_TLS_; // Vector of energy thresholds for the reaction
167  static G4ThreadLocal std::vector <G4double> *colCS_G4MT_TLS_ = 0 ; if (!colCS_G4MT_TLS_) colCS_G4MT_TLS_ = new std::vector <G4double> ; std::vector <G4double> &colCS = *colCS_G4MT_TLS_; // Vector of last cross sections for the reaction
168  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
169 
170  G4double pEn=pMom;
171  onlyCS=false;
172 
173  G4bool in=false; // By default the isotope must be found in the AMDB
174  lastP = 0.; // New momentum history (nothing to compare with)
175  lastN = tgN; // The last N of the calculated nucleus
176  lastZ = tgZ; // The last Z of the calculated nucleus
177  lastI = colN.size(); // Size of the Associative Memory DB in the heap
178  if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
179  { // The nucleus with projPDG is found in AMDB
180  if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
181  {
182  lastI=i;
183  lastTH =colTH[i]; // Last THreshold (A-dependent)
184  if(pEn<=lastTH)
185  {
186  return 0.; // Energy is below the Threshold value
187  }
188  lastP =colP [i]; // Last Momentum (A-dependent)
189  lastCS =colCS[i]; // Last CrossSect (A-dependent)
190  if(lastP == pMom) // Do not recalculate
191  {
192  CalculateCrossSection(onlyCS,-1,i,2212,lastZ,lastN,pMom); // Update param's only
193  return lastCS*millibarn; // Use theLastCS
194  }
195  in = true; // This is the case when the isotop is found in DB
196  // Momentum pMom is in IU ! @@ Units
197  lastCS=CalculateCrossSection(onlyCS,-1,i,2212,lastZ,lastN,pMom); // read & update
198  if(lastCS<=0. && pEn>lastTH) // Correct the threshold
199  {
200  lastTH=pEn;
201  }
202  break; // Go out of the LOOP with found lastI
203  }
204  } // End of attampt to find the nucleus in DB
205  if(!in) // This nucleus has not been calculated previously
206  {
207  //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
208  lastCS=CalculateCrossSection(onlyCS,0,lastI,2212,lastZ,lastN,pMom);//calculate&create
209  if(lastCS<=0.)
210  {
211  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
212  if(pEn>lastTH)
213  {
214  lastTH=pEn;
215  }
216  }
217  colN.push_back(tgN);
218  colZ.push_back(tgZ);
219  colP.push_back(pMom);
220  colTH.push_back(lastTH);
221  colCS.push_back(lastCS);
222  return lastCS*millibarn;
223  } // End of creation of the new set of parameters
224  else
225  {
226  colP[lastI]=pMom;
227  colCS[lastI]=lastCS;
228  }
229  return lastCS*millibarn;
230 }
#define G4ThreadLocal
Definition: tls.hh:52
int G4int
Definition: G4Types.hh:78
int millibarn
Definition: hepunit.py:40
bool G4bool
Definition: G4Types.hh:79
double G4double
Definition: G4Types.hh:76
G4double G4ChipsProtonElasticXS::GetExchangeT ( G4int  tZ,
G4int  tN,
G4int  pPDG 
)

Definition at line 615 of file G4ChipsProtonElasticXS.cc.

References G4cout, G4endl, G4UniformRand, and python.hepunit::gigaelectronvolt.

Referenced by G4QuasiElRatios::ChExer(), G4ChipsElasticModel::SampleInvariantT(), and G4QuasiElRatios::Scatter().

616 {
617  static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
618  static const G4double third=1./3.;
619  static const G4double fifth=1./5.;
620  static const G4double sevth=1./7.;
621  if(PDG!=2212) G4cout<<"**Warning*G4ChipsProtonElasticXS::GetExT:PDG="<<PDG<<G4endl;
622  if(onlyCS) G4cout<<"**Warning*G4ChipsProtonElasticXS::GetExchanT:onlyCS=1"<<G4endl;
623  if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
624  G4double q2=0.;
625  if(tgZ==1 && tgN==0) // ===> p+p=p+p
626  {
627  G4double E1=lastTM*theB1;
628  G4double R1=(1.-std::exp(-E1));
629  G4double E2=lastTM*theB2;
630  G4double R2=(1.-std::exp(-E2*E2*E2));
631  G4double E3=lastTM*theB3;
632  G4double R3=(1.-std::exp(-E3));
633  G4double I1=R1*theS1/theB1;
634  G4double I2=R2*theS2;
635  G4double I3=R3*theS3;
636  G4double I12=I1+I2;
637  G4double rand=(I12+I3)*G4UniformRand();
638  if (rand<I1 )
639  {
640  G4double ran=R1*G4UniformRand();
641  if(ran>1.) ran=1.;
642  q2=-std::log(1.-ran)/theB1;
643  }
644  else if(rand<I12)
645  {
646  G4double ran=R2*G4UniformRand();
647  if(ran>1.) ran=1.;
648  q2=-std::log(1.-ran);
649  if(q2<0.) q2=0.;
650  q2=std::pow(q2,third)/theB2;
651  }
652  else
653  {
654  G4double ran=R3*G4UniformRand();
655  if(ran>1.) ran=1.;
656  q2=-std::log(1.-ran)/theB3;
657  }
658  }
659  else
660  {
661  G4double a=tgZ+tgN;
662  G4double E1=lastTM*(theB1+lastTM*theSS);
663  G4double R1=(1.-std::exp(-E1));
664  G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
665  G4double tm2=lastTM*lastTM;
666  G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
667  if(a>6.5)E2*=tm2; // for heavy nuclei
668  G4double R2=(1.-std::exp(-E2));
669  G4double E3=lastTM*theB3;
670  if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
671  G4double R3=(1.-std::exp(-E3));
672  G4double E4=lastTM*theB4;
673  G4double R4=(1.-std::exp(-E4));
674  G4double I1=R1*theS1;
675  G4double I2=R2*theS2;
676  G4double I3=R3*theS3;
677  G4double I4=R4*theS4;
678  G4double I12=I1+I2;
679  G4double I13=I12+I3;
680  G4double rand=(I13+I4)*G4UniformRand();
681  if(rand<I1)
682  {
683  G4double ran=R1*G4UniformRand();
684  if(ran>1.) ran=1.;
685  q2=-std::log(1.-ran)/theB1;
686  if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
687  }
688  else if(rand<I12)
689  {
690  G4double ran=R2*G4UniformRand();
691  if(ran>1.) ran=1.;
692  q2=-std::log(1.-ran)/theB2;
693  if(q2<0.) q2=0.;
694  if(a<6.5) q2=std::pow(q2,third);
695  else q2=std::pow(q2,fifth);
696  }
697  else if(rand<I13)
698  {
699  G4double ran=R3*G4UniformRand();
700  if(ran>1.) ran=1.;
701  q2=-std::log(1.-ran)/theB3;
702  if(q2<0.) q2=0.;
703  if(a>6.5) q2=std::pow(q2,sevth);
704  }
705  else
706  {
707  G4double ran=R4*G4UniformRand();
708  if(ran>1.) ran=1.;
709  q2=-std::log(1.-ran)/theB4;
710  if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
711  }
712  }
713  if(q2<0.) q2=0.;
714  if(!(q2>=-1.||q2<=1.)) G4cout<<"*NAN*G4QElasticCrossSect::GetExchangeT: -t="<<q2<<G4endl;
715  if(q2>lastTM)
716  {
717  q2=lastTM;
718  }
719  return q2*GeVSQ;
720 }
int gigaelectronvolt
Definition: hepunit.py:110
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4ChipsProtonElasticXS::GetHMaxT ( )

Definition at line 745 of file G4ChipsProtonElasticXS.cc.

Referenced by G4QuasiElRatios::ChExer(), and G4QuasiElRatios::Scatter().

746 {
747  static const G4double HGeVSQ=gigaelectronvolt*gigaelectronvolt/2.;
748  return lastTM*HGeVSQ;
749 }
int gigaelectronvolt
Definition: hepunit.py:110
double G4double
Definition: G4Types.hh:76
G4double G4ChipsProtonElasticXS::GetIsoCrossSection ( const G4DynamicParticle Pt,
G4int  tgZ,
G4int  A,
const G4Isotope iso = 0,
const G4Element elm = 0,
const G4Material mat = 0 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 147 of file G4ChipsProtonElasticXS.cc.

References GetChipsCrossSection(), and G4DynamicParticle::GetTotalMomentum().

151 {
152  G4double pMom=Pt->GetTotalMomentum();
153  G4int tgN = A - tgZ;
154 
155  return GetChipsCrossSection(pMom, tgZ, tgN, 2212);
156 }
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
double G4double
Definition: G4Types.hh:76
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
G4bool G4ChipsProtonElasticXS::IsIsoApplicable ( const G4DynamicParticle Pt,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 137 of file G4ChipsProtonElasticXS.cc.

References G4DynamicParticle::GetDefinition(), and G4Proton::Proton().

140 {
141  G4ParticleDefinition* particle = Pt->GetDefinition();
142  if (particle == G4Proton::Proton() ) return true;
143  return false;
144 }
G4ParticleDefinition * GetDefinition() const
static G4Proton * Proton()
Definition: G4Proton.cc:93

The documentation for this class was generated from the following files: