Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4NeutronGeneralProcess Class Reference

#include <G4NeutronGeneralProcess.hh>

+ Inheritance diagram for G4NeutronGeneralProcess:

Public Member Functions

 G4NeutronGeneralProcess (const G4String &pname="NeutronGeneralProc")
 
 ~G4NeutronGeneralProcess () override
 
G4bool IsApplicable (const G4ParticleDefinition &) override
 
void ProcessDescription (std::ostream &outFile) const override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
G4bool StorePhysicsTable (const G4ParticleDefinition *part, const G4String &directory, G4bool ascii) override
 
void StartTracking (G4Track *) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
const G4VProcessGetCreatorProcess () const override
 
const G4StringGetSubProcessName () const
 
G4int GetSubProcessSubType () const
 
void SetInelasticProcess (G4HadronicProcess *)
 
void SetElasticProcess (G4HadronicProcess *)
 
void SetCaptureProcess (G4HadronicProcess *)
 
G4VCrossSectionDataSetGetXSection (G4int type)
 
G4HadronicProcessGetHadronicProcess (G4int type)
 
const G4VProcessGetSelectedProcess () const
 
void SetTimeLimit (G4double val)
 
void SetMinEnergyLimit (G4double val)
 
 G4NeutronGeneralProcess (G4NeutronGeneralProcess &)=delete
 
G4NeutronGeneralProcessoperator= (const G4NeutronGeneralProcess &right)=delete
 
- Public Member Functions inherited from G4HadronicProcess
 G4HadronicProcess (const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
 
 G4HadronicProcess (const G4String &processName, G4HadronicProcessType subType)
 
 ~G4HadronicProcess () override
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=nullptr)
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
 
