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

#include <G4CRCoalescence.hh>

+ Inheritance diagram for G4CRCoalescence:

Public Member Functions

 G4CRCoalescence ()
 
 ~G4CRCoalescence () override
 
 G4CRCoalescence (const G4CRCoalescence &right)=delete
 
const G4CRCoalescenceoperator= (const G4CRCoalescence &right)=delete
 
G4bool operator== (const G4CRCoalescence &right) const =delete
 
G4bool operator!= (const G4CRCoalescence &right) const =delete
 
void SetP0Coalescence (const G4HadProjectile &thePrimary, G4String)
 
void GenerateDeuterons (G4ReactionProductVector *result)
 
- 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 77 of file G4CRCoalescence.hh.

Constructor & Destructor Documentation

◆ G4CRCoalescence() [1/2]

G4CRCoalescence::G4CRCoalescence ( )
explicit

Definition at line 81 of file G4CRCoalescence.cc.

81 : G4HadronicInteraction("G4CRCoalescence" ),
82 fP0_d( 0.0 ), fP0_dbar( 0.0 ) {}

◆ ~G4CRCoalescence()

G4CRCoalescence::~G4CRCoalescence ( )
override

Definition at line 85 of file G4CRCoalescence.cc.

85{}

◆ G4CRCoalescence() [2/2]

G4CRCoalescence::G4CRCoalescence ( const G4CRCoalescence right)
delete

Member Function Documentation

◆ GenerateDeuterons()

void G4CRCoalescence::GenerateDeuterons ( G4ReactionProductVector result)

Definition at line 110 of file G4CRCoalescence.cc.

110 {
111 // Deuteron clusters are made with the first nucleon pair that fulfills
112 // the coalescence conditions, starting with the protons.
113 // A deuteron is a pair (i,j) where i is the proton and j the neutron in current event
114 // with the relative momentum less than p0 (i.e. within a sphere of radius p0).
115 // The same applies for antideuteron clusters, with antiprotons and antineutrons,
116 // instead of protons and neutrons, respectively.
117
118 // Vectors of index-position and 3-momentum pairs for, respectively:
119 // protons, neutrons, antiprotons and antineutrons
120 std::vector< std::pair< G4int, G4ThreeVector > > proton;
121 std::vector< std::pair< G4int, G4ThreeVector > > neutron;
122 std::vector< std::pair< G4int, G4ThreeVector > > antiproton;
123 std::vector< std::pair< G4int, G4ThreeVector > > antineutron;
124 for ( unsigned int i = 0; i < result->size(); ++i ) {
125 G4int pdgid = result->operator[](i)->GetDefinition()->GetPDGEncoding();
126 if ( pdgid == 2212 ) { // proton
127 proton.push_back( std::make_pair( i, result->operator[](i)->GetMomentum() ) );
128 result->erase( result->begin() + i );
129 }
130 }
131 for ( unsigned int i = 0; i < result->size(); ++i ) {
132 G4int pdgid = result->operator[](i)->GetDefinition()->GetPDGEncoding();
133 if ( pdgid == 2112 ) { // neutron
134 neutron.push_back( std::make_pair( i, result->operator[](i)->GetMomentum() ) );
135 result->erase( result->begin() + i );
136 }
137 }
138 for ( unsigned int i = 0; i < result->size(); ++i ) {
139 G4int pdgid = result->operator[](i)->GetDefinition()->GetPDGEncoding();
140 if ( pdgid == -2212 ) { // antiproton
141 antiproton.push_back( std::make_pair( i, result->operator[](i)->GetMomentum() ) );
142 result->erase( result->begin() + i );
143 }
144 }
145 for ( unsigned int i = 0; i < result->size(); ++i ) {
146 G4int pdgid = result->operator[](i)->GetDefinition()->GetPDGEncoding();
147 if ( pdgid == -2112 ) { // antineutron
148 antineutron.push_back( std::make_pair( i, result->operator[](i)->GetMomentum() ) );
149 result->erase( result->begin() + i );
150 }
151 }
152
153 for ( unsigned int i = 0; i < proton.size(); ++i ) { // loop over protons
154 if ( proton.at(i).first == -1 ) continue;
155 G4ThreeVector p1 = proton.at(i).second;
156 int partner1 = FindPartner( p1, G4Proton::Proton()->GetPDGMass(), neutron,
157 G4Neutron::Neutron()->GetPDGMass(), 1 );
158 if ( partner1 == -1 ) { // if no partner found, then the proton is a final-state secondary
161 finalp->SetDefinition( prt );
162 G4double massp = prt->GetPDGMass();
163 G4double totalEnergy = std::sqrt( p1.mag()*p1.mag() + massp*massp );
164 finalp->SetMomentum( p1 );
165 finalp->SetTotalEnergy( totalEnergy );
166 finalp->SetMass( massp );
167 result->push_back( finalp );
168 continue;
169 }
170 G4ThreeVector p2 = neutron.at(partner1).second;
171 PushDeuteron( p1, p2, 1, result );
172 neutron.at(partner1).first = -1; // tag the bound neutron
173 }
174
175 for ( unsigned int i = 0; i < neutron.size(); ++i ) { // loop over neutrons
176 if ( neutron.at(i).first == -1 ) continue; // Skip already bound neutron, else it is a final-state secondary
179 finaln->SetDefinition( nrt );
180 G4ThreeVector p2 = neutron.at(i).second;
181 G4double massn = nrt->GetPDGMass();
182 G4double totalEnergy = std::sqrt( p2.mag()*p2.mag() + massn*massn );
183 finaln->SetMomentum( p2 );
184 finaln->SetTotalEnergy( totalEnergy );
185 finaln->SetMass( massn );
186 result->push_back( finaln );
187 }
188
189 for ( unsigned int i = 0; i < antiproton.size(); ++i ) { // loop over antiprotons
190 if ( antiproton.at(i).first == -1 ) continue;
191 G4ThreeVector p1 = antiproton.at(i).second;
192 int partner1 = FindPartner( p1, G4Proton::Proton()->GetPDGMass(), antineutron,
193 G4Neutron::Neutron()->GetPDGMass(), -1 );
194 if ( partner1 == -1 ) { // if no partner found, then the antiproton is a final-state secondary
196 G4ReactionProduct* finalpbar = new G4ReactionProduct;
197 finalpbar->SetDefinition( pbar );
198 G4double massp = pbar->GetPDGMass();
199 G4double totalEnergy = std::sqrt( p1.mag()*p1.mag() + massp*massp );
200 finalpbar->SetMomentum( p1 );
201 finalpbar->SetTotalEnergy( totalEnergy );
202 finalpbar->SetMass( massp );
203 result->push_back( finalpbar );
204 continue;
205 }
206 G4ThreeVector p2 = antineutron.at(partner1).second;
207 PushDeuteron( p1, p2, -1, result );
208 antineutron.at(partner1).first = -1; // tag the bound antineutron
209 }
210
211 for ( unsigned int i = 0; i < antineutron.size(); ++i ) { // loop over antineutrons
212 if ( antineutron.at(i).first == -1 ) continue; // Skip already bound antineutron, else it is a final-state secondary
214 G4ReactionProduct* finalnbar = new G4ReactionProduct;
215 finalnbar->SetDefinition( nbar );
216 G4ThreeVector p2 = antineutron.at(i).second;
217 G4double massn = nbar->GetPDGMass();
218 G4double totalEnergy = std::sqrt( p2.mag()*p2.mag() + massn*massn );
219 finalnbar->SetMomentum( p2 );
220 finalnbar->SetTotalEnergy( totalEnergy );
221 finalnbar->SetMass( massn );
222 result->push_back( finalnbar );
223 }
224}
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
double mag() const
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
static G4ParticleTable * GetParticleTable()
G4ParticleDefinition * FindAntiParticle(G4int PDGEncoding)
static G4Proton * Proton()
Definition: G4Proton.cc:92
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetMass(const G4double mas)

