Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes
G4ParticleHPInelastic Class Reference

#include <G4ParticleHPInelastic.hh>

Inheritance diagram for G4ParticleHPInelastic:
G4HadronicInteraction

Public Member Functions

void ActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Material *aMaterial)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void DeActivateFor (const G4Element *anElement)
 
void DeActivateFor (const G4Material *aMaterial)
 
 G4ParticleHPInelastic (G4ParticleDefinition *projectile=G4Neutron::Neutron(), const char *name="NeutronHPInelastic")
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
const G4StringGetModelName () const
 
G4double GetRecoilEnergyThreshold () const
 
G4int GetVerboseLevel () const
 
virtual void InitialiseModel ()
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4bool IsBlocked (const G4Element *anElement) const
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
virtual void ModelDescription (std::ostream &outFile) const
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
void SetRecoilEnergyThreshold (G4double val)
 
void SetVerboseLevel (G4int)
 
 ~G4ParticleHPInelastic ()
 

Protected Member Functions

void Block ()
 
G4bool IsBlocked () const
 
void SetModelName (const G4String &nam)
 

Protected Attributes

G4String dataDirVariable
 
G4String dirName
 
G4bool isBlocked
 
G4int numEle
 
std::vector< G4ParticleHPChannelList * > * theInelastic
 
G4double theMaxEnergy
 
G4double theMinEnergy
 
G4HadFinalState theParticleChange
 
G4int verboseLevel
 

Private Attributes

std::pair< G4double, G4doubleepCheckLevels
 
G4double recoilEnergyThreshold
 
G4HadronicInteractionRegistryregistry
 
std::vector< const G4Material * > theBlockedList
 
std::vector< const G4Element * > theBlockedListElements
 
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements
 
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList
 
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
 
G4String theModelName
 
G4ParticleDefinitiontheProjectile
 

Detailed Description

Definition at line 89 of file G4ParticleHPInelastic.hh.

Constructor & Destructor Documentation

◆ G4ParticleHPInelastic()

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

References G4Alpha::Alpha(), dataDirVariable, G4Deuteron::Deuteron(), dirName, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), GetVerboseLevel(), G4He3::He3(), G4HadronicParameters::Instance(), MeV, G4Neutron::Neutron(), G4Proton::Proton(), G4HadronicInteraction::SetMaxEnergy(), G4HadronicInteraction::SetMinEnergy(), theProjectile, and G4Triton::Triton().

◆ ~G4ParticleHPInelastic()

G4ParticleHPInelastic::~G4ParticleHPInelastic ( )

Definition at line 241 of file G4ParticleHPInelastic.cc.

242 {
243// delete [] theInelastic;
244 //Vector is shared, only master deletes
246 if ( theInelastic != NULL ) {
247 for ( std::vector<G4ParticleHPChannelList*>::iterator
248 it = theInelastic->begin() ; it != theInelastic->end() ; it++ ) {
249 delete *it;
250 }
251 theInelastic->clear();
252 }
253 }
254 }
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

References G4Threading::IsWorkerThread(), and theInelastic.

Member Function Documentation

◆ ActivateFor() [1/2]

void G4HadronicInteraction::ActivateFor ( const G4Element anElement)
inlineinherited

◆ ActivateFor() [2/2]

void G4HadronicInteraction::ActivateFor ( const G4Material aMaterial)
inlineinherited

◆ ApplyYourself()

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

Reimplemented from G4HadronicInteraction.

Definition at line 258 of file G4ParticleHPInelastic.cc.