G4HadronicInteractionGetHadronicModel (const G4String &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
const G4NucleusGetTargetNucleus () const
 
G4NucleusGetTargetNucleusPointer ()
 
const G4IsotopeGetTargetIsotope ()
 
G4double ComputeCrossSection (const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
 
G4HadXSType CrossSectionType () const
 
void SetCrossSectionType (G4HadXSType val)
 
void BiasCrossSectionByFactor (G4double aScale)
 
void MultiplyCrossSectionBy (G4double factor)
 
G4double CrossSectionFactor () const
 
void SetIntegral (G4bool val)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
std::vector< G4TwoPeaksHadXS * > * TwoPeaksXS () const
 
std::vector< G4double > * EnergyOfCrossSectionMax () const
 
G4HadronicProcessoperator= (const G4HadronicProcess &right)=delete
 
 G4HadronicProcess (const G4HadronicProcess &)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4double ComputeGeneralLambda (size_t idxe, size_t idxt)
 
G4double GetProbability (size_t idxt)
 
void SelectedProcess (const G4Step &step, G4HadronicProcess *ptr, G4CrossSectionDataStore *)
 
- Protected Member Functions inherited from G4HadronicProcess
G4HadronicInteractionChooseHadronicInteraction (const G4HadProjectile &aHadProjectile, G4Nucleus &aTargetNucleus, const G4Material *aMaterial, const G4Element *anElement)
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
- Protected Member Functions inherited from G4VDiscreteProcess
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4HadronicProcess
G4HadProjectile thePro
 
G4ParticleChangetheTotalResult
 
G4CrossSectionDataStoretheCrossSectionDataStore
 
G4double fWeight = 1.0
 
G4double aScaleFactor = 1.0
 
G4double theLastCrossSection = 0.0
 
G4double mfpKinEnergy = DBL_MAX
 
G4int epReportLevel = 0
 
G4HadXSType fXSType = fHadNoIntegral
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 64 of file G4NeutronGeneralProcess.hh.

Constructor & Destructor Documentation

◆ G4NeutronGeneralProcess() [1/2]

G4NeutronGeneralProcess::G4NeutronGeneralProcess ( const G4String & pname = "NeutronGeneralProc")
explicit

Definition at line 80 of file G4NeutronGeneralProcess.cc.

81: G4HadronicProcess(pname),
82 fMinEnergy(1*CLHEP::keV),
83 fMiddleEnergy(20*CLHEP::MeV),
84 fMaxEnergy(100*CLHEP::TeV),
85 fTimeLimit(10*CLHEP::microsecond)
86{
89
90 fNeutron = G4Neutron::Neutron();
91
93 isMaster = false;
94 }
95}
@ fNeutronGeneral
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
void SetVerboseLevel(G4int value)
void SetProcessSubType(G4int)
G4bool IsWorkerThread()

◆ ~G4NeutronGeneralProcess()

G4NeutronGeneralProcess::~G4NeutronGeneralProcess ( )
override

Definition at line 99 of file G4NeutronGeneralProcess.cc.

100{
101 if(isMaster) {
102 delete theHandler;
103 theHandler = nullptr;
104 }
105}

◆ G4NeutronGeneralProcess() [2/2]

G4NeutronGeneralProcess::G4NeutronGeneralProcess ( G4NeutronGeneralProcess & )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronGeneralProcess::BuildPhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 271 of file G4NeutronGeneralProcess.cc.

272{
273 if(1 < verboseLevel) {
274 G4cout << "### G4NeutronGeneralProcess::BuildPhysicsTable() for "
275 << GetProcessName()
276 << " and particle " << part.GetParticleName()
277 << G4endl;
278 }
279 fElasticP->BuildPhysicsTable(part);
280 fInelasticP->BuildPhysicsTable(part);
281 fCaptureP->BuildPhysicsTable(part);
282
283 if(isMaster) {
284 std::size_t nmat = G4Material::GetNumberOfMaterials();
286
287 auto tables = theHandler->GetTables();
288
289 G4double sigEl(0.), sigInel(0.), sigCap(0.), val(0.), sum(0.);
290
291 for(std::size_t i=0; i<nmat; ++i) {
292 const G4Material* mat = (*matTable)[i];
293
294 // energy interval 0
295 std::size_t nn = (*(tables[0]))[i]->GetVectorLength();
296 if(1 < verboseLevel) {
297 G4cout << "======= Zone 0 ======= N= " << nn
298 << " for " << mat->GetName() << G4endl;
299 }
300 for(std::size_t j=0; j<nn; ++j) {
301 G4double e = (*(tables[0]))[i]->Energy(j);
302 G4double loge = G4Log(e);
303 sigEl = fXSFactorEl*ComputeCrossSection(fElasticXS, mat, e, loge);
304 sigInel = fXSFactorInel*ComputeCrossSection(fInelasticXS, mat, e, loge);
305 sigCap = ComputeCrossSection(fCaptureXS, mat, e, loge);
306 sum = sigEl + sigInel + sigCap;
307 if(1 < verboseLevel) {
308 G4cout << j << ". E= " << e << " xs=" << sum << " sigEl=" << sigEl
309 << " sigInel=" << sigInel << " sigCap=" << sigCap << G4endl;
310 }
311 (*(tables[0]))[i]->PutValue(j, sum);
312 val = sigEl/sum;
313 (*(tables[1]))[i]->PutValue(j, val);
314 val = (sigEl + sigInel)/sum;
315 (*(tables[2]))[i]->PutValue(j, val);
316 }
317
318 // energy interval 1
319 nn = (*(tables[3]))[0]->GetVectorLength();
320 if(1 < verboseLevel) {
321 G4cout << "======= Zone 1 ======= N= " << nn << G4endl;
322 }
323 for(std::size_t j=0; j<nn; ++j) {
324 G4double e = (*(tables[3]))[i]->Energy(j);
325 G4double loge = G4Log(e);
326 sigEl = fXSFactorEl*ComputeCrossSection(fElasticXS, mat, e, loge);
327 sigInel = fXSFactorInel*ComputeCrossSection(fInelasticXS, mat, e, loge);
328 sum = sigEl + sigInel;
329 if(1 < verboseLevel) {
330 G4cout << j << ". E= " << e << " xs=" << sum << " sigEl=" << sigEl
331 << " sigInel=" << sigInel << " factInel=" << fXSFactorInel
332 << G4endl;
333 }
334 (*(tables[3]))[i]->PutValue(j, sum);
335 val = sigInel/sum;
336 (*(tables[4]))[i]->PutValue(j, val);
337 }
338 }
339 }
340 if(1 < verboseLevel) {
341 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
342 << GetProcessName()
343 << " and particle " << part.GetParticleName()
344 << G4endl;
345 }
346}
G4double G4Log(G4double x)
Definition G4Log.hh:227
std::vector< G4Material * > G4MaterialTable
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
const std::vector< G4PhysicsTable * > & GetTables() const
void BuildPhysicsTable(const G4ParticleDefinition &) override
static std::size_t GetNumberOfMaterials()
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const
const G4String & GetParticleName() const
G4int verboseLevel
const G4String & GetProcessName() const

◆ ComputeGeneralLambda()

G4double G4NeutronGeneralProcess::ComputeGeneralLambda ( size_t idxe,
size_t idxt )
inlineprotected

Definition at line 193 of file G4NeutronGeneralProcess.hh.

194{
195 idxEnergy = idxe;
196 return theHandler->GetVector(idxt, matIndex)
197 ->LogVectorValue(fCurrE, fCurrLogE);
198}
const G4PhysicsVector * GetVector(std::size_t itable, std::size_t ivec) const
G4double LogVectorValue(const G4double energy, const G4double theLogEnergy) const

◆ GetCreatorProcess()

const G4VProcess * G4NeutronGeneralProcess::GetCreatorProcess ( ) const
overridevirtual

Reimplemented from G4VProcess.

Definition at line 521 of file G4NeutronGeneralProcess.cc.

522{
523 return fSelectedProc;
524}

◆ GetHadronicProcess()

G4HadronicProcess * G4NeutronGeneralProcess::GetHadronicProcess ( G4int type)

Definition at line 180 of file G4NeutronGeneralProcess.cc.

181{
182 G4HadronicProcess* ptr = nullptr;
184 if(type == fHadronElastic) { ptr = fElasticP; }
185 else if(type == fHadronInelastic) { ptr = fInelasticP; }
186 else if(type == fCapture) { ptr = fCaptureP; }
187 return ptr;
188}
G4HadronicProcessType
@ fHadronInelastic

◆ GetMeanFreePath()

G4double G4NeutronGeneralProcess::GetMeanFreePath ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overrideprotectedvirtual

Reimplemented from G4HadronicProcess.

Definition at line 484 of file G4NeutronGeneralProcess.cc.

487{
489 // recompute total cross section if needed
490 CurrentCrossSection(track);
492}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double currentInteractionLength

◆ GetProbability()

G4double G4NeutronGeneralProcess::GetProbability ( size_t idxt)
inlineprotected

Definition at line 202 of file G4NeutronGeneralProcess.hh.

203{
204 return theHandler->GetVector(idxt, matIndex)
205 ->LogVectorValue(fCurrE, fCurrLogE);
206}

Referenced by PostStepDoIt().

◆ GetSelectedProcess()

const G4VProcess * G4NeutronGeneralProcess::GetSelectedProcess ( ) const
inline

Definition at line 223 of file G4NeutronGeneralProcess.hh.

224{
225 return fSelectedProc;
226}

◆ GetSubProcessName()

const G4String & G4NeutronGeneralProcess::GetSubProcessName ( ) const

Definition at line 505 of file G4NeutronGeneralProcess.cc.

506{
507 return (nullptr != fSelectedProc) ? fSelectedProc->GetProcessName()
509}

◆ GetSubProcessSubType()

G4int G4NeutronGeneralProcess::GetSubProcessSubType ( ) const

Definition at line 513 of file G4NeutronGeneralProcess.cc.

514{
515 return (nullptr != fSelectedProc) ? fSelectedProc->GetProcessSubType()
517}
G4int GetProcessSubType() const

◆ GetXSection()

G4VCrossSectionDataSet * G4NeutronGeneralProcess::GetXSection ( G4int type)

Definition at line 168 of file G4NeutronGeneralProcess.cc.

169{
170 G4VCrossSectionDataSet* xs = nullptr;
172 if(type == fHadronElastic) { xs = fElasticXS; }
173 else if(type == fHadronInelastic) { xs = fInelasticXS; }
174 else if(type == fCapture) { xs = fCaptureXS; }
175 return xs;
176}

◆ IsApplicable()

G4bool G4NeutronGeneralProcess::IsApplicable ( const G4ParticleDefinition & )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 109 of file G4NeutronGeneralProcess.cc.

110{
111 return true;
112}

◆ operator=()

G4NeutronGeneralProcess & G4NeutronGeneralProcess::operator= ( const G4NeutronGeneralProcess & right)
delete

◆ PostStepDoIt()

G4VParticleChange * G4NeutronGeneralProcess::PostStepDoIt ( const G4Track & track,
const G4Step & step )
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 416 of file G4NeutronGeneralProcess.cc.

418{
419 fSelectedProc = this;
420 // time limit
421 if(0.0 == fLambda) {
424 return theTotalResult;
425 }
426 // In all cases clear number of interaction lengths
429 /*
430 G4cout << "PostStep: preStepLambda= " << fLambda << " idxE= " << idxEnergy
431 << " matIndex=" << matIndex << G4endl;
432 */
433 if (0 == idxEnergy) {
434 if(q <= GetProbability(1)) {
435 SelectedProcess(step, fElasticP, fXSSElastic);
436 } else if(q <= GetProbability(2)) {
437 SelectedProcess(step, fInelasticP, fXSSInelastic);
438 } else {
439 SelectedProcess(step, fCaptureP, fXSSCapture);
440 }
441 } else {
442 if(q <= GetProbability(4)) {
443 SelectedProcess(step, fInelasticP, fXSSInelastic);
444 } else {
445 SelectedProcess(step, fElasticP, fXSSElastic);
446 }
447 }
448 // total cross section is needed for selection of an element
449 if(fCurrMat->GetNumberOfElements() > 1) {
450 fCurrentXSS->ComputeCrossSection(track.GetDynamicParticle(), fCurrMat);
451 }
452 /*
453 G4cout << "## neutron E(MeV)=" << fCurrE << " inside " << fCurrMat->GetName()
454 << fSelectedProc->GetProcessName()
455 << " time(ns)=" << track.GetGlobalTime()/ns << G4endl;
456 */
457 // sample secondaries
458 return fSelectedProc->PostStepDoIt(track, step);
459}
@ fStopAndKill
#define G4UniformRand()
Definition Randomize.hh:52
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4ParticleChange * theTotalResult
std::size_t GetNumberOfElements() const
void SelectedProcess(const G4Step &step, G4HadronicProcess *ptr, G4CrossSectionDataStore *)
G4double GetProbability(size_t idxt)
void Initialize(const G4Track &) override
const G4DynamicParticle * GetDynamicParticle() const
void ProposeTrackStatus(G4TrackStatus status)
G4double theNumberOfInteractionLengthLeft

◆ PostStepGetPhysicalInteractionLength()

G4double G4NeutronGeneralProcess::PostStepGetPhysicalInteractionLength ( const G4Track & track,
G4double previousStepSize,
G4ForceCondition * condition )
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 375 of file G4NeutronGeneralProcess.cc.

379{
381
382 // time limit
383 if(track.GetGlobalTime() >= fTimeLimit) {
384 fLambda = 0.0;
385 return 0.0;
386 }
387
388 // recompute total cross section if needed
389 CurrentCrossSection(track);
390
392
393 // beggining of tracking (or just after DoIt of this process)
396
397 } else {
398
400 previousStepSize/currentInteractionLength;
403 }
404
406 /*
407 G4cout << "PostStepGetPhysicalInteractionLength: e= " << energy
408 << " idxe= " << idxEnergy << " xs= " << fLambda
409 << " x= " << x << G4endl;
410 */
411 return x;
412}
G4double GetGlobalTime() const
G4double theInitialNumberOfInteractionLength

