Geant4 11.4.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 82 of file G4CRCoalescence.cc.

82 : G4HadronicInteraction("G4CRCoalescence" ),
83 fP0_d( 0.0 ), fP0_dbar( 0.0 ), secID( -1 ) {
84 secID = G4PhysicsModelCatalog::GetModelID( "model_G4CRCoalescence" );
85}
G4HadronicInteraction(const G4String &modelName="HadronicModel")
static G4int GetModelID(const G4int modelIndex)

Referenced by G4CRCoalescence(), operator!=(), operator=(), and operator==().

◆ ~G4CRCoalescence()

G4CRCoalescence::~G4CRCoalescence ( )
override

Definition at line 88 of file G4CRCoalescence.cc.

88{}

◆ G4CRCoalescence() [2/2]

G4CRCoalescence::G4CRCoalescence ( const G4CRCoalescence & right)
delete

Member Function Documentation

◆ GenerateDeuterons()

void G4CRCoalescence::GenerateDeuterons ( G4ReactionProductVector * result)

Definition at line 113 of file G4CRCoalescence.cc.

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

◆ 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 91 of file G4CRCoalescence.cc.

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

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