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

#include <G4ParticleHPInelastic.hh>

Inheritance diagram for G4ParticleHPInelastic:
Collaboration diagram for G4ParticleHPInelastic:

Public Member Functions

 G4ParticleHPInelastic (G4ParticleDefinition *projectile=G4Neutron::Neutron(), const char *name="NeutronHPInelastic")
 
 ~G4ParticleHPInelastic ()
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
virtual const std::pair
< G4double, G4double
GetFatalEnergyCheckLevels () const
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void ModelDescription (std::ostream &outFile) const
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual std::pair< G4double,
G4double
GetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void InitialiseModel ()
 

Protected Attributes

std::vector
< G4ParticleHPChannelList * > * 
theInelastic
 
G4String dataDirVariable
 
G4String dirName
 
G4int numEle
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 

Detailed Description

Definition at line 89 of file G4ParticleHPInelastic.hh.

Constructor & Destructor Documentation

G4ParticleHPInelastic::G4ParticleHPInelastic ( G4ParticleDefinition projectile = G4Neutron::Neutron(),
const char *  name = "NeutronHPInelastic" 
)

Definition at line 47 of file G4ParticleHPInelastic.cc.

49  ,theInelastic(NULL)
50  ,numEle(0)
51  ,theProjectile(projectile)
52 {
53  G4String baseName;
54  if ( getenv("G4PARTICLEHPDATA") ) {
55  baseName = getenv( "G4PARTICLEHPDATA" );
56  }
57  //const char* dataDirVariable;
58  G4String particleName;
59  if( theProjectile == G4Neutron::Neutron() ) {
60  dataDirVariable = "G4NEUTRONHPDATA";
61  }else if( theProjectile == G4Proton::Proton() ) {
62  dataDirVariable = "G4PROTONHPDATA";
63  particleName = "Proton";
64  }else if( theProjectile == G4Deuteron::Deuteron() ) {
65  dataDirVariable = "G4DEUTERONHPDATA";
66  particleName = "Deuteron";
67  }else if( theProjectile == G4Triton::Triton() ) {
68  dataDirVariable = "G4TRITONHPDATA";
69  particleName = "Triton";
70  }else if( theProjectile == G4He3::He3() ) {
71  dataDirVariable = "G4HE3HPDATA";
72  particleName = "He3";
73  }else if( theProjectile == G4Alpha::Alpha() ) {
74  dataDirVariable = "G4ALPHAHPDATA";
75  particleName = "Alpha";
76  } else {
77  G4String message("G4ParticleHPInelastic may only be called for neutron, proton, deuteron, triton, He3 or alpha, while it is called for " + theProjectile->GetParticleName());
78  throw G4HadronicException(__FILE__, __LINE__,message.c_str());
79  }
80 
81  SetMinEnergy( 0.0 );
82  SetMaxEnergy( 20.*MeV );
83 
84 // G4cout << " entering G4ParticleHPInelastic constructor"<<G4endl;
85  if ( !getenv("G4PARTICLEHPDATA") && !getenv(dataDirVariable) ) {
86  G4String message( "Please set the environement variable " + G4String(dataDirVariable) + " to point to the " + theProjectile->GetParticleName() + " cross-section files." );
87  throw G4HadronicException(__FILE__, __LINE__,message.c_str());
88  }
89  if ( getenv(dataDirVariable) ) {
90  dirName = getenv(dataDirVariable);
91  } else {
92  dirName = baseName + "/" + particleName;
93  }
94 G4cout << dirName << G4endl;
95 
96  G4String tString = "/Inelastic";
97  dirName = dirName + tString;
98  //numEle = G4Element::GetNumberOfElements();
99 
100  G4cout << "@@@ G4ParticleHPInelastic instantiated for particle " << theProjectile->GetParticleName() << " data directory variable is " << dataDirVariable << " pointing to " << dirName << G4endl;
101 
102 /*
103  theInelastic = new G4ParticleHPChannelList[numEle];
104  for (G4int i=0; i<numEle; i++)
105  {
106  theInelastic[i].Init((*(G4Element::GetElementTable()))[i], dirName);
107  G4int itry = 0;
108  do
109  {
110  theInelastic[i].Register(&theNFS, "F01"); // has
111  theInelastic[i].Register(&theNXFS, "F02");
112  theInelastic[i].Register(&the2NDFS, "F03");
113  theInelastic[i].Register(&the2NFS, "F04"); // has, E Done
114  theInelastic[i].Register(&the3NFS, "F05"); // has, E Done
115  theInelastic[i].Register(&theNAFS, "F06");
116  theInelastic[i].Register(&theN3AFS, "F07");
117  theInelastic[i].Register(&the2NAFS, "F08");
118  theInelastic[i].Register(&the3NAFS, "F09");
119  theInelastic[i].Register(&theNPFS, "F10");
120  theInelastic[i].Register(&theN2AFS, "F11");
121  theInelastic[i].Register(&the2N2AFS, "F12");
122  theInelastic[i].Register(&theNDFS, "F13");
123  theInelastic[i].Register(&theNTFS, "F14");
124  theInelastic[i].Register(&theNHe3FS, "F15");
125  theInelastic[i].Register(&theND2AFS, "F16");
126  theInelastic[i].Register(&theNT2AFS, "F17");
127  theInelastic[i].Register(&the4NFS, "F18"); // has, E Done
128  theInelastic[i].Register(&the2NPFS, "F19");
129  theInelastic[i].Register(&the3NPFS, "F20");
130  theInelastic[i].Register(&theN2PFS, "F21");
131  theInelastic[i].Register(&theNPAFS, "F22");
132  theInelastic[i].Register(&thePFS, "F23");
133  theInelastic[i].Register(&theDFS, "F24");
134  theInelastic[i].Register(&theTFS, "F25");
135  theInelastic[i].Register(&theHe3FS, "F26");
136  theInelastic[i].Register(&theAFS, "F27");
137  theInelastic[i].Register(&the2AFS, "F28");
138  theInelastic[i].Register(&the3AFS, "F29");
139  theInelastic[i].Register(&the2PFS, "F30");
140  theInelastic[i].Register(&thePAFS, "F31");
141  theInelastic[i].Register(&theD2AFS, "F32");
142  theInelastic[i].Register(&theT2AFS, "F33");
143  theInelastic[i].Register(&thePDFS, "F34");
144  theInelastic[i].Register(&thePTFS, "F35");
145  theInelastic[i].Register(&theDAFS, "F36");
146  theInelastic[i].RestartRegistration();
147  itry++;
148  }
149  //while(!theInelastic[i].HasDataInAnyFinalState());
150  while( !theInelastic[i].HasDataInAnyFinalState() && itry < 6 );
151  // 6 is corresponding to the value(5) of G4ParticleHPChannel. TK
152 
153  if ( itry == 6 )
154  {
155  // No Final State at all.
156  G4bool exceptional = false;
157  if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
158  {
159  if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true; //1H
160  }
161  if ( !exceptional ) throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
162  }
163  }
164 */
165 /*
166  for (G4int i=0; i<numEle; i++)
167  {
168  theInelastic.push_back( new G4ParticleHPChannelList );
169  (*theInelastic[i]).Init((*(G4Element::GetElementTable()))[i], dirName, theProjectile);
170  G4int itry = 0;
171  do
172  {
173  (*theInelastic[i]).Register(&theNFS, "F01"); // has
174  (*theInelastic[i]).Register(&theNXFS, "F02");
175  (*theInelastic[i]).Register(&the2NDFS, "F03");
176  (*theInelastic[i]).Register(&the2NFS, "F04"); // has, E Done
177  (*theInelastic[i]).Register(&the3NFS, "F05"); // has, E Done
178  (*theInelastic[i]).Register(&theNAFS, "F06");
179  (*theInelastic[i]).Register(&theN3AFS, "F07");
180  (*theInelastic[i]).Register(&the2NAFS, "F08");
181  (*theInelastic[i]).Register(&the3NAFS, "F09");
182  (*theInelastic[i]).Register(&theNPFS, "F10");
183  (*theInelastic[i]).Register(&theN2AFS, "F11");
184  (*theInelastic[i]).Register(&the2N2AFS, "F12");
185  (*theInelastic[i]).Register(&theNDFS, "F13");
186  (*theInelastic[i]).Register(&theNTFS, "F14");
187  (*theInelastic[i]).Register(&theNHe3FS, "F15");
188  (*theInelastic[i]).Register(&theND2AFS, "F16");
189  (*theInelastic[i]).Register(&theNT2AFS, "F17");
190  (*theInelastic[i]).Register(&the4NFS, "F18"); // has, E Done
191  (*theInelastic[i]).Register(&the2NPFS, "F19");
192  (*theInelastic[i]).Register(&the3NPFS, "F20");
193  (*theInelastic[i]).Register(&theN2PFS, "F21");
194  (*theInelastic[i]).Register(&theNPAFS, "F22");
195  (*theInelastic[i]).Register(&thePFS, "F23");
196  (*theInelastic[i]).Register(&theDFS, "F24");
197  (*theInelastic[i]).Register(&theTFS, "F25");
198  (*theInelastic[i]).Register(&theHe3FS, "F26");
199  (*theInelastic[i]).Register(&theAFS, "F27");
200  (*theInelastic[i]).Register(&the2AFS, "F28");
201  (*theInelastic[i]).Register(&the3AFS, "F29");
202  (*theInelastic[i]).Register(&the2PFS, "F30");
203  (*theInelastic[i]).Register(&thePAFS, "F31");
204  (*theInelastic[i]).Register(&theD2AFS, "F32");
205  (*theInelastic[i]).Register(&theT2AFS, "F33");
206  (*theInelastic[i]).Register(&thePDFS, "F34");
207  (*theInelastic[i]).Register(&thePTFS, "F35");
208  (*theInelastic[i]).Register(&theDAFS, "F36");
209  (*theInelastic[i]).RestartRegistration();
210  itry++;
211  }
212  while( !(*theInelastic[i]).HasDataInAnyFinalState() && itry < 6 );
213  // 6 is corresponding to the value(5) of G4ParticleHPChannel. TK
214 
215  if ( itry == 6 )
216  {
217  // No Final State at all.
218  G4bool exceptional = false;
219  if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
220  {
221  if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true; //1H
222  }
223  if ( !exceptional ) {
224  G4cerr << " ELEMENT Z " << (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() << " N " << (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() << G4endl; //1H
225 throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
226  }
227  }
228 
229  }
230 */
231 
232  }
const XML_Char * name
Definition: expat.h:151
const G4String & GetParticleName() const
void SetMinEnergy(G4double anEnergy)
G4GLOB_DLL std::ostream G4cout
std::vector< G4ParticleHPChannelList * > * theInelastic
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
void SetMaxEnergy(const G4double anEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4He3 * He3()
Definition: G4He3.cc:94

Here is the call graph for this function:

G4ParticleHPInelastic::~G4ParticleHPInelastic ( )

Definition at line 234 of file G4ParticleHPInelastic.cc.

235  {
236 // delete [] theInelastic;
237  if ( theInelastic != NULL ) {
238  for ( std::vector<G4ParticleHPChannelList*>::iterator
239  it = theInelastic->begin() ; it != theInelastic->end() ; it++ ) {
240  delete *it;
241  }
242  theInelastic->clear();
243  }
244  }
std::vector< G4ParticleHPChannelList * > * theInelastic

Member Function Documentation

G4HadFinalState * G4ParticleHPInelastic::ApplyYourself ( const G4HadProjectile aTrack,
G4Nucleus aTargetNucleus 
)
virtual

Implements G4HadronicInteraction.

Definition at line 248 of file G4ParticleHPInelastic.cc.

249  {
251  const G4Material * theMaterial = aTrack.GetMaterial();
252  G4int n = theMaterial->GetNumberOfElements();
253  G4int index = theMaterial->GetElement(0)->GetIndex();
254  G4int it=0;
255  if(n!=1)
256  {
257  G4double* xSec = new G4double[n];
258  G4double sum=0;
259  G4int i;
260  const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
261  G4double rWeight;
262  G4ParticleHPThermalBoost aThermalE;
263  for (i=0; i<n; i++)
264  {
265  index = theMaterial->GetElement(i)->GetIndex();
266  rWeight = NumAtomsPerVolume[i];
267  if ( aTrack.GetDefinition() == G4Neutron::Neutron() ) {
268  xSec[i] = ((*theInelastic)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
269  theMaterial->GetElement(i),
270  theMaterial->GetTemperature()));
271  } else {
272  xSec[i] = ((*theInelastic)[index])->GetXsec(aTrack.GetKineticEnergy());
273  }
274  xSec[i] *= rWeight;
275  sum+=xSec[i];
276 #ifdef G4PHPDEBUG
277  if( getenv("G4ParticleHPDebug") ) G4cout << " G4ParticleHPInelastic XSEC ELEM " << i << " = " << xSec[i] << G4endl;
278 #endif
279 
280  }
281  G4double random = G4UniformRand();
282  G4double running = 0;
283  for (i=0; i<n; i++)
284  {
285  running += xSec[i];
286  index = theMaterial->GetElement(i)->GetIndex();
287  it = i;
288  //if(random<=running/sum) break;
289  if( sum == 0 || random<=running/sum) break;
290  }
291  delete [] xSec;
292  }
293 
294 #ifdef G4PHPDEBUG
295  if( getenv("G4ParticleHPDebug") ) G4cout << " G4ParticleHPInelastic SELECTED ELEM " << it << " = " << theMaterial->GetElement(it)->GetName() << " FROM MATERIAL " << theMaterial->GetName() << G4endl;
296 #endif
297  //return theInelastic[index].ApplyYourself(theMaterial->GetElement(it), aTrack);
298  G4HadFinalState* result = ((*theInelastic)[index])->ApplyYourself(theMaterial->GetElement(it), aTrack);
299  //
300  aNucleus.SetParameters(G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA(),G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargZ());
301  const G4Element* target_element = (*G4Element::GetElementTable())[index];
302  const G4Isotope* target_isotope=NULL;
303  G4int iele = target_element->GetNumberOfIsotopes();
304  for ( G4int j = 0 ; j != iele ; j++ ) {
305  target_isotope=target_element->GetIsotope( j );
306  if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
307  }
308  //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
309  //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
310  //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
311  aNucleus.SetIsotope( target_isotope );
312 
314 
315  //GDEB
316  if( getenv("G4PHPTEST") ) {
317  G4HadSecondary* seco = result->GetSecondary(0);
318  if(seco) {
319  G4ThreeVector secoMom = seco->GetParticle()->GetMomentum();
320  G4cout << " G4ParticleHPinelastic COS THETA " << std::cos(secoMom.theta()) <<" " << secoMom << G4endl;
321  }
322  }
323 
324  return result;
325  }
G4double G4ParticleHPJENDLHEData::G4double result
static G4ParticleHPManager * GetInstance()
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
G4HadSecondary * GetSecondary(size_t i)
const G4String & GetName() const
Definition: G4Material.hh:178
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
int G4int
Definition: G4Types.hh:78
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:206
G4int GetN() const
Definition: G4Isotope.hh:94
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
size_t GetIndex() const
Definition: G4Element.hh:182
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
double theta() const
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)
G4DynamicParticle * GetParticle()
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
G4double GetTemperature() const
Definition: G4Material.hh:182
#define G4endl
Definition: G4ios.hh:61
const G4Material * GetMaterial() const
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
double G4double
Definition: G4Types.hh:76
const G4String & GetName() const
Definition: G4Element.hh:127
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
G4ThreeVector GetMomentum() const

