Geant4-11
Public Member Functions | Protected Member Functions | Protected Attributes
G4ParticleHP2NDInelasticFS Class Reference

#include <G4ParticleHP2NDInelasticFS.hh>

Inheritance diagram for G4ParticleHP2NDInelasticFS:
G4ParticleHPInelasticBaseFS G4ParticleHPFinalState

Public Member Functions

G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack)
 
void BaseApply (const G4HadProjectile &theTrack, G4ParticleDefinition **theDefs, G4int nDef)
 
 G4ParticleHP2NDInelasticFS ()
 
G4double GetA ()
 
G4int GetM ()
 
G4double GetN ()
 
virtual G4ParticleHPVectorGetXsec ()
 
virtual G4double GetXsec (G4double anEnergy)
 
G4double GetZ ()
 
G4bool HasAnyData ()
 
G4bool HasFSData ()
 
G4bool HasXsec ()
 
void Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *)
 
void Init (G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile)
 
void InitGammas (G4double AR, G4double ZR)
 
G4ParticleHPFinalStateNew ()
 
void SetA_Z (G4double anA, G4double aZ, G4int aM=0)
 
void SetAZMs (G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
 
void SetProjectile (G4ParticleDefinition *projectile)
 
 ~G4ParticleHP2NDInelasticFS ()
 

Protected Member Functions

void adjust_final_state (G4LorentzVector)
 

Protected Attributes

G4String gammaPath
 
G4bool hasAnyData
 
G4bool hasFSData
 
G4bool hasXsec
 
G4double Qvalue
 
G4int secID
 
G4ParticleHPAngulartheAngularDistribution
 
G4double theBaseA
 
G4int theBaseM
 
G4double theBaseZ
 
G4ParticleHPEnAngCorrelationtheEnergyAngData
 
G4ParticleHPEnergyDistributiontheEnergyDistribution
 
G4ParticleHPPhotonDisttheFinalStatePhotons
 
G4ParticleHPDeExGammas theGammas
 
G4ParticleHPNames theNames
 
G4int theNDLDataA
 
G4int theNDLDataM
 
G4int theNDLDataZ
 
G4double theNuclearMassDifference
 
G4ParticleDefinitiontheProjectile
 
G4Cache< G4HadFinalState * > theResult
 
G4ParticleHPVectortheXsection
 

Detailed Description

Definition at line 41 of file G4ParticleHP2NDInelasticFS.hh.

Constructor & Destructor Documentation

◆ G4ParticleHP2NDInelasticFS()

G4ParticleHP2NDInelasticFS::G4ParticleHP2NDInelasticFS ( )

Definition at line 37 of file G4ParticleHP2NDInelasticFS.cc.

38{
39 secID = G4PhysicsModelCatalog::GetModelID( "model_G4ParticleHP2NDInelasticFS_F03" );
40}
static G4int GetModelID(const G4int modelIndex)

References G4PhysicsModelCatalog::GetModelID(), and G4ParticleHPFinalState::secID.

Referenced by New().

◆ ~G4ParticleHP2NDInelasticFS()

G4ParticleHP2NDInelasticFS::~G4ParticleHP2NDInelasticFS ( )
inline

Definition at line 46 of file G4ParticleHP2NDInelasticFS.hh.

46{}

Member Function Documentation

◆ adjust_final_state()

void G4ParticleHPFinalState::adjust_final_state ( G4LorentzVector  init_4p_lab)
protectedinherited

Definition at line 47 of file G4ParticleHPFinalState.cc.

48{
49
50 G4double minimum_energy = 1*keV;
51
52 if ( G4ParticleHPManager::GetInstance()->GetDoNotAdjustFinalState() ) return;
53
54 G4int nSecondaries = theResult.Get()->GetNumberOfSecondaries();
55
56 G4int sum_Z = 0;
57 G4int sum_A = 0;
58 G4int max_SecZ = 0;
59 G4int max_SecA = 0;
60 G4int imaxA = -1;
61 for ( int i = 0 ; i < nSecondaries ; i++ )
62 {
63 //G4cout << "G4ParticleHPFinalState::adjust_final_state theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetParticleName() = " << theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetParticleName() << G4endl;
65 max_SecZ = std::max ( max_SecZ , theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicNumber() );
67 max_SecA = std::max ( max_SecA , theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() );
68 if ( theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() == max_SecA ) imaxA = i;
69#ifdef G4PHPDEBUG
70 if( std::getenv("G4ParticleHPDebug")) G4cout << "G4ParticleHPFinalState::adjust_final_stat SECO " << i << " " <<theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetParticleName() << G4endl;
71#endif
72
73 }
74
75 G4ParticleDefinition* resi_pd = 0;
76
77 G4double baseZNew = theBaseZ;
78 G4double baseANew = theBaseA;
80 baseANew ++;
81 } else if( theProjectile == G4Proton::Proton() ) {
82 baseZNew ++;
83 baseANew ++;
84 } else if( theProjectile == G4Deuteron::Deuteron() ) {
85 baseZNew ++;
86 baseANew += 2;
87 } else if( theProjectile == G4Triton::Triton() ) {
88 baseZNew ++;
89 baseANew += 3;
90 } else if( theProjectile == G4He3::He3() ) {
91 baseZNew += 2;
92 baseANew += 3;
93 } else if( theProjectile == G4Alpha::Alpha() ) {
94 baseZNew += 2;
95 baseANew += 4;
96 }
97
98#ifdef G4PHPDEBUG
99 if( std::getenv("G4ParticleHPDebug")) G4cout << "G4ParticleHPFinalState::adjust_final_stat BaseZ " << baseZNew << " BaseA " << baseANew << " sum_Z " << sum_Z << " sum_A " << sum_A << G4endl;
100#endif
101
102 G4bool needOneMoreSec = false;
103 G4ParticleDefinition* oneMoreSec_pd = 0;
104 if ( (int)(baseZNew - sum_Z) == 0 && (int)(baseANew - sum_A) == 0 )
105 {
106 //All secondaries are already created;
107 resi_pd = theResult.Get()->GetSecondary( imaxA )->GetParticle()->GetDefinition();
108 }
109 else
110 {
111 if ( max_SecA >= int(baseANew - sum_A) )
112 {
113 //Most heavy secondary is interpreted as residual
114 resi_pd = theResult.Get()->GetSecondary( imaxA )->GetParticle()->GetDefinition();
115 needOneMoreSec = true;
116 }
117 else
118 {
119 //creation of residual is requierd
120 resi_pd = G4IonTable::GetIonTable()->GetIon ( int(baseZNew - sum_Z) , (int)(baseANew - sum_A) , 0.0 );
121 }
122
123 if ( needOneMoreSec )
124 {
125 if ( int(baseZNew - sum_Z) == 0 && (int)(baseANew - sum_A) > 0 )
126 {
127 //In this case, one neutron is added to secondaries
128 if ( int(baseANew - sum_A) > 1 ) G4cout << "More than one neutron is required for the balance of baryon number!" << G4endl;
129 oneMoreSec_pd = G4Neutron::Neutron();
130 }
131 else
132 {
133#ifdef G4PHPDEBUG
134 if( std::getenv("G4ParticleHPDebug")) G4cout << this << "G4ParticleHPFinalState oneMoreSec_pd Z " << baseZNew << " - " << sum_Z << " A " << baseANew << " - " << sum_A << " projectile " << theProjectile->GetParticleName() << G4endl;
135#endif
136 oneMoreSec_pd = G4IonTable::GetIonTable()->GetIon ( int(baseZNew - sum_Z) , (int)(baseANew - sum_A) , 0.0 );
137 if( !oneMoreSec_pd ) {
138 G4cerr << this << "G4ParticleHPFinalState oneMoreSec_pd Z " << baseZNew << " - " << sum_Z << " A " << baseANew << " - " << sum_A << " projectile " << theProjectile->GetParticleName() << G4endl;
139 G4Exception("G4ParticleHPFinalState:adjust_final_state",
140 "Warning",
142 "No adjustment will be done!");
143 return;
144 }
145 }
146 }
147
148 if ( resi_pd == 0 )
149 {
150 // theNDLDataZ,A has the Z and A of used NDL file
151 G4double ndlZNew = theNDLDataZ;
152 G4double ndlANew = theNDLDataA;
154 ndlANew ++;
155 } else if( theProjectile == G4Proton::Proton() ) {
156 ndlZNew ++;
157 ndlANew ++;
158 } else if( theProjectile == G4Deuteron::Deuteron() ) {
159 ndlZNew ++;
160 ndlANew += 2;
161 } else if( theProjectile == G4Triton::Triton() ) {
162 ndlZNew ++;
163 ndlANew += 3;
164 } else if( theProjectile == G4He3::He3() ) {
165 ndlZNew += 2;
166 ndlANew += 3;
167 } else if( theProjectile == G4Alpha::Alpha() ) {
168 ndlZNew += 2;
169 ndlANew += 4;
170 }
171 // theNDLDataZ,A has the Z and A of used NDL file
172 if ( (int)(ndlZNew - sum_Z) == 0 && (int)(ndlANew - sum_A) == 0 )
173 {
174 G4int dif_Z = ( int ) ( theNDLDataZ - theBaseZ );
175 G4int dif_A = ( int ) ( theNDLDataA - theBaseA );
176 resi_pd = G4IonTable::GetIonTable()->GetIon ( max_SecZ - dif_Z , max_SecA - dif_A , 0.0 );
177 if( !resi_pd ) {
178 G4cerr << "G4ParticleHPFinalState resi_pd Z " << max_SecZ << " - " << dif_Z << " A " << max_SecA << " - " << dif_A << " projectile " << theProjectile->GetParticleName() << G4endl;
179 G4Exception("G4ParticleHPFinalState:adjust_final_state",
180 "Warning",
182 "No adjustment will be done!");
183 return;
184 }
185
186 for ( int i = 0 ; i < nSecondaries ; i++ )
187 {
188 if ( theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicNumber() == max_SecZ
189 && theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition()->GetAtomicMass() == max_SecA )
190 {
192 p = p * resi_pd->GetPDGMass()/ G4IonTable::GetIonTable()->GetIon ( max_SecZ , max_SecA , 0.0 )->GetPDGMass();
193 theResult.Get()->GetSecondary( i )->GetParticle()->SetDefinition( resi_pd );
195 }
196 }
197 }
198 }
199 }
200
201
202 G4LorentzVector secs_4p_lab( 0.0 );
203
205 G4double fast = 0;
206 G4double slow = 1;
207 G4int ifast = 0;
208 G4int islow = 0;
209 G4int ires = -1;
210
211 for ( G4int i = 0 ; i < n_sec ; i++ )
212 {
213
214 //G4cout << "HP_DB " << i
215 // << " " << theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetParticleName()
216 // << " 4p " << theResult.GetSecondary( i )->GetParticle()->Get4Momentum()
217 // << " ke " << theResult.GetSecondary( i )->GetParticle()->Get4Momentum().e() - theResult.GetSecondary( i )->GetParticle()->GetDefinition()->GetPDGMass()
218 // << G4endl;
219
220 secs_4p_lab += theResult.Get()->GetSecondary( i )->GetParticle()->Get4Momentum();
221
222 G4double beta = 0;
224 {
226 }
227 else
228 {
229 beta = 1;
230 }
231
232 if ( theResult.Get()->GetSecondary( i )->GetParticle()->GetDefinition() == resi_pd ) ires = i;
233
234 if ( slow > beta && beta != 0 )
235 {
236 slow = beta;
237 islow = i;
238 }
239
240 if ( fast <= beta )
241 {
242 if ( fast != 1 )
243 {
244 fast = beta;
245 ifast = i;
246 }
247 else
248 {
249// fast is already photon then check E
251 if ( e > theResult.Get()->GetSecondary( ifast )->GetParticle()->Get4Momentum().e() )
252 {
253// among photons, the highest E becomes the fastest
254 ifast = i;
255 }
256 }
257 }
258 }
259
260
261 G4LorentzVector dif_4p = init_4p_lab - secs_4p_lab;
262
263 //G4cout << "HP_DB dif_4p " << init_4p_lab - secs_4p_lab << G4endl;
264 //G4cout << "HP_DB dif_3p mag " << ( dif_4p.v() ).mag() << G4endl;
265 //G4cout << "HP_DB dif_e " << dif_4p.e() - ( dif_4p.v() ).mag()<< G4endl;
266
267 G4LorentzVector p4(0);
268 if ( ires == -1 )
269 {
270// Create and Add Residual Nucleus
271 ires = nSecondaries;
272 nSecondaries += 1;
273
274 G4DynamicParticle* res = new G4DynamicParticle ( resi_pd , dif_4p.v() );
275 theResult.Get()->AddSecondary ( res, secID );
276
277 p4 = res->Get4Momentum();
278 if ( slow > p4.beta() )
279 {
280 slow = p4.beta();
281 islow = ires;
282 }
283 dif_4p = init_4p_lab - ( secs_4p_lab + p4 );
284 }
285
286 if ( needOneMoreSec && oneMoreSec_pd)
287 //
288 // fca: this is not a fix, this is a crash avoidance...
289 // fca: the baryon number is still wrong, most probably because it
290 // fca: should have been decreased, but since we could not create a particle
291 // fca: we just do not add it
292 //
293 {
294 nSecondaries += 1;
295 G4DynamicParticle* one = new G4DynamicParticle ( oneMoreSec_pd , dif_4p.v() );
296 theResult.Get()->AddSecondary ( one, secID );
297 p4 = one->Get4Momentum();
298 if ( slow > p4.beta() )
299 {
300 slow = p4.beta();
301 islow = nSecondaries-1; //Because the first is 0th, so the last becomes "nSecondaries-1"
302 }
303 dif_4p = init_4p_lab - ( secs_4p_lab + p4 );
304 }
305
306 //Which is bigger dif_p or dif_e
307
308 if ( dif_4p.v().mag() < std::abs( dif_4p.e() ) )
309 {
310
311 // Adjust p
312 //if ( dif_4p.v().mag() < 1*MeV )
313 if ( minimum_energy < dif_4p.v().mag() && dif_4p.v().mag() < 1*MeV )
314 {
315
316 nSecondaries += 1;
317 theResult.Get()->AddSecondary ( new G4DynamicParticle ( G4Gamma::Gamma() , dif_4p.v() ), secID );
318
319 }
320 else
321 {
322 //G4cout << "HP_DB Difference in dif_p is too large (>1MeV) or too small(<1keV) to adjust, so that give up tuning" << G4endl;
323 }
324
325 }
326 else
327 {
328
329 // dif_p > dif_e
330 // at first momentum
331 // Move residual momentum
332
333 p4 = theResult.Get()->GetSecondary( ires )->GetParticle()->Get4Momentum();
334 theResult.Get()->GetSecondary( ires )->GetParticle()->SetMomentum( p4.v() + dif_4p.v() );
335 dif_4p = init_4p_lab - ( secs_4p_lab - p4 + theResult.Get()->GetSecondary( ires )->GetParticle()->Get4Momentum() );
336
337 //G4cout << "HP_DB new residual kinetic energy " << theResult.GetSecondary( ires )->GetParticle()->GetKineticEnergy() << G4endl;
338
339 }
340
341 G4double dif_e = dif_4p.e() - ( dif_4p.v() ).mag();
342 //G4cout << "HP_DB dif_e " << dif_e << G4endl;
343
344 if ( dif_e > 0 )
345 {
346
347// create 2 gamma
348
349 nSecondaries += 2;
350 G4double e1 = ( dif_4p.e() -dif_4p.v().mag() ) / 2;
351
352 if ( minimum_energy < e1 )
353 {
354 G4double costh = 2.*G4UniformRand()-1.;
356 G4ThreeVector dir( std::sin(std::acos(costh))*std::cos(phi),
357 std::sin(std::acos(costh))*std::sin(phi),
358 costh);
361 }
362 else
363 {
364 //G4cout << "HP_DB Difference is too small(<1keV) to adjust, so that neglect it" << G4endl;
365 }
366
367 }
368 else //dif_e < 0
369 {
370
371// At first reduce KE of the fastest secondary;
374 G4ThreeVector dir = ( theResult.Get()->GetSecondary( ifast )->GetParticle()->GetMomentum() ).unit();
375
376 //G4cout << "HP_DB ifast " << ifast << " ke0 " << ke0 << G4endl;
377
378 if ( ke0 + dif_e > 0 )
379 {
380 theResult.Get()->GetSecondary( ifast )->GetParticle()->SetKineticEnergy( ke0 + dif_e );
381 G4ThreeVector dp = p0 - theResult.Get()->GetSecondary( ifast )->GetParticle()->GetMomentum();
382
384 //theResult.GetSecondary( islow )->GetParticle()->SetMomentum( p - dif_e*dir );
385 theResult.Get()->GetSecondary( islow )->GetParticle()->SetMomentum( p + dp );
386 }
387 else
388 {
389 //G4cout << "HP_DB Difference in dif_e too large ( <0MeV ) to adjust, so that give up tuning" << G4endl;
390 }
391
392 }
393
394}
static const G4double e1[44]
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double keV
Definition: G4SIunits.hh:202
static constexpr double MeV
Definition: G4SIunits.hh:200
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
double mag() const
Hep3Vector v() const
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
value_type & Get() const
Definition: G4Cache.hh:315
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4ParticleDefinition * GetDefinition() const
G4LorentzVector Get4Momentum() const
G4double GetKineticEnergy() const
void SetMomentum(const G4ThreeVector &momentum)
G4ThreeVector GetMomentum() const
void SetKineticEnergy(G4double aEnergy)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
std::size_t GetNumberOfSecondaries() const
G4HadSecondary * GetSecondary(size_t i)
G4DynamicParticle * GetParticle()
static G4He3 * He3()
Definition: G4He3.cc:93
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4IonTable * GetIonTable()
Definition: G4IonTable.cc:170
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
const G4String & GetParticleName() const
G4ParticleDefinition * theProjectile
G4Cache< G4HadFinalState * > theResult
static G4ParticleHPManager * GetInstance()
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4Triton * Triton()
Definition: G4Triton.cc:93
T max(const T t1, const T t2)
brief Return the largest of the two arguments