259 {
261 const G4Material * theMaterial = aTrack.GetMaterial();
262 G4int n = theMaterial->GetNumberOfElements();
263 G4int index = theMaterial->GetElement(0)->GetIndex();
264 G4int it=0;
265 if(n!=1)
266 {
267 G4double* xSec = new G4double[n];
268 G4double sum=0;
269 G4int i;
270 const G4double * NumAtomsPerVolume = theMaterial->GetVecNbOfAtomsPerVolume();
271 G4double rWeight;
272 G4ParticleHPThermalBoost aThermalE;
273 for (i=0; i<n; i++)
274 {
275 index = theMaterial->GetElement(i)->GetIndex();
276 rWeight = NumAtomsPerVolume[i];
277 if ( aTrack.GetDefinition() == G4Neutron::Neutron() ) {
278 xSec[i] = ((*theInelastic)[index])->GetXsec(aThermalE.GetThermalEnergy(aTrack,
279 theMaterial->GetElement(i),
280 theMaterial->GetTemperature()));
281 } else {
282 xSec[i] = ((*theInelastic)[index])->GetXsec(aTrack.GetKineticEnergy());
283 }
284 xSec[i] *= rWeight;
285 sum+=xSec[i];
286#ifdef G4PHPDEBUG
287 #ifdef G4VERBOSE
288 if( std::getenv("G4ParticleHPDebug") && G4HadronicParameters::Instance()->GetVerboseLevel() > 0 )
289 G4cout << " G4ParticleHPInelastic XSEC ELEM " << i << " = " << xSec[i] << G4endl;
290 #endif
291#endif
292
293 }
294 G4double random = G4UniformRand();
295 G4double running = 0;
296 for (i=0; i<n; i++)
297 {
298 running += xSec[i];
299 index = theMaterial->GetElement(i)->GetIndex();
300 it = i;
301 //if(random<=running/sum) break;
302 if( sum == 0 || random<=running/sum) break;
303 }
304 delete [] xSec;
305 }
306
307#ifdef G4PHPDEBUG
308 #ifdef G4VERBOSE
309 if( std::getenv("G4ParticleHPDebug") && G4HadronicParameters::Instance()->GetVerboseLevel() > 0 )
310 G4cout << " G4ParticleHPInelastic SELECTED ELEM " << it << " = " << theMaterial->GetElement(it)->GetName() << " FROM MATERIAL " << theMaterial->GetName() << G4endl;
311 #endif
312#endif
313 //return theInelastic[index].ApplyYourself(theMaterial->GetElement(it), aTrack);
314 G4HadFinalState* result = ((*theInelastic)[index])->ApplyYourself(theMaterial->GetElement(it), aTrack);
315 //
317 const G4Element* target_element = (*G4Element::GetElementTable())[index];
318 const G4Isotope* target_isotope=NULL;
319 G4int iele = target_element->GetNumberOfIsotopes();
320 for ( G4int j = 0 ; j != iele ; j++ ) {
321 target_isotope=target_element->GetIsotope( j );
322 if ( target_isotope->GetN() == G4ParticleHPManager::GetInstance()->GetReactionWhiteBoard()->GetTargA() ) break;
323 }
324 //G4cout << "Target Material of this reaction is " << theMaterial->GetName() << G4endl;
325 //G4cout << "Target Element of this reaction is " << target_element->GetName() << G4endl;
326 //G4cout << "Target Isotope of this reaction is " << target_isotope->GetName() << G4endl;
327 aNucleus.SetIsotope( target_isotope );
328
330
331 //GDEB
332 if( std::getenv("G4PHPTEST") ) {
333 G4HadSecondary* seco = result->GetSecondary(0);
334 if(seco) {
335 G4ThreeVector secoMom = seco->GetParticle()->GetMomentum();
336 #ifdef G4VERBOSE
338 G4cout << " G4ParticleHPinelastic COS THETA " << std::cos(secoMom.theta()) <<" " << secoMom << G4endl;
339 #endif
340 }
341 }
342
343 return result;
344 }
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
#define G4UniformRand()
Definition: Randomize.hh:52
double theta() const
G4ThreeVector GetMomentum() const
static G4ElementTable * GetElementTable()
Definition: G4Element.cc:397
const G4Isotope * GetIsotope(G4int iso) const
Definition: G4Element.hh:170
size_t GetIndex() const
Definition: G4Element.hh:182
size_t GetNumberOfIsotopes() const
Definition: G4Element.hh:159
const G4String & GetName() const
Definition: G4Element.hh:127
G4HadSecondary * GetSecondary(size_t i)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
G4DynamicParticle * GetParticle()
G4int GetN() const
Definition: G4Isotope.hh:93
G4double GetTemperature() const
Definition: G4Material.hh:178
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:198
size_t GetNumberOfElements() const
Definition: G4Material.hh:182
const G4double * GetVecNbOfAtomsPerVolume() const
Definition: G4Material.hh:202
const G4String & GetName() const
Definition: G4Material.hh:173
G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &aTargetNucleus)
static G4ParticleHPManager * GetInstance()
G4ParticleHPReactionWhiteBoard * GetReactionWhiteBoard()
G4double GetThermalEnergy(const G4HadProjectile &aP, const G4Element *anE, G4double aT)

References ApplyYourself(), G4ParticleHPManager::CloseReactionWhiteBoard(), G4cout, G4endl, G4UniformRand, G4HadProjectile::GetDefinition(), G4Material::GetElement(), G4Element::GetElementTable(), G4Element::GetIndex(), G4ParticleHPManager::GetInstance(), G4Element::GetIsotope(), G4HadProjectile::GetKineticEnergy(), G4HadProjectile::GetMaterial(), G4DynamicParticle::GetMomentum(), G4Isotope::GetN(), G4Element::GetName(), G4Material::GetName(), G4Material::GetNumberOfElements(), G4Element::GetNumberOfIsotopes(), G4HadSecondary::GetParticle(), G4ParticleHPManager::GetReactionWhiteBoard(), G4HadFinalState::GetSecondary(), G4ParticleHPReactionWhiteBoard::GetTargA(), G4ParticleHPReactionWhiteBoard::GetTargZ(), G4Material::GetTemperature(), G4ParticleHPThermalBoost::GetThermalEnergy(), G4Material::GetVecNbOfAtomsPerVolume(), GetVerboseLevel(), G4HadronicParameters::GetVerboseLevel(), G4HadronicParameters::Instance(), CLHEP::detail::n, G4Neutron::Neutron(), G4ParticleHPManager::OpenReactionWhiteBoard(), G4Nucleus::SetIsotope(), G4Nucleus::SetParameters(), and CLHEP::Hep3Vector::theta().

