Geant4 11.1.1
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 *)
 
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 StartTracking (G4Track *track) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double, G4ForceCondition *) override
 
G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep) override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
 
G4HadronicInteractionGetHadronicModel (const G4String &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *) override
 
const G4NucleusGetTargetNucleus () const
 
G4NucleusGetTargetNucleusPointer ()
 
const G4IsotopeGetTargetIsotope ()
 
G4double ComputeCrossSection (const G4ParticleDefinition *, const G4Material *, const G4double kinEnergy)
 
G4HadXSType CrossSectionType () const
 
void SetCrossSectionType (G4HadXSType val)
 
void ProcessDescription (std::ostream &outFile) const override
 
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 PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
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 const G4VProcessGetCreatorProcess () const
 
virtual void StartTracking (G4Track *)
 
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
 
virtual void ProcessDescription (std::ostream &outfile) 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 &)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- 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
 
G4long 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 79 of file G4NeutronGeneralProcess.cc.

80: G4HadronicProcess(pname),
81 fMinEnergy(1*CLHEP::keV),
82 fMiddleEnergy(20*CLHEP::MeV),
83 fMaxEnergy(100*CLHEP::TeV),
84 fTimeLimit(10*CLHEP::microsecond)
85{
88
89 fNeutron = G4Neutron::Neutron();
90
92 isMaster = false;
93 }
94}
@ fNeutronGeneral
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
void SetVerboseLevel(G4int value)
Definition: G4VProcess.hh:416
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410
G4bool IsWorkerThread()
Definition: G4Threading.cc:123

◆ ~G4NeutronGeneralProcess()

G4NeutronGeneralProcess::~G4NeutronGeneralProcess ( )
override

Definition at line 98 of file G4NeutronGeneralProcess.cc.

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

◆ G4NeutronGeneralProcess() [2/2]

G4NeutronGeneralProcess::G4NeutronGeneralProcess ( G4NeutronGeneralProcess )
delete

Member Function Documentation

◆ BuildPhysicsTable()

void G4NeutronGeneralProcess::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 246 of file G4NeutronGeneralProcess.cc.

247{
248 if(1 < verboseLevel) {
249 G4cout << "### G4NeutronGeneralProcess::BuildPhysicsTable() for "
250 << GetProcessName()
251 << " and particle " << part.GetParticleName()
252 << G4endl;
253 }
254 fElastic->BuildPhysicsTable(part);
255 fInelastic->BuildPhysicsTable(part);
256 fCapture->BuildPhysicsTable(part);
257
258 if(isMaster) {
259 std::size_t nmat = G4Material::GetNumberOfMaterials();
261
262 auto tables = theHandler->GetTables();
263
264 G4double sigEl(0.), sigInel(0.), sigCap(0.), val(0.), sum(0.);
265
266 for(std::size_t i=0; i<nmat; ++i) {
267 const G4Material* mat = (*matTable)[i];
268
269 // energy interval 0
270 std::size_t nn = (*(tables[0]))[i]->GetVectorLength();
271 if(1 < verboseLevel) {
272 G4cout << "======= Zone 0 ======= N= " << nn
273 << " for " << mat->GetName() << G4endl;
274 }
275 for(std::size_t j=0; j<nn; ++j) {
276 G4double e = (*(tables[0]))[i]->Energy(j);
277 G4double loge = G4Log(e);
278 sigEl = fXSFactorEl*ComputeCrossSection(fElasticXS, mat, e, loge);
279 sigInel = fXSFactorInel*ComputeCrossSection(fInelasticXS, mat, e, loge);
280 sigCap = ComputeCrossSection(fCaptureXS, mat, e, loge);
281 sum = sigEl + sigInel + sigCap;
282 if(1 < verboseLevel) {
283 G4cout << j << ". E= " << e << " xs=" << sum << " sigEl=" << sigEl
284 << " sigInel=" << sigInel << " sigCap=" << sigCap << G4endl;
285 }
286 (*(tables[0]))[i]->PutValue(j, sum);
287 val = sigEl/sum;
288 (*(tables[1]))[i]->PutValue(j, val);
289 val = (sigEl + sigInel)/sum;
290 (*(tables[2]))[i]->PutValue(j, val);
291 }
292
293 // energy interval 1
294 nn = (*(tables[3]))[0]->GetVectorLength();
295 if(1 < verboseLevel) {
296 G4cout << "======= Zone 1 ======= N= " << nn << G4endl;
297 }
298 for(std::size_t j=0; j<nn; ++j) {
299 G4double e = (*(tables[3]))[i]->Energy(j);
300 G4double loge = G4Log(e);
301 sigEl = fXSFactorEl*ComputeCrossSection(fElasticXS, mat, e, loge);
302 sigInel = fXSFactorInel*ComputeCrossSection(fInelasticXS, mat, e, loge);
303 sum = sigEl + sigInel;
304 if(1 < verboseLevel) {
305 G4cout << j << ". E= " << e << " xs=" << sum << " sigEl=" << sigEl
306 << " sigInel=" << sigInel << " factInel=" << fXSFactorInel
307 << G4endl;
308 }
309 (*(tables[3]))[i]->PutValue(j, sum);
310 val = sigInel/sum;
311 (*(tables[4]))[i]->PutValue(j, val);
312 }
313 }
314 }
315 if(1 < verboseLevel) {
316 G4cout << "### G4VEmProcess::BuildPhysicsTable() done for "
317 << GetProcessName()
318 << " and particle " << part.GetParticleName()
319 << G4endl;
320 }
321}
G4double G4Log(G4double x)
Definition: G4Log.hh:227
std::vector< G4Material * > G4MaterialTable
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const std::vector< G4PhysicsTable * > & GetTables() const
void BuildPhysicsTable(const G4ParticleDefinition &) override
static size_t GetNumberOfMaterials()
Definition: G4Material.cc:684
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
const G4String & GetParticleName() const
G4int verboseLevel
Definition: G4VProcess.hh:360
const G4String & GetProcessName() const
Definition: G4VProcess.hh:386