References G4HadFinalState::AddSecondary(), G4Alpha::Alpha(), CLHEP::HepLorentzVector::beta(), anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, G4Deuteron::Deuteron(), CLHEP::HepLorentzVector::e(), e1, G4cerr, G4cout, G4endl, G4Exception(), G4UniformRand, G4Gamma::Gamma(), G4Cache< VALTYPE >::Get(), G4DynamicParticle::Get4Momentum(), G4ParticleDefinition::GetAtomicMass(), G4ParticleDefinition::GetAtomicNumber(), G4DynamicParticle::GetDefinition(), G4ParticleHPManager::GetInstance(), G4IonTable::GetIon(), G4IonTable::GetIonTable(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentum(), G4HadFinalState::GetNumberOfSecondaries(), G4HadSecondary::GetParticle(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGMass(), G4HadFinalState::GetSecondary(), G4He3::He3(), JustWarning, keV, CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), MeV, G4Neutron::Neutron(), G4Proton::Proton(), G4ParticleHPFinalState::secID, G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetKineticEnergy(), G4DynamicParticle::SetMomentum(), G4ParticleHPFinalState::theBaseA, G4ParticleHPFinalState::theBaseZ, G4ParticleHPFinalState::theNDLDataA, G4ParticleHPFinalState::theNDLDataZ, G4ParticleHPFinalState::theProjectile, G4ParticleHPFinalState::theResult, G4Triton::Triton(), twopi, and CLHEP::HepLorentzVector::v().

Referenced by G4ParticleHPInelasticBaseFS::BaseApply(), and G4ParticleHPInelasticCompFS::CompositeApply().

◆ ApplyYourself()

G4HadFinalState * G4ParticleHP2NDInelasticFS::ApplyYourself ( const G4HadProjectile theTrack)
virtual

Implements G4ParticleHPInelasticBaseFS.

Definition at line 42 of file G4ParticleHP2NDInelasticFS.cc.

43{
44// these are the particle types in the final state
45
46 G4ParticleDefinition * theDefs[3];
47 theDefs[0] = G4Neutron::Neutron();
48 theDefs[1] = G4Neutron::Neutron();
49 theDefs[2] = G4Deuteron::Deuteron();
50
51// fill the final state
52 G4ParticleHPInelasticBaseFS::BaseApply(theTrack, theDefs, 3);
53
54// return the result
55 return theResult.Get();
56}
void BaseApply(const G4HadProjectile &theTrack, G4ParticleDefinition **theDefs, G4int nDef)

References G4ParticleHPInelasticBaseFS::BaseApply(), G4Deuteron::Deuteron(), G4Cache< VALTYPE >::Get(), G4Neutron::Neutron(), and G4ParticleHPFinalState::theResult.

◆ BaseApply()

void G4ParticleHPInelasticBaseFS::BaseApply ( const G4HadProjectile theTrack,
G4ParticleDefinition **  theDefs,
G4int  nDef 
)
inherited

Definition at line 187 of file G4ParticleHPInelasticBaseFS.cc.

190{
191
192// prepare neutron
193 if ( theResult.Get() == NULL ) theResult.Put( new G4HadFinalState );
194 theResult.Get()->Clear();
195 G4double eKinetic = theTrack.GetKineticEnergy();
196 const G4HadProjectile *hadProjectile = &theTrack;
197 G4ReactionProduct incidReactionProduct( const_cast<G4ParticleDefinition *>(hadProjectile->GetDefinition()) );
198 incidReactionProduct.SetMomentum( hadProjectile->Get4Momentum().vect() );
199 incidReactionProduct.SetKineticEnergy( eKinetic );
200
201// prepare target
202 G4double targetMass;
203 G4double eps = 0.0001;
204 targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) /
205 //theProjectile->GetPDGMass();
207
208 //give priority to ENDF vales for target mass
209 if(theEnergyAngData!=0)
210 { targetMass = theEnergyAngData->GetTargetMass(); }
212 { targetMass = theAngularDistribution->GetTargetMass(); }
213
214//080731a
215//110512 TKDB ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded.
216//if ( targetMass == 0 ) G4cout << "080731a It looks like something wrong value in G4NDL, please update the latest version. If you use the latest, then please report this problem to Geant4 Hyper news." << G4endl;
217 if ( targetMass == 0 )
218 {
219 //G4cout << "TKDB targetMass = 0; ENDF-VII.0 21Sc45 has trouble in MF4MT22 (n,np) targetMass is not properly recorded. This could be a similar situation." << G4endl;
220 //targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / theProjectile->GetPDGMass();
221 targetMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps))) / G4Neutron::Neutron()->GetPDGMass();
222 }
223
224 G4Nucleus aNucleus;
225 G4ReactionProduct theTarget;
226 //G4ThreeVector neuVelo = (1./hadProjectile->GetDefinition()->GetPDGMass())*incidReactionProduct.GetMomentum();
227 //theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
228 //G4Nucleus::GetBiasedThermalNucleus requests normalization of mass and velocity in neutron mass
229 G4ThreeVector neuVelo = (1./G4Neutron::Neutron()->GetPDGMass())*incidReactionProduct.GetMomentum();
230 theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
231
233
234// prepare energy in target rest frame
235 G4ReactionProduct boosted;
236 boosted.Lorentz(incidReactionProduct, theTarget);
237 eKinetic = boosted.GetKineticEnergy();
238 G4double orgMomentum = boosted.GetMomentum().mag();
239
240// Take N-body phase-space distribution, if no other data present.
241 if(!HasFSData()) // adding the residual is trivial here @@@
242 {
243 G4ParticleHPNBodyPhaseSpace thePhaseSpaceDistribution;
244 G4double aPhaseMass=0;
245 G4int ii;
246 for(ii=0; ii<nDef; ii++)
247 {
248 aPhaseMass+=theDefs[ii]->GetPDGMass();
249 }
250
251 //----------------------------------------------------------------------------
252 if(Qvalue<1.*CLHEP::keV && Qvalue>-1.*CLHEP::keV){ //Not in the G4NDL lib or not calculated yet:
253 //Calculate residual:
254 G4int ResidualA=theBaseA;
255 G4int ResidualZ=theBaseZ;
256 for (ii = 0; ii < nDef; ii++) {
257 ResidualZ -= theDefs[ii]->GetAtomicNumber();
258 ResidualA -= theDefs[ii]->GetBaryonNumber();
259 }
260
261 if (ResidualA > 0 && ResidualZ > 0) {
262 G4ParticleDefinition* resid = G4IonTable::GetIonTable()->GetIon(ResidualZ,ResidualA);
263 Qvalue = incidReactionProduct.GetMass()+theTarget.GetMass()-aPhaseMass-resid->GetPDGMass();
264 }
265
266 if (Qvalue > 400*CLHEP::MeV || Qvalue < -400*CLHEP::MeV) {
267 //Then Q value is probably too large ...
268 Qvalue = 1.1*CLHEP::keV;
269 }
270 }
271 //----------------------------------------------------------------------------
272
273 thePhaseSpaceDistribution.Init(aPhaseMass, nDef);
274 thePhaseSpaceDistribution.SetProjectileRP(&incidReactionProduct);
275 thePhaseSpaceDistribution.SetTarget(&theTarget);
276 thePhaseSpaceDistribution.SetQValue(Qvalue);
277
278 for(ii=0; ii<nDef; ii++)
279 {
280 G4double massCode = 1000.*std::abs(theDefs[ii]->GetPDGCharge());
281 massCode += theDefs[ii]->GetBaryonNumber();
282 G4double dummy = 0;
283 G4ReactionProduct * aSec = thePhaseSpaceDistribution.Sample(eKinetic, massCode, dummy);
284 aSec->Lorentz(*aSec, -1.*theTarget);
285 G4DynamicParticle * aPart = new G4DynamicParticle();
286 aPart->SetDefinition(aSec->GetDefinition());
287 aPart->SetMomentum(aSec->GetMomentum());
288 delete aSec;
289 theResult.Get()->AddSecondary(aPart, secID);
290#ifdef G4PHPDEBUG
291 if( std::getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply NoFSData add secondary " << aPart->GetParticleDefinition()->GetParticleName() << " E= " << aPart->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
292#endif
293 }
295 //TK120607
296 //Final momentum check should be done before return
298 G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
299 G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
300 G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
301 adjust_final_state ( init_4p_lab );
302
303 return;
304 }
305
306// set target and neutron in the relevant exit channel
308 {
310 theAngularDistribution->SetProjectileRP(incidReactionProduct);
311 }
312 else if(theEnergyAngData!=0)
313 {
314 theEnergyAngData->SetTarget(theTarget);
315 theEnergyAngData->SetProjectileRP(incidReactionProduct);
316 }
317
318 G4ReactionProductVector * tmpHadrons = 0;
319#ifdef G4PHPDEBUG
320 //To avoid compilation error around line 532.
321 G4int ii(0);
322#endif
323 G4int dummy;
324 unsigned int i;
325
326 if(theEnergyAngData != 0)
327 {
328 tmpHadrons = theEnergyAngData->Sample(eKinetic);
329
330 if ( ! G4ParticleHPManager::GetInstance()->GetDoNotAdjustFinalState() ) {
331 //141017 Fix BEGIN
332 //Adjust A and Z in the case of miss much between selected data and target nucleus
333 if ( tmpHadrons != NULL ) {
334 G4int sumA = 0;
335 G4int sumZ = 0;
336 G4int maxA = 0;
337 G4int jAtMaxA = 0;
338 for ( G4int j = 0 ; j != (G4int)tmpHadrons->size() ; j++ ) {
339 //G4cout << __FILE__ << " " << __LINE__ << "th line: tmpHadrons->at(j)->GetDefinition()->GetParticleName() = " << tmpHadrons->at(j)->GetDefinition()->GetParticleName() << G4endl;
340 if ( tmpHadrons->at(j)->GetDefinition()->GetBaryonNumber() > maxA ) {
341 maxA = tmpHadrons->at(j)->GetDefinition()->GetBaryonNumber();
342 jAtMaxA = j;
343 }
344 sumA += tmpHadrons->at(j)->GetDefinition()->GetBaryonNumber();
345 sumZ += G4int( tmpHadrons->at(j)->GetDefinition()->GetPDGCharge() + eps );
346 }
347 G4int dA = (G4int)theBaseA + hadProjectile->GetDefinition()->GetBaryonNumber() - sumA;
348 G4int dZ = (G4int)theBaseZ + G4int( hadProjectile->GetDefinition()->GetPDGCharge() + eps ) - sumZ;
349 if ( dA < 0 || dZ < 0 ) {
350 G4int newA = tmpHadrons->at(jAtMaxA)->GetDefinition()->GetBaryonNumber() + dA ;
351 G4int newZ = G4int( tmpHadrons->at(jAtMaxA)->GetDefinition()->GetPDGCharge() + eps ) + dZ;
352 if(newA>newZ && newZ>0){
353 G4ParticleDefinition* pd = G4IonTable::GetIonTable()->GetIon ( newZ , newA );
354 tmpHadrons->at( jAtMaxA )->SetDefinition( pd );
355 }
356 }
357 }
358 //141017 Fix END
359 }
360 }
361 else if(theAngularDistribution!= 0)
362 {
363 G4bool * Done = new G4bool[nDef];
364 G4int i0;
365 for(i0=0; i0<nDef; i0++) Done[i0] = false;
366
367//Following lines are commented out to fix coverity defeat 58622
368// if(tmpHadrons == 0)
369// {
370 tmpHadrons = new G4ReactionProductVector;
371// }
372// else
373// {
374// for(i=0; i<tmpHadrons->size(); i++)
375// {
376// for(ii=0; ii<nDef; ii++)
377// if(!Done[ii] && tmpHadrons->operator[](i)->GetDefinition() == theDefs[ii])
378// Done[ii] = true;
379// }
380#ifdef G4PHPDEBUG
381// if( getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply secondary previously added " << tmpHadrons->operator[](i)->GetDefinition()->GetParticleName() << " E= " << tmpHadrons->operator[](i)->GetKineticEnergy() << G4endl;
382#endif
383// }
384 G4ReactionProduct * aHadron;
385 G4double localMass = ( G4NucleiProperties::GetNuclearMass(static_cast<G4int>(theBaseA+eps), static_cast<G4int>(theBaseZ+eps)));
386 G4ThreeVector bufferedDirection(0,0,0);
387 for(i0=0; i0<nDef; i0++)
388 {
389 if(!Done[i0])
390 {
391 aHadron = new G4ReactionProduct;
393 {
394 aHadron->SetDefinition(theDefs[i0]);
395 aHadron->SetKineticEnergy(theEnergyDistribution->Sample(eKinetic, dummy));
396 }
397 else if(nDef == 1)
398 {
399 aHadron->SetDefinition(theDefs[i0]);
400 aHadron->SetKineticEnergy(eKinetic);
401 }
402 else if(nDef == 2)
403 {
404 aHadron->SetDefinition(theDefs[i0]);
405 aHadron->SetKineticEnergy(50*CLHEP::MeV);
406 }
407 else
408 {
409 throw G4HadronicException(__FILE__, __LINE__, "No energy distribution to sample from in InelasticBaseFS::BaseApply");
410 }
412 if(theEnergyDistribution==0 && nDef == 2)
413 {
414 if(i0==0)
415 {
416 G4double mass1 = theDefs[0]->GetPDGMass();
417 G4double mass2 = theDefs[1]->GetPDGMass();
419 G4int z1 = static_cast<G4int>(theBaseZ+eps-theDefs[0]->GetPDGCharge()-theDefs[1]->GetPDGCharge());
420 G4int a1 = static_cast<G4int>(theBaseA+eps)-theDefs[0]->GetBaryonNumber()-theDefs[1]->GetBaryonNumber();
421 G4double concreteMass = G4NucleiProperties::GetNuclearMass(a1, z1);
422 G4double availableEnergy = eKinetic+massn+localMass-mass1-mass2-concreteMass;
423 // available kinetic energy in CMS (non relativistic)
424 G4double emin = availableEnergy+mass1+mass2 - std::sqrt((mass1+mass2)*(mass1+mass2)+orgMomentum*orgMomentum);
425 G4double p1=std::sqrt(2.*mass2*emin);
426 bufferedDirection = p1*aHadron->GetMomentum().unit();
427#ifdef G4PHPDEBUG
428 if(std::getenv("G4ParticleHPDebug")) // @@@@@ verify the nucleon counting...
429 {
430 G4cout << "G4ParticleHPInelasticBaseFS "<<z1<<" "<<theBaseZ<<" "<<a1<<" "<<theBaseA<<" "<<availableEnergy<<" "
431 << emin<<G4endl;
432 }
433#endif
434 }
435 else
436 {
437 bufferedDirection = -bufferedDirection;
438 }
439 // boost from cms to lab
440#ifdef G4PHPDEBUG
441 if(std::getenv("G4ParticleHPDebug"))
442 {
443 G4cout << " G4ParticleHPInelasticBaseFS "<<bufferedDirection.mag2()<<G4endl;
444 }
445#endif
446 aHadron->SetTotalEnergy( std::sqrt(aHadron->GetMass()*aHadron->GetMass()
447 +bufferedDirection.mag2()) );
448 aHadron->SetMomentum(bufferedDirection);
449 aHadron->Lorentz(*aHadron, -1.*(theTarget+incidReactionProduct));
450#ifdef G4PHPDEBUG
451 if(std::getenv("G4ParticleHPDebug"))
452 {
453 G4cout << " G4ParticleHPInelasticBaseFS "<<aHadron->GetTotalEnergy()<<" "<<aHadron->GetMomentum()<<G4endl;
454 }
455#endif
456 }
457 tmpHadrons->push_back(aHadron);
458#ifdef G4PHPDEBUG
459 if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply FSData add secondary " << aHadron->GetDefinition()->GetParticleName() << " E= " << aHadron->GetKineticEnergy() << G4endl;
460#endif
461 }
462 }
463 delete [] Done;
464 }
465 else
466 {
467 throw G4HadronicException(__FILE__, __LINE__, "No data to create the neutrons in NInelasticFS");
468 }
469
470 G4ReactionProductVector * thePhotons = 0;
472 {
473 // the photon distributions are in the Nucleus rest frame.
474 G4ReactionProduct boosted_tmp;
475 boosted_tmp.Lorentz(incidReactionProduct, theTarget);
476 G4double anEnergy = boosted_tmp.GetKineticEnergy();
477 thePhotons = theFinalStatePhotons->GetPhotons(anEnergy);
478 if(thePhotons!=0)
479 {
480 for(i=0; i<thePhotons->size(); i++)
481 {
482 // back to lab
483 thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), -1.*theTarget);
484 }
485 }
486 }
487 else if(theEnergyAngData!=0)
488 {
489
490 // PA130927: do not create photons to adjust binding energy
491 G4bool bAdjustPhotons = true;
492#ifdef PHP_AS_HP
493 bAdjustPhotons = true;
494#else
495 if ( G4ParticleHPManager::GetInstance()->GetDoNotAdjustFinalState() ) bAdjustPhotons = false;
496#endif
497
498 if( bAdjustPhotons ) {
499 G4double theGammaEnergy = theEnergyAngData->GetTotalMeanEnergy();
500 G4double anEnergy = boosted.GetKineticEnergy();
501 theGammaEnergy = anEnergy-theGammaEnergy;
502 theGammaEnergy += theNuclearMassDifference;
503 G4double eBindProducts = 0;
504 G4double eBindN = 0;
505 G4double eBindP = 0;
510 G4int ia=0;
511 for(i=0; i<tmpHadrons->size(); i++)
512 {
513 if(tmpHadrons->operator[](i)->GetDefinition() == G4Neutron::Neutron())
514 {
515 eBindProducts+=eBindN;
516 }
517 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Proton::Proton())
518 {
519 eBindProducts+=eBindP;
520 }
521 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Deuteron::Deuteron())
522 {
523 eBindProducts+=eBindD;
524 }
525 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Triton::Triton())
526 {
527 eBindProducts+=eBindT;
528 }
529 else if(tmpHadrons->operator[](i)->GetDefinition() == G4He3::He3())
530 {
531 eBindProducts+=eBindHe3;
532 }
533 else if(tmpHadrons->operator[](i)->GetDefinition() == G4Alpha::Alpha())
534 {
535 eBindProducts+=eBindA;
536 ia++;
537 }
538 }
539
540 theGammaEnergy += eBindProducts;
541
542#ifdef G4PHPDEBUG
543 if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply gamma Energy " << theGammaEnergy << " eBindProducts " << eBindProducts << G4endl;
544#endif
545
546 //101111
547 //Special treatment for Be9 + n -> 2n + Be8 -> 2n + a + a
548 if ( (G4int)(theBaseZ+eps) == 4 && (G4int)(theBaseA+eps) == 9 )
549 {
550 // This only valid for G4NDL3.13,,,
551 if ( std::abs( theNuclearMassDifference -
554 && ia == 2 )
555 {
556 theGammaEnergy -= (2*eBindA);
557 }
558 }
559
560 G4ReactionProductVector * theOtherPhotons = 0;
561 G4int iLevel;
562 while(theGammaEnergy>=theGammas.GetLevelEnergy(0)) // Loop checking, 11.05.2015, T. Koi
563 {
564 for(iLevel=theGammas.GetNumberOfLevels()-1; iLevel>=0; iLevel--)
565 {
566 if(theGammas.GetLevelEnergy(iLevel)<theGammaEnergy) break;
567 }
568 if(iLevel==0||iLevel==theGammas.GetNumberOfLevels()-1)
569 {
570 theOtherPhotons = theGammas.GetDecayGammas(iLevel);
571#ifdef G4PHPDEBUG
572 if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply adding gamma from level " << iLevel << theOtherPhotons->operator[](ii)->GetKineticEnergy() << G4endl;
573#endif
574 }
575 else
576 {
577 G4double random = G4UniformRand();
578 G4double eLow = theGammas.GetLevelEnergy(iLevel);
579 G4double eHigh = theGammas.GetLevelEnergy(iLevel+1);
580 if(random > (eHigh-eLow)/(theGammaEnergy-eLow)) iLevel++;
581 theOtherPhotons = theGammas.GetDecayGammas(iLevel);
582 }
583 if(thePhotons==0) thePhotons = new G4ReactionProductVector;
584 if(theOtherPhotons != 0)
585 {
586 for(unsigned int iii=0; iii<theOtherPhotons->size(); iii++)
587 {
588 thePhotons->push_back(theOtherPhotons->operator[](iii));
589#ifdef G4PHPDEBUG
590 if( std::getenv("G4ParticleHPDebug"))
591 G4cout << iii << " G4ParticleHPInelasticBaseFS::BaseApply adding gamma " << theOtherPhotons->operator[](iii)->GetKineticEnergy() << G4endl;
592#endif
593 }
594 delete theOtherPhotons;
595 }
596 theGammaEnergy -= theGammas.GetLevelEnergy(iLevel);
597 if(iLevel == -1) break;
598 }
599 }
600 }
601
602 // fill the result
603 unsigned int nSecondaries = tmpHadrons->size();
604 unsigned int nPhotons = 0;
605 if(thePhotons!=0) { nPhotons = thePhotons->size(); }
606 nSecondaries += nPhotons;
607 G4DynamicParticle * theSec;
608#ifdef G4PHPDEBUG
609 if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply N hadrons " << nSecondaries-nPhotons << G4endl;
610#endif
611
612 for(i=0; i<nSecondaries-nPhotons; i++)
613 {
614 theSec = new G4DynamicParticle;
615 theSec->SetDefinition(tmpHadrons->operator[](i)->GetDefinition());
616 theSec->SetMomentum(tmpHadrons->operator[](i)->GetMomentum());
617 theResult.Get()->AddSecondary(theSec, secID);
618#ifdef G4PHPDEBUG
619 if( std::getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply add secondary2 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
620#endif
621 if( std::getenv("G4PHPTEST") ) G4cout << " InelasticBaseFS COS THETA " << std::cos(theSec->GetMomentum().theta()) << " " << (theSec->GetMomentum().theta()) << " " << theSec->GetMomentum() << " E "<< theSec->GetKineticEnergy() << " " << theSec->GetDefinition()->GetParticleName() << G4endl; //GDEB
622 delete tmpHadrons->operator[](i);
623 }
624#ifdef G4PHPDEBUG
625 if( std::getenv("G4ParticleHPDebug")) G4cout << " G4ParticleHPInelasticBaseFS::BaseApply N photons " << nPhotons << G4endl;
626#endif
627 if(thePhotons != 0)
628 {
629 for(i=0; i<nPhotons; i++)
630 {
631 theSec = new G4DynamicParticle;
632 theSec->SetDefinition(thePhotons->operator[](i)->GetDefinition());
633 theSec->SetMomentum(thePhotons->operator[](i)->GetMomentum());
634 theResult.Get()->AddSecondary(theSec, secID);
635#ifdef G4PHPDEBUG
636 if( std::getenv("G4ParticleHPDebug")) G4cout << this << " G4ParticleHPInelasticBaseFS::BaseApply add secondary3 " << theSec->GetParticleDefinition()->GetParticleName() << " E= " << theSec->GetKineticEnergy() << " NSECO " << theResult.Get()->GetNumberOfSecondaries() << G4endl;
637#endif
638 delete thePhotons->operator[](i);
639 }
640 }
641
642// some garbage collection
643 delete thePhotons;
644 delete tmpHadrons;
645
646//080721
648 G4LorentzVector targ_4p_lab ( theTarget.GetMomentum() , std::sqrt( targ_pd->GetPDGMass()*targ_pd->GetPDGMass() + theTarget.GetMomentum().mag2() ) );
649 G4LorentzVector proj_4p_lab = theTrack.Get4Momentum();
650 G4LorentzVector init_4p_lab = proj_4p_lab + targ_4p_lab;
651
652 //if data in MF=6 format (no correlated particle emission), then adjust_final_state can give severe errors:
653 if(theEnergyAngData==0){adjust_final_state ( init_4p_lab );}
654
655// clean up the primary neutron
657}
static const G4double eps
static const G4int maxA
@ stopAndKill
std::vector< G4ReactionProduct * > G4ReactionProductVector
Hep3Vector unit() const
double theta() const
double mag2() const
Hep3Vector vect() const
void Put(const value_type &val) const
Definition: G4Cache.hh:321
const G4ParticleDefinition * GetParticleDefinition() const
void SetStatusChange(G4HadFinalStateStatus aS)
const G4Material * GetMaterial() const
const G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4LorentzVector & Get4Momentum() const
G4double GetTemperature() const
Definition: G4Material.hh:178
static G4double GetBindingEnergy(const G4int A, const G4int Z)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:118
G4double GetPDGCharge() const
void SetTarget(const G4ReactionProduct &aTarget)
void SetProjectileRP(const G4ReactionProduct &anIncidentParticleRP)
void SampleAndUpdate(G4ReactionProduct &anIncidentParticle)
G4ReactionProductVector * GetDecayGammas(G4int aLevel)
G4double GetLevelEnergy(G4int aLevel)
void SetProjectileRP(G4ReactionProduct &aIncidentPart)
void SetTarget(G4ReactionProduct &aTarget)
G4ReactionProductVector * Sample(G4double anEnergy)
G4double Sample(G4double anEnergy, G4int &it)
void adjust_final_state(G4LorentzVector)
G4ParticleHPEnergyDistribution * theEnergyDistribution
G4ParticleHPPhotonDist * theFinalStatePhotons
G4ParticleHPAngular * theAngularDistribution
G4ParticleHPEnAngCorrelation * theEnergyAngData
G4ReactionProduct * Sample(G4double anEnergy, G4double massCode, G4double mass)
void Init(G4double aMass, G4int aCount)
G4ReactionProductVector * GetPhotons(G4double anEnergy)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
G4double GetKineticEnergy() const
const G4ParticleDefinition * GetDefinition() const
G4double GetTotalEnergy() const
G4ThreeVector GetMomentum() const
void Lorentz(const G4ReactionProduct &p1, const G4ReactionProduct &p2)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
void SetKineticEnergy(const G4double en)
G4double GetMass() const
void SetProjectileRP(G4ReactionProduct *aIncidentParticleRP)
void SetTarget(G4ReactionProduct *aTarget)
static constexpr double keV
static constexpr double MeV

