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

#include <G4ChipsHyperonElasticXS.hh>

Inheritance diagram for G4ChipsHyperonElasticXS:
G4VCrossSectionDataSet

Public Member Functions

 G4ChipsHyperonElasticXS ()
 
 ~G4ChipsHyperonElasticXS ()
 
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)
 
- 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 46 of file G4ChipsHyperonElasticXS.hh.

Constructor & Destructor Documentation

G4ChipsHyperonElasticXS::G4ChipsHyperonElasticXS ( )

Definition at line 61 of file G4ChipsHyperonElasticXS.cc.

61  :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
62 {
63  lPMin=-8.; //Min tabulatedLogarithmMomentum(D)
64  lPMax= 8.; //Max tabulatedLogarithmMomentum(D)
65  dlnP=(lPMax-lPMin)/nLast;// LogStep inTable (D)
66  onlyCS=true;//Flag toCalculOnlyCS(not Si/Bi)(L)
67  lastSIG=0.; //Last calculated cross section (L)
68  lastLP=-10.;//LastLog(mom_of IncidentHadron)(L)
69  lastTM=0.; //Last t_maximum (L)
70  theSS=0.; //TheLastSqSlope of 1st difr.Max(L)
71  theS1=0.; //TheLastMantissa of 1st difrMax(L)
72  theB1=0.; //TheLastSlope of 1st difructMax(L)
73  theS2=0.; //TheLastMantissa of 2nd difrMax(L)
74  theB2=0.; //TheLastSlope of 2nd difructMax(L)
75  theS3=0.; //TheLastMantissa of 3d difr.Max(L)
76  theB3=0.; //TheLastSlope of 3d difruct.Max(L)
77  theS4=0.; //TheLastMantissa of 4th difrMax(L)
78  theB4=0.; //TheLastSlope of 4th difructMax(L)
79  lastTZ=0; // Last atomic number of the target
80  lastTN=0; // Last # of neutrons in the target
81  lastPIN=0.; // Last initialized max momentum
82  lastCST=0; // Elastic cross-section table
83  lastPAR=0; // ParametersForFunctionCalculation
84  lastSST=0; // E-dep ofSqardSlope of 1st difMax
85  lastS1T=0; // E-dep of mantissa of 1st dif.Max
86  lastB1T=0; // E-dep of the slope of 1st difMax
87  lastS2T=0; // E-dep of mantissa of 2nd difrMax
88  lastB2T=0; // E-dep of the slope of 2nd difMax
89  lastS3T=0; // E-dep of mantissa of 3d difr.Max
90  lastB3T=0; // E-dep of the slope of 3d difrMax
91  lastS4T=0; // E-dep of mantissa of 4th difrMax
92  lastB4T=0; // E-dep of the slope of 4th difMax
93  lastN=0; // The last N of calculated nucleus
94  lastZ=0; // The last Z of calculated nucleus
95  lastP=0.; // LastUsed inCrossSection Momentum
96  lastTH=0.; // Last threshold momentum
97  lastCS=0.; // Last value of the Cross Section
98  lastI=0; // The last position in the DAMDB
99 }
G4VCrossSectionDataSet(const G4String &nam="")
static const char * Default_Name()
G4ChipsHyperonElasticXS::~G4ChipsHyperonElasticXS ( )

Definition at line 101 of file G4ChipsHyperonElasticXS.cc.

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

Member Function Documentation

static const char* G4ChipsHyperonElasticXS::Default_Name ( )
inlinestatic

Definition at line 54 of file G4ChipsHyperonElasticXS.hh.

Referenced by G4ChipsComponentXS::G4ChipsComponentXS().

54 {return "ChipsHyperonElasticXS";}
G4double G4ChipsHyperonElasticXS::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 189 of file G4ChipsHyperonElasticXS.cc.

References G4ThreadLocal, and python.hepunit::millibarn.

Referenced by GetIsoCrossSection().

