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

#include <G4QuasiElasticChannel.hh>

+ Inheritance diagram for G4QuasiElasticChannel:

Public Member Functions

 G4QuasiElasticChannel ()
 
 ~G4QuasiElasticChannel () override
 
G4double GetFraction (G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 
G4KineticTrackVectorScatter (G4Nucleus &theNucleus, const G4DynamicParticle &thePrimary)
 
 G4QuasiElasticChannel (const G4QuasiElasticChannel &)=delete
 
const G4QuasiElasticChanneloperator= (const G4QuasiElasticChannel &)=delete
 
G4bool operator== (const G4QuasiElasticChannel &) const =delete
 
G4bool operator!= (const G4QuasiElasticChannel &) const =delete
 
- Public Member Functions inherited from G4HadronicInteraction
 G4HadronicInteraction (const G4String &modelName="HadronicModel")
 
virtual ~G4HadronicInteraction ()
 
virtual G4HadFinalStateApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
virtual G4double SampleInvariantT (const G4ParticleDefinition *p, G4double plab, G4int Z, G4int A)
 
virtual G4bool IsApplicable (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)
 
G4double GetMinEnergy () const
 
G4double GetMinEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMinEnergy (G4double anEnergy)
 
void SetMinEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMinEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4double GetMaxEnergy () const
 
G4double GetMaxEnergy (const G4Material *aMaterial, const G4Element *anElement) const
 
void SetMaxEnergy (const G4double anEnergy)
 
void SetMaxEnergy (G4double anEnergy, const G4Element *anElement)
 
void SetMaxEnergy (G4double anEnergy, const G4Material *aMaterial)
 
G4int GetVerboseLevel () const
 
void SetVerboseLevel (G4int value)
 
const G4StringGetModelName () const
 
void DeActivateFor (const G4Material *aMaterial)
 
void ActivateFor (const G4Material *aMaterial)
 
void DeActivateFor (const G4Element *anElement)
 
void ActivateFor (const G4Element *anElement)
 
G4bool IsBlocked (const G4Material *aMaterial) const
 
G4bool IsBlocked (const G4Element *anElement) const
 
void SetRecoilEnergyThreshold (G4double val)
 
G4double GetRecoilEnergyThreshold () const
 
virtual const std::pair< G4double, G4doubleGetFatalEnergyCheckLevels () const
 
virtual std::pair< G4double, G4doubleGetEnergyMomentumCheckLevels () const
 
void SetEnergyMomentumCheckLevels (G4double relativeLevel, G4double absoluteLevel)
 
virtual void ModelDescription (std::ostream &outFile) const
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void InitialiseModel ()
 
 G4HadronicInteraction (const G4HadronicInteraction &right)=delete
 
const G4HadronicInteractionoperator= (const G4HadronicInteraction &right)=delete
 
G4bool operator== (const G4HadronicInteraction &right) const =delete
 
G4bool operator!= (const G4HadronicInteraction &right) const =delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4HadronicInteraction
void SetModelName (const G4String &nam)
 
G4bool IsBlocked () const
 
void Block ()
 
- Protected Attributes inherited from G4HadronicInteraction
G4HadFinalState theParticleChange
 
G4int verboseLevel
 
G4double theMinEnergy
 
G4double theMaxEnergy
 
G4bool isBlocked
 

Detailed Description

Definition at line 51 of file G4QuasiElasticChannel.hh.

Constructor & Destructor Documentation

◆ G4QuasiElasticChannel() [1/2]

G4QuasiElasticChannel::G4QuasiElasticChannel ( )
explicit

Definition at line 62 of file G4QuasiElasticChannel.cc.

63 : G4HadronicInteraction("QuasiElastic"),
64 theQuasiElastic(new G4QuasiElRatios),
65 the3DNucleus(new G4Fancy3DNucleus),
66 secID(-1) {
67 secID = G4PhysicsModelCatalog::GetModelID( "model_QuasiElastic" );
68}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4int GetModelID(const G4int modelIndex)

◆ ~G4QuasiElasticChannel()

G4QuasiElasticChannel::~G4QuasiElasticChannel ( )
override

Definition at line 70 of file G4QuasiElasticChannel.cc.

71{
72 delete the3DNucleus;
73 delete theQuasiElastic;
74}

◆ G4QuasiElasticChannel() [2/2]

G4QuasiElasticChannel::G4QuasiElasticChannel ( const G4QuasiElasticChannel & )
delete

Member Function Documentation

◆ GetFraction()

G4double G4QuasiElasticChannel::GetFraction ( G4Nucleus & theNucleus,
const G4DynamicParticle & thePrimary )

Definition at line 76 of file G4QuasiElasticChannel.cc.

78{
79 #ifdef debug_scatter
80 G4cout << "G4QuasiElasticChannel:: P=" << thePrimary.GetTotalMomentum()
81 << ", pPDG=" << thePrimary.GetDefinition()->GetPDGEncoding()
82 << ", Z = " << theNucleus.GetZ_asInt())
83 << ", N = " << theNucleus.GetN_asInt())
84 << ", A = " << theNucleus.GetA_asInt() << G4endl;
85 #endif
86
87 std::pair<G4double,G4double> ratios;
88 ratios=theQuasiElastic->GetRatios(thePrimary.GetTotalMomentum(),
89 thePrimary.GetDefinition()->GetPDGEncoding(),
90 theNucleus.GetZ_asInt(),
91 theNucleus.GetN_asInt());
92 #ifdef debug_scatter
93 G4cout << "G4QuasiElasticChannel::ratios " << ratios.first << " x " <<ratios.second
94 << " = " << ratios.first*ratios.second << G4endl;
95 #endif
96
97 return ratios.first*ratios.second;
98}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4ParticleDefinition * GetDefinition() const
G4double GetTotalMomentum() const
G4int GetA_asInt() const
Definition G4Nucleus.hh:99
G4int GetZ_asInt() const
Definition G4Nucleus.hh:105
G4int GetN_asInt() const
Definition G4Nucleus.hh:102
std::pair< G4double, G4double > GetRatios(G4double pIU, G4int prPDG, G4int tgZ, G4int tgN)