References G4HadFinalState::AddSecondary(), G4ParticleHPFinalState::adjust_final_state(), G4Alpha::Alpha(), G4HadFinalState::Clear(), G4Deuteron::Deuteron(), eps, G4cout, G4endl, G4UniformRand, G4Cache< VALTYPE >::Get(), G4HadProjectile::Get4Momentum(), G4ParticleDefinition::GetAtomicNumber(), G4ParticleDefinition::GetBaryonNumber(), G4Nucleus::GetBiasedThermalNucleus(), G4NucleiProperties::GetBindingEnergy(), G4ParticleHPDeExGammas::GetDecayGammas(), G4DynamicParticle::GetDefinition(), G4HadProjectile::GetDefinition(), G4ReactionProduct::GetDefinition(), G4ParticleHPManager::GetDoNotAdjustFinalState(), G4ParticleHPManager::GetInstance(), G4IonTable::GetIon(), G4IonTable::GetIonTable(), G4DynamicParticle::GetKineticEnergy(), G4HadProjectile::GetKineticEnergy(), G4ReactionProduct::GetKineticEnergy(), G4ParticleHPDeExGammas::GetLevelEnergy(), G4ReactionProduct::GetMass(), G4HadProjectile::GetMaterial(), G4DynamicParticle::GetMomentum(), G4ReactionProduct::GetMomentum(), G4NucleiProperties::GetNuclearMass(), G4ParticleHPDeExGammas::GetNumberOfLevels(), G4HadFinalState::GetNumberOfSecondaries(), G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4ParticleHPPhotonDist::GetPhotons(), G4ParticleHPAngular::GetTargetMass(), G4ParticleHPEnAngCorrelation::GetTargetMass(), G4Material::GetTemperature(), G4ReactionProduct::GetTotalEnergy(), G4ParticleHPEnAngCorrelation::GetTotalMeanEnergy(), G4ParticleHPFinalState::HasFSData(), G4He3::He3(), G4ParticleHPNBodyPhaseSpace::Init(), CLHEP::keV, G4ReactionProduct::Lorentz(), CLHEP::Hep3Vector::mag(), CLHEP::Hep3Vector::mag2(), maxA, CLHEP::MeV, G4Neutron::Neutron(), G4Proton::Proton(), G4Cache< VALTYPE >::Put(), G4ParticleHPInelasticBaseFS::Qvalue, G4ParticleHPEnAngCorrelation::Sample(), G4ParticleHPNBodyPhaseSpace::Sample(), G4ParticleHPEnergyDistribution::Sample(), G4ParticleHPAngular::SampleAndUpdate(), G4ParticleHPFinalState::secID, G4DynamicParticle::SetDefinition(), G4ReactionProduct::SetDefinition(), G4ReactionProduct::SetKineticEnergy(), G4ReactionProduct::SetMomentum(), G4DynamicParticle::SetMomentum(), G4ParticleHPAngular::SetProjectileRP(), G4ParticleHPEnAngCorrelation::SetProjectileRP(), G4VParticleHPEnergyAngular::SetProjectileRP(), G4VParticleHPEnergyAngular::SetQValue(), G4HadFinalState::SetStatusChange(), G4ParticleHPAngular::SetTarget(), G4ParticleHPEnAngCorrelation::SetTarget(), G4VParticleHPEnergyAngular::SetTarget(), G4ReactionProduct::SetTotalEnergy(), stopAndKill, G4ParticleHPInelasticBaseFS::theAngularDistribution, G4ParticleHPFinalState::theBaseA, G4ParticleHPFinalState::theBaseZ, G4ParticleHPInelasticBaseFS::theEnergyAngData, G4ParticleHPInelasticBaseFS::theEnergyDistribution, G4ParticleHPInelasticBaseFS::theFinalStatePhotons, G4ParticleHPInelasticBaseFS::theGammas, G4ParticleHPInelasticBaseFS::theNuclearMassDifference, G4ParticleHPFinalState::theProjectile, G4ParticleHPFinalState::theResult, CLHEP::Hep3Vector::theta(), G4Triton::Triton(), CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