190 {
191  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)
192  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)
193  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
194  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
195  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
196  // ***---*** End of the mandatory Static Definitions of the Associative Memory ***---***
197 
198  G4bool fCS = false;
199  G4double pEn=pMom;
200 
201  onlyCS=fCS;
202 
203  G4bool in=false; // By default the isotope must be found in the AMDB
204  lastP = 0.; // New momentum history (nothing to compare with)
205  lastN = tgN; // The last N of the calculated nucleus
206  lastZ = tgZ; // The last Z of the calculated nucleus
207  lastI = colN.size(); // Size of the Associative Memory DB in the heap
208  if(lastI) for(G4int i=0; i<lastI; i++) // Loop over proj/tgZ/tgN lines of DB
209  { // The nucleus with projPDG is found in AMDB
210  if(colN[i]==tgN && colZ[i]==tgZ) // Isotope is foind in AMDB
211  {
212  lastI=i;
213  lastTH =colTH[i]; // Last THreshold (A-dependent)
214  if(pEn<=lastTH)
215  {
216  return 0.; // Energy is below the Threshold value
217  }
218  lastP =colP [i]; // Last Momentum (A-dependent)
219  lastCS =colCS[i]; // Last CrossSect (A-dependent)
220  // if(std::fabs(lastP/pMom-1.)<tolerance) //VI (do not use tolerance)
221  if(lastP == pMom) // Do not recalculate
222  {
223  CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // Update param's only
224  return lastCS*millibarn; // Use theLastCS
225  }
226  in = true; // This is the case when the isotop is found in DB
227  // Momentum pMom is in IU ! @@ Units
228  lastCS=CalculateCrossSection(fCS,-1,i,pPDG,lastZ,lastN,pMom); // read & update
229  if(lastCS<=0. && pEn>lastTH) // Correct the threshold
230  {
231  lastTH=pEn;
232  }
233  break; // Go out of the LOOP with found lastI
234  }
235  } // End of attampt to find the nucleus in DB
236  if(!in) // This nucleus has not been calculated previously
237  {
238  //!!The slave functions must provide cross-sections in millibarns (mb) !! (not in IU)
239  lastCS=CalculateCrossSection(fCS,0,lastI,pPDG,lastZ,lastN,pMom);//calculate&create
240  if(lastCS<=0.)
241  {
242  lastTH = 0; //ThresholdEnergy(tgZ, tgN); // The Threshold Energy which is now the last
243  if(pEn>lastTH)
244  {
245  lastTH=pEn;
246  }
247  }
248  colN.push_back(tgN);
249  colZ.push_back(tgZ);
250  colP.push_back(pMom);
251  colTH.push_back(lastTH);
252  colCS.push_back(lastCS);
253  return lastCS*millibarn;
254  } // End of creation of the new set of parameters
255  else
256  {
257  colP[lastI]=pMom;
258  colCS[lastI]=lastCS;
259  }
260  return lastCS*millibarn;
261 }
#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 G4ChipsHyperonElasticXS::GetExchangeT ( G4int  tZ,
G4int  tN,
G4int  pPDG 
)

Definition at line 635 of file G4ChipsHyperonElasticXS.cc.

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