Referenced by ApplyYourself().

◆ Block()

void G4HadronicInteraction::Block ( )
inlineprotectedinherited

◆ BuildPhysicsTable()

void G4ParticleHPInelastic::BuildPhysicsTable ( const G4ParticleDefinition projectile)
virtual

Reimplemented from G4HadronicInteraction.

Definition at line 468 of file G4ParticleHPInelastic.cc.

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

References dataDirVariable, dirName, G4ParticleHPManager::DumpSetting(), G4cout, G4endl, G4Element::GetElementTable(), G4ParticleHPManager::GetInelasticFinalStates(), G4ParticleHPManager::GetInstance(), G4Element::GetNumberOfElements(), G4ParticleDefinition::GetParticleName(), G4ParticleHPManager::GetVerboseLevel(), GetVerboseLevel(), G4HadronicParameters::GetVerboseLevel(), lightions::Init(), G4HadronicParameters::Instance(), G4Threading::IsMasterThread(), numEle, G4AutoDelete::Register(), G4ParticleHPManager::RegisterInelasticFinalStates(), and theInelastic.

◆ DeActivateFor() [1/2]

void G4HadronicInteraction::DeActivateFor ( const G4Element anElement)
inherited

Definition at line 186 of file G4HadronicInteraction.cc.

187{
188 Block();
189 theBlockedListElements.push_back(anElement);
190}
std::vector< const G4Element * > theBlockedListElements

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedListElements.

◆ DeActivateFor() [2/2]

void G4HadronicInteraction::DeActivateFor ( const G4Material aMaterial)
inherited

Definition at line 180 of file G4HadronicInteraction.cc.

181{
182 Block();
183 theBlockedList.push_back(aMaterial);
184}
std::vector< const G4Material * > theBlockedList

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theBlockedList.

Referenced by G4HadronHElasticPhysics::ConstructProcess().

◆ GetEnergyMomentumCheckLevels()

std::pair< G4double, G4double > G4HadronicInteraction::GetEnergyMomentumCheckLevels ( ) const
virtualinherited

◆ GetFatalEnergyCheckLevels()

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

Reimplemented from G4HadronicInteraction.

Definition at line 346 of file G4ParticleHPInelastic.cc.

347{
348 // max energy non-conservation is mass of heavy nucleus
349 // This should be same to the hadron default value
350 return std::pair<G4double, G4double>(10*perCent,DBL_MAX);
351}
static constexpr double perCent
Definition: G4SIunits.hh:325
#define DBL_MAX
Definition: templates.hh:62

References DBL_MAX, and perCent.

◆ GetMaxEnergy() [1/2]

G4double G4HadronicInteraction::GetMaxEnergy ( ) const
inlineinherited

◆ GetMaxEnergy() [2/2]

G4double G4HadronicInteraction::GetMaxEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 131 of file G4HadronicInteraction.cc.

133{
134 if(!IsBlocked()) { return theMaxEnergy; }
135 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return 0.0; }
136 if(!theMaxEnergyListElements.empty()) {
137 for(auto const& elmlist : theMaxEnergyListElements) {
138 if( anElement == elmlist.second )
139 { return elmlist.first; }
140 }
141 }
142 if(!theMaxEnergyList.empty()) {
143 for(auto const& matlist : theMaxEnergyList) {
144 if( aMaterial == matlist.second )
145 { return matlist.first; }
146 }
147 }
148 return theMaxEnergy;
149}
std::vector< std::pair< G4double, const G4Material * > > theMaxEnergyList
std::vector< std::pair< G4double, const G4Element * > > theMaxEnergyListElements

References G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMaxEnergy, G4HadronicInteraction::theMaxEnergyList, and G4HadronicInteraction::theMaxEnergyListElements.

◆ GetMinEnergy() [1/2]

G4double G4HadronicInteraction::GetMinEnergy ( ) const
inlineinherited

◆ GetMinEnergy() [2/2]

G4double G4HadronicInteraction::GetMinEnergy ( const G4Material aMaterial,
const G4Element anElement 
) const
inherited

Definition at line 81 of file G4HadronicInteraction.cc.