◆ PreparePhysicsTable()

void G4NeutronGeneralProcess::PreparePhysicsTable ( const G4ParticleDefinition & part)
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 192 of file G4NeutronGeneralProcess.cc.

193{
194 if(1 < verboseLevel) {
195 G4cout << "G4NeutronGeneralProcess::PreparePhysicsTable() for "
196 << GetProcessName()
197 << " and particle " << part.GetParticleName()
198 << " isMaster: " << isMaster << G4endl;
199 }
200 G4bool noEl = (nullptr == fElasticP);
201 G4bool noInel = (nullptr == fInelasticP);
202 G4bool noCap = (nullptr == fCaptureP);
203 if(noEl || noInel || noCap) {
205 ed << "Incomplete configuration of the neutron general process." << G4endl;
206 if(noEl) { ed << "Neutron elastic process is not defined" << G4endl; }
207 if(noInel) { ed << "Neutron inelastic process is not defined" << G4endl; }
208 if(noCap) { ed << "Neutron capture process is not defined" << G4endl; }
209 G4Exception ("G4NeutronGeneralProcess::PreparePhysicsTable(..)", "had001",
210 FatalException, ed, "");
211 return;
212 }
213
215
217 fMaxEnergy = std::max(100*MeV, param->GetMaxEnergy());
218 if(param->ApplyFactorXS()) {
219 fXSFactorEl = param->XSFactorNucleonElastic();
220 fXSFactorInel = param->XSFactorNucleonInelastic();
221 }
222
223 fElasticP->PreparePhysicsTable(part);
224 fInelasticP->PreparePhysicsTable(part);
225 fCaptureP->PreparePhysicsTable(part);
226
227 std::size_t nmat = G4Material::GetNumberOfMaterials();
229
230 std::size_t nmax = 0;
231 for(std::size_t i=0; i<nmat; ++i) {
232 std::size_t nelm = (*matTable)[i]->GetNumberOfElements();
233 nmax = std::max(nmax, nelm);
234 }
235 fXsec.resize(nmax);
236
237 if(isMaster) {
238 if(nullptr == theHandler) {
239 theHandler = new G4HadDataHandler(nTables);
240 }
241
242 fMaxEnergy = std::max(fMaxEnergy, param->GetMaxEnergy());
243 nLowE *= G4lrint(std::log10(fMiddleEnergy/fMinEnergy));
244 nHighE *= G4lrint(std::log10(fMaxEnergy/fMiddleEnergy));
245
246 G4PhysicsVector* vec = nullptr;
247 G4PhysicsLogVector aVector(fMinEnergy, fMiddleEnergy, nLowE, false);
248 G4PhysicsLogVector bVector(fMiddleEnergy, fMaxEnergy, nHighE, false);
249
250 for(std::size_t i=0; i<nTables; ++i) {
251 G4PhysicsTable* table = new G4PhysicsTable();
252 theHandler->UpdateTable(table, i);
253 table->resize(nmat);
254 for(std::size_t j=0; j<nmat; ++j) {
255 vec = (*table)[j];
256 if (nullptr == vec) {
257 if(i <= 2) {
258 vec = new G4PhysicsVector(aVector);
259 } else {
260 vec = new G4PhysicsVector(bVector);
261 }
263 }
264 }
265 }
266 }
267}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
bool G4bool
Definition G4Types.hh:86
void UpdateTable(G4PhysicsTable *, std::size_t idx)
static G4HadronicParameters * Instance()
G4double XSFactorNucleonElastic() const
G4double XSFactorNucleonInelastic() const
void PreparePhysicsTable(const G4ParticleDefinition &) override
static void SetPhysicsVector(G4PhysicsTable *physTable, std::size_t idx, G4PhysicsVector *vec)
void resize(std::size_t, G4PhysicsVector *vec=nullptr)
int G4lrint(double ad)
Definition templates.hh:134