Referenced by G4TheoFSGenerator::ApplyYourself().

◆ operator!=()

G4bool G4CRCoalescence::operator!= ( const G4CRCoalescence right) const
delete

◆ operator=()

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

◆ operator==()

G4bool G4CRCoalescence::operator== ( const G4CRCoalescence right) const
delete

◆ SetP0Coalescence()

void G4CRCoalescence::SetP0Coalescence ( const G4HadProjectile thePrimary,
G4String   
)

Definition at line 88 of file G4CRCoalescence.cc.

88 {
89 // Note by A.R. : in the present version, the coalescence parameters are set only for
90 // proton projectile. If we want to extend this coalescence algorithm
91 // for other applications, besides cosmic rays, we need to set these
92 // coalescence parameters also for all projectiles.
93 // (Note that the method "GenerateDeuterons", instead, can be already used
94 // as it is for all projectiles.)
95 fP0_dbar = 0.0;
96 fP0_d = 0.0;
97 if ( thePrimary.GetDefinition()->GetPDGEncoding() == 2212 ) { // proton
98 G4double mproj = thePrimary.GetDefinition()->GetPDGMass();
99 G4double pz = thePrimary.Get4Momentum().z();
100 G4double ekin = std::sqrt( pz*pz + mproj*mproj ) - mproj;
101 if ( ekin > 10.0 ) {
102 fP0_dbar = 130.0 / ( 1.0 + std::exp( 21.6 - std::log( 0.001*ekin )/0.089 ) ); // set p0 for antideuteron
103 fP0_d = 118.1 * ( 1.0 + std::exp( 5.53 - std::log( 0.001*ekin )/0.43 ) ); // set p0 for deuteron
104 }
105 }
106 //G4cout << "Coalescence parameter p0 deuteron / antideuteron: " << fP0_d << " / " << fP0_dbar << G4endl;
107}
const G4ParticleDefinition * GetDefinition() const
const G4LorentzVector & Get4Momentum() const

Referenced by G4TheoFSGenerator::ApplyYourself().


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