Geant4  10.02.p03
G4ChargeExchangeProcess Class Reference

#include <G4ChargeExchangeProcess.hh>

Inheritance diagram for G4ChargeExchangeProcess:
Collaboration 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)
 
virtual ~G4HadronicProcess ()
 
void RegisterMe (G4HadronicInteraction *a)
 
G4double GetElementCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
G4double GetMicroscopicCrossSection (const G4DynamicParticle *part, const G4Element *elm, const G4Material *mat=0)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
void DumpPhysicsTable (const G4ParticleDefinition &p)
 
void AddDataSet (G4VCrossSectionDataSet *aDataSet)
 
std::vector< G4HadronicInteraction * > & GetHadronicInteractionList ()
 
G4double GetMeanFreePath (const G4Track &aTrack, G4double, G4ForceCondition *)
 
const G4NucleusGetTargetNucleus () const
 
const G4IsotopeGetTargetIsotope ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
void BiasCrossSectionByFactor (G4double aScale)
 
void SetEpReportLevel (G4int level)
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
G4CrossSectionDataStoreGetCrossSectionDataStore ()
 
void MultiplyCrossSectionBy (G4double factor)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChange * AlongStepDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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 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
 
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 &)
 

Private Attributes

const G4ParticleDefinitiontheParticle
 
const G4ParticleDefinitiontheProton
 
const G4ParticleDefinitiontheNeutron
 
const G4ParticleDefinitiontheAProton
 
const G4ParticleDefinitiontheANeutron
 
const G4ParticleDefinitionthePiPlus
 
const G4ParticleDefinitionthePiMinus
 
const G4ParticleDefinitionthePiZero
 
const G4ParticleDefinitiontheKPlus
 
const G4ParticleDefinitiontheKMinus
 
const G4ParticleDefinitiontheK0S
 
const G4ParticleDefinitiontheK0L
 
const G4ParticleDefinitiontheL
 
const G4ParticleDefinitiontheAntiL
 
const G4ParticleDefinitiontheSPlus
 
const G4ParticleDefinitiontheASPlus
 
const G4ParticleDefinitiontheSMinus
 
const G4ParticleDefinitiontheASMinus
 
const G4ParticleDefinitiontheS0
 
const G4ParticleDefinitiontheAS0
 
const G4ParticleDefinitiontheXiMinus
 
const G4ParticleDefinitiontheXi0
 
const G4ParticleDefinitiontheAXiMinus
 
const G4ParticleDefinitiontheAXi0
 
const G4ParticleDefinitiontheOmega
 
const G4ParticleDefinitiontheAOmega
 
const G4ParticleDefinitiontheD
 
const G4ParticleDefinitiontheT
 
const G4ParticleDefinitiontheA
 
const G4ParticleDefinitiontheHe3
 
G4CrossSectionDataStorestore
 
G4PhysicsLinearVectorfactors
 
G4double thEnergy
 
G4int pPDG
 
G4bool first
 

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, G4Material *aMaterial, G4Element *anElement)
 
G4NucleusGetTargetNucleusPointer ()
 
void DumpState (const G4Track &, const G4String &, G4ExceptionDescription &)
 
G4HadronicInteractionGetHadronicInteraction () const
 
G4double GetLastCrossSection ()
 
void FillResult (G4HadFinalState *aR, const G4Track &aT)
 
G4HadFinalStateCheckResult (const G4HadProjectile &thePro, const G4Nucleus &targetNucleus, G4HadFinalState *result)
 
void CheckEnergyMomentumConservation (const G4Track &, const G4Nucleus &)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4HadronicProcess
G4HadProjectile thePro
 
G4ParticleChange * theTotalResult
 
G4int epReportLevel
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 

Detailed Description