◆ ProcessDescription()

void G4NeutronGeneralProcess::ProcessDescription ( std::ostream & outFile) const
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 496 of file G4NeutronGeneralProcess.cc.

497{
498 fElasticP->ProcessDescription(out);
499 fInelasticP->ProcessDescription(out);
500 fCaptureP->ProcessDescription(out);
501}
void ProcessDescription(std::ostream &outFile) const override

◆ SelectedProcess()

void G4NeutronGeneralProcess::SelectedProcess ( const G4Step & step,
G4HadronicProcess * ptr,
G4CrossSectionDataStore * xs )
inlineprotected

Definition at line 211 of file G4NeutronGeneralProcess.hh.

215{
216 fSelectedProc = ptr;
217 fCurrentXSS = xs;
219}
void SetProcessDefinedStep(const G4VProcess *aValue)
G4StepPoint * GetPostStepPoint() const

Referenced by PostStepDoIt().

◆ SetCaptureProcess()

void G4NeutronGeneralProcess::SetCaptureProcess ( G4HadronicProcess * ptr)

Definition at line 142 of file G4NeutronGeneralProcess.cc.

143{
144 fCaptureP = ptr;
145 fXSSCapture = ptr->GetCrossSectionDataStore();
146 fCaptureXS = InitialisationXS(ptr);
147 if(nullptr == fCaptureXS) {
148 fCaptureXS = new G4NeutronCaptureXS();
149 ptr->AddDataSet(fCaptureXS);
150 }
151}
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
G4CrossSectionDataStore * GetCrossSectionDataStore()

