source: trunk/source/processes/hadronic/models/neutron_hp/src/G4NeutronHPCaptureFS.cc @ 1196

Last change on this file since 1196 was 962, checked in by garnier, 15 years ago

update processes

File size: 7.8 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26// neutron_hp -- source file
27// J.P. Wellisch, Nov-1996
28// A prototype of the low energy neutron transport model.
29//
30// 12-April-06 Enable IC electron emissions T. Koi
31// 26-January-07 Add G4NEUTRONHP_USE_ONLY_PHOTONEVAPORATION flag
32// 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
33//
34#include "G4NeutronHPCaptureFS.hh"
35#include "G4Gamma.hh"
36#include "G4ReactionProduct.hh"
37#include "G4Nucleus.hh"
38#include "G4PhotonEvaporation.hh"
39#include "G4Fragment.hh"
40#include "G4ParticleTable.hh"
41#include "G4NeutronHPDataUsed.hh"
42
43  G4HadFinalState * G4NeutronHPCaptureFS::ApplyYourself(const G4HadProjectile & theTrack)
44  {
45    G4int i;
46    theResult.Clear();
47// prepare neutron
48    G4double eKinetic = theTrack.GetKineticEnergy();
49    const G4HadProjectile *incidentParticle = &theTrack;
50    G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
51    theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
52    theNeutron.SetKineticEnergy( eKinetic );
53
54// prepare target
55    G4ReactionProduct theTarget; 
56    G4Nucleus aNucleus;
57    G4double eps = 0.0001;
58    if(targetMass<500*MeV)
59      targetMass = ( G4NucleiProperties::GetNuclearMass( static_cast<G4int>(theBaseA+eps) , static_cast<G4int>(theBaseZ+eps) )) /
60                     G4Neutron::Neutron()->GetPDGMass();
61    G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
62    G4double temperature = theTrack.GetMaterial()->GetTemperature();
63    theTarget = aNucleus.GetBiasedThermalNucleus(targetMass, neutronVelocity, temperature);
64
65// go to nucleus rest system
66    theNeutron.Lorentz(theNeutron, -1*theTarget);
67    eKinetic = theNeutron.GetKineticEnergy();
68
69// dice the photons
70
71    G4ReactionProductVector * thePhotons = 0;
72    if ( HasFSData() && !getenv ( "G4NEUTRONHP_USE_ONLY_PHOTONEVAPORATION" ) ) 
73    { 
74      thePhotons = theFinalStatePhotons.GetPhotons(eKinetic);
75    }
76    else
77    {
78      G4ThreeVector aCMSMomentum = theNeutron.GetMomentum()+theTarget.GetMomentum();
79      G4LorentzVector p4(aCMSMomentum, theTarget.GetTotalEnergy() + theNeutron.GetTotalEnergy());
80      G4Fragment nucleus(static_cast<G4int>(theBaseA+1), static_cast<G4int>(theBaseZ) ,p4);
81      G4PhotonEvaporation photonEvaporation;
82      // T. K. add
83      photonEvaporation.SetICM( TRUE );
84      G4FragmentVector* products = photonEvaporation.BreakItUp(nucleus);
85      G4FragmentVector::iterator i;
86      thePhotons = new G4ReactionProductVector;
87      for(i=products->begin(); i!=products->end(); i++)
88      {
89        G4ReactionProduct * theOne = new G4ReactionProduct;
90        // T. K. add
91        if ( (*i)->GetParticleDefinition() != 0 ) 
92           theOne->SetDefinition( (*i)->GetParticleDefinition() );
93        else
94           theOne->SetDefinition( G4Gamma::Gamma() ); // this definiion will be over writen
95       
96        // T. K. comment out below line
97        //theOne->SetDefinition( G4Gamma::Gamma() );
98        G4ParticleTable* theTable = G4ParticleTable::GetParticleTable();
99        if((*i)->GetMomentum().mag() > 10*MeV) 
100                 theOne->SetDefinition( 
101                 theTable->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ)) );
102        theOne->SetMomentum( (*i)->GetMomentum().vect() ) ;
103        theOne->SetTotalEnergy( (*i)->GetMomentum().t() );
104        thePhotons->push_back(theOne);
105        delete *i;
106      } 
107      delete products;
108    }
109
110// add them to the final state
111
112    G4int nPhotons = 0;
113    if(thePhotons!=0) nPhotons=thePhotons->size();
114    G4int nParticles = nPhotons;
115    if(1==nPhotons) nParticles = 2;
116
117    // back to lab system
118    for(i=0; i<nPhotons; i++)
119    {
120      thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), theTarget);
121    }
122   
123    // Recoil, if only one gamma
124    if (1==nPhotons)
125    {
126       G4DynamicParticle * theOne = new G4DynamicParticle;
127       G4ParticleDefinition * aRecoil = G4ParticleTable::GetParticleTable()
128                                        ->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ));
129       theOne->SetDefinition(aRecoil);
130       // Now energy;
131       // Can be done slightly better @
132       G4ThreeVector aMomentum =  theTrack.Get4Momentum().vect()
133                                 +theTarget.GetMomentum()
134                                 -thePhotons->operator[](0)->GetMomentum();
135
136       G4ThreeVector theMomUnit = aMomentum.unit();
137       G4double aKinEnergy =  theTrack.GetKineticEnergy()
138                             +theTarget.GetKineticEnergy(); // gammas come from Q-value
139       G4double theResMass = aRecoil->GetPDGMass();
140       G4double theResE = aRecoil->GetPDGMass()+aKinEnergy;
141       G4double theAbsMom = std::sqrt(theResE*theResE - theResMass*theResMass);
142       G4ThreeVector theMomentum = theAbsMom*theMomUnit;
143       theOne->SetMomentum(theMomentum);
144       theResult.AddSecondary(theOne);
145    }
146
147    // Now fill in the gammas.
148    for(i=0; i<nPhotons; i++)
149    {
150      // back to lab system
151      G4DynamicParticle * theOne = new G4DynamicParticle;
152      theOne->SetDefinition(thePhotons->operator[](i)->GetDefinition());
153      theOne->SetMomentum(thePhotons->operator[](i)->GetMomentum());
154      theResult.AddSecondary(theOne);
155      delete thePhotons->operator[](i);
156    }
157    delete thePhotons; 
158// clean up the primary neutron
159    theResult.SetStatusChange(stopAndKill);
160    return &theResult;
161  }
162
163  void G4NeutronHPCaptureFS::Init (G4double A, G4double Z, G4String & dirName, G4String & )
164  {
165    G4String tString = "/FS/";
166    G4bool dbool;
167    G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), dirName, tString, dbool);
168    G4String filename = aFile.GetName();
169    theBaseA = A;
170    theBaseZ = G4int(Z+.5);
171    if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
172    {
173      hasAnyData = false;
174      hasFSData = false; 
175      hasXsec = false;
176      return;
177    }
178    std::ifstream theData(filename, std::ios::in);
179   
180    hasFSData = theFinalStatePhotons.InitMean(theData); 
181    if(hasFSData)
182    {
183      targetMass = theFinalStatePhotons.GetTargetMass();
184      theFinalStatePhotons.InitAngular(theData); 
185      theFinalStatePhotons.InitEnergies(theData); 
186    }
187    theData.close();
188  }
Note: See TracBrowser for help on using the repository browser.