Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4ParticleHPCaptureFS Class Reference

#include <G4ParticleHPCaptureFS.hh>

Inheritance diagram for G4ParticleHPCaptureFS:
Collaboration diagram for G4ParticleHPCaptureFS:

Public Member Functions

 G4ParticleHPCaptureFS ()
 
 ~G4ParticleHPCaptureFS ()
 
void Init (G4double A, G4double Z, G4int M, G4String &dirName, G4String &aFSType, G4ParticleDefinition *)
 
G4HadFinalStateApplyYourself (const G4HadProjectile &theTrack)
 
G4ParticleHPFinalStateNew ()
 
- Public Member Functions inherited from G4ParticleHPFinalState
 G4ParticleHPFinalState ()
 
virtual ~G4ParticleHPFinalState ()
 
void Init (G4double A, G4double Z, G4String &dirName, G4String &aFSType, G4ParticleDefinition *projectile)
 
G4bool HasXsec ()
 
G4bool HasFSData ()
 
G4bool HasAnyData ()
 
virtual G4double GetXsec (G4double)
 
virtual G4ParticleHPVectorGetXsec ()
 
void SetA_Z (G4double anA, G4double aZ, G4int aM=0)
 
G4double GetZ ()
 
G4double GetN ()
 
G4double GetA ()
 
G4int GetM ()
 
