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

#include <G4ChipsKaonMinusElasticXS.hh>

Inheritance diagram for G4ChipsKaonMinusElasticXS:
G4VCrossSectionDataSet

Public Member Functions

 G4ChipsKaonMinusElasticXS ()
 
 ~G4ChipsKaonMinusElasticXS ()
 
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 G4ChipsKaonMinusElasticXS.hh.

Constructor & Destructor Documentation

G4ChipsKaonMinusElasticXS::G4ChipsKaonMinusElasticXS ( )

Definition at line 76 of file G4ChipsKaonMinusElasticXS.cc.

References G4ParticleDefinition::GetPDGMass(), G4KaonMinus::KaonMinus(), and G4TemplateAutoLock< M, L, U >::unlock().

76  :G4VCrossSectionDataSet(Default_Name()), nPoints(128), nLast(nPoints-1)
77 {
78  G4AutoLock l(&initM);
79  mK = G4KaonMinus::KaonMinus()->GetPDGMass()*.001;
80  mK2 = mK*mK;
81  l.unlock();
82  lPMin=-8.; //Min tabulatedLogarithmMomentum/D
83  lPMax= 8.; //Max tabulatedLogarithmMomentum/D
84  dlnP=(lPMax-lPMin)/nLast;// LogStep inTable /D
85  onlyCS=true;//Flag toCalculOnlyCS(not Si/Bi)/L
86  lastSIG=0.; //Last calculated cross section /L
87  lastLP=-10.;//LastLog(mom_of IncidentHadron)/L
88  lastTM=0.; //Last t_maximum /L
89  theSS=0.; //TheLastSqSlope of 1st difr.Max/L
90  theS1=0.; //TheLastMantissa of 1st difrMax/L
91  theB1=0.; //TheLastSlope of 1st difructMax/L
92  theS2=0.; //TheLastMantissa of 2nd difrMax/L
93  theB2=0.; //TheLastSlope of 2nd difructMax/L
94  theS3=0.; //TheLastMantissa of 3d difr.Max/L
95  theB3=0.; //TheLastSlope of 3d difruct.Max/L
96  theS4=0.; //TheLastMantissa of 4th difrMax/L
97  theB4=0.; //TheLastSlope of 4th difructMax/L
98  lastTZ=0; // Last atomic number of theTarget
99  lastTN=0; // Last # of neutrons in theTarget
100  lastPIN=0.;// Last initialized max momentum
101  lastCST=0; // Elastic cross-section table
102  lastPAR=0; // ParametersForFunctionCalculation
103  lastSST=0; // E-dep ofSqardSlope of 1st difMax
104  lastS1T=0; // E-dep of mantissa of 1st dif.Max
105  lastB1T=0; // E-dep of the slope of 1st difMax
106  lastS2T=0; // E-dep of mantissa of 2nd difrMax
107  lastB2T=0; // E-dep of the slope of 2nd difMax
108  lastS3T=0; // E-dep of mantissa of 3d difr.Max
109  lastB3T=0; // E-dep of the slope of 3d difrMax
110  lastS4T=0; // E-dep of mantissa of 4th difrMax
111  lastB4T=0; // E-dep of the slope of 4th difMax
112  lastN=0; // The last N of calculated nucleus
113  lastZ=0; // The last Z of calculated nucleus
114  lastP=0.; // LastUsed inCrossSection Momentum
115  lastTH=0.; // Last threshold momentum
116  lastCS=0.; // Last value of the Cross Section
117  lastI=0; // The last position in the DAMDB
118 }
G4VCrossSectionDataSet(const G4String &nam="")
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
G4double GetPDGMass() const
G4ChipsKaonMinusElasticXS::~G4ChipsKaonMinusElasticXS ( )

Definition at line 120 of file G4ChipsKaonMinusElasticXS.cc.

121 {
122  std::vector<G4double*>::iterator pos;
123  for (pos=CST.begin(); pos<CST.end(); pos++)
124  { delete [] *pos; }
125  CST.clear();
126  for (pos=PAR.begin(); pos<PAR.end(); pos++)
127  { delete [] *pos; }
128  PAR.clear();
129  for (pos=SST.begin(); pos<SST.end(); pos++)
130  { delete [] *pos; }
131  SST.clear();
132  for (pos=S1T.begin(); pos<S1T.end(); pos++)
133  { delete [] *pos; }
134  S1T.clear();
135  for (pos=B1T.begin(); pos<B1T.end(); pos++)
136  { delete [] *pos; }
137  B1T.clear();
138  for (pos=S2T.begin(); pos<S2T.end(); pos++)
139  { delete [] *pos; }
140  S2T.clear();
141  for (pos=B2T.begin(); pos<B2T.end(); pos++)
142  { delete [] *pos; }
143  B2T.clear();
144  for (pos=S3T.begin(); pos<S3T.end(); pos++)
145  { delete [] *pos; }
146  S3T.clear();
147  for (pos=B3T.begin(); pos<B3T.end(); pos++)
148  { delete [] *pos; }
149  B3T.clear();
150  for (pos=S4T.begin(); pos<S4T.end(); pos++)
151  { delete [] *pos; }
152  S4T.clear();
153  for (pos=B4T.begin(); pos<B4T.end(); pos++)
154  { delete [] *pos; }
155  B4T.clear();
156 }

Member Function Documentation

static const char* G4ChipsKaonMinusElasticXS::Default_Name ( )
inlinestatic
G4double G4ChipsKaonMinusElasticXS::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 180 of file G4ChipsKaonMinusElasticXS.cc.

References G4ThreadLocal, and python.hepunit::millibarn.

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

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

Definition at line 619 of file G4ChipsKaonMinusElasticXS.cc.

References G4cout, G4endl, and G4UniformRand.

Referenced by G4ChipsElasticModel::SampleInvariantT().

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

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

173 {
174  G4double pMom=Pt->GetTotalMomentum();
175  G4int tgN = A - tgZ;
176 
177  return GetChipsCrossSection(pMom, tgZ, tgN, -321);
178 }
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 G4ChipsKaonMinusElasticXS::IsIsoApplicable ( const G4DynamicParticle Pt,
G4int  Z,
G4int  A,
const G4Element elm,
const G4Material mat 
)
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 158 of file G4ChipsKaonMinusElasticXS.cc.

References G4DynamicParticle::GetDefinition(), and G4KaonMinus::KaonMinus().

161 {
162  G4ParticleDefinition* particle = Pt->GetDefinition();
163  if (particle == G4KaonMinus::KaonMinus() ) return true;
164  return false;
165 }
G4ParticleDefinition * GetDefinition() const
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113

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