◆ ComputeGeneralLambda()

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

Definition at line 189 of file G4NeutronGeneralProcess.hh.

190{
191 idxEnergy = idxe;
192 return theHandler->GetVector(idxt, matIndex)
193 ->LogVectorValue(fCurrE, fCurrLogE);
194}
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 496 of file G4NeutronGeneralProcess.cc.

497{
498 return fSelectedProc;
499}

◆ GetMeanFreePath()

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

Reimplemented from G4HadronicProcess.

Definition at line 459 of file G4NeutronGeneralProcess.cc.

462{
464 // recompute total cross section if needed
465 CurrentCrossSection(track);
467}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced
G4double currentInteractionLength
Definition: G4VProcess.hh:339

◆ GetProbability()

G4double G4NeutronGeneralProcess::GetProbability ( size_t  idxt)
inlineprotected

Definition at line 198 of file G4NeutronGeneralProcess.hh.

199{
200 return theHandler->GetVector(idxt, matIndex)
201 ->LogVectorValue(fCurrE, fCurrLogE);
202}

Referenced by PostStepDoIt().

◆ GetSelectedProcess()

const G4VProcess * G4NeutronGeneralProcess::GetSelectedProcess ( ) const
inline

Definition at line 219 of file G4NeutronGeneralProcess.hh.

220{
221 return fSelectedProc;
222}

◆ GetSubProcessName()

const G4String & G4NeutronGeneralProcess::GetSubProcessName ( ) const

Definition at line 480 of file G4NeutronGeneralProcess.cc.

481{
482 return (nullptr != fSelectedProc) ? fSelectedProc->GetProcessName()
484}

◆ GetSubProcessSubType()

G4int G4NeutronGeneralProcess::GetSubProcessSubType ( ) const

Definition at line 488 of file G4NeutronGeneralProcess.cc.

489{
490 return (nullptr != fSelectedProc) ? fSelectedProc->GetProcessSubType()
492}
G4int GetProcessSubType() const
Definition: G4VProcess.hh:404

◆ IsApplicable()

G4bool G4NeutronGeneralProcess::IsApplicable ( const G4ParticleDefinition )
overridevirtual

Reimplemented from G4VProcess.

Definition at line 108 of file G4NeutronGeneralProcess.cc.

109{
110 return true;
111}

◆ 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 391 of file G4NeutronGeneralProcess.cc.