void SetAZMs (G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
 
void SetProjectile (G4ParticleDefinition *projectile)
 

Additional Inherited Members

- Protected Member Functions inherited from G4ParticleHPFinalState
void adjust_final_state (G4LorentzVector)
 
G4bool DoNotAdjustFinalState ()
 
- Protected Attributes inherited from G4ParticleHPFinalState
G4bool hasXsec
 
G4bool hasFSData
 
G4bool hasAnyData
 
G4ParticleHPNames theNames
 
G4Cache< G4HadFinalState * > theResult
 
G4double theBaseA
 
G4double theBaseZ
 
G4int theBaseM
 
G4int theNDLDataZ
 
G4int theNDLDataA
 
G4int theNDLDataM
 
G4ParticleDefinitiontheProjectile
 

Detailed Description

Definition at line 42 of file G4ParticleHPCaptureFS.hh.

Constructor & Destructor Documentation

G4ParticleHPCaptureFS::G4ParticleHPCaptureFS ( )
inline

Definition at line 46 of file G4ParticleHPCaptureFS.hh.

47  {
48  hasXsec = false;
49  hasExactMF6 = false;
50  targetMass = 0;
51  }

Here is the caller graph for this function:

G4ParticleHPCaptureFS::~G4ParticleHPCaptureFS ( )
inline

Definition at line 53 of file G4ParticleHPCaptureFS.hh.

54  {
55  }

Member Function Documentation

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

Reimplemented from G4ParticleHPFinalState.

Definition at line 50 of file G4ParticleHPCaptureFS.cc.

51  {
52 
53  if ( theResult.Get() == NULL ) theResult.Put( new G4HadFinalState );
54  theResult.Get()->Clear();
55 
56  G4int i;
57 
58 // prepare neutron
59  G4double eKinetic = theTrack.GetKineticEnergy();
60  const G4HadProjectile *incidentParticle = &theTrack;
61  G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition() ) );
62  theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
63  theNeutron.SetKineticEnergy( eKinetic );
64 
65 // prepare target
66  G4ReactionProduct theTarget;
67  G4Nucleus aNucleus;
68  G4double eps = 0.0001;
69  if(targetMass<500*MeV)
70  targetMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theBaseA+eps) , static_cast<G4int>(theBaseZ+eps) )) /
72  G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
73  G4double temperature = theTrack.GetMaterial()->GetTemperature();
74  theTarget = aNucleus.GetBiasedThermalNucleus(targetMass, neutronVelocity, temperature);
75  theTarget.SetDefinition( G4IonTable::GetIonTable()->GetIon( G4int(theBaseZ), G4int(theBaseA) , 0.0 ) ); //TESTPHP
76 
77 // go to nucleus rest system
78  theNeutron.Lorentz(theNeutron, -1*theTarget);
79  eKinetic = theNeutron.GetKineticEnergy();
80 
81 // dice the photons
82 
83  G4ReactionProductVector * thePhotons = 0;
84  if ( HasFSData() && !G4ParticleHPManager::GetInstance()->GetUseOnlyPhotoEvaporation() )
85  {
86  //NDL has final state data
87  if ( hasExactMF6 ) {
88  theMF6FinalState.SetTarget(theTarget);
89  theMF6FinalState.SetProjectileRP(theNeutron);
90  thePhotons = theMF6FinalState.Sample( eKinetic );
91  } else {
92  thePhotons = theFinalStatePhotons.GetPhotons(eKinetic);
93  }
94  if ( thePhotons == NULL ) {
95  throw G4HadronicException(__FILE__, __LINE__, "Final state data for photon is not properly allocated");
96  }
97  }
98  else
99  {
100  //NDL does not have final state data or forced to use PhotoEvaporation model
101  G4ThreeVector aCMSMomentum = theNeutron.GetMomentum()+theTarget.GetMomentum();
102  G4LorentzVector p4(aCMSMomentum, theTarget.GetTotalEnergy() + theNeutron.GetTotalEnergy());
103  G4Fragment nucleus(static_cast<G4int>(theBaseA+1), static_cast<G4int>(theBaseZ) ,p4);
104  G4PhotonEvaporation photonEvaporation;
105  // T. K. add
106  photonEvaporation.SetICM( TRUE );
107  G4FragmentVector* products = photonEvaporation.BreakItUp(nucleus);
108  G4FragmentVector::iterator it;
109  thePhotons = new G4ReactionProductVector;
110  for(it=products->begin(); it!=products->end(); it++)
111  {
112  G4ReactionProduct * theOne = new G4ReactionProduct;
113  // T. K. add
114  if ( (*it)->GetParticleDefinition() != 0 )
115  theOne->SetDefinition( (*it)->GetParticleDefinition() );
116  else
117  theOne->SetDefinition( G4Gamma::Gamma() ); // this definiion will be over writen
118 
119  // T. K. comment out below line
120  //theOne->SetDefinition( G4Gamma::Gamma() );
121  G4IonTable* theTable = G4IonTable::GetIonTable();
122  if( (*it)->GetMomentum().mag() > 10*MeV ) theOne->SetDefinition( theTable->GetIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0 ) );
123 
124  //if ( (*i)->GetExcitationEnergy() > 0 )
125  if ( (*it)->GetExcitationEnergy() > 1.0e-2*eV )
126  {
127  G4double ex = (*it)->GetExcitationEnergy();
128  G4ReactionProduct* aPhoton = new G4ReactionProduct;
129  aPhoton->SetDefinition( G4Gamma::Gamma() );
130  aPhoton->SetMomentum( (*it)->GetMomentum().vect().unit() * ex );
131  //aPhoton->SetTotalEnergy( ex ); //will be calculated from momentum
132  thePhotons->push_back(aPhoton);
133  }
134 
135  theOne->SetMomentum( (*it)->GetMomentum().vect() * ( (*it)->GetMomentum().t() - (*it)->GetExcitationEnergy() ) / (*it)->GetMomentum().t() ) ;
136  //theOne->SetTotalEnergy( (*i)->GetMomentum().t() - (*i)->GetExcitationEnergy() ); //will be calculated from momentum
137  thePhotons->push_back(theOne);
138  delete *it;
139  }
140  delete products;
141  }
142 
143 // add them to the final state
144 
145  G4int nPhotons = 0;
146  nPhotons=thePhotons->size();
147 
149  if ( DoNotAdjustFinalState() ) {
150 //Make at least one photon
151 //101203 TK
152  if ( nPhotons == 0 )
153  {
154  G4ReactionProduct * theOne = new G4ReactionProduct;
155  theOne->SetDefinition( G4Gamma::Gamma() );
156  G4double theta = pi*G4UniformRand();
157  G4double phi = twopi*G4UniformRand();
158  G4double sinth = std::sin(theta);
159  G4ThreeVector direction( sinth*std::cos(phi), sinth*std::sin(phi), std::cos(theta) );
160  theOne->SetMomentum( direction ) ;
161  thePhotons->push_back(theOne);
162  nPhotons++; // 0 -> 1
163  }
164 //One photon case: energy set to Q-value
165 //101203 TK
166  //if ( nPhotons == 1 )
167  if ( nPhotons == 1 && thePhotons->operator[](0)->GetDefinition()->GetBaryonNumber() == 0 )
168  {
169  G4ThreeVector direction = thePhotons->operator[](0)->GetMomentum().unit();
170 
171  G4double Q = G4IonTable::GetIonTable()->GetIonMass(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA), 0) + G4Neutron::Neutron()->GetPDGMass()
172  - G4IonTable::GetIonTable()->GetIonMass(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0);
173 
174  thePhotons->operator[](0)->SetMomentum( Q*direction );
175  }
176 //
177  }
178 
179  // back to lab system
180  for(i=0; i<nPhotons; i++)
181  {
182  thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), theTarget);
183  }
184 
185  // Recoil, if only one gamma
186  //if (1==nPhotons)
187  if ( nPhotons == 1 && thePhotons->operator[](0)->GetDefinition()->GetBaryonNumber() == 0 )
188  {
189  G4DynamicParticle * theOne = new G4DynamicParticle;
191  ->GetIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0);
192  theOne->SetDefinition(aRecoil);
193  // Now energy;
194  // Can be done slightly better @
195  G4ThreeVector aMomentum = theTrack.Get4Momentum().vect()
196  +theTarget.GetMomentum()
197  -thePhotons->operator[](0)->GetMomentum();
198 
199  //TKDB 140520
200  //G4ThreeVector theMomUnit = aMomentum.unit();
201  //G4double aKinEnergy = theTrack.GetKineticEnergy()
202  // +theTarget.GetKineticEnergy(); // gammas come from Q-value
203  //G4double theResMass = aRecoil->GetPDGMass();
204  //G4double theResE = aRecoil->GetPDGMass()+aKinEnergy;
205  //G4double theAbsMom = std::sqrt(theResE*theResE - theResMass*theResMass);
206  //G4ThreeVector theMomentum = theAbsMom*theMomUnit;
207  //theOne->SetMomentum(theMomentum);
208 
209  theOne->SetMomentum(aMomentum);
210  theResult.Get()->AddSecondary(theOne);
211  }
212 
213  // Now fill in the gammas.
214  for(i=0; i<nPhotons; i++)
215  {
216  // back to lab system
217  G4DynamicParticle * theOne = new G4DynamicParticle;
218  theOne->SetDefinition(thePhotons->operator[](i)->GetDefinition());
219  theOne->SetMomentum(thePhotons->operator[](i)->GetMomentum());
220  theResult.Get()->AddSecondary(theOne);
221  delete thePhotons->operator[](i);
222  }
223  delete thePhotons;
224 
225 //101203TK
226  G4bool residual = false;
228  ->GetIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0);
229  for ( G4int j = 0 ; j != theResult.Get()->GetNumberOfSecondaries() ; j++ )
230  {
231  if ( theResult.Get()->GetSecondary(j)->GetParticle()->GetDefinition() == aRecoil ) residual = true;
232  }
233 
234  if ( residual == false )
235  {
236  G4int nNonZero = 0;
237  G4LorentzVector p_photons(0,0,0,0);
238  for ( G4int j = 0 ; j != theResult.Get()->GetNumberOfSecondaries() ; j++ )
239  {
240  p_photons += theResult.Get()->GetSecondary(j)->GetParticle()->Get4Momentum();
241  // To many 0 momentum photons -> Check PhotonDist
242  if ( theResult.Get()->GetSecondary(j)->GetParticle()->Get4Momentum().e() > 0 ) nNonZero++;
243  }
244 
245  // Can we include kinetic energy here?
246  G4double deltaE = ( theTrack.Get4Momentum().e() + theTarget.GetTotalEnergy() )
247  - ( p_photons.e() + aRecoil->GetPDGMass() );
248 
249 //Add photons
250  if ( nPhotons - nNonZero > 0 )
251  {
252  //G4cout << "TKDB G4ParticleHPCaptureFS::ApplyYourself we will create additional " << nPhotons - nNonZero << " photons" << G4endl;
253  std::vector<G4double> vRand;
254  vRand.push_back( 0.0 );
255  for ( G4int j = 0 ; j != nPhotons - nNonZero - 1 ; j++ )
256  {
257  vRand.push_back( G4UniformRand() );
258  }
259  vRand.push_back( 1.0 );
260  std::sort( vRand.begin(), vRand.end() );
261 
262  std::vector<G4double> vEPhoton;
263  for ( G4int j = 0 ; j < (G4int)vRand.size() - 1 ; j++ )
264  {
265  vEPhoton.push_back( deltaE * ( vRand[j+1] - vRand[j] ) );
266  }
267  std::sort( vEPhoton.begin(), vEPhoton.end() );
268 
269  for ( G4int j = 0 ; j < nPhotons - nNonZero - 1 ; j++ )
270  {
271  //Isotopic in LAB OK?
272  G4double theta = pi*G4UniformRand();
273  G4double phi = twopi*G4UniformRand();
274  G4double sinth = std::sin(theta);
275  G4double en = vEPhoton[j];
276  G4ThreeVector tempVector(en*sinth*std::cos(phi), en*sinth*std::sin(phi), en*std::cos(theta) );
277 
278  p_photons += G4LorentzVector ( tempVector, tempVector.mag() );
279  G4DynamicParticle * theOne = new G4DynamicParticle;
280  theOne->SetDefinition( G4Gamma::Gamma() );
281  theOne->SetMomentum( tempVector );
282  theResult.Get()->AddSecondary(theOne);
283  }
284 
285 // Add last photon
286  G4DynamicParticle * theOne = new G4DynamicParticle;
287  theOne->SetDefinition( G4Gamma::Gamma() );
288 // For better momentum conservation
289  G4ThreeVector lastPhoton = -p_photons.vect().unit()*vEPhoton.back();
290  p_photons += G4LorentzVector( lastPhoton , lastPhoton.mag() );
291  theOne->SetMomentum( lastPhoton );
292  theResult.Get()->AddSecondary(theOne);
293  }
294 
295 //Add residual
296  G4DynamicParticle * theOne = new G4DynamicParticle;
297  G4ThreeVector aMomentum = theTrack.Get4Momentum().vect() + theTarget.GetMomentum()
298  - p_photons.vect();
299  theOne->SetDefinition(aRecoil);
300  theOne->SetMomentum( aMomentum );
301  theResult.Get()->AddSecondary(theOne);
302 
303  }
304 //101203TK END
305 
306 // clean up the primary neutron
308  return theResult.Get();
309  }
virtual void SetICM(G4bool)
static G4ParticleHPManager * GetInstance()
G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)
static G4double GetNuclearMass(const G4double A, const G4double Z)
G4Cache< G4HadFinalState * > theResult
void SetMomentum(const G4ThreeVector &momentum)
G4HadSecondary * GetSecondary(size_t i)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetProjectileRP(G4ReactionProduct &aIncidentPart)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:503
value_type & Get() const
Definition: G4Cache.hh:282
static double Q[]
static const G4double eps
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
G4ReactionProductVector * Sample(G4double anEnergy)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
static constexpr double twopi
Definition: G4SIunits.hh:76
void SetStatusChange(G4HadFinalStateStatus aS)
std::vector< G4ReactionProduct * > G4ReactionProductVector
G4ReactionProduct GetBiasedThermalNucleus(G4double aMass, G4ThreeVector aVelocity, G4double temp=-1) const
Definition: G4Nucleus.cc:113
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ParticleDefinition * GetDefinition() const
bool G4bool
Definition: G4Types.hh:79
G4double GetIonMass(G4int Z, G4int A, G4int L=0, G4int lvl=0) const
Definition: G4IonTable.cc:1335
G4double GetKineticEnergy() const
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:63
static constexpr double eV
Definition: G4SIunits.hh:215
#define TRUE
Definition: globals.hh:55
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4LorentzVector & Get4Momentum() const
G4LorentzVector Get4Momentum() const
static G4IonTable * GetIonTable()
Definition: G4IonTable.hh:78
G4double GetTotalEnergy() const
G4double GetPDGMass() const
Hep3Vector unit() const
G4DynamicParticle * GetParticle()
void SetTarget(G4ReactionProduct &aTarget)
G4ReactionProductVector * GetPhotons(G4double anEnergy)
G4ThreeVector GetMomentum() const
G4double GetTemperature() const
Definition: G4Material.hh:182
static constexpr double MeV
Definition: G4SIunits.hh:214
const G4Material * GetMaterial() const
static constexpr double pi
Definition: G4SIunits.hh:75
void AddSecondary(G4DynamicParticle *aP, G4int mod=-1)
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
double mag() const
void Put(const value_type &val) const
Definition: G4Cache.hh:286
G4int GetNumberOfSecondaries() const
CLHEP::HepLorentzVector G4LorentzVector