Here is the call graph for this function:

void G4ParticleHPInelastic::BuildPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 451 of file G4ParticleHPInelastic.cc.

451  {
452 
454 
455  theInelastic = hpmanager->GetInelasticFinalStates( &projectile );
456 
457  if ( G4Threading::IsMasterThread() ) {
458 
459  if ( theInelastic == NULL ) theInelastic = new std::vector<G4ParticleHPChannelList*>;
460 
461  if ( numEle == (G4int)G4Element::GetNumberOfElements() ) return;
462 
463  if ( theInelastic->size() == G4Element::GetNumberOfElements() ) {
465  return;
466  }
467 /*
468  const char* dataDirVariable;
469  if( &projectile == G4Neutron::Neutron() ) {
470  dataDirVariable = "G4NEUTRONHPDATA";
471  }else if( &projectile == G4Proton::Proton() ) {
472  dataDirVariable = "G4PROTONHPDATA";
473  }else if( &projectile == G4Deuteron::Deuteron() ) {
474  dataDirVariable = "G4DEUTERONHPDATA";
475  }else if( &projectile == G4Triton::Triton() ) {
476  dataDirVariable = "G4TRITONHPDATA";
477  }else if( &projectile == G4He3::He3() ) {
478  dataDirVariable = "G4HE3HPDATA";
479  }else if( &projectile == G4Alpha::Alpha() ) {
480  dataDirVariable = "G4ALPHAHPDATA";
481  } else {
482  G4String message("G4ParticleHPInelastic may only be called for neutron, proton, deuteron, triton, He3 or alpha, while it is called for " + projectile.GetParticleName());
483  throw G4HadronicException(__FILE__, __LINE__,message.c_str());
484  }
485  if(!getenv(dataDirVariable)){
486  G4String message("Please set the environement variable " + G4String(dataDirVariable) + " to point to the " + projectile.GetParticleName() + " cross-section files.");
487  throw G4HadronicException(__FILE__, __LINE__,message.c_str());
488  }
489  dirName = getenv(dataDirVariable);
490  G4cout << dirName << G4endl;
491 
492  G4String tString = "/Inelastic";
493  dirName = dirName + tString;
494 
495 */
496  G4cout << "@@@ G4ParticleHPInelastic instantiated for particle " << projectile.GetParticleName() << " data directory variable is " << dataDirVariable << " pointing to " << dirName << G4endl;
497  for (G4int i = numEle ; i < (G4int)G4Element::GetNumberOfElements(); i++)
498  {
499  theInelastic->push_back( new G4ParticleHPChannelList );
500  ((*theInelastic)[i])->Init((*(G4Element::GetElementTable()))[i], dirName, const_cast<G4ParticleDefinition*>(&projectile));
501  G4int itry = 0;
502  do
503  {
504  ((*theInelastic)[i])->Register( new G4ParticleHPNInelasticFS , "F01"); // has
505  ((*theInelastic)[i])->Register( new G4ParticleHPNXInelasticFS , "F02");
506  ((*theInelastic)[i])->Register( new G4ParticleHP2NDInelasticFS , "F03");
507  ((*theInelastic)[i])->Register( new G4ParticleHP2NInelasticFS , "F04"); // has, E Done
508  ((*theInelastic)[i])->Register( new G4ParticleHP3NInelasticFS , "F05"); // has, E Done
509  ((*theInelastic)[i])->Register( new G4ParticleHPNAInelasticFS , "F06");
510  ((*theInelastic)[i])->Register( new G4ParticleHPN3AInelasticFS , "F07");
511  ((*theInelastic)[i])->Register( new G4ParticleHP2NAInelasticFS , "F08");
512  ((*theInelastic)[i])->Register( new G4ParticleHP3NAInelasticFS , "F09");
513  ((*theInelastic)[i])->Register( new G4ParticleHPNPInelasticFS , "F10");
514  ((*theInelastic)[i])->Register( new G4ParticleHPN2AInelasticFS , "F11");
515  ((*theInelastic)[i])->Register( new G4ParticleHP2N2AInelasticFS , "F12");
516  ((*theInelastic)[i])->Register( new G4ParticleHPNDInelasticFS , "F13");
517  ((*theInelastic)[i])->Register( new G4ParticleHPNTInelasticFS , "F14");
518  ((*theInelastic)[i])->Register( new G4ParticleHPNHe3InelasticFS , "F15");
519  ((*theInelastic)[i])->Register( new G4ParticleHPND2AInelasticFS , "F16");
520  ((*theInelastic)[i])->Register( new G4ParticleHPNT2AInelasticFS , "F17");
521  ((*theInelastic)[i])->Register( new G4ParticleHP4NInelasticFS , "F18"); // has, E Done
522  ((*theInelastic)[i])->Register( new G4ParticleHP2NPInelasticFS , "F19");
523  ((*theInelastic)[i])->Register( new G4ParticleHP3NPInelasticFS , "F20");
524  ((*theInelastic)[i])->Register( new G4ParticleHPN2PInelasticFS , "F21");
525  ((*theInelastic)[i])->Register( new G4ParticleHPNPAInelasticFS , "F22");
526  ((*theInelastic)[i])->Register( new G4ParticleHPPInelasticFS , "F23");
527  ((*theInelastic)[i])->Register( new G4ParticleHPDInelasticFS , "F24");
528  ((*theInelastic)[i])->Register( new G4ParticleHPTInelasticFS , "F25");
529  ((*theInelastic)[i])->Register( new G4ParticleHPHe3InelasticFS , "F26");
530  ((*theInelastic)[i])->Register( new G4ParticleHPAInelasticFS , "F27");
531  ((*theInelastic)[i])->Register( new G4ParticleHP2AInelasticFS , "F28");
532  ((*theInelastic)[i])->Register( new G4ParticleHP3AInelasticFS , "F29");
533  ((*theInelastic)[i])->Register( new G4ParticleHP2PInelasticFS , "F30");
534  ((*theInelastic)[i])->Register( new G4ParticleHPPAInelasticFS , "F31");
535  ((*theInelastic)[i])->Register( new G4ParticleHPD2AInelasticFS , "F32");
536  ((*theInelastic)[i])->Register( new G4ParticleHPT2AInelasticFS , "F33");
537  ((*theInelastic)[i])->Register( new G4ParticleHPPDInelasticFS , "F34");
538  ((*theInelastic)[i])->Register( new G4ParticleHPPTInelasticFS , "F35");
539  ((*theInelastic)[i])->Register( new G4ParticleHPDAInelasticFS , "F36");
540  ((*theInelastic)[i])->RestartRegistration();
541  itry++;
542  }
543  while( !((*theInelastic)[i])->HasDataInAnyFinalState() && itry < 6 ); // Loop checking, 11.05.2015, T. Koi
544  // 6 is corresponding to the value(5) of G4ParticleHPChannel. TK
545 
546  if ( itry == 6 ) {
547  // No Final State at all.
548  /*
549  G4bool exceptional = false;
550  if ( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() == 1 )
551  {
552  if ( (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() == 1 && (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() == 1 ) exceptional = true; //1H
553  }
554  if ( !exceptional ) {
555  G4cerr << " ELEMENT Z " << (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetZ() << " N " << (*(G4Element::GetElementTable()))[i]->GetIsotope( 0 )->GetN() << G4endl; //1H
556  throw G4HadronicException(__FILE__, __LINE__, "Channel: Do not know what to do with this element");
557  }
558  */
560  G4cout << "ParticleHP::Inelastic for " << projectile.GetParticleName() << ". Do not know what to do with element of \"" << (*(G4Element::GetElementTable()))[i]->GetName() << "\"." << G4endl;
561  G4cout << "The components of the element are" << G4endl;
562  //G4cout << "TKDB dataDirVariable = " << dataDirVariable << G4endl;
563  for ( G4int ii = 0 ; ii < (G4int)( (*(G4Element::GetElementTable()))[i]->GetNumberOfIsotopes() ) ; ii++ ) {
564  G4cout << " Z: " << (*(G4Element::GetElementTable()))[i]->GetIsotope( ii )->GetZ() << ", A: " << (*(G4Element::GetElementTable()))[i]->GetIsotope( ii )->GetN() << G4endl;
565  }
566  G4cout << "No possible final state data of the element is found in " << dataDirVariable << "." << G4endl;
567  }
568  }
569  }
570  hpmanager->RegisterInelasticFinalStates( &projectile , theInelastic );
571  }
573 }
static G4ParticleHPManager * GetInstance()
void Init()
Definition: G4IonTable.cc:90
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
void Register(T *inst)
Definition: G4AutoDelete.hh:65
G4GLOB_DLL std::ostream G4cout
static size_t GetNumberOfElements()
Definition: G4Element.cc:405
std::vector< G4ParticleHPChannelList * > * theInelastic
std::vector< G4ParticleHPChannelList * > * GetInelasticFinalStates(const G4ParticleDefinition *)
#define G4endl
Definition: G4ios.hh:61
G4bool IsMasterThread()
Definition: G4Threading.cc:146
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:398
void RegisterInelasticFinalStates(const G4ParticleDefinition *, std::vector< G4ParticleHPChannelList * > *)