Referenced by G4ParticleHP2AInelasticFS::ApplyYourself(), G4ParticleHP2N2AInelasticFS::ApplyYourself(), G4ParticleHP2NAInelasticFS::ApplyYourself(), ApplyYourself(), G4ParticleHP2NInelasticFS::ApplyYourself(), G4ParticleHP2NPInelasticFS::ApplyYourself(), G4ParticleHP2PInelasticFS::ApplyYourself(), G4ParticleHP3AInelasticFS::ApplyYourself(), G4ParticleHP3NAInelasticFS::ApplyYourself(), G4ParticleHP3NInelasticFS::ApplyYourself(), G4ParticleHP3NPInelasticFS::ApplyYourself(), G4ParticleHP4NInelasticFS::ApplyYourself(), G4ParticleHPD2AInelasticFS::ApplyYourself(), G4ParticleHPDAInelasticFS::ApplyYourself(), G4ParticleHPN2AInelasticFS::ApplyYourself(), G4ParticleHPN2PInelasticFS::ApplyYourself(), G4ParticleHPN3AInelasticFS::ApplyYourself(), G4ParticleHPNAInelasticFS::ApplyYourself(), G4ParticleHPND2AInelasticFS::ApplyYourself(), G4ParticleHPNDInelasticFS::ApplyYourself(), G4ParticleHPNHe3InelasticFS::ApplyYourself(), G4ParticleHPNPAInelasticFS::ApplyYourself(), G4ParticleHPNPInelasticFS::ApplyYourself(), G4ParticleHPNT2AInelasticFS::ApplyYourself(), G4ParticleHPNTInelasticFS::ApplyYourself(), G4ParticleHPNXInelasticFS::ApplyYourself(), G4ParticleHPPAInelasticFS::ApplyYourself(), G4ParticleHPPDInelasticFS::ApplyYourself(), G4ParticleHPPTInelasticFS::ApplyYourself(), and G4ParticleHPT2AInelasticFS::ApplyYourself().