393{
394 fSelectedProc = this;
395 // time limit
396 if(0.0 == fLambda) {
399 return theTotalResult;
400 }
401 // In all cases clear number of interaction lengths
404 /*
405 G4cout << "PostStep: preStepLambda= " << fLambda << " idxE= " << idxEnergy
406 << " matIndex=" << matIndex << G4endl;
407 */
408 if (0 == idxEnergy) {
409 if(q <= GetProbability(1)) {
410 SelectedProcess(step, fElastic, fXSSElastic);
411 } else if(q <= GetProbability(2)) {
412 SelectedProcess(step, fInelastic, fXSSInelastic);
413 } else {
414 SelectedProcess(step, fCapture, fXSSCapture);
415 }
416 } else {
417 if(q <= GetProbability(4)) {
418 SelectedProcess(step, fInelastic, fXSSInelastic);
419 } else {
420 SelectedProcess(step, fElastic, fXSSElastic);
421 }
422 }
423 // total cross section is needed for selection of an element
424 if(fCurrMat->GetNumberOfElements() > 1) {
425 fCurrentXSS->ComputeCrossSection(track.GetDynamicParticle(), fCurrMat);
426 }
427 /*
428 G4cout << "## neutron E(MeV)=" << fCurrE << " inside " << fCurrMat->GetName()
429 << fSelectedProc->GetProcessName()
430 << " time(ns)=" << track.GetGlobalTime()/ns << G4endl;
431 */
432 // sample secondaries
433 return fSelectedProc->PostStepDoIt(track, step);
434}
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:52
G4double ComputeCrossSection(const G4DynamicParticle *, const G4Material *)
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
G4ParticleChange * theTotalResult
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
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
Definition: G4VProcess.hh:335

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4HadronicProcess.

Definition at line 350 of file G4NeutronGeneralProcess.cc.

354{
356
357 // time limit
358 if(track.GetGlobalTime() >= fTimeLimit) {
359 fLambda = 0.0;
360 return 0.0;
361 }
362
363 // recompute total cross section if needed
364 CurrentCrossSection(track);
365
367
368 // beggining of tracking (or just after DoIt of this process)
371
372 } else {
373
375 previousStepSize/currentInteractionLength;
378 }
379
381 /*
382 G4cout << "PostStepGetPhysicalInteractionLength: e= " << energy
383 << " idxe= " << idxEnergy << " xs= " << fLambda
384 << " x= " << x << G4endl;
385 */
386 return x;
387}
G4double GetGlobalTime() const
G4double theInitialNumberOfInteractionLength
Definition: G4VProcess.hh:342

◆ PreparePhysicsTable()

void G4NeutronGeneralProcess::PreparePhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 167 of file G4NeutronGeneralProcess.cc.

168{
169 if(1 < verboseLevel) {
170 G4cout << "G4NeutronGeneralProcess::PreparePhysicsTable() for "
171 << GetProcessName()
172 << " and particle " << part.GetParticleName()
173 << " isMaster: " << isMaster << G4endl;
174 }
175 G4bool noEl = (nullptr == fElastic);
176 G4bool noInel = (nullptr == fInelastic);
177 G4bool noCap = (nullptr == fCapture);
178 if(noEl || noInel || noCap) {
180 ed << "Incomplete configuration of the neutron general process." << G4endl;
181 if(noEl) { ed << "Neutron elastic process is not defined" << G4endl; }
182 if(noInel) { ed << "Neutron inelastic process is not defined" << G4endl; }
183 if(noCap) { ed << "Neutron capture process is not defined" << G4endl; }
184 G4Exception ("G4NeutronGeneralProcess::PreparePhysicsTable(..)", "had001",
185 FatalException, ed, "");
186 return;
187 }
188
190
192 fMaxEnergy = std::max(100*MeV, param->GetMaxEnergy());
193 if(param->ApplyFactorXS()) {
194 fXSFactorEl = param->XSFactorNucleonElastic();
195 fXSFactorInel = param->XSFactorNucleonInelastic();
196 }
197
198 fElastic->PreparePhysicsTable(part);
199 fInelastic->PreparePhysicsTable(part);
200 fCapture->PreparePhysicsTable(part);
201
202 std::size_t nmat = G4Material::GetNumberOfMaterials();
204
205 std::size_t nmax = 0;
206 for(std::size_t i=0; i<nmat; ++i) {
207 std::size_t nelm = (*matTable)[i]->GetNumberOfElements();
208 nmax = std::max(nmax, nelm);
209 }
210 fXsec.resize(nmax);
211
212 if(isMaster) {
213 if(nullptr == theHandler) {
214 theHandler = new G4HadDataHandler(nTables);
215 }
216
217 fMaxEnergy = std::max(fMaxEnergy, param->GetMaxEnergy());
218 nLowE *= G4lrint(std::log10(fMiddleEnergy/fMinEnergy));
219 nHighE *= G4lrint(std::log10(fMaxEnergy/fMiddleEnergy));
220
221 G4PhysicsVector* vec = nullptr;
222 G4PhysicsLogVector aVector(fMinEnergy, fMiddleEnergy, nLowE, false);
223 G4PhysicsLogVector bVector(fMiddleEnergy, fMaxEnergy, nHighE, false);
224
225 for(std::size_t i=0; i<nTables; ++i) {
226 G4PhysicsTable* table = new G4PhysicsTable();
227 theHandler->UpdateTable(table, i);
228 table->resize(nmat);
229 for(std::size_t j=0; j<nmat; ++j) {
230 vec = (*table)[j];
231 if (nullptr == vec) {
232 if(i <= 2) {
233 vec = new G4PhysicsVector(aVector);
234 } else {
235 vec = new G4PhysicsVector(bVector);
236 }
238 }
239 }
240 }
241 }
242}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
bool G4bool
Definition: G4Types.hh:86
void UpdateTable(G4PhysicsTable *, std::size_t idx)
static G4HadronicParameters * Instance()
G4double XSFactorNucleonElastic() const
G4double GetMaxEnergy() 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 471 of file G4NeutronGeneralProcess.cc.

