Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ChipsKaonPlusElasticXS Class Reference

#include <G4ChipsKaonPlusElasticXS.hh>

Inheritance diagram for G4ChipsKaonPlusElasticXS:
Collaboration diagram for G4ChipsKaonPlusElasticXS:

Public Member Functions

 G4ChipsKaonPlusElasticXS ()
 
 ~G4ChipsKaonPlusElasticXS ()
 
virtual void CrossSectionDescription (std::ostream &) const
 
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 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 G4ChipsKaonPlusElasticXS.hh.

Constructor & Destructor Documentation

G4ChipsKaonPlusElasticXS::G4ChipsKaonPlusElasticXS ( )

Definition at line 77 of file G4ChipsKaonPlusElasticXS.cc.

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

Here is the call graph for this function:

G4ChipsKaonPlusElasticXS::~G4ChipsKaonPlusElasticXS ( )

Definition at line 121 of file G4ChipsKaonPlusElasticXS.cc.

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

Member Function Documentation

void G4ChipsKaonPlusElasticXS::CrossSectionDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4VCrossSectionDataSet.

Definition at line 160 of file G4ChipsKaonPlusElasticXS.cc.

161 {
162  outFile << "G4ChipsKaonPlusElasticXS provides the elastic cross\n"
163  << "section for K+ nucleus scattering as a function of incident\n"
164  << "momentum. The cross section is calculated using M. Kossov's\n"
165  << "CHIPS parameterization of cross section data.\n";
166 }
static const char* G4ChipsKaonPlusElasticXS::Default_Name ( )
inlinestatic

Definition at line 54 of file G4ChipsKaonPlusElasticXS.hh.

54 {return "ChipsKaonPlusElasticXS";}

Here is the caller graph for this function:

G4double G4ChipsKaonPlusElasticXS::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 G4ChipsKaonPlusElasticXS.cc.

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

Here is the caller graph for this function:

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

Definition at line 617 of file G4ChipsKaonPlusElasticXS.cc.

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

Here is the caller graph for this function:

G4double G4ChipsKaonPlusElasticXS::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 178 of file G4ChipsKaonPlusElasticXS.cc.

182 {
183  G4double pMom=Pt->GetTotalMomentum();
184  G4int tgN = A - tgZ;
185 
186  return GetChipsCrossSection(pMom, tgZ, tgN, 321);
187 }
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
double A(double temperature)
virtual G4double GetChipsCrossSection(G4double momentum, G4int Z, G4int N, G4int pdg)
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Reimplemented from G4VCrossSectionDataSet.

Definition at line 169 of file G4ChipsKaonPlusElasticXS.cc.

172 {
173  return true;
174 }

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