◆ GetA()

G4double G4ParticleHPFinalState::GetA ( void  )
inlineinherited

Definition at line 104 of file G4ParticleHPFinalState.hh.

104{ return theBaseA; }

References G4ParticleHPFinalState::theBaseA.

◆ GetM()

G4int G4ParticleHPFinalState::GetM ( )
inlineinherited

◆ GetN()

G4double G4ParticleHPFinalState::GetN ( )
inlineinherited

◆ GetXsec() [1/2]

virtual G4ParticleHPVector * G4ParticleHPInelasticBaseFS::GetXsec ( )
inlinevirtualinherited

Reimplemented from G4ParticleHPFinalState.

Definition at line 87 of file G4ParticleHPInelasticBaseFS.hh.

87{return theXsection;}

References G4ParticleHPInelasticBaseFS::theXsection.

◆ GetXsec() [2/2]

virtual G4double G4ParticleHPInelasticBaseFS::GetXsec ( G4double  anEnergy)
inlinevirtualinherited

Reimplemented from G4ParticleHPFinalState.

Definition at line 82 of file G4ParticleHPInelasticBaseFS.hh.

83 {
84 return std::max(0., theXsection->GetY(anEnergy));
85 }
G4double GetY(G4double x)

References G4ParticleHPVector::GetY(), G4INCL::Math::max(), and G4ParticleHPInelasticBaseFS::theXsection.

