source: trunk/source/processes/hadronic/models/lll_fission/src/G4FissionLibrary.cc@ 1350

Last change on this file since 1350 was 819, checked in by garnier, 17 years ago

import all except CVS

  • Property svn:executable set to *
File size: 10.1 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//
27// This software was developed by Lawrence Livermore National Laboratory.
28//
29// Redistribution and use in source and binary forms, with or without
30// modification, are permitted provided that the following conditions are met:
31//
32// 1. Redistributions of source code must retain the above copyright notice,
33// this list of conditions and the following disclaimer.
34// 2. Redistributions in binary form must reproduce the above copyright notice,
35// this list of conditions and the following disclaimer in the documentation
36// and/or other materials provided with the distribution.
37// 3. The name of the author may not be used to endorse or promote products
38// derived from this software without specific prior written permission.
39//
40// THIS SOFTWARE IS PROVIDED BY THE AUTHOR ``AS IS'' AND ANY EXPRESS OR IMPLIED
41// WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED WARRANTIES OF
42// MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO
43// EVENT SHALL THE AUTHOR BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL,
44// SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
45// PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR PROFITS;
46// OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY,
47// WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR
48// OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF
49// ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
50//
51// Copyright (c) 2006 The Regents of the University of California.
52// All rights reserved.
53// UCRL-CODE-224807
54//
55//
56// $Id: G4FissionLibrary.cc,v 1.4 2007/06/01 14:02:08 gcosmo Exp $
57//
58// neutron_hp -- source file
59// J.M. Verbeke, Jan-2007
60// A low energy neutron-induced fission model.
61//
62
63#include "G4FissionLibrary.hh"
64
65G4FissionLibrary::G4FissionLibrary()
66 : G4NeutronHPFinalState()
67{
68 hasXsec = false;
69}
70
71G4FissionLibrary::~G4FissionLibrary()
72{
73}
74
75G4NeutronHPFinalState * G4FissionLibrary::New()
76{
77 G4FissionLibrary * theNew = new G4FissionLibrary;
78 return theNew;
79}
80
81void G4FissionLibrary::Init (G4double A, G4double Z, G4String & dirName, G4String &)
82{
83 G4String tString = "/FS/";
84 G4bool dbool;
85 theIsotope = static_cast<G4int>(1000*Z+A);
86 G4NeutronHPDataUsed aFile = theNames.GetName(static_cast<G4int>(A), static_cast<G4int>(Z), dirName, tString, dbool);
87 G4String filename = aFile.GetName();
88
89 if(!dbool)
90 {
91 hasAnyData = false;
92 hasFSData = false;
93 hasXsec = false;
94 return;
95 }
96 std::ifstream theData(filename, std::ios::in);
97
98 // here it comes
99 G4int infoType, dataType;
100 hasFSData = false;
101 while (theData >> infoType)
102 {
103 hasFSData = true;
104 theData >> dataType;
105 switch(infoType)
106 {
107 case 1:
108 if(dataType==4) theNeutronAngularDis.Init(theData);
109 if(dataType==5) thePromptNeutronEnDis.Init(theData);
110 if(dataType==12) theFinalStatePhotons.InitMean(theData);
111 if(dataType==14) theFinalStatePhotons.InitAngular(theData);
112 if(dataType==15) theFinalStatePhotons.InitEnergies(theData);
113 break;
114 case 2:
115 if(dataType==1) theFinalStateNeutrons.InitMean(theData);
116 break;
117 case 3:
118 if(dataType==1) theFinalStateNeutrons.InitDelayed(theData);
119 if(dataType==5) theDelayedNeutronEnDis.Init(theData);
120 break;
121 case 4:
122 if(dataType==1) theFinalStateNeutrons.InitPrompt(theData);
123 break;
124 case 5:
125 if(dataType==1) theEnergyRelease.Init(theData);
126 break;
127 default:
128 G4cout << "G4FissionLibrary::Init: unknown data type"<<dataType<<G4endl;
129 throw G4HadronicException(__FILE__, __LINE__, "G4FissionLibrary::Init: unknown data type");
130 break;
131 }
132 }
133 targetMass = theFinalStateNeutrons.GetTargetMass();
134 theData.close();
135}
136
137G4HadFinalState * G4FissionLibrary::ApplyYourself(const G4HadProjectile & theTrack)
138{
139 theResult.Clear();
140
141// prepare neutron
142 G4double eKinetic = theTrack.GetKineticEnergy();
143 const G4HadProjectile *incidentParticle = &theTrack;
144 G4ReactionProduct theNeutron( const_cast<G4ParticleDefinition *>(incidentParticle->GetDefinition()) );
145 theNeutron.SetMomentum( incidentParticle->Get4Momentum().vect() );
146 theNeutron.SetKineticEnergy( eKinetic );
147
148// prepare target
149 G4Nucleus aNucleus;
150 G4ReactionProduct theTarget;
151 G4ThreeVector neuVelo = (1./incidentParticle->GetDefinition()->GetPDGMass())*theNeutron.GetMomentum();
152 theTarget = aNucleus.GetBiasedThermalNucleus( targetMass, neuVelo, theTrack.GetMaterial()->GetTemperature());
153
154// set neutron and target in the FS classes
155 theNeutronAngularDis.SetNeutron(theNeutron);
156 theNeutronAngularDis.SetTarget(theTarget);
157
158// boost to target rest system
159 theNeutron.Lorentz(theNeutron, -1*theTarget);
160
161 eKinetic = theNeutron.GetKineticEnergy();
162
163// dice neutron and gamma multiplicities, energies and momenta in Lab. @@
164// no energy conservation on an event-to-event basis. we rely on the data to be ok. @@
165// also for mean, we rely on the consistency of the data. @@
166
167 G4int nPrompt=0, gPrompt=0;
168 SampleMult(theTrack, &nPrompt, &gPrompt, eKinetic);
169
170// Build neutrons and add them to dynamic particle vector
171 G4double momentum;
172 for(G4int i=0; i<nPrompt; i++)
173 {
174 G4DynamicParticle * it = new G4DynamicParticle;
175 it->SetDefinition(G4Neutron::Neutron());
176 it->SetKineticEnergy(getneng_(&i)*MeV);
177 momentum = it->GetTotalMomentum();
178 G4ThreeVector temp(momentum*getndircosu_(&i),
179 momentum*getndircosv_(&i),
180 momentum*getndircosw_(&i));
181 it->SetMomentum( temp );
182// it->SetGlobalTime(getnage_(&i)*second);
183 theResult.AddSecondary(it);
184// G4cout <<"G4FissionLibrary::ApplyYourself: energy of prompt neutron " << i << " = " << it->GetKineticEnergy()<<G4endl;
185 }
186
187// Build gammas, lorentz transform them, and add them to dynamic particle vector
188 for(G4int i=0; i<gPrompt; i++)
189 {
190 G4ReactionProduct * thePhoton = new G4ReactionProduct;
191 thePhoton->SetDefinition(G4Gamma::Gamma());
192 thePhoton->SetKineticEnergy(getpeng_(&i)*MeV);
193 momentum = thePhoton->GetTotalMomentum();
194 G4ThreeVector temp(momentum*getpdircosu_(&i),
195 momentum*getpdircosv_(&i),
196 momentum*getpdircosw_(&i));
197 thePhoton->SetMomentum( temp );
198 thePhoton->Lorentz(*thePhoton, -1.*theTarget);
199
200 G4DynamicParticle * it = new G4DynamicParticle;
201 it->SetDefinition(thePhoton->GetDefinition());
202 it->SetMomentum(thePhoton->GetMomentum());
203// it->SetGlobalTime(getpage_(&i)*second);
204// G4cout <<"G4FissionLibrary::ApplyYourself: energy of prompt photon " << i << " = " << it->GetKineticEnergy()<<G4endl;
205 theResult.AddSecondary(it);
206 delete thePhoton;
207 }
208// G4cout <<"G4FissionLibrary::ApplyYourself: Number of secondaries = "<<theResult.GetNumberOfSecondaries()<< G4endl;
209// G4cout <<"G4FissionLibrary::ApplyYourself: Number of induced prompt neutron = "<<nPrompt<<G4endl;
210// G4cout <<"G4FissionLibrary::ApplyYourself: Number of induced prompt photons = "<<gPrompt<<G4endl;
211
212// finally deal with local energy depositions.
213 G4double eDepByFragments = theEnergyRelease.GetFragmentKinetic();
214 theResult.SetLocalEnergyDeposit(eDepByFragments);
215// G4cout << "G4FissionLibrary::local energy deposit" << eDepByFragments<<G4endl;
216// clean up the primary neutron
217 theResult.SetStatusChange(stopAndKill);
218 return &theResult;
219}
220
221void G4FissionLibrary::SampleMult(const G4HadProjectile & theTrack, G4int* nPrompt,
222 G4int* gPrompt, G4double eKinetic)
223{
224 G4double promptNeutronMulti = 0;
225 promptNeutronMulti = theFinalStateNeutrons.GetPrompt(eKinetic); // prompt nubar from Geant
226 G4double delayedNeutronMulti = 0;
227 delayedNeutronMulti = theFinalStateNeutrons.GetDelayed(eKinetic); // delayed nubar from Geant
228
229 G4double time = theTrack.GetGlobalTime()/second;
230 if(delayedNeutronMulti==0&&promptNeutronMulti==0) {
231 // no data for prompt and delayed neutrons in Geant
232 // but there is perhaps data for the total neutron multiplicity, in which case
233 // we use it for prompt neutron emission
234 G4double totalNeutronMulti = theFinalStateNeutrons.GetMean(eKinetic);
235 genfissevt_(&theIsotope, &time, &totalNeutronMulti, &eKinetic);
236 } else {
237 // prompt nubar != 0 || delayed nubar != 0
238 genfissevt_(&theIsotope, &time, &promptNeutronMulti, &eKinetic);
239 }
240 *nPrompt = getnnu_();
241 if (*nPrompt == -1) *nPrompt = 0; // the fission library libFission.a has no data for neutrons
242 *gPrompt = getpnu_();
243 if (*gPrompt == -1) *gPrompt = 0; // the fission library libFission.a has no data for gammas
244}
245
Note: See TracBrowser for help on using the repository browser.