Geant4  10.02.p03
G4ParticleHPChannelList Class Reference

#include <G4ParticleHPChannelList.hh>

Collaboration diagram for G4ParticleHPChannelList:

Public Member Functions

 G4ParticleHPChannelList (G4int n, G4ParticleDefinition *projectile)
 
 G4ParticleHPChannelList ()
 
void Init (G4int n)
 
 ~G4ParticleHPChannelList ()
 
G4HadFinalStateApplyYourself (const G4Element *theElement, const G4HadProjectile &aTrack)
 
void Init (G4Element *anElement, const G4String &dirName, G4ParticleDefinition *projectile)
 
void Register (G4ParticleHPFinalState *theFS, const G4String &aName)
 
G4double GetXsec (G4double anEnergy)
 
G4int GetNumberOfChannels ()
 
G4bool HasDataInAnyFinalState ()
 
void RestartRegistration ()
 
void DumpInfo ()
 

Private Attributes

G4ParticleHPChannel ** theChannels
 
G4int nChannels
 
G4String theDir
 
G4ElementtheElement
 
G4bool allChannelsCreated
 
G4int theInitCount
 
G4StableIsotopes theStableOnes
 
G4ParticleDefinitiontheProjectile
 
G4HadFinalState unChanged
 

Static Private Attributes

static G4ThreadLocal G4int trycounter = 0
 

Detailed Description

Definition at line 47 of file G4ParticleHPChannelList.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPChannelList() [1/2]

G4ParticleHPChannelList::G4ParticleHPChannelList ( G4int  n,
G4ParticleDefinition projectile 
)

Definition at line 42 of file G4ParticleHPChannelList.cc.

◆ G4ParticleHPChannelList() [2/2]

G4ParticleHPChannelList::G4ParticleHPChannelList ( )

Definition at line 51 of file G4ParticleHPChannelList.cc.

◆ ~G4ParticleHPChannelList()

G4ParticleHPChannelList::~G4ParticleHPChannelList ( )

Definition at line 59 of file G4ParticleHPChannelList.cc.

60  {
61  if(theChannels!=0)
62  {
63  for(G4int i=0;i<nChannels; i++)
64  {
65  delete theChannels[i];
66  }
67  delete [] theChannels;
68  }
69  }
int G4int
Definition: G4Types.hh:78
G4ParticleHPChannel ** theChannels

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * G4ParticleHPChannelList::ApplyYourself ( const G4Element theElement,
const G4HadProjectile aTrack 
)

Definition at line 73 of file G4ParticleHPChannelList.cc.