◆ GetZ()

G4double G4ParticleHPFinalState::GetZ ( void  )
inlineinherited

◆ HasAnyData()

G4bool G4ParticleHPFinalState::HasAnyData ( )
inlineinherited

◆ HasFSData()

G4bool G4ParticleHPFinalState::HasFSData ( )
inlineinherited

◆ HasXsec()

G4bool G4ParticleHPFinalState::HasXsec ( )
inlineinherited

◆ Init() [1/2]

void G4ParticleHP2NDInelasticFS::Init ( G4double  A,
G4double  Z,
G4int  M,
G4String dirName,
G4String aFSType,
G4ParticleDefinition projectile 
)
virtual

Implements G4ParticleHPFinalState.

Definition at line 58 of file G4ParticleHP2NDInelasticFS.cc.

60{
61 G4ParticleHPInelasticBaseFS::Init(A, Z, M, dirName, aFSType, projectile);
62 G4double ResidualA = 0;
63 G4double ResidualZ = 0;
64 if( projectile == G4Neutron::Neutron() ) {
65 ResidualA = A-3;
66 ResidualZ = Z-1;
67 } else if( projectile == G4Proton::Proton() ) {
68 ResidualA = A-3;
69 ResidualZ = Z;
70 } else if( projectile == G4Deuteron::Deuteron() ) {
71 ResidualA = A-2;
72 ResidualZ = Z;
73 } else if( projectile == G4Triton::Triton() ) {
74 ResidualA = A-1;
75 ResidualZ = Z;
76 } else if( projectile == G4He3::He3() ) {
77 ResidualA = A-1;
78 ResidualZ = Z+1;
79 } else if( projectile == G4Alpha::Alpha() ) {
80 ResidualA = A;
81 ResidualZ = Z+1;
82 }
83
84 G4ParticleHPInelasticBaseFS::InitGammas(ResidualA, ResidualZ);
85}
#define M(row, col)
const G4int Z[17]
const G4double A[17]
void InitGammas(G4double AR, G4double ZR)
void Init(G4double A, G4double Z, G4int M, G4String &dirName, G4String &bit, G4ParticleDefinition *)

