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

Last change on this file since 847 was 819, checked in by garnier, 16 years ago

import all except CVS

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//
33#include "G4NeutronHPCaptureFS.hh"
34#include "G4Gamma.hh"
35#include "G4ReactionProduct.hh"
36#include "G4Nucleus.hh"
37#include "G4PhotonEvaporation.hh"
38#include "G4Fragment.hh"
39#include "G4ParticleTable.hh"
40#include "G4NeutronHPDataUsed.hh"
41
42  G4HadFinalState * G4NeutronHPCaptureFS::ApplyYourself(const G4HadProjectile & theTrack)
43  {
44    G4int i;
45    theResult.Clear();
46// prepare neutron
47    G4double eKinetic = theTrack.GetKineticEnergy();
48    const G4HadProjectile *incidentParticle = &theTrack;
49    G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
50    theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
51    theNeutron.SetKineticEnergy( eKinetic );
52
53// prepare target
54    G4ReactionProduct theTarget; 
55    G4Nucleus aNucleus;
56    G4double eps = 0.0001;
57    if(targetMass<500*MeV)
58      targetMass = ( G4NucleiPropertiesTable::GetNuclearMass(static_cast<G4int>(theBaseZ+eps), static_cast<G4int>(theBaseA+eps))) /
59                     G4Neutron::Neutron()->GetPDGMass();
60    G4ThreeVector neutronVelocity = 1./G4Neutron::Neutron()->GetPDGMass()*theNeutron.GetMomentum();
61    G4double temperature = theTrack.GetMaterial()->GetTemperature();
62    theTarget = aNucleus.GetBiasedThermalNucleus(targetMass, neutronVelocity, temperature);
63
64// go to nucleus rest system
65    theNeutron.Lorentz(theNeutron, -1*theTarget);
66    eKinetic = theNeutron.GetKineticEnergy();
67
68// dice the photons
69
70    G4ReactionProductVector * thePhotons = 0;
71    if ( HasFSData() && !getenv ( "G4NEUTRONHP_USE_ONLY_PHOTONEVAPORATION" ) ) 
72    { 
73      thePhotons = theFinalStatePhotons.GetPhotons(eKinetic);
74    }
75    else
76    {
77      G4ThreeVector aCMSMomentum = theNeutron.GetMomentum()+theTarget.GetMomentum();
78      G4LorentzVector p4(aCMSMomentum, theTarget.GetTotalEnergy() + theNeutron.GetTotalEnergy());
79      G4Fragment nucleus(static_cast<G4int>(theBaseA+1), static_cast<G4int>(theBaseZ) ,p4);
80      G4PhotonEvaporation photonEvaporation;
81      // T. K. add
82      photonEvaporation.SetICM( TRUE );
83      G4FragmentVector* products = photonEvaporation.BreakItUp(nucleus);
84      G4FragmentVector::iterator i;
85      thePhotons = new G4ReactionProductVector;
86      for(i=products->begin(); i!=products->end(); i++)
87      {
88        G4ReactionProduct * theOne = new G4ReactionProduct;
89        // T. K. add
90        if ( (*i)->GetParticleDefinition() != 0 ) 
91           theOne->SetDefinition( (*i)->GetParticleDefinition() );
92        else
93           theOne->SetDefinition( G4Gamma::Gamma() ); // this definiion will be over writen
94       
95        // T. K. comment out below line
96        //theOne->SetDefinition( G4Gamma::Gamma() );
97        G4ParticleTable* theTable = G4ParticleTable::GetParticleTable();
98        if((*i)->GetMomentum().mag() > 10*MeV) 
99                 theOne->SetDefinition( 
100                 theTable->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ)) );
101        theOne->SetMomentum( (*i)->GetMomentum().vect() ) ;
102        theOne->SetTotalEnergy( (*i)->GetMomentum().t() );
103        thePhotons->push_back(theOne);
104        delete *i;
105      } 
106      delete products;
107    }
108
109// add them to the final state
110
111    G4int nPhotons = 0;
112    if(thePhotons!=0) nPhotons=thePhotons->size();
113    G4int nParticles = nPhotons;
114    if(1==nPhotons) nParticles = 2;
115
116    // back to lab system
117    for(i=0; i<nPhotons; i++)
118    {
119      thePhotons->operator[](i)->Lorentz(*(thePhotons->operator[](i)), theTarget);
120    }
121   
122    // Recoil, if only one gamma
123    if (1==nPhotons)
124    {
125       G4DynamicParticle * theOne = new G4DynamicParticle;
126       G4ParticleDefinition * aRecoil = G4ParticleTable::GetParticleTable()
127                                        ->FindIon(static_cast<G4int>(theBaseZ), static_cast<G4int>(theBaseA+1), 0, static_cast<G4int>(theBaseZ));
128       theOne->SetDefinition(aRecoil);
129       // Now energy;
130       // Can be done slightly better @
131       G4ThreeVector aMomentum =  theTrack.Get4Momentum().vect()
132                                 +theTarget.GetMomentum()
133                                 -thePhotons->operator[](0)->GetMomentum();
134
135       G4ThreeVector theMomUnit = aMomentum.unit();
136       G4double aKinEnergy =  theTrack.GetKineticEnergy()
137                             +theTarget.GetKineticEnergy(); // gammas come from Q-value
138       G4double theResMass = aRecoil->GetPDGMass();
139       G4double theResE = aRecoil->GetPDGMass()+aKinEnergy;
140       G4double theAbsMom = std::sqrt(theResE*theResE - theResMass*theResMass);
141       G4ThreeVector theMomentum = theAbsMom*theMomUnit;
142       theOne->SetMomentum(theMomentum);
143       theResult.AddSecondary(theOne);
144    }
145
146    // Now fill in the gammas.
147    for(i=0; i<nPhotons; i++)
148    {
149      // back to lab system
150      G4DynamicParticle * theOne = new G4DynamicParticle;
151      theOne->SetDefinition(thePhotons->operator[](i)->GetDefinition());
152      theOne->SetMomentum(thePhotons->operator[](i)->GetMomentum());
153      theResult.AddSecondary(theOne);
154      delete thePhotons->operator[](i);
155    }
156    delete thePhotons; 
157// clean up the primary neutron
158    theResult.SetStatusChange(stopAndKill);
159    return &theResult;
160  }
161
162  void G4NeutronHPCaptureFS::Init (G4double A, G4double Z, G4String & dirName, G4String & )
163  {
164    G4String tString = "/FS/";
165    G4bool dbool;
166    G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), dirName, tString, dbool);
167    G4String filename = aFile.GetName();
168    theBaseA = A;
169    theBaseZ = G4int(Z+.5);
170    if(!dbool || ( Z<2.5 && ( std::abs(theBaseZ - Z)>0.0001 || std::abs(theBaseA - A)>0.0001)))
171    {
172      hasAnyData = false;
173      hasFSData = false; 
174      hasXsec = false;
175      return;
176    }
177    std::ifstream theData(filename, std::ios::in);
178   
179    hasFSData = theFinalStatePhotons.InitMean(theData); 
180    if(hasFSData)
181    {
182      targetMass = theFinalStatePhotons.GetTargetMass();
183      theFinalStatePhotons.InitAngular(theData); 
184      theFinalStatePhotons.InitEnergies(theData); 
185    }
186    theData.close();
187  }
Note: See TracBrowser for help on using the repository browser.