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

#include <G4ChargeExchangeProcess.hh>

+ Inheritance diagram for G4ChargeExchangeProcess:

Public Member Functions

 G4ChargeExchangeProcess (const G4String &procName="chargeExchange")
 
virtual ~G4ChargeExchangeProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &aParticleType)
 
virtual G4double GetElementCrossSection (const G4DynamicParticle *aParticle, const G4Element *anElement, const G4Material *mat=0)
 
- 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 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 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 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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- 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 ()
 
- 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 52 of file G4ChargeExchangeProcess.hh.

Constructor & Destructor Documentation

◆ G4ChargeExchangeProcess()

G4ChargeExchangeProcess::G4ChargeExchangeProcess ( const G4String & procName = "chargeExchange")

Definition at line 58 of file G4ChargeExchangeProcess.cc.

59 : G4HadronicProcess(procName,fChargeExchange), first(true)
60{
61 thEnergy = 20.*MeV;
62 pPDG = 0;
63 verboseLevel= 1;
65 theProton = G4Proton::Proton();
66 theNeutron = G4Neutron::Neutron();
67 theAProton = G4AntiProton::AntiProton();
68 theANeutron = G4AntiNeutron::AntiNeutron();
69 thePiPlus = G4PionPlus::PionPlus();
70 thePiMinus = G4PionMinus::PionMinus();
71 thePiZero = G4PionZero::PionZero();
72 theKPlus = G4KaonPlus::KaonPlus();
73 theKMinus = G4KaonMinus::KaonMinus();
76 theL = G4Lambda::Lambda();
77 theAntiL = G4AntiLambda::AntiLambda();
78 theSPlus = G4SigmaPlus::SigmaPlus();
80 theSMinus = G4SigmaMinus::SigmaMinus();
82 theS0 = G4SigmaZero::SigmaZero();
84 theXiMinus = G4XiMinus::XiMinus();
85 theXi0 = G4XiZero::XiZero();
86 theAXiMinus = G4AntiXiMinus::AntiXiMinus();
87 theAXi0 = G4AntiXiZero::AntiXiZero();
88 theOmega = G4OmegaMinus::OmegaMinus();
90 theD = G4Deuteron::Deuteron();
91 theT = G4Triton::Triton();
92 theA = G4Alpha::Alpha();
93 theHe3 = G4He3::He3();
94}
@ fChargeExchange
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4AntiLambda * AntiLambda()
static G4AntiNeutron * AntiNeutron()
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4AntiProton * AntiProton()
static G4AntiSigmaMinus * AntiSigmaMinus()
static G4AntiSigmaPlus * AntiSigmaPlus()
static G4AntiSigmaZero * AntiSigmaZero()
static G4AntiXiMinus * AntiXiMinus()
static G4AntiXiZero * AntiXiZero()
static G4Deuteron * Deuteron()
Definition G4Deuteron.cc:90
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
static G4He3 * He3()
Definition G4He3.cc:90
static G4KaonMinus * KaonMinus()
static G4KaonPlus * KaonPlus()
static G4KaonZeroLong * KaonZeroLong()
static G4KaonZeroShort * KaonZeroShort()
static G4Lambda * Lambda()
Definition G4Lambda.cc:105
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
static G4OmegaMinus * OmegaMinus()
static G4PionMinus * PionMinus()
static G4PionPlus * PionPlus()
Definition G4PionPlus.cc:93
static G4PionZero * PionZero()
static G4Proton * Proton()
Definition G4Proton.cc:90
static G4SigmaMinus * SigmaMinus()
static G4SigmaPlus * SigmaPlus()
static G4SigmaZero * SigmaZero()
static G4Triton * Triton()
Definition G4Triton.cc:90
G4int verboseLevel
static G4XiMinus * XiMinus()
Definition G4XiMinus.cc:103
static G4XiZero * XiZero()
Definition G4XiZero.cc:103

◆ ~G4ChargeExchangeProcess()

G4ChargeExchangeProcess::~G4ChargeExchangeProcess ( )
virtual

Definition at line 96 of file G4ChargeExchangeProcess.cc.

97{
98 if (factors) delete factors;
99}

Member Function Documentation

◆ BuildPhysicsTable()