Definition at line 53 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;
92  theA = G4Alpha::Alpha();
93  theHe3 = G4He3::He3();
94 }
const G4ParticleDefinition * theSMinus
static const double MeV
Definition: G4SIunits.hh:211
G4int verboseLevel
Definition: G4VProcess.hh:368
static G4AntiOmegaMinus * AntiOmegaMinus()
static G4OmegaMinus * OmegaMinus()
const G4ParticleDefinition * thePiZero
static G4KaonZeroLong * KaonZeroLong()
const G4ParticleDefinition * theHe3
const G4ParticleDefinition * theD
const G4ParticleDefinition * thePiPlus
const G4ParticleDefinition * theOmega
static G4AntiSigmaPlus * AntiSigmaPlus()
const G4ParticleDefinition * theAS0
const G4ParticleDefinition * theASMinus
const G4ParticleDefinition * theL
static G4SigmaZero * SigmaZero()
Definition: G4SigmaZero.cc:102
static G4KaonMinus * KaonMinus()
Definition: G4KaonMinus.cc:113
const G4ParticleDefinition * theAntiL
static G4AntiSigmaMinus * AntiSigmaMinus()
const G4ParticleDefinition * theSPlus
static G4XiZero * XiZero()
Definition: G4XiZero.cc:106
void AddDataSet(G4VCrossSectionDataSet *aDataSet)
const G4ParticleDefinition * theAXi0
static G4KaonZeroShort * KaonZeroShort()
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
const G4ParticleDefinition * theXi0
const G4ParticleDefinition * theAProton
const G4ParticleDefinition * theASPlus
static G4XiMinus * XiMinus()
Definition: G4XiMinus.cc:106
static G4AntiXiMinus * AntiXiMinus()
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4PionPlus * PionPlus()
Definition: G4PionPlus.cc:98
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4PionZero * PionZero()
Definition: G4PionZero.cc:108
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
const G4ParticleDefinition * thePiMinus
static G4SigmaMinus * SigmaMinus()
const G4ParticleDefinition * theKMinus
const G4ParticleDefinition * theKPlus
const G4ParticleDefinition * theK0S
static G4AntiLambda * AntiLambda()
const G4ParticleDefinition * theXiMinus
static G4PionMinus * PionMinus()
Definition: G4PionMinus.cc:98
static G4AntiSigmaZero * AntiSigmaZero()
const G4ParticleDefinition * theAOmega
const G4ParticleDefinition * theProton
static G4AntiXiZero * AntiXiZero()
G4HadronicProcess(const G4String &processName="Hadronic", G4ProcessType procType=fHadronic)
const G4ParticleDefinition * theNeutron
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
const G4ParticleDefinition * theA
static G4SigmaPlus * SigmaPlus()
Definition: G4SigmaPlus.cc:108
static G4Lambda * Lambda()
Definition: G4Lambda.cc:108
static G4KaonPlus * KaonPlus()
Definition: G4KaonPlus.cc:113
const G4ParticleDefinition * theAXiMinus
static G4He3 * He3()
Definition: G4He3.cc:94
static G4AntiNeutron * AntiNeutron()
const G4ParticleDefinition * theK0L
const G4ParticleDefinition * theANeutron
const G4ParticleDefinition * theT
const G4ParticleDefinition * theS0
Here is the call graph for this function:

◆ ~G4ChargeExchangeProcess()

G4ChargeExchangeProcess::~G4ChargeExchangeProcess ( )
virtual

Definition at line 96 of file G4ChargeExchangeProcess.cc.

97 {
98  if (factors) delete factors;
99 }
G4PhysicsLinearVector * factors
Here is the call graph for this function:

Member Function Documentation

◆ BuildPhysicsTable()

void G4ChargeExchangeProcess::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4HadronicProcess.

Definition at line 102 of file G4ChargeExchangeProcess.cc.

103 {
104  if(first) {
105  first = false;
106  theParticle = &aParticleType;
108 
110 
111  const size_t n = 10;
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  //factors->SetSpline(true);
127 
128  if(verboseLevel>1)
129  G4cout << "G4ChargeExchangeProcess for "
131  << G4endl;
132  }
134 }
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4ParticleDefinition * thePiPlus
const G4ParticleDefinition * theParticle
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
Char_t n[5]
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
G4CrossSectionDataStore * GetCrossSectionDataStore()
void PutValue(size_t index, G4double theValue)
static const double GeV
Definition: G4SIunits.hh:214
const G4ParticleDefinition * thePiMinus
const G4ParticleDefinition * theKMinus
const G4ParticleDefinition * theKPlus
const G4ParticleDefinition * theK0S
#define G4endl
Definition: G4ios.hh:61
G4CrossSectionDataStore * store
double G4double
Definition: G4Types.hh:76
G4PhysicsLinearVector * factors
const G4ParticleDefinition * theK0L
Here is the call graph for this function:
Here is the caller graph for this function:

◆ DumpPhysicsTable()

void G4ChargeExchangeProcess::DumpPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Definition at line 196 of file G4ChargeExchangeProcess.cc.

197 {
198  store->DumpPhysicsTable(aParticleType);
199 }
void DumpPhysicsTable(const G4ParticleDefinition &)
G4CrossSectionDataStore * store
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetElementCrossSection()

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

Definition at line 136 of file G4ChargeExchangeProcess.cc.