References A, G4Alpha::Alpha(), G4Deuteron::Deuteron(), G4He3::He3(), G4ParticleHPInelasticBaseFS::Init(), G4ParticleHPInelasticBaseFS::InitGammas(), M, G4Neutron::Neutron(), G4Proton::Proton(), G4Triton::Triton(), and Z.

◆ Init() [2/2]

void G4ParticleHPFinalState::Init ( G4double  A,
G4double  Z,
G4String dirName,
G4String aFSType,
G4ParticleDefinition projectile 
)
inlineinherited

Definition at line 76 of file G4ParticleHPFinalState.hh.

78 {
79 G4int M = 0;
80 Init ( A, Z, M, dirName, aFSType,const_cast<G4ParticleDefinition*>(projectile));
81 }
void Init(G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile)

References A, G4ParticleHPFinalState::Init(), M, and Z.

Referenced by G4ParticleHPFinalState::Init(), and G4ParticleHPChannel::UpdateData().

◆ InitGammas()

void G4ParticleHPInelasticBaseFS::InitGammas ( G4double  AR,
G4double  ZR 
)
inherited

Definition at line 51 of file G4ParticleHPInelasticBaseFS.cc.

52{
53 // char the[100] = {""};
54 // std::ostrstream ost(the, 100, std::ios::out);
55 // ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
56 // G4String * aName = new G4String(the);
57 // std::ifstream from(*aName, std::ios::in);
58
59 std::ostringstream ost;
60 ost <<gammaPath<<"z"<<ZR<<".a"<<AR;
61 G4String aName = ost.str();
62 std::ifstream from(aName, std::ios::in);
63
64 if(!from) return; // no data found for this isotope
65 // std::ifstream theGammaData(*aName, std::ios::in);
66 std::ifstream theGammaData(aName, std::ios::in);
67
68 G4double eps = 0.001;
70 G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(AR+eps),static_cast<G4int>(ZR+eps)) -
72 theGammas.Init(theGammaData);
73 // delete aName;
74
75}
void Init(std::istream &aDataFile)

References eps, G4ParticleHPInelasticBaseFS::gammaPath, G4NucleiProperties::GetBindingEnergy(), G4ParticleHPDeExGammas::Init(), G4ParticleHPFinalState::theBaseA, G4ParticleHPFinalState::theBaseZ, G4ParticleHPInelasticBaseFS::theGammas, and G4ParticleHPInelasticBaseFS::theNuclearMassDifference.

Referenced by G4ParticleHP2AInelasticFS::Init(), G4ParticleHP2N2AInelasticFS::Init(), G4ParticleHP2NAInelasticFS::Init(), Init(), G4ParticleHP2NInelasticFS::Init(), G4ParticleHP2NPInelasticFS::Init(), G4ParticleHP2PInelasticFS::Init(), G4ParticleHP3AInelasticFS::Init(), G4ParticleHP3NAInelasticFS::Init(), G4ParticleHP3NInelasticFS::Init(), G4ParticleHP3NPInelasticFS::Init(), G4ParticleHP4NInelasticFS::Init(), G4ParticleHPD2AInelasticFS::Init(), G4ParticleHPDAInelasticFS::Init(), G4ParticleHPN2AInelasticFS::Init(), G4ParticleHPN2PInelasticFS::Init(), G4ParticleHPN3AInelasticFS::Init(), G4ParticleHPNAInelasticFS::Init(), G4ParticleHPND2AInelasticFS::Init(), G4ParticleHPNDInelasticFS::Init(), G4ParticleHPNHe3InelasticFS::Init(), G4ParticleHPNPAInelasticFS::Init(), G4ParticleHPNPInelasticFS::Init(), G4ParticleHPNT2AInelasticFS::Init(), G4ParticleHPNTInelasticFS::Init(), G4ParticleHPNXInelasticFS::Init(), G4ParticleHPPAInelasticFS::Init(), G4ParticleHPPDInelasticFS::Init(), G4ParticleHPPTInelasticFS::Init(), and G4ParticleHPT2AInelasticFS::Init().

◆ New()

G4ParticleHPFinalState * G4ParticleHP2NDInelasticFS::New ( )
inlinevirtual

◆ SetA_Z()

void G4ParticleHPFinalState::SetA_Z ( G4double  anA,
G4double  aZ,
G4int  aM = 0 
)
inlineinherited

◆ SetAZMs()

void G4ParticleHPFinalState::SetAZMs ( G4double  anA,
G4double  aZ,
G4int  aM,
G4ParticleHPDataUsed  used 
)
inlineinherited

◆ SetProjectile()

void G4ParticleHPFinalState::SetProjectile ( G4ParticleDefinition projectile)
inlineinherited

Definition at line 115 of file G4ParticleHPFinalState.hh.

116 {
117 theProjectile = projectile;
118 }

References G4ParticleHPFinalState::theProjectile.

Referenced by G4ParticleHPChannel::Register().

Field Documentation

◆ gammaPath

G4String G4ParticleHPInelasticBaseFS::gammaPath
protectedinherited

◆ hasAnyData

G4bool G4ParticleHPFinalState::hasAnyData
protectedinherited

◆ hasFSData

G4bool G4ParticleHPFinalState::hasFSData
protectedinherited

◆ hasXsec

G4bool G4ParticleHPFinalState::hasXsec
protectedinherited

◆ Qvalue

G4double G4ParticleHPInelasticBaseFS::Qvalue
protectedinherited

◆ secID

G4int G4ParticleHPFinalState::secID
protectedinherited

Definition at line 140 of file G4ParticleHPFinalState.hh.