void G4ChargeExchangeProcess::BuildPhysicsTable ( const G4ParticleDefinition & aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 101 of file G4ChargeExchangeProcess.cc.

103{
104 if(first) {
105 first = false;
106 theParticle = &aParticleType;
107 pPDG = theParticle->GetPDGEncoding();
108
110
111 const size_t n = 10;
112 if(theParticle == thePiPlus || theParticle == thePiMinus ||
113 theParticle == theKPlus || theParticle == theKMinus ||
114 theParticle == theK0S || theParticle == theK0L) {
115
116 G4double F[n] = {0.33,0.27,0.29,0.31,0.27,0.18,0.13,0.1,0.09,0.07};
117 factors = new G4PhysicsLinearVector(0.0,2.0*GeV,n);
118 for(size_t i=0; i<n; i++) {factors->PutValue(i,F[i]);}
119
120 } else {
121
122 G4double F[n] = {0.50,0.45,0.40,0.35,0.30,0.25,0.06,0.04,0.005,0.0};
123 factors = new G4PhysicsLinearVector(0.0,4.0*GeV,n);
124 for(size_t i=0; i<n; i++) {factors->PutValue(i,F[i]);}
125 }
126
127 if(verboseLevel>1)
128 G4cout << "G4ChargeExchangeProcess for "
129 << theParticle->GetParticleName()
130 << G4endl;
131 }
133}
double G4double
Definition G4Types.hh:83
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void BuildPhysicsTable(const G4ParticleDefinition &) override
G4CrossSectionDataStore * GetCrossSectionDataStore()
const G4String & GetParticleName() const
void PutValue(const std::size_t index, const G4double value)

◆ DumpPhysicsTable()

void G4ChargeExchangeProcess::DumpPhysicsTable ( const G4ParticleDefinition & aParticleType)
virtual

Definition at line 194 of file G4ChargeExchangeProcess.cc.

196{
197 store->DumpPhysicsTable(aParticleType);
198}
void DumpPhysicsTable(const G4ParticleDefinition &)

◆ GetElementCrossSection()

G4double G4ChargeExchangeProcess::GetElementCrossSection ( const G4DynamicParticle * aParticle,
const G4Element * anElement,
const G4Material * mat = 0 )
virtual

Definition at line 135 of file G4ChargeExchangeProcess.cc.

139{
140 // gives the microscopic cross section in GEANT4 internal units
141 G4double Z = elm->GetZ();
142 G4int iz = G4int(Z);
143 G4double x = 0.0;
144
145 // The process is effective only above the threshold
146 if(iz == 1 || dp->GetKineticEnergy() < thEnergy) return x;
147
148 if(verboseLevel>1)
149 G4cout << "G4ChargeExchangeProcess compute GHAD CS for element "
150 << elm->GetName()
151 << G4endl;
152 x = store->GetCrossSection(dp, elm, mat);
153
154 if(verboseLevel>1)
155 G4cout << "G4ChargeExchangeProcess cross(mb)= " << x/millibarn
156 << " E(MeV)= " << dp->GetKineticEnergy()
157 << " " << theParticle->GetParticleName()
158 << " in Z= " << iz
159 << G4endl;
160 G4bool b;
161 G4double A = elm->GetN();
162 G4double ptot = dp->GetTotalMomentum();
163 x *= factors->GetValue(ptot, b)/G4Pow::GetInstance()->powA(A, 0.42);
164 if(theParticle == thePiPlus || theParticle == theProton ||
165 theParticle == theKPlus || theParticle == theANeutron)
166 { x *= (1.0 - Z/A); }
167
168 else if(theParticle == thePiMinus || theParticle == theNeutron ||
169 theParticle == theKMinus || theParticle == theAProton)
170 { x *= Z/A; }
171
172 if(theParticle->GetPDGMass() < GeV) {
173 if(ptot > 2.*GeV) x *= 4.0*GeV*GeV/(ptot*ptot);
174 }
175
176 if(verboseLevel>1)
177 G4cout << "Corrected cross(mb)= " << x/millibarn << G4endl;
178
179 return x;
180}
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
const G4double A[17]
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
G4double GetValue(const G4double energy, G4bool &isOutRange) const
static G4Pow * GetInstance()
Definition G4Pow.cc:41
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230

◆ IsApplicable()

G4bool G4ChargeExchangeProcess::IsApplicable ( const G4ParticleDefinition & aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 182 of file G4ChargeExchangeProcess.cc.

184{
185 const G4ParticleDefinition* p = &aParticleType;
186 return (p == thePiPlus || p == thePiMinus ||
187 p == theProton || p == theNeutron ||
188 p == theAProton|| p == theANeutron||
189 p == theKPlus || p == theKMinus ||
190 p == theK0S || p == theK0L ||
191 p == theL);
192}

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