636 {
637  static const G4double GeVSQ=gigaelectronvolt*gigaelectronvolt;
638  static const G4double third=1./3.;
639  static const G4double fifth=1./5.;
640  static const G4double sevth=1./7.;
641  if(PDG==3222 || PDG<3000 || PDG>3334)G4cout<<"*Warning*G4QHyElCS::GET:PDG="<<PDG<<G4endl;
642  if(onlyCS)G4cout<<"*Warning*G4ChipsHyperonElasticXS::GetExchanT: onlyCS=1"<<G4endl;
643  if(lastLP<-4.3) return lastTM*GeVSQ*G4UniformRand();// S-wave for p<14 MeV/c (kinE<.1MeV)
644  G4double q2=0.;
645  if(tgZ==1 && tgN==0) // ===> p+p=p+p
646  {
647  G4double E1=lastTM*theB1;
648  G4double R1=(1.-std::exp(-E1));
649  G4double E2=lastTM*theB2;
650  G4double R2=(1.-std::exp(-E2*E2*E2));
651  G4double E3=lastTM*theB3;
652  G4double R3=(1.-std::exp(-E3));
653  G4double I1=R1*theS1/theB1;
654  G4double I2=R2*theS2;
655  G4double I3=R3*theS3;
656  G4double I12=I1+I2;
657  G4double rand=(I12+I3)*G4UniformRand();
658  if (rand<I1 )
659  {
660  G4double ran=R1*G4UniformRand();
661  if(ran>1.) ran=1.;
662  q2=-std::log(1.-ran)/theB1;
663  }
664  else if(rand<I12)
665  {
666  G4double ran=R2*G4UniformRand();
667  if(ran>1.) ran=1.;
668  q2=-std::log(1.-ran);
669  if(q2<0.) q2=0.;
670  q2=std::pow(q2,third)/theB2;
671  }
672  else
673  {
674  G4double ran=R3*G4UniformRand();
675  if(ran>1.) ran=1.;
676  q2=-std::log(1.-ran)/theB3;
677  }
678  }
679  else
680  {
681  G4double a=tgZ+tgN;
682  G4double E1=lastTM*(theB1+lastTM*theSS);
683  G4double R1=(1.-std::exp(-E1));
684  G4double tss=theSS+theSS; // for future solution of quadratic equation (imediate check)
685  G4double tm2=lastTM*lastTM;
686  G4double E2=lastTM*tm2*theB2; // power 3 for lowA, 5 for HighA (1st)
687  if(a>6.5)E2*=tm2; // for heavy nuclei
688  G4double R2=(1.-std::exp(-E2));
689  G4double E3=lastTM*theB3;
690  if(a>6.5)E3*=tm2*tm2*tm2; // power 1 for lowA, 7 (2nd) for HighA
691  G4double R3=(1.-std::exp(-E3));
692  G4double E4=lastTM*theB4;
693  G4double R4=(1.-std::exp(-E4));
694  G4double I1=R1*theS1;
695  G4double I2=R2*theS2;
696  G4double I3=R3*theS3;
697  G4double I4=R4*theS4;
698  G4double I12=I1+I2;
699  G4double I13=I12+I3;
700  G4double rand=(I13+I4)*G4UniformRand();
701  if(rand<I1)
702  {
703  G4double ran=R1*G4UniformRand();
704  if(ran>1.) ran=1.;
705  q2=-std::log(1.-ran)/theB1;
706  if(std::fabs(tss)>1.e-7) q2=(std::sqrt(theB1*(theB1+(tss+tss)*q2))-theB1)/tss;
707  }
708  else if(rand<I12)
709  {
710  G4double ran=R2*G4UniformRand();
711  if(ran>1.) ran=1.;
712  q2=-std::log(1.-ran)/theB2;
713  if(q2<0.) q2=0.;
714  if(a<6.5) q2=std::pow(q2,third);
715  else q2=std::pow(q2,fifth);
716  }
717  else if(rand<I13)
718  {
719  G4double ran=R3*G4UniformRand();
720  if(ran>1.) ran=1.;
721  q2=-std::log(1.-ran)/theB3;
722  if(q2<0.) q2=0.;
723  if(a>6.5) q2=std::pow(q2,sevth);
724  }
725  else
726  {
727  G4double ran=R4*G4UniformRand();
728  if(ran>1.) ran=1.;
729  q2=-std::log(1.-ran)/theB4;
730  if(a<6.5) q2=lastTM-q2; // u reduced for lightA (starts from 0)
731  }
732  }
733  if(q2<0.) q2=0.;
734  if(!(q2>=-1.||q2<=1.))G4cout<<"*NAN*G4QHyElasticCrossSect::GetExchangeT:-t="<<q2<<G4endl;
735  if(q2>lastTM)
736  {
737  q2=lastTM;
738  }
739  return q2*GeVSQ;
740 }
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 G4ChipsHyperonElasticXS::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 177 of file G4ChipsHyperonElasticXS.cc.

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

181 {
182  G4double pMom=Pt->GetTotalMomentum();
183  G4int tgN = A - tgZ;
184  G4int pdg = Pt->GetDefinition()->GetPDGEncoding();
185 
186  return GetChipsCrossSection(pMom, tgZ, tgN, pdg);
187 }
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
double G4double
Definition: G4Types.hh:76
G4bool G4ChipsHyperonElasticXS::IsIsoApplicable ( const G4DynamicParticle Pt,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 139 of file G4ChipsHyperonElasticXS.cc.

References G4DynamicParticle::GetDefinition(), G4Lambda::Lambda(), G4OmegaMinus::OmegaMinus(), G4SigmaMinus::SigmaMinus(), G4SigmaPlus::SigmaPlus(), G4SigmaZero::SigmaZero(), G4XiMinus::XiMinus(), and G4XiZero::XiZero().

142 {
143  G4ParticleDefinition* particle = Pt->GetDefinition();
144  if (particle == G4Lambda::Lambda())
145  {
146  return true;
147  }
148  else if(particle == G4SigmaPlus::SigmaPlus())
149  {
150  return true;
151  }
152  else if(particle == G4SigmaMinus::SigmaMinus())
153  {
154  return true;
155  }
156  else if(particle == G4SigmaZero::SigmaZero())
157  {
158  return true;
159  }
160  else if(particle == G4XiMinus::XiMinus())
161  {
162  return true;
163  }
164  else if(particle == G4XiZero::XiZero())
165  {
166  return true;
167  }
168  else if(particle == G4OmegaMinus::OmegaMinus())
169  {
170  return true;
171  }
172  return false;
173 }
static G4OmegaMinus * OmegaMinus()
G4ParticleDefinition * GetDefinition() const
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:99
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108

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