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

#include <G4NeutronHPChannelList.hh>

Public Member Functions

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

Detailed Description

Definition at line 44 of file G4NeutronHPChannelList.hh.

Constructor & Destructor Documentation

G4NeutronHPChannelList::G4NeutronHPChannelList ( G4int  n)

Definition at line 41 of file G4NeutronHPChannelList.cc.

References n.

42  {
43  nChannels = n;
44  theChannels = new G4NeutronHPChannel * [n];
45  allChannelsCreated = false;
46  theInitCount = 0;
47  }
const G4int n
G4NeutronHPChannelList::G4NeutronHPChannelList ( )

Definition at line 49 of file G4NeutronHPChannelList.cc.

50  {
51  nChannels = 0;
52  theChannels = 0;
53  allChannelsCreated = false;
54  theInitCount = 0;
55  }
G4NeutronHPChannelList::~G4NeutronHPChannelList ( )

Definition at line 57 of file G4NeutronHPChannelList.cc.

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

Member Function Documentation

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

Definition at line 71 of file G4NeutronHPChannelList.cc.

References G4HadFinalState::AddSecondary(), G4NeutronHPChannel::ApplyYourself(), G4HadFinalState::Clear(), G4cout, G4endl, G4UniformRand, G4HadProjectile::Get4Momentum(), G4NeutronHPChannel::GetFSCrossSection(), G4NeutronHPManager::GetInstance(), G4IonTable::GetIon(), G4IonTable::GetIonTable(), G4HadProjectile::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4NeutronHPChannel::GetN(), G4NeutronHPChannel::GetNiso(), G4NeutronHPManager::GetReactionWhiteBoard(), G4Material::GetTemperature(), G4NeutronHPThermalBoost::GetThermalEnergy(), G4NeutronHPChannel::GetWeightedXsec(), G4NeutronHPChannel::GetZ(), G4HadFinalState::SetEnergyChange(), G4HadFinalState::SetMomentumChange(), G4NeutronHPReactionWhiteBoard::SetTargA(), G4NeutronHPReactionWhiteBoard::SetTargZ(), and CLHEP::HepLorentzVector::vect().

72  {
73  G4NeutronHPThermalBoost aThermalE;
74  G4int i, ii;
75  // decide on the isotope
76  G4int numberOfIsos(0);
77  for(ii=0; ii<nChannels; ii++)
78  {
79  numberOfIsos = theChannels[ii]->GetNiso();
80  if(numberOfIsos!=0) break;
81  }
82  G4double * running= new G4double [numberOfIsos];
83  running[0] = 0;
84  for(i=0;i<numberOfIsos; i++)
85  {
86  if(i!=0) running[i] = running[i-1];
87  for(ii=0; ii<nChannels; ii++)
88  {
89  if(theChannels[ii]->HasAnyData(i))
90  {
91  running[i] +=theChannels[ii]->GetWeightedXsec(aThermalE.GetThermalEnergy(aTrack,
92  theChannels[ii]->GetN(i),
93  theChannels[ii]->GetZ(i),
94  aTrack.GetMaterial()->GetTemperature()),
95  i);
96  }
97  }
98  }
99  G4int isotope=nChannels-1;
100  G4double random=G4UniformRand();
101  for(i=0;i<numberOfIsos; i++)
102  {
103  isotope = i;
104  //if(random<running[i]/running[numberOfIsos-1]) break;
105  if(running[numberOfIsos-1] != 0) if(random<running[i]/running[numberOfIsos-1]) break;
106  }
107  delete [] running;
108 
109  // decide on the channel
110  running = new G4double[nChannels];
111  running[0]=0;
112  G4int targA=-1; // For production of unChanged
113  G4int targZ=-1;
114  for(i=0; i<nChannels; i++)
115  {
116  if(i!=0) running[i] = running[i-1];
117  if(theChannels[i]->HasAnyData(isotope))
118  {
119  running[i] += theChannels[i]->GetFSCrossSection(aThermalE.GetThermalEnergy(aTrack,
120  theChannels[i]->GetN(isotope),
121  theChannels[i]->GetZ(isotope),
122  aTrack.GetMaterial()->GetTemperature()),
123  isotope);
124  targA=(G4int)theChannels[i]->GetN(isotope); //Will be simply used the last valid value
125  targZ=(G4int)theChannels[i]->GetZ(isotope);
126  }
127  }
128 
129  //TK120607
130  if ( running[nChannels-1] == 0 )
131  {
132  //This happened usually by the miss match between the cross section data and model
133  if ( targA == -1 && targZ == -1 ) {
134  throw G4HadronicException(__FILE__, __LINE__, "NeutronHP model encounter lethal discrepancy with cross section data");
135  }
136 
137  //TK121106
138  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;
139  unChanged.Clear();
140 
141  //For Ep Check create unchanged final state including rest target
142  G4ParticleDefinition* targ_pd = G4IonTable::GetIonTable()->GetIon ( targZ , targA , 0.0 );
143  G4DynamicParticle* targ_dp = new G4DynamicParticle( targ_pd , G4ThreeVector(1,0,0), 0.0 );
144  unChanged.SetEnergyChange(aTrack.GetKineticEnergy());
145  unChanged.SetMomentumChange(aTrack.Get4Momentum().vect() );
146  unChanged.AddSecondary(targ_dp);
147  //TK121106
150  return &unChanged;
151  }
152  //TK120607
153 
154 
155  G4int lChan=0;
156  random=G4UniformRand();
157  for(i=0; i<nChannels; i++)
158  {
159  lChan = i;
160  if(running[nChannels-1] != 0) if(random<running[i]/running[nChannels-1]) break;
161  }
162  delete [] running;
163  return theChannels[lChan]->ApplyYourself(aTrack, isotope);
164  }
G4double GetWeightedXsec(G4double energy, G4int isoNumber)
CLHEP::Hep3Vector G4ThreeVector
G4double GetZ(G4int i)
static G4NeutronHPManager * GetInstance()
G4NeutronHPReactionWhiteBoard * GetReactionWhiteBoard()
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:449
int G4int
Definition: G4Types.hh:78
G4HadFinalState * ApplyYourself(const G4HadProjectile &theTrack, G4int isoNumber=-1)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:80
void SetEnergyChange(G4double anEnergy)
G4double GetFSCrossSection(G4double energy, G4int isoNumber)
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4double GetTemperature() const
Definition: G4Material.hh:180
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
double G4double
Definition: G4Types.hh:76
void SetMomentumChange(const G4ThreeVector &aV)
void AddSecondary(G4DynamicParticle *aP)
G4int G4NeutronHPChannelList::GetNumberOfChannels ( )
inline