74  {
75  G4ParticleHPThermalBoost aThermalE;
76  G4int i, ii;
77  // decide on the isotope
78  G4int numberOfIsos(0);
79  for(ii=0; ii<nChannels; ii++)
80  {
81  numberOfIsos = theChannels[ii]->GetNiso();
82  if(numberOfIsos!=0) break;
83  }
84  G4double * running= new G4double [numberOfIsos];
85  running[0] = 0;
86  for(i=0;i<numberOfIsos; i++)
87  {
88  if(i!=0) running[i] = running[i-1];
89  for(ii=0; ii<nChannels; ii++)
90  {
91  if(theChannels[ii]->HasAnyData(i))
92  {
93  running[i] +=theChannels[ii]->GetWeightedXsec(aThermalE.GetThermalEnergy(aTrack,
94  theChannels[ii]->GetN(i),
95  theChannels[ii]->GetZ(i),
96  aTrack.GetMaterial()->GetTemperature()),
97  i);
98  }
99  }
100  }
101  G4int isotope=nChannels-1;
102  G4double random=G4UniformRand();
103  for(i=0;i<numberOfIsos; i++)
104  {
105  isotope = i;
106  //if(random<running[i]/running[numberOfIsos-1]) break;
107  if(running[numberOfIsos-1] != 0) if(random<running[i]/running[numberOfIsos-1]) break;
108  }
109  delete [] running;
110 
111  // decide on the channel
112  running = new G4double[nChannels];
113  running[0]=0;
114  G4int targA=-1; // For production of unChanged
115  G4int targZ=-1;
116  for(i=0; i<nChannels; i++)
117  {
118  if(i!=0) running[i] = running[i-1];
119  if(theChannels[i]->HasAnyData(isotope))
120  {
121  targA=(G4int)theChannels[i]->GetN(isotope); //Will be simply used the last valid value
122  targZ=(G4int)theChannels[i]->GetZ(isotope);
123  running[i] += theChannels[i]->GetFSCrossSection(aThermalE.GetThermalEnergy(aTrack,
124  targA,
125  targZ,
126  aTrack.GetMaterial()->GetTemperature()),
127  isotope);
128  targA=(G4int)theChannels[i]->GetN(isotope); //Will be simply used the last valid value
129  targZ=(G4int)theChannels[i]->GetZ(isotope);
130  // G4cout << " G4ParticleHPChannelList " << i << " targA " << targA << " targZ " << targZ << " running " << running[i] << G4endl;
131  }
132  }
133 
134  //TK120607
135  if ( running[nChannels-1] == 0 )
136  {
137  //This happened usually by the miss match between the cross section data and model
138  if ( targA == -1 && targZ == -1 ) {
139  throw G4HadronicException(__FILE__, __LINE__, "ParticleHP model encounter lethal discrepancy with cross section data");
140  }
141 
142  //TK121106
143  G4cout << "Warning from NeutronHP: could not find proper reaction channel. This may cause by inconsistency between cross section and model. Unchanged final states are returned." << G4endl;
144  unChanged.Clear();
145 
146  //For Ep Check create unchanged final state including rest target
147  G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon ( targZ , targA , 0.0 );
148  G4DynamicParticle* targ_dp = new G4DynamicParticle( targ_pd , G4ThreeVector(1,0,0), 0.0 );
151  unChanged.AddSecondary(targ_dp);
152  //TK121106
155  return &unChanged;
156  }
157  //TK120607
158 
159 
160  G4int lChan=0;
161  random=G4UniformRand();
162  for(i=0; i<nChannels; i++)
163  {
164  lChan = i;
165  if(running[nChannels-1] != 0) if(random<running[i]/running[nChannels-1]) break;
166  }
167  delete [] running;
168 #ifdef G4PHPDEBUG
169  if( getenv("G4ParticleHPDebug") ) G4cout << " G4ParticleHPChannelList SELECTED ISOTOPE " << isotope << " SELECTED CHANNEL " << lChan << G4endl;
170 #endif
171  return theChannels[lChan]->ApplyYourself(aTrack, isotope);
172  }
static G4ParticleHPManager * GetInstance()
const G4Material * GetMaterial() const
CLHEP::Hep3Vector G4ThreeVector
const G4LorentzVector & Get4Momentum() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:491
G4double GetN(G4int i)
int G4int
Definition: G4Types.hh:78
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
G4double GetZ(G4int i)
G4double GetTemperature() const
Definition: G4Material.hh:182
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
void SetEnergyChange(G4double anEnergy)
G4double GetKineticEnergy() const
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
#define G4endl
Definition: G4ios.hh:61
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4ParticleHPChannel ** theChannels
void SetMomentumChange(const G4ThreeVector &aV)
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
Here is the call graph for this function:

◆ DumpInfo()

void G4ParticleHPChannelList::DumpInfo ( )

Definition at line 218 of file G4ParticleHPChannelList.cc.

218  {
219 
220  G4cout<<"================================================================"<<G4endl;
221  G4cout<<" Element: "<<theElement->GetName()<<G4endl;
222  G4cout<<" Number of channels: "<<nChannels<<G4endl;
223  G4cout<<" Projectile: "<<theProjectile->GetParticleName()<<G4endl;
224  G4cout<<" Directory name: "<<theDir<<G4endl;
225  for(int i=0;i<nChannels;i++){
227  G4cout<<"----------------------------------------------------------------"<<G4endl;
228  theChannels[i]->DumpInfo();
229  G4cout<<"----------------------------------------------------------------"<<G4endl;
230  }
231  }
232  G4cout<<"================================================================"<<G4endl;
233 
234 }
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
Definition: G4Element.hh:127
#define G4endl
Definition: G4ios.hh:61
G4ParticleHPChannel ** theChannels
G4ParticleDefinition * theProjectile
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetNumberOfChannels()

G4int G4ParticleHPChannelList::GetNumberOfChannels ( )
inline

Definition at line 77 of file G4ParticleHPChannelList.hh.

◆ GetXsec()

G4double G4ParticleHPChannelList::GetXsec ( G4double  anEnergy)
inline

Definition at line 66 of file G4ParticleHPChannelList.hh.

67  {
68  G4double result=0;
69  G4int i;
70  for(i=0; i<nChannels; i++)
71  {
72  result+=std::max(0., theChannels[i]->GetXsec(anEnergy));
73  }
74  return result;
75  }
int G4int
Definition: G4Types.hh:78
G4double GetXsec(G4double anEnergy)
double G4double
Definition: G4Types.hh:76
G4ParticleHPChannel ** theChannels