◆ SetElasticProcess()

void G4NeutronGeneralProcess::SetElasticProcess ( G4HadronicProcess * ptr)

Definition at line 129 of file G4NeutronGeneralProcess.cc.

130{
131 fElasticP = ptr;
132 fXSSElastic = ptr->GetCrossSectionDataStore();
133 fElasticXS = InitialisationXS(ptr);
134 if(nullptr == fElasticXS) {
135 fElasticXS = new G4NeutronElasticXS();
136 ptr->AddDataSet(fElasticXS);
137 }
138}

◆ SetInelasticProcess()

void G4NeutronGeneralProcess::SetInelasticProcess ( G4HadronicProcess * ptr)

Definition at line 116 of file G4NeutronGeneralProcess.cc.

117{
118 fInelasticP = ptr;
119 fXSSInelastic = ptr->GetCrossSectionDataStore();
120 fInelasticXS = InitialisationXS(ptr);
121 if(nullptr == fInelasticXS) {
122 fInelasticXS = new G4NeutronInelasticXS();
123 ptr->AddDataSet(fInelasticXS);
124 }
125}

◆ SetMinEnergyLimit()

void G4NeutronGeneralProcess::SetMinEnergyLimit ( G4double val)
inline

Definition at line 254 of file G4NeutronGeneralProcess.hh.