Here is the call graph for this function:

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

Implements G4ParticleHPFinalState.

Definition at line 312 of file G4ParticleHPCaptureFS.cc.

313  {
314 
315  //TK110430 BEGIN
316  std::stringstream ss;
317  ss << static_cast<G4int>(Z);
318  G4String sZ;
319  ss >> sZ;
320  ss.clear();
321  ss << static_cast<G4int>(A);
322  G4String sA;
323  ss >> sA;
324 
325  ss.clear();
326  G4String sM;
327  if ( M > 0 )
328  {
329  ss << "m";
330  ss << M;
331  ss >> sM;
332  ss.clear();
333  }
334 
335  G4String element_name = theNames.GetName( static_cast<G4int>(Z)-1 );
336  G4String filenameMF6 = dirName+"/FSMF6/"+sZ+"_"+sA+sM+"_"+element_name;
337  //std::ifstream dummyIFS(filenameMF6, std::ios::in);
338  //if ( dummyIFS.good() == true ) hasExactMF6=true;
339  std::istringstream theData(std::ios::in);
340  G4ParticleHPManager::GetInstance()->GetDataStream(filenameMF6,theData);
341 
342  //TK110430 Only use MF6MT102 which has exactly same A and Z
343  //Even _nat_ do not select and there is no _nat_ case in ENDF-VII.0
344  if ( theData.good() == true ) {
345  hasExactMF6=true;
346  theMF6FinalState.Init(theData);
347  //theData.close();
348  return;
349  }
350  //TK110430 END
351 
352 
353  G4String tString = "/FS";
354  G4bool dbool;
355  G4ParticleHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), M, dirName, tString, dbool);
356 
357  G4String filename = aFile.GetName();
358  SetAZMs( A, Z, M, aFile );
359  //theBaseA = A;
360  //theBaseZ = G4int(Z+.5);
361  if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
362  {
363  hasAnyData = false;
364  hasFSData = false;
365  hasXsec = false;
366  return;
367  }
368  //std::ifstream theData(filename, std::ios::in);
369  //std::istringstream theData(std::ios::in);
370  theData.clear();
371  G4ParticleHPManager::GetInstance()->GetDataStream(filename,theData);
372  hasFSData = theFinalStatePhotons.InitMean(theData);
373  if(hasFSData)
374  {
375  targetMass = theFinalStatePhotons.GetTargetMass();
376  theFinalStatePhotons.InitAngular(theData);
377  theFinalStatePhotons.InitEnergies(theData);
378  }
379  //theData.close();
380  }
static G4ParticleHPManager * GetInstance()
void SetAZMs(G4double anA, G4double aZ, G4int aM, G4ParticleHPDataUsed used)
void GetDataStream(G4String, std::istringstream &iss)
double A(double temperature)
G4bool InitMean(std::istream &aDataFile)
bool G4bool
Definition: G4Types.hh:79
G4ParticleHPDataUsed GetName(G4int A, G4int Z, G4String base, G4String rest, G4bool &active)
void InitEnergies(std::istream &aDataFile)
void InitAngular(std::istream &aDataFile)
void Init(std::istream &aDataFile)

Here is the call graph for this function:

G4ParticleHPFinalState* G4ParticleHPCaptureFS::New ( )
inlinevirtual

Implements G4ParticleHPFinalState.

Definition at line 59 of file G4ParticleHPCaptureFS.hh.

Here is the call graph for this function:


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