Referenced by G4ParticleHPFinalState::adjust_final_state(), G4ParticleHPCaptureFS::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4ParticleHPFissionFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4ParticleHPInelasticCompFS::CompositeApply(), G4ParticleHP2AInelasticFS::G4ParticleHP2AInelasticFS(), G4ParticleHP2N2AInelasticFS::G4ParticleHP2N2AInelasticFS(), G4ParticleHP2NAInelasticFS::G4ParticleHP2NAInelasticFS(), G4ParticleHP2NDInelasticFS(), G4ParticleHP2NInelasticFS::G4ParticleHP2NInelasticFS(), G4ParticleHP2NPInelasticFS::G4ParticleHP2NPInelasticFS(), G4ParticleHP2PInelasticFS::G4ParticleHP2PInelasticFS(), G4ParticleHP3AInelasticFS::G4ParticleHP3AInelasticFS(), G4ParticleHP3NAInelasticFS::G4ParticleHP3NAInelasticFS(), G4ParticleHP3NInelasticFS::G4ParticleHP3NInelasticFS(), G4ParticleHP3NPInelasticFS::G4ParticleHP3NPInelasticFS(), G4ParticleHP4NInelasticFS::G4ParticleHP4NInelasticFS(), G4ParticleHPAInelasticFS::G4ParticleHPAInelasticFS(), G4ParticleHPCaptureFS::G4ParticleHPCaptureFS(), G4ParticleHPD2AInelasticFS::G4ParticleHPD2AInelasticFS(), G4ParticleHPDAInelasticFS::G4ParticleHPDAInelasticFS(), G4ParticleHPDInelasticFS::G4ParticleHPDInelasticFS(), G4ParticleHPElasticFS::G4ParticleHPElasticFS(), G4ParticleHPFinalState::G4ParticleHPFinalState(), G4ParticleHPFissionFS::G4ParticleHPFissionFS(), G4ParticleHPHe3InelasticFS::G4ParticleHPHe3InelasticFS(), G4ParticleHPN2AInelasticFS::G4ParticleHPN2AInelasticFS(), G4ParticleHPN2PInelasticFS::G4ParticleHPN2PInelasticFS(), G4ParticleHPN3AInelasticFS::G4ParticleHPN3AInelasticFS(), G4ParticleHPNAInelasticFS::G4ParticleHPNAInelasticFS(), G4ParticleHPND2AInelasticFS::G4ParticleHPND2AInelasticFS(), G4ParticleHPNDInelasticFS::G4ParticleHPNDInelasticFS(), G4ParticleHPNHe3InelasticFS::G4ParticleHPNHe3InelasticFS(), G4ParticleHPNInelasticFS::G4ParticleHPNInelasticFS(), G4ParticleHPNPAInelasticFS::G4ParticleHPNPAInelasticFS(), G4ParticleHPNPInelasticFS::G4ParticleHPNPInelasticFS(), G4ParticleHPNT2AInelasticFS::G4ParticleHPNT2AInelasticFS(), G4ParticleHPNTInelasticFS::G4ParticleHPNTInelasticFS(), G4ParticleHPNXInelasticFS::G4ParticleHPNXInelasticFS(), G4ParticleHPPAInelasticFS::G4ParticleHPPAInelasticFS(), G4ParticleHPPDInelasticFS::G4ParticleHPPDInelasticFS(), G4ParticleHPPInelasticFS::G4ParticleHPPInelasticFS(), G4ParticleHPPTInelasticFS::G4ParticleHPPTInelasticFS(), G4ParticleHPT2AInelasticFS::G4ParticleHPT2AInelasticFS(), G4ParticleHPTInelasticFS::G4ParticleHPTInelasticFS(), and G4ParticleHPInelasticCompFS::use_nresp71_model().

◆ theAngularDistribution

G4ParticleHPAngular* G4ParticleHPInelasticBaseFS::theAngularDistribution
protectedinherited

◆ theBaseA

G4double G4ParticleHPFinalState::theBaseA
protectedinherited

◆ theBaseM

G4int G4ParticleHPFinalState::theBaseM
protectedinherited

◆ theBaseZ

G4double G4ParticleHPFinalState::theBaseZ
protectedinherited

◆ theEnergyAngData

G4ParticleHPEnAngCorrelation* G4ParticleHPInelasticBaseFS::theEnergyAngData
protectedinherited

◆ theEnergyDistribution

G4ParticleHPEnergyDistribution* G4ParticleHPInelasticBaseFS::theEnergyDistribution
protectedinherited

◆ theFinalStatePhotons

G4ParticleHPPhotonDist* G4ParticleHPInelasticBaseFS::theFinalStatePhotons
protectedinherited

◆ theGammas

G4ParticleHPDeExGammas G4ParticleHPInelasticBaseFS::theGammas
protectedinherited

◆ theNames

G4ParticleHPNames G4ParticleHPFinalState::theNames
protectedinherited

◆ theNDLDataA

G4int G4ParticleHPFinalState::theNDLDataA
protectedinherited

◆ theNDLDataM

G4int G4ParticleHPFinalState::theNDLDataM
protectedinherited

◆ theNDLDataZ

G4int G4ParticleHPFinalState::theNDLDataZ
protectedinherited

◆ theNuclearMassDifference

G4double G4ParticleHPInelasticBaseFS::theNuclearMassDifference
protectedinherited

◆ theProjectile

G4ParticleDefinition* G4ParticleHPFinalState::theProjectile
protectedinherited

◆ theResult

G4Cache< G4HadFinalState* > G4ParticleHPFinalState::theResult
protectedinherited

Definition at line 129 of file G4ParticleHPFinalState.hh.

Referenced by G4ParticleHPFinalState::adjust_final_state(), G4FissionLibrary::ApplyYourself(), G4ParticleHP2AInelasticFS::ApplyYourself(), G4ParticleHP2N2AInelasticFS::ApplyYourself(), G4ParticleHP2NAInelasticFS::ApplyYourself(), ApplyYourself(), G4ParticleHP2NInelasticFS::ApplyYourself(), G4ParticleHP2NPInelasticFS::ApplyYourself(), G4ParticleHP2PInelasticFS::ApplyYourself(), G4ParticleHP3AInelasticFS::ApplyYourself(), G4ParticleHP3NAInelasticFS::ApplyYourself(), G4ParticleHP3NInelasticFS::ApplyYourself(), G4ParticleHP3NPInelasticFS::ApplyYourself(), G4ParticleHP4NInelasticFS::ApplyYourself(), G4ParticleHPAInelasticFS::ApplyYourself(), G4ParticleHPCaptureFS::ApplyYourself(), G4ParticleHPD2AInelasticFS::ApplyYourself(), G4ParticleHPDAInelasticFS::ApplyYourself(), G4ParticleHPDInelasticFS::ApplyYourself(), G4ParticleHPElasticFS::ApplyYourself(), G4ParticleHPFissionFS::ApplyYourself(), G4ParticleHPHe3InelasticFS::ApplyYourself(), G4ParticleHPN2AInelasticFS::ApplyYourself(), G4ParticleHPN2PInelasticFS::ApplyYourself(), G4ParticleHPN3AInelasticFS::ApplyYourself(), G4ParticleHPNAInelasticFS::ApplyYourself(), G4ParticleHPND2AInelasticFS::ApplyYourself(), G4ParticleHPNDInelasticFS::ApplyYourself(), G4ParticleHPNHe3InelasticFS::ApplyYourself(), G4ParticleHPNInelasticFS::ApplyYourself(), G4ParticleHPNPAInelasticFS::ApplyYourself(), G4ParticleHPNPInelasticFS::ApplyYourself(), G4ParticleHPNT2AInelasticFS::ApplyYourself(), G4ParticleHPNTInelasticFS::ApplyYourself(), G4ParticleHPNXInelasticFS::ApplyYourself(), G4ParticleHPPAInelasticFS::ApplyYourself(), G4ParticleHPPDInelasticFS::ApplyYourself(), G4ParticleHPPInelasticFS::ApplyYourself(), G4ParticleHPPTInelasticFS::ApplyYourself(), G4ParticleHPT2AInelasticFS::ApplyYourself(), G4ParticleHPTInelasticFS::ApplyYourself(), G4ParticleHPInelasticBaseFS::BaseApply(), G4ParticleHPInelasticCompFS::CompositeApply(), G4ParticleHPFinalState::G4ParticleHPFinalState(), G4ParticleHPInelasticCompFS::use_nresp71_model(), and G4ParticleHPFinalState::~G4ParticleHPFinalState().

◆ theXsection

G4ParticleHPVector* G4ParticleHPInelasticBaseFS::theXsection
protectedinherited

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