Referenced by G4TheoFSGenerator::ApplyYourself().

◆ operator!=()

G4bool G4QuasiElasticChannel::operator!= ( const G4QuasiElasticChannel & ) const
delete

◆ operator=()

const G4QuasiElasticChannel & G4QuasiElasticChannel::operator= ( const G4QuasiElasticChannel & )
delete

◆ operator==()

G4bool G4QuasiElasticChannel::operator== ( const G4QuasiElasticChannel & ) const
delete

◆ Scatter()

G4KineticTrackVector * G4QuasiElasticChannel::Scatter ( G4Nucleus & theNucleus,
const G4DynamicParticle & thePrimary )

Definition at line 100 of file G4QuasiElasticChannel.cc.

102{
103 G4int A=theNucleus.GetA_asInt();
104 G4int Z=theNucleus.GetZ_asInt();
105 // build Nucleus and choose random nucleon to scatter with
106 the3DNucleus->Init(theNucleus.GetA_asInt(),theNucleus.GetZ_asInt());
107 const std::vector<G4Nucleon>& nucleons=the3DNucleus->GetNucleons();
108 G4double targetNucleusMass=the3DNucleus->GetMass();
109 G4LorentzVector targetNucleus4Mom(0.,0.,0.,targetNucleusMass);
110 G4int index;
111 do {
112 index=G4lrint((A-1)*G4UniformRand());
113 } while (index < 0 || index >= static_cast<G4int>(nucleons.size())); /* Loop checking, 07.08.2015, A.Ribon */
114
115 const G4ParticleDefinition * pDef= nucleons[index].GetDefinition();
116
117 G4int resA=A - 1;
118 G4int resZ=Z - static_cast<int>(pDef->GetPDGCharge());
119 const G4ParticleDefinition* resDef;
120 G4double residualNucleusMass;
121 if(resZ)
122 {
123 resDef=G4ParticleTable::GetParticleTable()->GetIonTable()->GetIon(resZ,resA,0);
124 residualNucleusMass=resDef->GetPDGMass();
125 }
126 else {
127 resDef=G4Neutron::Neutron();
128 residualNucleusMass=resA * G4Neutron::Neutron()->GetPDGMass();
129 }
130 #ifdef debug_scatter
131 G4cout<<"G4QElChan::Scatter: neutron - proton? A ="<<A<<", Z="<<Z<<", projName="
132 <<pDef->GetParticleName()<<G4endl;
133 #endif
134
135 G4LorentzVector pNucleon=nucleons[index].Get4Momentum();
136 G4double residualNucleusEnergy=std::sqrt(sqr(residualNucleusMass) +
137 pNucleon.vect().mag2());
138 pNucleon.setE(targetNucleusMass-residualNucleusEnergy);
139 G4LorentzVector residualNucleus4Mom=targetNucleus4Mom-pNucleon;
140
141 std::pair<G4LorentzVector,G4LorentzVector> result;
142
143 result=theQuasiElastic->Scatter(pDef->GetPDGEncoding(),pNucleon,
144 thePrimary.GetDefinition()->GetPDGEncoding(),
145 thePrimary.Get4Momentum());
146 G4LorentzVector scatteredHadron4Mom;
147 if (result.first.e() > 0.)
148 scatteredHadron4Mom=result.second;
149 else { //scatter failed
150 //G4cout << "Warning - G4QuasiElasticChannel::Scatter no scattering" << G4endl;
151 //return 0; //no scatter
152 scatteredHadron4Mom=thePrimary.Get4Momentum();
153 residualNucleus4Mom=G4LorentzVector(0.,0.,0.,targetNucleusMass);
155 }
156
157#ifdef debug_scatter
158 G4LorentzVector EpConservation=pNucleon+thePrimary.Get4Momentum()
159 - result.first - result.second;
160 if ( (EpConservation.vect().mag2() > .01*MeV*MeV )
161 || (std::abs(EpConservation.e()) > 0.1 * MeV ) )
162 {
163 G4cout << "Warning - G4QuasiElasticChannel::Scatter E-p non conservation : "
164 << EpConservation << G4endl;
165 }
166#endif
167
169 G4KineticTrack * sPrim=new G4KineticTrack(thePrimary.GetDefinition(),
170 0.,G4ThreeVector(0), scatteredHadron4Mom);
171 sPrim->SetCreatorModelID( secID );
172 ktv->push_back(sPrim);
173 if (result.first.e() > 0.)
174 {
175 G4KineticTrack * sNuc=new G4KineticTrack(pDef, 0.,G4ThreeVector(0), result.first);
176 sNuc->SetCreatorModelID( secID );
177 ktv->push_back(sNuc);
178 }
179 if(resZ || resA==1) // For the only neutron or for tnuclei with Z>0
180 {
181 G4KineticTrack * rNuc=new G4KineticTrack(resDef,
182 0.,G4ThreeVector(0), residualNucleus4Mom);
183 rNuc->SetCreatorModelID( secID );
184 ktv->push_back(rNuc);
185 }
186 else // The residual nucleus consists of only neutrons
187 {
188 residualNucleus4Mom/=resA; // Split 4-mom of A*n system equally
189 for(G4int in=0; in<resA; in++) // Loop over neutrons in A*n system.
190 {
191 G4KineticTrack* rNuc=new G4KineticTrack(resDef,
192 0.,G4ThreeVector(0), residualNucleus4Mom);
193 rNuc->SetCreatorModelID( secID );
194 ktv->push_back(rNuc);
195 }
196 }
197#ifdef debug_scatter
198 G4cout<<"G4QElC::Scat: Nucleon: "<<result.first <<" mass "<<result.first.mag() << G4endl;
199 G4cout<<"G4QElC::Scat: Project: "<<result.second<<" mass "<<result.second.mag()<< G4endl;
200#endif
201 return ktv;
202}
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
const G4double A[17]
#define G4UniformRand()
Definition Randomize.hh:52
double mag2() const
Hep3Vector vect() const
G4LorentzVector Get4Momentum() const
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
void SetCreatorModelID(G4int id)
static G4Neutron * Neutron()
Definition G4Neutron.cc:101
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
std::pair< G4LorentzVector, G4LorentzVector > Scatter(G4int NPDG, G4LorentzVector N4M, G4int pPDG, G4LorentzVector p4M)
virtual G4double GetMass()=0
virtual void Init(G4int theA, G4int theZ, G4int numberOfLambdas=0)=0
virtual const std::vector< G4Nucleon > & GetNucleons()=0
int G4lrint(double ad)
Definition templates.hh:134
T sqr(const T &x)
Definition templates.hh:128

Referenced by G4TheoFSGenerator::ApplyYourself().


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