◆ HasDataInAnyFinalState()

G4bool G4ParticleHPChannelList::HasDataInAnyFinalState ( )
inline

Definition at line 79 of file G4ParticleHPChannelList.hh.

80  {
81  G4bool result = false;
82  G4int i;
83  for(i=0; i<nChannels; i++)
84  {
85  if(theChannels[i]->HasDataInAnyFinalState()) result = true;
86  }
87  return result;
88  }
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
G4ParticleHPChannel ** theChannels
Here is the caller graph for this function:

◆ Init() [1/2]

void G4ParticleHPChannelList::Init ( G4int  n)

◆ Init() [2/2]

void G4ParticleHPChannelList::Init ( G4Element anElement,
const G4String dirName,
G4ParticleDefinition projectile 
)

Definition at line 174 of file G4ParticleHPChannelList.cc.

175  {
176  theDir = dirName;
177 // G4cout << theDir << G4endl;
178  theElement = anElement;
179 // G4cout << theElement << G4endl;
180  theProjectile = projectile;
181  }
G4ParticleDefinition * theProjectile

◆ Register()

void G4ParticleHPChannelList::Register ( G4ParticleHPFinalState theFS,
const G4String aName 
)

Definition at line 183 of file G4ParticleHPChannelList.cc.

185  {
186  if(!allChannelsCreated)
187  {
188  if(nChannels!=0)
189  {
190  G4ParticleHPChannel ** theBuffer = new G4ParticleHPChannel * [nChannels+1];
191  G4int i;
192  for(i=0; i<nChannels; i++)
193  {
194  theBuffer[i] = theChannels[i];
195  }
196  delete [] theChannels;
197  theChannels = theBuffer;
198  }
199  else
200  {
202  }
203  G4String name;
204  name = aName+"/";
207  // theChannels[nChannels]->SetProjectile( theProjectile );
208  nChannels++;
209  }
210 
211  //110527TKDB Unnessary codes, Detected by gcc4.6 compiler
212  //G4bool result;
213  //result = theChannels[theInitCount]->Register(theFS);
215  theInitCount++;
216  }
G4String name
Definition: TRTMaterials.hh:40
int G4int
Definition: G4Types.hh:78
G4bool Register(G4ParticleHPFinalState *theFS)
void Init(G4Element *theElement, const G4String dirName)
G4ParticleHPChannel ** theChannels
G4ParticleDefinition * theProjectile
Here is the call graph for this function:

◆ RestartRegistration()

void G4ParticleHPChannelList::RestartRegistration ( )
inline

Definition at line 89 of file G4ParticleHPChannelList.hh.

Here is the call graph for this function:

Member Data Documentation

◆ allChannelsCreated

G4bool G4ParticleHPChannelList::allChannelsCreated
private

Definition at line 106 of file G4ParticleHPChannelList.hh.

◆ nChannels

G4int G4ParticleHPChannelList::nChannels
private

Definition at line 102 of file G4ParticleHPChannelList.hh.

◆ theChannels

G4ParticleHPChannel** G4ParticleHPChannelList::theChannels
private

Definition at line 101 of file G4ParticleHPChannelList.hh.

◆ theDir

G4String G4ParticleHPChannelList::theDir
private

Definition at line 103 of file G4ParticleHPChannelList.hh.

◆ theElement

G4Element* G4ParticleHPChannelList::theElement
private

Definition at line 104 of file G4ParticleHPChannelList.hh.

◆ theInitCount

G4int G4ParticleHPChannelList::theInitCount
private

Definition at line 107 of file G4ParticleHPChannelList.hh.

◆ theProjectile

G4ParticleDefinition* G4ParticleHPChannelList::theProjectile
private

Definition at line 111 of file G4ParticleHPChannelList.hh.

◆ theStableOnes

G4StableIsotopes G4ParticleHPChannelList::theStableOnes
private

Definition at line 109 of file G4ParticleHPChannelList.hh.

◆ trycounter

G4ThreadLocal G4int G4ParticleHPChannelList::trycounter = 0
staticprivate

Definition at line 100 of file G4ParticleHPChannelList.hh.

◆ unChanged

G4HadFinalState G4ParticleHPChannelList::unChanged
private

Definition at line 113 of file G4ParticleHPChannelList.hh.


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