140 {
141  // gives the microscopic cross section in GEANT4 internal units
142  G4double Z = elm->GetZ();
143  G4int iz = G4int(Z);
144  G4double x = 0.0;
145 
146  // The process is effective only above the threshold
147  if(iz == 1 || dp->GetKineticEnergy() < thEnergy) return x;
148 
149  if(verboseLevel>1)
150  G4cout << "G4ChargeExchangeProcess compute GHAD CS for element "
151  << elm->GetName()
152  << G4endl;
153  x = store->GetCrossSection(dp, elm, mat);
154 
155  if(verboseLevel>1)
156  G4cout << "G4ChargeExchangeProcess cross(mb)= " << x/millibarn
157  << " E(MeV)= " << dp->GetKineticEnergy()
158  << " " << theParticle->GetParticleName()
159  << " in Z= " << iz
160  << G4endl;
161  G4bool b;
162  G4double A = elm->GetN();
163  G4double ptot = dp->GetTotalMomentum();
164  x *= factors->GetValue(ptot, b)/G4Pow::GetInstance()->powA(A, 0.42);
167  { x *= (1.0 - Z/A); }
168 
169  else if(theParticle == thePiMinus || theParticle == theNeutron ||
171  { x *= Z/A; }
172 
173  if(theParticle->GetPDGMass() < GeV) {
174  if(ptot > 2.*GeV) x *= 4.0*GeV*GeV/(ptot*ptot);
175  }
176 
177  if(verboseLevel>1)
178  G4cout << "Corrected cross(mb)= " << x/millibarn << G4endl;
179 
180  return x;
181 }
static G4Pow * GetInstance()
Definition: G4Pow.cc:55
G4int verboseLevel
Definition: G4VProcess.hh:368
const G4ParticleDefinition * thePiPlus
const G4ParticleDefinition * theParticle
int G4int
Definition: G4Types.hh:78
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
double A(double temperature)
Float_t Z
const G4ParticleDefinition * theAProton
bool G4bool
Definition: G4Types.hh:79
G4double iz
Definition: TRTMaterials.hh:39
G4double GetCrossSection(const G4DynamicParticle *, const G4Material *)
static const double GeV
Definition: G4SIunits.hh:214
const G4ParticleDefinition * thePiMinus
const G4ParticleDefinition * theKMinus
G4double GetValue(G4double theEnergy, G4bool &isOutRange) const
const G4ParticleDefinition * theKPlus
static const double millibarn
Definition: G4SIunits.hh:105
const G4ParticleDefinition * theProton
#define G4endl
Definition: G4ios.hh:61
G4CrossSectionDataStore * store
const G4ParticleDefinition * theNeutron
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:259
double G4double
Definition: G4Types.hh:76
G4PhysicsLinearVector * factors
const G4ParticleDefinition * theANeutron
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4ChargeExchangeProcess::IsApplicable ( const G4ParticleDefinition aParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 184 of file G4ChargeExchangeProcess.cc.

185 {
186  const G4ParticleDefinition* p = &aParticleType;
187  return (p == thePiPlus || p == thePiMinus ||
188  p == theProton || p == theNeutron ||
189  p == theAProton|| p == theANeutron||
190  p == theKPlus || p == theKMinus ||
191  p == theK0S || p == theK0L ||
192  p == theL);
193 }
const G4ParticleDefinition * thePiPlus
const G4ParticleDefinition * theL
const G4ParticleDefinition * theAProton
const G4ParticleDefinition * thePiMinus
const G4ParticleDefinition * theKMinus
const G4ParticleDefinition * theKPlus
const G4ParticleDefinition * theK0S
const G4ParticleDefinition * theProton
const G4ParticleDefinition * theNeutron
const G4ParticleDefinition * theK0L
const G4ParticleDefinition * theANeutron
Here is the call graph for this function:
Here is the caller graph for this function:

Member Data Documentation

◆ factors

G4PhysicsLinearVector* G4ChargeExchangeProcess::factors
private

Definition at line 106 of file G4ChargeExchangeProcess.hh.

◆ first

G4bool G4ChargeExchangeProcess::first
private

Definition at line 111 of file G4ChargeExchangeProcess.hh.

◆ pPDG

G4int G4ChargeExchangeProcess::pPDG
private

Definition at line 110 of file G4ChargeExchangeProcess.hh.

◆ store

G4CrossSectionDataStore* G4ChargeExchangeProcess::store
private

Definition at line 105 of file G4ChargeExchangeProcess.hh.

◆ theA

const G4ParticleDefinition* G4ChargeExchangeProcess::theA
private

Definition at line 102 of file G4ChargeExchangeProcess.hh.

◆ theANeutron

const G4ParticleDefinition* G4ChargeExchangeProcess::theANeutron
private

Definition at line 78 of file G4ChargeExchangeProcess.hh.

◆ theAntiL

const G4ParticleDefinition* G4ChargeExchangeProcess::theAntiL
private

Definition at line 87 of file G4ChargeExchangeProcess.hh.

◆ theAOmega

const G4ParticleDefinition* G4ChargeExchangeProcess::theAOmega
private

Definition at line 99 of file G4ChargeExchangeProcess.hh.

◆ theAProton

const G4ParticleDefinition* G4ChargeExchangeProcess::theAProton
private

Definition at line 77 of file G4ChargeExchangeProcess.hh.

◆ theAS0

const G4ParticleDefinition* G4ChargeExchangeProcess::theAS0
private

Definition at line 93 of file G4ChargeExchangeProcess.hh.

◆ theASMinus

const G4ParticleDefinition* G4ChargeExchangeProcess::theASMinus
private

Definition at line 91 of file G4ChargeExchangeProcess.hh.

◆ theASPlus

const G4ParticleDefinition* G4ChargeExchangeProcess::theASPlus
private

Definition at line 89 of file G4ChargeExchangeProcess.hh.

◆ theAXi0

const G4ParticleDefinition* G4ChargeExchangeProcess::theAXi0
private

Definition at line 97 of file G4ChargeExchangeProcess.hh.

◆ theAXiMinus

const G4ParticleDefinition* G4ChargeExchangeProcess::theAXiMinus
private

Definition at line 96 of file G4ChargeExchangeProcess.hh.

◆ theD

const G4ParticleDefinition* G4ChargeExchangeProcess::theD
private

Definition at line 100 of file G4ChargeExchangeProcess.hh.

◆ theHe3

const G4ParticleDefinition* G4ChargeExchangeProcess::theHe3
private

Definition at line 103 of file G4ChargeExchangeProcess.hh.

◆ theK0L

const G4ParticleDefinition* G4ChargeExchangeProcess::theK0L
private

Definition at line 85 of file G4ChargeExchangeProcess.hh.

◆ theK0S

const G4ParticleDefinition* G4ChargeExchangeProcess::theK0S
private

Definition at line 84 of file G4ChargeExchangeProcess.hh.

◆ theKMinus

const G4ParticleDefinition* G4ChargeExchangeProcess::theKMinus
private

Definition at line 83 of file G4ChargeExchangeProcess.hh.

◆ theKPlus

const G4ParticleDefinition* G4ChargeExchangeProcess::theKPlus
private

Definition at line 82 of file G4ChargeExchangeProcess.hh.

◆ theL

const G4ParticleDefinition* G4ChargeExchangeProcess::theL
private

Definition at line 86 of file G4ChargeExchangeProcess.hh.

◆ thEnergy

G4double G4ChargeExchangeProcess::thEnergy
private

Definition at line 108 of file G4ChargeExchangeProcess.hh.

◆ theNeutron

const G4ParticleDefinition* G4ChargeExchangeProcess::theNeutron
private

Definition at line 76 of file G4ChargeExchangeProcess.hh.

◆ theOmega

const G4ParticleDefinition* G4ChargeExchangeProcess::theOmega
private

Definition at line 98 of file G4ChargeExchangeProcess.hh.

◆ theParticle

const G4ParticleDefinition* G4ChargeExchangeProcess::theParticle
private

Definition at line 73 of file G4ChargeExchangeProcess.hh.

◆ thePiMinus

const G4ParticleDefinition* G4ChargeExchangeProcess::thePiMinus
private

Definition at line 80 of file G4ChargeExchangeProcess.hh.

◆ thePiPlus

const G4ParticleDefinition* G4ChargeExchangeProcess::thePiPlus
private

Definition at line 79 of file G4ChargeExchangeProcess.hh.

◆ thePiZero

const G4ParticleDefinition* G4ChargeExchangeProcess::thePiZero
private

Definition at line 81 of file G4ChargeExchangeProcess.hh.

◆ theProton

const G4ParticleDefinition* G4ChargeExchangeProcess::theProton
private

Definition at line 75 of file G4ChargeExchangeProcess.hh.

◆ theS0

const G4ParticleDefinition* G4ChargeExchangeProcess::theS0
private

Definition at line 92 of file G4ChargeExchangeProcess.hh.

◆ theSMinus

const G4ParticleDefinition* G4ChargeExchangeProcess::theSMinus
private

Definition at line 90 of file G4ChargeExchangeProcess.hh.

◆ theSPlus

const G4ParticleDefinition* G4ChargeExchangeProcess::theSPlus
private

Definition at line 88 of file G4ChargeExchangeProcess.hh.

◆ theT

const G4ParticleDefinition* G4ChargeExchangeProcess::theT
private

Definition at line 101 of file G4ChargeExchangeProcess.hh.

◆ theXi0

const G4ParticleDefinition* G4ChargeExchangeProcess::theXi0
private

Definition at line 95 of file G4ChargeExchangeProcess.hh.

◆ theXiMinus

const G4ParticleDefinition* G4ChargeExchangeProcess::theXiMinus
private

Definition at line 94 of file G4ChargeExchangeProcess.hh.


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