83{
84 if(!IsBlocked()) { return theMinEnergy; }
85 if( IsBlocked(aMaterial) || IsBlocked(anElement) ) { return DBL_MAX; }
86 if(!theMinEnergyListElements.empty()) {
87 for(auto const& elmlist : theMinEnergyListElements) {
88 if( anElement == elmlist.second )
89 { return elmlist.first; }
90 }
91 }
92 if(!theMinEnergyList.empty()) {
93 for(auto const & matlist : theMinEnergyList) {
94 if( aMaterial == matlist.second )
95 { return matlist.first; }
96 }
97 }
98 return theMinEnergy;
99}
std::vector< std::pair< G4double, const G4Element * > > theMinEnergyListElements
std::vector< std::pair< G4double, const G4Material * > > theMinEnergyList

References DBL_MAX, G4HadronicInteraction::IsBlocked(), G4HadronicInteraction::theMinEnergy, G4HadronicInteraction::theMinEnergyList, and G4HadronicInteraction::theMinEnergyListElements.

◆ GetModelName()

const G4String & G4HadronicInteraction::GetModelName ( ) const
inlineinherited

Definition at line 115 of file G4HadronicInteraction.hh.

116 { return theModelName; }

References G4HadronicInteraction::theModelName.

Referenced by G4MuMinusCapturePrecompound::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4INCLXXInterface::ApplyYourself(), G4TheoFSGenerator::ApplyYourself(), G4HadronStoppingProcess::AtRestDoIt(), G4VHadronPhysics::BuildModel(), G4HadronicProcess::CheckEnergyMomentumConservation(), G4HadronicProcess::CheckResult(), G4ChargeExchangePhysics::ConstructProcess(), G4MuonicAtomDecay::DecayIt(), G4LENDModel::DumpLENDTargetInfo(), G4AblaInterface::G4AblaInterface(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4ExcitedStringDecay::G4ExcitedStringDecay(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LENDorBERTModel::G4LENDorBERTModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4LowEIonFragmentation::G4LowEIonFragmentation(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4INCLXXInterface::GetDeExcitationModelName(), G4EnergyRangeManager::GetHadronicInteraction(), G4VHighEnergyGenerator::GetProjectileNucleus(), G4NeutronRadCapture::InitialiseModel(), G4BinaryCascade::ModelDescription(), G4LMsdGenerator::ModelDescription(), G4VPartonStringModel::ModelDescription(), G4TheoFSGenerator::ModelDescription(), G4VHadronPhysics::NewModel(), G4NeutrinoElectronProcess::PostStepDoIt(), G4HadronicProcess::PostStepDoIt(), G4ElNeutrinoNucleusProcess::PostStepDoIt(), G4HadronElasticProcess::PostStepDoIt(), G4MuNeutrinoNucleusProcess::PostStepDoIt(), G4HadronicProcessStore::PrintModelHtml(), G4BinaryCascade::PropagateModelDescription(), G4HadronicProcessStore::RegisterInteraction(), and G4LENDModel::returnUnchanged().

◆ GetRecoilEnergyThreshold()

G4double G4HadronicInteraction::GetRecoilEnergyThreshold ( ) const
inlineinherited

◆ GetVerboseLevel()

G4int G4ParticleHPInelastic::GetVerboseLevel ( ) const

◆ InitialiseModel()

void G4HadronicInteraction::InitialiseModel ( )
virtualinherited

◆ IsApplicable()

G4bool G4HadronicInteraction::IsApplicable ( const G4HadProjectile aTrack,
G4Nucleus targetNucleus 
)
virtualinherited

◆ IsBlocked() [1/3]

G4bool G4HadronicInteraction::IsBlocked ( ) const
inlineprotectedinherited

◆ IsBlocked() [2/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Element anElement) const
inherited

Definition at line 202 of file G4HadronicInteraction.cc.

203{
204 for (auto const& elm : theBlockedListElements) {
205 if (anElement == elm) return true;
206 }
207 return false;
208}

References G4HadronicInteraction::theBlockedListElements.

◆ IsBlocked() [3/3]

G4bool G4HadronicInteraction::IsBlocked ( const G4Material aMaterial) const
inherited

Definition at line 193 of file G4HadronicInteraction.cc.

194{
195 for (auto const& mat : theBlockedList) {
196 if (aMaterial == mat) return true;
197 }
198 return false;
199}

References G4HadronicInteraction::theBlockedList.

◆ ModelDescription()

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

Reimplemented from G4HadronicInteraction.

Definition at line 600 of file G4ParticleHPInelastic.cc.

601{
602 outFile << "Extension of High Precision model for inelastic reaction of proton, deuteron, triton, He3 and alpha below 20MeV\n";
603}

◆ operator!=()

G4bool G4HadronicInteraction::operator!= ( const G4HadronicInteraction right) const
deleteinherited

◆ operator==()

G4bool G4HadronicInteraction::operator== ( const G4HadronicInteraction right) const
deleteinherited

◆ SampleInvariantT()

G4double G4HadronicInteraction::SampleInvariantT ( const G4ParticleDefinition p,
G4double  plab,
G4int  Z,
G4int  A 
)
virtualinherited

◆ SetEnergyMomentumCheckLevels()

void G4HadronicInteraction::SetEnergyMomentumCheckLevels ( G4double  relativeLevel,
G4double  absoluteLevel 
)
inlineinherited

Definition at line 149 of file G4HadronicInteraction.hh.

150 { epCheckLevels.first = relativeLevel;
151 epCheckLevels.second = absoluteLevel; }

References G4HadronicInteraction::epCheckLevels.

Referenced by G4BinaryCascade::G4BinaryCascade(), G4CascadeInterface::G4CascadeInterface(), and G4FTFModel::G4FTFModel().

◆ SetMaxEnergy() [1/3]

void G4HadronicInteraction::SetMaxEnergy ( const G4double  anEnergy)
inlineinherited

Definition at line 102 of file G4HadronicInteraction.hh.

103 { theMaxEnergy = anEnergy; }

References G4HadronicInteraction::theMaxEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4IonINCLXXPhysics::AddProcess(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4HadronicBuilder::BuildFTFP_BERT(), G4HadronicBuilder::BuildFTFQGSP_BERT(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4HadronicBuilder::BuildQGSP_FTFP_BERT(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMaxEnergy() [2/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 151 of file G4HadronicInteraction.cc.

153{
154 Block();
155 if(!theMaxEnergyListElements.empty()) {
156 for(auto & elmlist : theMaxEnergyListElements) {
157 if( anElement == elmlist.second ) {
158 elmlist.first = anEnergy;
159 return;
160 }
161 }
162 }
163 theMaxEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
164}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyListElements.

◆ SetMaxEnergy() [3/3]

void G4HadronicInteraction::SetMaxEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 166 of file G4HadronicInteraction.cc.

167{
168 Block();
169 if(!theMaxEnergyList.empty()) {
170 for(auto & matlist: theMaxEnergyList) {
171 if( aMaterial == matlist.second ) {
172 matlist.first = anEnergy;
173 return;
174 }
175 }
176 }
177 theMaxEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
178}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMaxEnergyList.

◆ SetMinEnergy() [1/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy)
inlineinherited

Definition at line 89 of file G4HadronicInteraction.hh.

90 { theMinEnergy = anEnergy; }

References G4HadronicInteraction::theMinEnergy.

Referenced by G4HadronicInteraction::ActivateFor(), G4BertiniElectroNuclearBuilder::Build(), G4LENDBertiniGammaElectroNuclearBuilder::Build(), G4NeutronLENDBuilder::Build(), G4NeutronPHPBuilder::Build(), G4AlphaPHPBuilder::Build(), G4BertiniKaonBuilder::Build(), G4BertiniNeutronBuilder::Build(), G4BertiniPiKBuilder::Build(), G4BertiniPionBuilder::Build(), G4BertiniProtonBuilder::Build(), G4BinaryAlphaBuilder::Build(), G4BinaryDeuteronBuilder::Build(), G4BinaryHe3Builder::Build(), G4BinaryNeutronBuilder::Build(), G4BinaryPiKBuilder::Build(), G4BinaryPionBuilder::Build(), G4BinaryProtonBuilder::Build(), G4BinaryTritonBuilder::Build(), G4DeuteronPHPBuilder::Build(), G4FTFBinaryKaonBuilder::Build(), G4FTFBinaryNeutronBuilder::Build(), G4FTFBinaryPiKBuilder::Build(), G4FTFBinaryPionBuilder::Build(), G4FTFBinaryProtonBuilder::Build(), G4FTFPAntiBarionBuilder::Build(), G4FTFPKaonBuilder::Build(), G4FTFPNeutronBuilder::Build(), G4FTFPPiKBuilder::Build(), G4FTFPPionBuilder::Build(), G4FTFPProtonBuilder::Build(), G4He3PHPBuilder::Build(), G4HyperonFTFPBuilder::Build(), G4HyperonQGSPBuilder::Build(), G4INCLXXNeutronBuilder::Build(), G4INCLXXPionBuilder::Build(), G4INCLXXProtonBuilder::Build(), G4PrecoNeutronBuilder::Build(), G4PrecoProtonBuilder::Build(), G4ProtonPHPBuilder::Build(), G4QGSBinaryKaonBuilder::Build(), G4QGSBinaryNeutronBuilder::Build(), G4QGSBinaryPiKBuilder::Build(), G4QGSBinaryPionBuilder::Build(), G4QGSBinaryProtonBuilder::Build(), G4QGSPAntiBarionBuilder::Build(), G4QGSPKaonBuilder::Build(), G4QGSPLundStrFragmProtonBuilder::Build(), G4QGSPNeutronBuilder::Build(), G4QGSPPiKBuilder::Build(), G4QGSPPionBuilder::Build(), G4TritonPHPBuilder::Build(), G4QGSPProtonBuilder::Build(), G4QGSBuilder::BuildModel(), G4VHadronPhysics::BuildModel(), G4EmExtraPhysics::ConstructGammaElectroNuclear(), LBE::ConstructHad(), G4EmExtraPhysics::ConstructLENDGammaNuclear(), G4HadronElasticPhysicsHP::ConstructProcess(), G4HadronElasticPhysicsLEND::ConstructProcess(), G4HadronElasticPhysicsPHP::ConstructProcess(), G4HadronDElasticPhysics::ConstructProcess(), G4HadronElasticPhysics::ConstructProcess(), G4HadronHElasticPhysics::ConstructProcess(), G4IonElasticPhysics::ConstructProcess(), G4IonINCLXXPhysics::ConstructProcess(), G4IonPhysics::ConstructProcess(), G4IonPhysicsPHP::ConstructProcess(), G4IonQMDPhysics::ConstructProcess(), G4ANuElNucleusNcModel::G4ANuElNucleusNcModel(), G4ANuMuNucleusNcModel::G4ANuMuNucleusNcModel(), G4BertiniKaonBuilder::G4BertiniKaonBuilder(), G4BertiniPiKBuilder::G4BertiniPiKBuilder(), G4BertiniPionBuilder::G4BertiniPionBuilder(), G4BinaryCascade::G4BinaryCascade(), G4BinaryPiKBuilder::G4BinaryPiKBuilder(), G4BinaryPionBuilder::G4BinaryPionBuilder(), G4ChargeExchange::G4ChargeExchange(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElectroVDNuclearModel::G4ElectroVDNuclearModel(), G4EMDissociation::G4EMDissociation(), G4FissLib::G4FissLib(), G4FTFBinaryKaonBuilder::G4FTFBinaryKaonBuilder(), G4FTFBinaryNeutronBuilder::G4FTFBinaryNeutronBuilder(), G4FTFBinaryPiKBuilder::G4FTFBinaryPiKBuilder(), G4FTFBinaryPionBuilder::G4FTFBinaryPionBuilder(), G4FTFBinaryProtonBuilder::G4FTFBinaryProtonBuilder(), G4FTFPAntiBarionBuilder::G4FTFPAntiBarionBuilder(), G4FTFPKaonBuilder::G4FTFPKaonBuilder(), G4FTFPNeutronBuilder::G4FTFPNeutronBuilder(), G4FTFPPiKBuilder::G4FTFPPiKBuilder(), G4FTFPPionBuilder::G4FTFPPionBuilder(), G4FTFPProtonBuilder::G4FTFPProtonBuilder(), G4HadronElastic::G4HadronElastic(), G4HadronicAbsorptionBertini::G4HadronicAbsorptionBertini(), G4HadronicAbsorptionFritiof::G4HadronicAbsorptionFritiof(), G4HadronicAbsorptionFritiofWithBinaryCascade::G4HadronicAbsorptionFritiofWithBinaryCascade(), G4hhElastic::G4hhElastic(), G4HyperonFTFPBuilder::G4HyperonFTFPBuilder(), G4HyperonQGSPBuilder::G4HyperonQGSPBuilder(), G4INCLXXPionBuilder::G4INCLXXPionBuilder(), G4LEHadronProtonElastic::G4LEHadronProtonElastic(), G4LENDModel::G4LENDModel(), G4LEnp::G4LEnp(), G4LEpp::G4LEpp(), G4LFission::G4LFission(), G4LowEGammaNuclearModel::G4LowEGammaNuclearModel(), G4MuonVDNuclearModel::G4MuonVDNuclearModel(), G4NeutrinoElectronCcModel::G4NeutrinoElectronCcModel(), G4NeutrinoElectronNcModel::G4NeutrinoElectronNcModel(), G4NeutrinoNucleusModel::G4NeutrinoNucleusModel(), G4NeutronElectronElModel::G4NeutronElectronElModel(), G4NeutronRadCapture::G4NeutronRadCapture(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4NuElNucleusNcModel::G4NuElNucleusNcModel(), G4NuMuNucleusNcModel::G4NuMuNucleusNcModel(), G4ParticleHPCapture::G4ParticleHPCapture(), G4ParticleHPElastic::G4ParticleHPElastic(), G4ParticleHPFission::G4ParticleHPFission(), G4ParticleHPInelastic(), G4ParticleHPThermalScattering::G4ParticleHPThermalScattering(), G4QGSPAntiBarionBuilder::G4QGSPAntiBarionBuilder(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4NeutrinoElectronCcModel::IsApplicable(), G4HadronPhysicsFTFP_BERT_HP::Neutron(), G4HadronPhysicsINCLXX::Neutron(), G4HadronPhysicsQGSP_BERT_HP::Neutron(), G4HadronPhysicsQGSP_BIC_HP::Neutron(), G4HadronPhysicsShielding::Neutron(), and G4VHadronPhysics::NewModel().

◆ SetMinEnergy() [2/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Element anElement 
)
inherited

Definition at line 101 of file G4HadronicInteraction.cc.

103{
104 Block();
105 if(!theMinEnergyListElements.empty()) {
106 for(auto & elmlist : theMinEnergyListElements) {
107 if( anElement == elmlist.second ) {
108 elmlist.first = anEnergy;
109 return;
110 }
111 }
112 }
113 theMinEnergyListElements.push_back(std::pair<G4double, const G4Element *>(anEnergy, anElement));
114}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyListElements.

◆ SetMinEnergy() [3/3]

void G4HadronicInteraction::SetMinEnergy ( G4double  anEnergy,
const G4Material aMaterial 
)
inherited

Definition at line 116 of file G4HadronicInteraction.cc.

118{
119 Block();
120 if(!theMinEnergyList.empty()) {
121 for(auto & matlist : theMinEnergyList) {
122 if( aMaterial == matlist.second ) {
123 matlist.first = anEnergy;
124 return;
125 }
126 }
127 }
128 theMinEnergyList.push_back(std::pair<G4double, const G4Material *>(anEnergy, aMaterial));
129}

References G4HadronicInteraction::Block(), and G4HadronicInteraction::theMinEnergyList.

◆ SetModelName()

void G4HadronicInteraction::SetModelName ( const G4String nam)
inlineprotectedinherited

◆ SetRecoilEnergyThreshold()

void G4HadronicInteraction::SetRecoilEnergyThreshold ( G4double  val)
inlineinherited

◆ SetVerboseLevel()

void G4ParticleHPInelastic::SetVerboseLevel ( G4int  newValue)

Field Documentation

◆ dataDirVariable

G4String G4ParticleHPInelastic::dataDirVariable
protected

Definition at line 110 of file G4ParticleHPInelastic.hh.

Referenced by BuildPhysicsTable(), and G4ParticleHPInelastic().

◆ dirName

G4String G4ParticleHPInelastic::dirName
protected

Definition at line 111 of file G4ParticleHPInelastic.hh.

Referenced by BuildPhysicsTable(), and G4ParticleHPInelastic().

◆ epCheckLevels

std::pair<G4double, G4double> G4HadronicInteraction::epCheckLevels
privateinherited

◆ isBlocked

G4bool G4HadronicInteraction::isBlocked
protectedinherited

◆ numEle

G4int G4ParticleHPInelastic::numEle
protected

Definition at line 112 of file G4ParticleHPInelastic.hh.

Referenced by BuildPhysicsTable().

◆ recoilEnergyThreshold

G4double G4HadronicInteraction::recoilEnergyThreshold
privateinherited

◆ registry

G4HadronicInteractionRegistry* G4HadronicInteraction::registry
privateinherited

◆ theBlockedList

std::vector<const G4Material *> G4HadronicInteraction::theBlockedList
privateinherited

◆ theBlockedListElements

std::vector<const G4Element *> G4HadronicInteraction::theBlockedListElements
privateinherited

◆ theInelastic

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

Definition at line 109 of file G4ParticleHPInelastic.hh.

Referenced by BuildPhysicsTable(), and ~G4ParticleHPInelastic().

◆ theMaxEnergy

G4double G4HadronicInteraction::theMaxEnergy
protectedinherited

◆ theMaxEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMaxEnergyList
privateinherited

◆ theMaxEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMaxEnergyListElements
privateinherited

◆ theMinEnergy

G4double G4HadronicInteraction::theMinEnergy
protectedinherited

◆ theMinEnergyList

std::vector<std::pair<G4double, const G4Material *> > G4HadronicInteraction::theMinEnergyList
privateinherited

◆ theMinEnergyListElements

std::vector<std::pair<G4double, const G4Element *> > G4HadronicInteraction::theMinEnergyListElements
privateinherited

◆ theModelName

G4String G4HadronicInteraction::theModelName
privateinherited

◆ theParticleChange

G4HadFinalState G4HadronicInteraction::theParticleChange
protectedinherited

Definition at line 172 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LENDCapture::ApplyYourself(), G4LENDElastic::ApplyYourself(), G4LENDFission::ApplyYourself(), G4LENDInelastic::ApplyYourself(), G4ElectroVDNuclearModel::ApplyYourself(), G4ParticleHPThermalScattering::ApplyYourself(), G4NeutrinoElectronNcModel::ApplyYourself(), G4NeutronElectronElModel::ApplyYourself(), G4LFission::ApplyYourself(), G4ANuElNucleusCcModel::ApplyYourself(), G4ANuElNucleusNcModel::ApplyYourself(), G4ANuMuNucleusCcModel::ApplyYourself(), G4ANuMuNucleusNcModel::ApplyYourself(), G4MuonVDNuclearModel::ApplyYourself(), G4NeutrinoElectronCcModel::ApplyYourself(), G4NuElNucleusCcModel::ApplyYourself(), G4NuElNucleusNcModel::ApplyYourself(), G4NuMuNucleusCcModel::ApplyYourself(), G4NuMuNucleusNcModel::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4BinaryCascade::ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4LMsdGenerator::ApplyYourself(), G4ElectroVDNuclearModel::CalculateEMVertex(), G4MuonVDNuclearModel::CalculateEMVertex(), G4ElectroVDNuclearModel::CalculateHadronicVertex(), G4MuonVDNuclearModel::CalculateHadronicVertex(), G4NeutrinoNucleusModel::CoherentPion(), G4CascadeInterface::copyOutputToHadronicResult(), G4BinaryCascade::DebugEpConservation(), G4BinaryCascade::DebugFinalEpConservation(), G4NeutrinoNucleusModel::FinalBarion(), G4NeutrinoNucleusModel::FinalMeson(), G4WilsonAbrasionModel::GetAbradedNucleons(), G4CascadeInterface::NoInteraction(), G4CascadeInterface::Propagate(), G4NeutrinoNucleusModel::RecoilDeexcitation(), G4LEHadronProtonElastic::~G4LEHadronProtonElastic(), G4LEnp::~G4LEnp(), and G4LFission::~G4LFission().

◆ theProjectile

G4ParticleDefinition* G4ParticleHPInelastic::theProjectile
private

Definition at line 154 of file G4ParticleHPInelastic.hh.

Referenced by G4ParticleHPInelastic().

◆ verboseLevel

G4int G4HadronicInteraction::verboseLevel
protectedinherited

Definition at line 177 of file G4HadronicInteraction.hh.

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4LFission::ApplyYourself(), G4MuMinusCapturePrecompound::ApplyYourself(), G4NeutronRadCapture::ApplyYourself(), G4LowEGammaNuclearModel::ApplyYourself(), G4ChargeExchange::ApplyYourself(), G4HadronElastic::ApplyYourself(), G4LEHadronProtonElastic::ApplyYourself(), G4LEnp::ApplyYourself(), G4LEpp::ApplyYourself(), G4CascadeInterface::ApplyYourself(), G4CascadeInterface::checkFinalResult(), G4CascadeInterface::copyOutputToHadronicResult(), G4CascadeInterface::copyOutputToReactionProducts(), G4LENDModel::create_used_target_map(), G4CascadeInterface::createBullet(), G4CascadeInterface::createTarget(), G4ElasticHadrNucleusHE::DefineHadronValues(), G4ElasticHadrNucleusHE::FillData(), G4ElasticHadrNucleusHE::FillFq2(), G4DiffuseElastic::G4DiffuseElastic(), G4DiffuseElasticV2::G4DiffuseElasticV2(), G4ElasticHadrNucleusHE::G4ElasticHadrNucleusHE(), G4EMDissociation::G4EMDissociation(), G4hhElastic::G4hhElastic(), G4NuclNuclDiffuseElastic::G4NuclNuclDiffuseElastic(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), G4ElasticHadrNucleusHE::GetFt(), G4ElasticHadrNucleusHE::GetLightFq2(), G4ElasticHadrNucleusHE::GetQ2_2(), G4HadronicInteraction::GetVerboseLevel(), G4ElasticHadrNucleusHE::HadronNucleusQ2_2(), G4ElasticHadrNucleusHE::HadronProtonQ2(), G4LFission::init(), G4DiffuseElastic::Initialise(), G4DiffuseElasticV2::Initialise(), G4NuclNuclDiffuseElastic::Initialise(), G4DiffuseElastic::InitialiseOnFly(), G4DiffuseElasticV2::InitialiseOnFly(), G4NuclNuclDiffuseElastic::InitialiseOnFly(), G4CascadeInterface::makeDynamicParticle(), G4CascadeInterface::NoInteraction(), G4CascadeInterface::Propagate(), G4ElasticHadrNucleusHE::SampleInvariantT(), G4AntiNuclElastic::SampleThetaCMS(), G4DiffuseElastic::SampleThetaLab(), G4NuclNuclDiffuseElastic::SampleThetaLab(), G4AntiNuclElastic::SampleThetaLab(), G4WilsonAbrasionModel::SetUseAblation(), G4HadronicInteraction::SetVerboseLevel(), G4WilsonAbrasionModel::SetVerboseLevel(), G4DiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElasticV2::ThetaCMStoThetaLab(), G4NuclNuclDiffuseElastic::ThetaCMStoThetaLab(), G4DiffuseElastic::ThetaLabToThetaCMS(), G4DiffuseElasticV2::ThetaLabToThetaCMS(), and G4NuclNuclDiffuseElastic::ThetaLabToThetaCMS().


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