Here is the call graph for this function:

const std::pair< G4double, G4double > G4ParticleHPInelastic::GetFatalEnergyCheckLevels ( ) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 327 of file G4ParticleHPInelastic.cc.

328 {
329  // max energy non-conservation is mass of heavy nucleus
330 // if ( getenv("G4PHP_DO_NOT_ADJUST_FINAL_STATE") ) return std::pair<G4double, G4double>(5*perCent,250*GeV);
331  // This should be same to the hadron default value
332 // return std::pair<G4double, G4double>(10*perCent,10*GeV);
333  return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
334 }
static constexpr double perCent
Definition: G4SIunits.hh:332
#define DBL_MAX
Definition: templates.hh:83
G4int G4ParticleHPInelastic::GetVerboseLevel ( ) const

Definition at line 405 of file G4ParticleHPInelastic.cc.

406 {
408 }
static G4ParticleHPManager * GetInstance()

Here is the call graph for this function:

void G4ParticleHPInelastic::ModelDescription ( std::ostream &  outFile) const
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 575 of file G4ParticleHPInelastic.cc.

576 {
577  outFile << "Extension of High Precision model for inelastic reaction of proton, deuteron, triton, He3 and alpha below 20MeV\n";
578 }
void G4ParticleHPInelastic::SetVerboseLevel ( G4int  newValue)

Definition at line 409 of file G4ParticleHPInelastic.cc.

410 {
412 }
static G4ParticleHPManager * GetInstance()

Here is the call graph for this function:

Member Data Documentation

G4String G4ParticleHPInelastic::dataDirVariable
protected

Definition at line 110 of file G4ParticleHPInelastic.hh.

G4String G4ParticleHPInelastic::dirName
protected

Definition at line 111 of file G4ParticleHPInelastic.hh.

G4int G4ParticleHPInelastic::numEle
protected

Definition at line 112 of file G4ParticleHPInelastic.hh.

std::vector<G4ParticleHPChannelList*>* G4ParticleHPInelastic::theInelastic
protected

Definition at line 109 of file G4ParticleHPInelastic.hh.


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