472{
473 fElastic->ProcessDescription(out);
474 fInelastic->ProcessDescription(out);
475 fCapture->ProcessDescription(out);
476}
void ProcessDescription(std::ostream &outFile) const override

◆ SelectedProcess()

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

Definition at line 207 of file G4NeutronGeneralProcess.hh.

211{
212 fSelectedProc = ptr;
213 fCurrentXSS = xs;
215}
void SetProcessDefinedStep(const G4VProcess *aValue)
G4StepPoint * GetPostStepPoint() const

Referenced by PostStepDoIt().

◆ SetCaptureProcess()

void G4NeutronGeneralProcess::SetCaptureProcess ( G4HadronicProcess ptr)

Definition at line 141 of file G4NeutronGeneralProcess.cc.

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

◆ SetElasticProcess()

void G4NeutronGeneralProcess::SetElasticProcess ( G4HadronicProcess ptr)

Definition at line 128 of file G4NeutronGeneralProcess.cc.

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

◆ SetInelasticProcess()

void G4NeutronGeneralProcess::SetInelasticProcess ( G4HadronicProcess ptr)

Definition at line 115 of file G4NeutronGeneralProcess.cc.

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

◆ SetMinEnergyLimit()

void G4NeutronGeneralProcess::SetMinEnergyLimit ( G4double  val)
inline

Definition at line 250 of file G4NeutronGeneralProcess.hh.

251{
252 fMinEnergy = val;
253}

◆ SetTimeLimit()

void G4NeutronGeneralProcess::SetTimeLimit ( G4double  val)
inline

Definition at line 243 of file G4NeutronGeneralProcess.hh.

244{
245 fTimeLimit = val;
246}

◆ StartTracking()

void G4NeutronGeneralProcess::StartTracking ( G4Track )
overridevirtual

Reimplemented from G4HadronicProcess.

Definition at line 342 of file G4NeutronGeneralProcess.cc.

343{
345 fCurrMat = nullptr;
346}

◆ StorePhysicsTable()

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

Reimplemented from G4VProcess.

Definition at line 439 of file G4NeutronGeneralProcess.cc.

442{
443 G4bool yes = true;
444 if(!isMaster) { return yes; }
445 for(std::size_t i=0; i<nTables; ++i) {
446 G4String nam = (0==i || 3==i)
447 ? "LambdaNeutronGeneral" + nameT[i] : "ProbNeutronGeneral" + nameT[i];
448 G4String fnam = GetPhysicsTableFileName(part, directory, nam, ascii);
449 auto table = theHandler->Table(i);
450 if(nullptr == table || !table->StorePhysicsTable(fnam, ascii)) {
451 yes = false;
452 }
453 }
454 return yes;
455}
G4PhysicsTable * Table(std::size_t idx) const
const G4String & GetPhysicsTableFileName(const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
Definition: G4VProcess.cc:188

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