Definition at line 74 of file G4NeutronHPChannelList.hh.

74 { return nChannels; }
G4double G4NeutronHPChannelList::GetXsec ( G4double  anEnergy)
inline

Definition at line 63 of file G4NeutronHPChannelList.hh.

References G4INCL::Math::max().

64  {
65  G4double result=0;
66  G4int i;
67  for(i=0; i<nChannels; i++)
68  {
69  result+=std::max(0., theChannels[i]->GetXsec(anEnergy));
70  }
71  return result;
72  }
int G4int
Definition: G4Types.hh:78
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4double GetXsec(G4double anEnergy)
double G4double
Definition: G4Types.hh:76
G4bool G4NeutronHPChannelList::HasDataInAnyFinalState ( )
inline

Definition at line 76 of file G4NeutronHPChannelList.hh.

77  {
78  G4bool result = false;
79  G4int i;
80  for(i=0; i<nChannels; i++)
81  {
82  if(theChannels[i]->HasDataInAnyFinalState()) result = true;
83  }
84  return result;
85  }
int G4int
Definition: G4Types.hh:78
bool G4bool
Definition: G4Types.hh:79
void G4NeutronHPChannelList::Init ( G4int  n)
void G4NeutronHPChannelList::Init ( G4Element anElement,
const G4String dirName 
)

Definition at line 166 of file G4NeutronHPChannelList.cc.

167  {
168  theDir = dirName;
169 // G4cout << theDir << G4endl;
170  theElement = anElement;
171 // G4cout << theElement << G4endl;
172  ;
173  }
void G4NeutronHPChannelList::Register ( G4NeutronHPFinalState theFS,
const G4String aName 
)

Definition at line 175 of file G4NeutronHPChannelList.cc.

References G4NeutronHPChannel::Init(), and G4NeutronHPChannel::Register().

177  {
178  if(!allChannelsCreated)
179  {
180  if(nChannels!=0)
181  {
182  G4NeutronHPChannel ** theBuffer = new G4NeutronHPChannel * [nChannels+1];
183  G4int i;
184  for(i=0; i<nChannels; i++)
185  {
186  theBuffer[i] = theChannels[i];
187  }
188  delete [] theChannels;
189  theChannels = theBuffer;
190  }
191  else
192  {
193  theChannels = new G4NeutronHPChannel * [nChannels+1];
194  }
195  G4String name;
196  name = aName+"/";
197  theChannels[nChannels] = new G4NeutronHPChannel;
198  theChannels[nChannels]->Init(theElement, theDir, name);
199  nChannels++;
200  }
201 
202  //110527TKDB Unnessary codes, Detected by gcc4.6 compiler
203  //G4bool result;
204  //result = theChannels[theInitCount]->Register(theFS);
205  theChannels[theInitCount]->Register(theFS);
206  theInitCount++;
207  }
G4bool Register(G4NeutronHPFinalState *theFS)
void Init(G4Element *theElement, const G4String dirName)
const XML_Char * name
int G4int
Definition: G4Types.hh:78
void G4NeutronHPChannelList::RestartRegistration ( )
inline

Definition at line 86 of file G4NeutronHPChannelList.hh.

87  {
88  allChannelsCreated = true;
89  theInitCount = 0;
90  }

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