255{
256 fMinEnergy = val;
257}

◆ SetTimeLimit()

void G4NeutronGeneralProcess::SetTimeLimit ( G4double val)
inline

Definition at line 247 of file G4NeutronGeneralProcess.hh.

248{
249 fTimeLimit = val;
250}

◆ StartTracking()

void G4NeutronGeneralProcess::StartTracking ( G4Track * )
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 367 of file G4NeutronGeneralProcess.cc.

368{
370 fCurrMat = nullptr;
371}

◆ StorePhysicsTable()

G4bool G4NeutronGeneralProcess::StorePhysicsTable ( const G4ParticleDefinition * part,
const G4String & directory,
G4bool ascii )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 464 of file G4NeutronGeneralProcess.cc.

467{
468 G4bool yes = true;
469 if(!isMaster) { return yes; }
470 for(std::size_t i=0; i<nTables; ++i) {
471 G4String nam = (0==i || 3==i)
472 ? "LambdaNeutronGeneral" + nameT[i] : "ProbNeutronGeneral" + nameT[i];
473 G4String fnam = GetPhysicsTableFileName(part, directory, nam, ascii);
474 auto table = theHandler->Table(i);
475 if(nullptr == table || !table->StorePhysicsTable(fnam, ascii)) {
476 yes = false;
477 }
478 }
479 return yes;
480}
G4PhysicsTable * Table(std::size_t idx) const
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)

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