source: trunk/source/processes/hadronic/models/de_excitation/fission/src/G4CompetitiveFission.cc @ 962

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

update processes

File size: 19.0 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// $Id: G4CompetitiveFission.cc,v 1.9 2008/11/20 13:46:27 dennis Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30// Hadronic Process: Nuclear De-excitations
31// by V. Lara (Oct 1998)
32
33#include "G4CompetitiveFission.hh"
34#include "G4PairingCorrection.hh"
35
36G4CompetitiveFission::G4CompetitiveFission() : G4VEvaporationChannel("fission")
37{
38    theFissionBarrierPtr = new G4FissionBarrier;
39    MyOwnFissionBarrier = true;
40
41    theFissionProbabilityPtr = new G4FissionProbability;
42    MyOwnFissionProbability = true;
43 
44    theLevelDensityPtr = new G4FissionLevelDensityParameter;
45    MyOwnLevelDensity = true;
46
47    MaximalKineticEnergy = -1000.0*MeV;
48    FissionBarrier = 0.0;
49    FissionProbability = 0.0;
50    LevelDensityParameter = 0.0;
51}
52
53G4CompetitiveFission::G4CompetitiveFission(const G4CompetitiveFission &) : G4VEvaporationChannel()
54{
55}
56
57G4CompetitiveFission::~G4CompetitiveFission()
58{
59    if (MyOwnFissionBarrier) delete theFissionBarrierPtr;
60
61    if (MyOwnFissionProbability) delete theFissionProbabilityPtr;
62
63    if (MyOwnLevelDensity) delete theLevelDensityPtr;
64}
65
66const G4CompetitiveFission & G4CompetitiveFission::operator=(const G4CompetitiveFission &)
67{
68    throw G4HadronicException(__FILE__, __LINE__, "G4CompetitiveFission::operator= meant to not be accessable");
69    return *this;
70}
71
72G4bool G4CompetitiveFission::operator==(const G4CompetitiveFission &right) const
73{
74    return (this == (G4CompetitiveFission *) &right);
75}
76
77G4bool G4CompetitiveFission::operator!=(const G4CompetitiveFission &right) const
78{
79    return (this != (G4CompetitiveFission *) &right);
80}
81
82
83
84
85void G4CompetitiveFission::Initialize(const G4Fragment & fragment)
86{
87    G4int anA = static_cast<G4int>(fragment.GetA());
88    G4int aZ = static_cast<G4int>(fragment.GetZ());
89    G4double ExEnergy = fragment.GetExcitationEnergy() - 
90           G4PairingCorrection::GetInstance()->GetFissionPairingCorrection(anA,aZ);
91 
92
93    // Saddle point excitation energy ---> A = 65
94    // Fission is excluded for A < 65
95    if (anA >= 65 && ExEnergy > 0.0) {
96        FissionBarrier = theFissionBarrierPtr->FissionBarrier(anA,aZ,ExEnergy);
97        MaximalKineticEnergy = ExEnergy - FissionBarrier;
98        LevelDensityParameter = theLevelDensityPtr->LevelDensityParameter(anA,aZ,ExEnergy);
99        FissionProbability = theFissionProbabilityPtr->EmissionProbability(fragment,MaximalKineticEnergy);
100    }
101    else {
102        MaximalKineticEnergy = -1000.0*MeV;
103        LevelDensityParameter = 0.0;
104        FissionProbability = 0.0;
105    }
106
107    return;
108}
109
110
111
112G4FragmentVector * G4CompetitiveFission::BreakUp(const G4Fragment & theNucleus)
113{
114    // Nucleus data
115    // Atomic number of nucleus
116    G4int A = static_cast<G4int>(theNucleus.GetA());
117    // Charge of nucleus
118    G4int Z = static_cast<G4int>(theNucleus.GetZ());
119    //   Excitation energy (in MeV)
120    G4double U = theNucleus.GetExcitationEnergy() - 
121        G4PairingCorrection::GetInstance()->GetFissionPairingCorrection(A,Z);
122    // Check that U > 0
123    if (U <= 0.0) {
124        G4FragmentVector * theResult = new  G4FragmentVector;
125        theResult->push_back(new G4Fragment(theNucleus));
126        return theResult;
127    }
128
129    // Atomic Mass of Nucleus (in MeV)
130    G4double M = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(Z,A)/MeV;
131    // Nucleus Momentum
132    G4LorentzVector theNucleusMomentum = theNucleus.GetMomentum();
133
134    // Calculate fission parameters
135    G4FissionParameters theParameters(A,Z,U,FissionBarrier);
136 
137    // First fragment
138    G4int A1 = 0;
139    G4int Z1 = 0;
140    G4double M1 = 0.0;
141
142    // Second fragment
143    G4int A2 = 0;
144    G4int Z2 = 0;
145    G4double M2 = 0.0;
146
147    G4double FragmentsExcitationEnergy = 0.0;
148    G4double FragmentsKineticEnergy = 0.0;
149
150    G4int Trials = 0;
151    do {
152
153        // First fragment
154        A1 = FissionAtomicNumber(A,theParameters);
155        Z1 = FissionCharge(A,Z,A1);
156        M1 = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(Z1,A1);
157
158
159        // Second Fragment
160        A2 = A - A1;
161        Z2 = Z - Z1;
162        if (A2 < 1 || Z2 < 0) 
163            throw G4HadronicException(__FILE__, __LINE__, "G4CompetitiveFission::BreakUp: Can't define second fragment! ");
164        M2 = G4ParticleTable::GetParticleTable()->GetIonTable()->GetIonMass(Z2,A2)/MeV;
165
166        // Check that fragment masses are less or equal than total energy
167        //  if (M1 + M2 > theNucleusMomentum.mag()/MeV)
168        if (M1 + M2 > theNucleusMomentum.e()/MeV)
169            throw G4HadronicException(__FILE__, __LINE__, "G4CompetitiveFission::BreakUp: Fragments Mass > Total Energy");
170
171        // Maximal Kinetic Energy (available energy for fragments)
172        //  G4double Tmax = theNucleusMomentum.mag()/MeV - M1 - M2;
173        G4double Tmax = M + U - M1 - M2;
174
175        FragmentsKineticEnergy = FissionKineticEnergy( A , Z,
176                                                       A1, Z1,
177                                                       A2, Z2,
178                                                       U , Tmax,
179                                                       theParameters);
180   
181        // Excitation Energy
182        FragmentsExcitationEnergy = Tmax - FragmentsKineticEnergy;
183   
184    } while (FragmentsExcitationEnergy < 0.0 && Trials++ < 100);
185 
186 
187
188    if (FragmentsExcitationEnergy <= 0.0) 
189        throw G4HadronicException(__FILE__, __LINE__, "G4CompetitiveFission::BreakItUp: Excitation energy for fragments < 0.0!");
190 
191
192  // while (FragmentsExcitationEnergy < 0 && Trials < 100);
193 
194  // Fragment 1
195    G4double U1 = FragmentsExcitationEnergy * (static_cast<G4double>(A1)/static_cast<G4double>(A));
196    // Fragment 2
197    G4double U2 = FragmentsExcitationEnergy * (static_cast<G4double>(A2)/static_cast<G4double>(A));
198
199
200    G4double Pmax = std::sqrt( 2 * ( ( (M1+U1)*(M2+U2) ) /
201                                ( (M1+U1)+(M2+U2) ) ) * FragmentsKineticEnergy);
202
203    G4ParticleMomentum momentum1 = IsotropicVector( Pmax );
204    G4ParticleMomentum momentum2( -momentum1 );
205
206    // Perform a Galileo boost for fragments
207    momentum1 += (theNucleusMomentum.boostVector() * (M1+U1));
208    momentum2 += (theNucleusMomentum.boostVector() * (M2+U2));
209
210
211    // Create 4-momentum for first fragment
212    // Warning!! Energy conservation is broken
213    G4LorentzVector FourMomentum1( momentum1 , std::sqrt(momentum1.mag2() + (M1+U1)*(M1+U1)));
214
215    // Create 4-momentum for second fragment
216    // Warning!! Energy conservation is broken
217    G4LorentzVector FourMomentum2( momentum2 , std::sqrt(momentum2.mag2() + (M2+U2)*(M2+U2)));
218
219    // Create Fragments
220    G4Fragment * Fragment1 = new G4Fragment( A1, Z1, FourMomentum1);
221    if (!Fragment1) throw G4HadronicException(__FILE__, __LINE__, "G4CompetitiveFission::BreakItUp: Can't create Fragment1! ");
222    G4Fragment * Fragment2 = new G4Fragment( A2, Z2, FourMomentum2);
223    if (!Fragment2) throw G4HadronicException(__FILE__, __LINE__, "G4CompetitiveFission::BreakItUp: Can't create Fragment2! ");
224
225#ifdef PRECOMPOUND_TEST
226    Fragment1->SetCreatorModel(G4String("G4CompetitiveFission"));
227    Fragment2->SetCreatorModel(G4String("G4CompetitiveFission"));
228#endif
229  // Create Fragment Vector
230    G4FragmentVector * theResult = new G4FragmentVector;
231
232    theResult->push_back(Fragment1);
233    theResult->push_back(Fragment2);
234
235#ifdef debug
236    CheckConservation(theNucleus,theResult);
237#endif
238
239    return theResult;
240}
241
242
243
244G4int G4CompetitiveFission::FissionAtomicNumber(const G4int A, const G4FissionParameters & theParam)
245    // Calculates the atomic number of a fission product
246{
247
248    // For Simplicity reading code
249    const G4double A1 = theParam.GetA1();
250    const G4double A2 = theParam.GetA2();
251    const G4double As = theParam.GetAs();
252//    const G4double Sigma1 = theParam.GetSigma1();
253    const G4double Sigma2 = theParam.GetSigma2();
254    const G4double SigmaS = theParam.GetSigmaS();
255    const G4double w = theParam.GetW();
256
257 
258//    G4double FasymAsym = 2.0*std::exp(-((A2-As)*(A2-As))/(2.0*Sigma2*Sigma2)) +
259//      std::exp(-((A1-As)*(A1-As))/(2.0*Sigma1*Sigma1));
260
261//    G4double FsymA1A2 = std::exp(-((As-(A1+A2))*(As-(A1+A2)))/(2.0*SigmaS*SigmaS));
262
263
264    G4double C2A = A2 + 3.72*Sigma2;
265    G4double C2S = As + 3.72*SigmaS;
266 
267    G4double C2 = 0.0;
268    if (w > 1000.0 ) C2 = C2S;
269    else if (w < 0.001) C2 = C2A;
270    else C2 =  std::max(C2A,C2S);
271
272    G4double C1 = A-C2;
273    if (C1 < 30.0) {
274        C2 = A-30.0;
275        C1 = 30.0;
276    }
277
278    G4double Am1 = (As + A1)/2.0;
279    G4double Am2 = (A1 + A2)/2.0;
280
281    // Get Mass distributions as sum of symmetric and asymmetric Gasussians
282    G4double Mass1 = MassDistribution(As,A,theParam); 
283    G4double Mass2 = MassDistribution(Am1,A,theParam); 
284    G4double Mass3 = MassDistribution(A1,A,theParam); 
285    G4double Mass4 = MassDistribution(Am2,A,theParam); 
286    G4double Mass5 = MassDistribution(A2,A,theParam); 
287    // get maximal value among Mass1,...,Mass5
288    G4double MassMax = Mass1;
289    if (Mass2 > MassMax) MassMax = Mass2;
290    if (Mass3 > MassMax) MassMax = Mass3;
291    if (Mass4 > MassMax) MassMax = Mass4;
292    if (Mass5 > MassMax) MassMax = Mass5;
293
294    // Sample a fragment mass number, which lies between C1 and C2
295    G4double m;
296    G4double Pm;
297    do {
298        m = C1+G4UniformRand()*(C2-C1);
299        Pm = MassDistribution(m,A,theParam); 
300    } while (G4UniformRand() > Pm/MassMax);
301
302    return static_cast<G4int>(m+0.5);
303}
304
305
306
307
308
309
310G4double G4CompetitiveFission::MassDistribution(const G4double x, const G4double A, 
311                                                const G4FissionParameters & theParam)
312    // This method gives mass distribution F(x) = F_{asym}(x)+w*F_{sym}(x)
313    // which consist of symmetric and asymmetric sum of gaussians components.
314{
315    G4double Xsym = std::exp(-0.5*(x-theParam.GetAs())*(x-theParam.GetAs())/
316                        (theParam.GetSigmaS()*theParam.GetSigmaS()));
317
318    G4double Xasym = std::exp(-0.5*(x-theParam.GetA2())*(x-theParam.GetA2())/
319                         (theParam.GetSigma2()*theParam.GetSigma2())) + 
320        std::exp(-0.5*(x-(A-theParam.GetA2()))*(x-(A-theParam.GetA2()))/
321            (theParam.GetSigma2()*theParam.GetSigma2())) +
322        0.5*std::exp(-0.5*(x-theParam.GetA1())*(x-theParam.GetA1())/
323                (theParam.GetSigma1()*theParam.GetSigma1())) +
324        0.5*std::exp(-0.5*(x-(A-theParam.GetA1()))*(x-(A-theParam.GetA1()))/
325                (theParam.GetSigma1()*theParam.GetSigma1()));
326
327    if (theParam.GetW() > 1000) return Xsym;
328    else if (theParam.GetW() < 0.001) return Xasym;
329    else return theParam.GetW()*Xsym+Xasym;
330}
331
332
333
334
335G4int G4CompetitiveFission::FissionCharge(const G4double A,
336                                          const G4double Z,
337                                          const G4double Af)
338    // Calculates the charge of a fission product for a given atomic number Af
339{
340    const G4double sigma = 0.6;
341    G4double DeltaZ = 0.0;
342    if (Af >= 134.0) DeltaZ = -0.45;                    //                      134 <= Af
343    else if (Af <= (A-134.0)) DeltaZ = 0.45;             // Af <= (A-134)
344    else DeltaZ = -0.45*(Af-(A/2.0))/(134.0-(A/2.0));   //       (A-134) < Af < 134
345
346    G4double Zmean = (Af/A)*Z + DeltaZ;
347 
348    G4double theZ;
349    do {
350        theZ = G4RandGauss::shoot(Zmean,sigma);
351    } while (theZ  < 1.0 || theZ > (Z-1.0) || theZ > Af);
352    //  return static_cast<G4int>(theZ+0.5);
353    return static_cast<G4int>(theZ+0.5);
354}
355
356
357
358
359G4double G4CompetitiveFission::FissionKineticEnergy(const G4double A, const G4double Z,
360                                                    const G4double Af1, const G4double /*Zf1*/,
361                                                    const G4double Af2, const G4double /*Zf2*/,
362                                                    const G4double /*U*/, const G4double Tmax,
363                                                    const G4FissionParameters & theParam)
364    // Gives the kinetic energy of fission products
365{
366    // Find maximal value of A for fragments
367    G4double AfMax = std::max(Af1,Af2);
368    if (AfMax < (A/2.0)) AfMax = A - AfMax;
369
370    // Weights for symmetric and asymmetric components
371    G4double Pas;
372    if (theParam.GetW() > 1000) Pas = 0.0;
373    else {
374        G4double P1 = 0.5*std::exp(-0.5*(AfMax-theParam.GetA1())*(AfMax-theParam.GetA1())/
375                              (theParam.GetSigma1()*theParam.GetSigma1()));
376
377        G4double P2 = std::exp(-0.5*(AfMax-theParam.GetA2())*(AfMax-theParam.GetA2())/
378                          (theParam.GetSigma2()*theParam.GetSigma2()));
379
380        Pas = P1+P2;
381    }
382
383
384    G4double Ps;
385    if (theParam.GetW() < 0.001) Ps = 0.0;
386    else 
387        Ps = theParam.GetW()*std::exp(-0.5*(AfMax-theParam.GetAs())*(AfMax-theParam.GetAs())/
388                                 (theParam.GetSigmaS()*theParam.GetSigmaS()));
389 
390
391
392    G4double Psy = Ps/(Pas+Ps);
393
394
395    // Fission fractions Xsy and Xas formed in symmetric and asymmetric modes
396    G4double PPas = theParam.GetSigma1() + 2.0 * theParam.GetSigma2();
397    G4double PPsy = theParam.GetW() * theParam.GetSigmaS();
398    G4double Xas = PPas / (PPas+PPsy);
399    G4double Xsy = PPsy / (PPas+PPsy);
400
401
402    // Average kinetic energy for symmetric and asymmetric components
403    G4double Eaverage = 0.1071*MeV*(Z*Z)/std::pow(A,1.0/3.0) + 22.2*MeV;
404
405
406    // Compute maximal average kinetic energy of fragments and Energy Dispersion (sqrt)
407    G4double TaverageAfMax;
408    G4double ESigma;
409    // Select randomly fission mode (symmetric or asymmetric)
410    if (G4UniformRand() > Psy) { // Asymmetric Mode
411        G4double A11 = theParam.GetA1()-0.7979*theParam.GetSigma1();
412        G4double A12 = theParam.GetA1()+0.7979*theParam.GetSigma1();
413        G4double A21 = theParam.GetA2()-0.7979*theParam.GetSigma2();
414        G4double A22 = theParam.GetA2()+0.7979*theParam.GetSigma2();
415        // scale factor
416        G4double ScaleFactor = 0.5*theParam.GetSigma1()*(AsymmetricRatio(A,A11)+AsymmetricRatio(A,A12))+
417            theParam.GetSigma2()*(AsymmetricRatio(A,A21)+AsymmetricRatio(A,A22));
418        // Compute average kinetic energy for fragment with AfMax
419        TaverageAfMax = (Eaverage + 12.5 * Xsy) * (PPas/ScaleFactor) * AsymmetricRatio(A,AfMax);
420        ESigma = 10.0*MeV; // MeV
421
422    } else { // Symmetric Mode
423        G4double As0 = theParam.GetAs() + 0.7979*theParam.GetSigmaS();
424        // scale factor
425        G4double ScaleFactor = theParam.GetW()*theParam.GetSigmaS()*SymmetricRatio(A,As0);
426        // Compute average kinetic energy for fragment with AfMax
427        TaverageAfMax = (Eaverage - 12.5*MeV*Xas) * (PPsy/ScaleFactor) * SymmetricRatio(A,AfMax);
428        ESigma = 8.0*MeV;
429    }
430
431
432    // Select randomly, in accordance with Gaussian distribution, fragment kinetic energy
433    G4double KineticEnergy;
434    G4int i = 0;
435    do {
436        KineticEnergy = G4RandGauss::shoot(TaverageAfMax,ESigma);
437        if (i++ > 100) return Eaverage;
438    } while (KineticEnergy < Eaverage-3.72*ESigma || 
439             KineticEnergy > Eaverage+3.72*ESigma ||
440             KineticEnergy > Tmax);
441
442    return KineticEnergy;
443}
444
445
446
447
448
449G4double G4CompetitiveFission::AsymmetricRatio(const G4double A,const G4double A11)
450{
451    const G4double B1 = 23.5;
452    const G4double A00 = 134.0;
453    return Ratio(A,A11,B1,A00);
454}
455
456G4double G4CompetitiveFission::SymmetricRatio(const G4double A,const G4double A11)
457{
458    const G4double B1 = 5.32;
459    const G4double A00 = A/2.0;
460    return Ratio(A,A11,B1,A00);
461}
462
463G4double G4CompetitiveFission::Ratio(const G4double A,const G4double A11,
464                                     const G4double B1,const G4double A00) 
465{
466    if (A == 0) throw G4HadronicException(__FILE__, __LINE__, "G4CompetitiveFission::Ratio: A == 0!");
467    if (A11 >= A/2.0 && A11 <= (A00+10.0)) return 1.0-B1*((A11-A00)/A)*((A11-A00)/A);
468    else return 1.0-B1*(10.0/A)*(10.0/A)-2.0*(10.0/A)*B1*((A11-A00-10.0)/A);
469}
470
471
472
473
474
475G4ThreeVector G4CompetitiveFission::IsotropicVector(const G4double Magnitude)
476    // Samples a isotropic random vectorwith a magnitud given by Magnitude.
477    // By default Magnitude = 1.0
478{
479    G4double CosTheta = 1.0 - 2.0*G4UniformRand();
480    G4double SinTheta = std::sqrt(1.0 - CosTheta*CosTheta);
481    G4double Phi = twopi*G4UniformRand();
482    G4ThreeVector Vector(Magnitude*std::cos(Phi)*SinTheta,
483                         Magnitude*std::sin(Phi)*SinTheta,
484                         Magnitude*CosTheta);
485    return Vector;
486}
487
488
489#ifdef debug
490void G4CompetitiveFission::CheckConservation(const G4Fragment & theInitialState,
491                                             G4FragmentVector * Result) const
492{
493    G4double ProductsEnergy =0;
494    G4ThreeVector ProductsMomentum;
495    G4int ProductsA = 0;
496    G4int ProductsZ = 0;
497    G4FragmentVector::iterator h;
498    for (h = Result->begin(); h != Result->end(); h++) {
499        G4LorentzVector tmp = (*h)->GetMomentum();
500        ProductsEnergy += tmp.e();
501        ProductsMomentum += tmp.vect();
502        ProductsA += static_cast<G4int>((*h)->GetA());
503        ProductsZ += static_cast<G4int>((*h)->GetZ());
504    }
505
506    if (ProductsA != theInitialState.GetA()) {
507        G4cout << "!!!!!!!!!! Baryonic Number Conservation Violation !!!!!!!!!!" << G4endl;
508        G4cout << "G4CompetitiveFission.cc: Barionic Number Conservation test for fission fragments" 
509               << G4endl; 
510        G4cout << "Initial A = " << theInitialState.GetA() 
511               << "   Fragments A = " << ProductsA << "   Diference --> " 
512               << theInitialState.GetA() - ProductsA << G4endl;
513    }
514    if (ProductsZ != theInitialState.GetZ()) {
515        G4cout << "!!!!!!!!!! Charge Conservation Violation !!!!!!!!!!" << G4endl;
516        G4cout << "G4CompetitiveFission.cc: Charge Conservation test for fission fragments" 
517               << G4endl; 
518        G4cout << "Initial Z = " << theInitialState.GetZ() 
519               << "   Fragments Z = " << ProductsZ << "   Diference --> " 
520               << theInitialState.GetZ() - ProductsZ << G4endl;
521    }
522    if (std::abs(ProductsEnergy-theInitialState.GetMomentum().e()) > 1.0*keV) {
523        G4cout << "!!!!!!!!!! Energy Conservation Violation !!!!!!!!!!" << G4endl;
524        G4cout << "G4CompetitiveFission.cc: Energy Conservation test for fission fragments" 
525               << G4endl; 
526        G4cout << "Initial E = " << theInitialState.GetMomentum().e()/MeV << " MeV"
527               << "   Fragments E = " << ProductsEnergy/MeV  << " MeV   Diference --> " 
528               << (theInitialState.GetMomentum().e() - ProductsEnergy)/MeV << " MeV" << G4endl;
529    } 
530    if (std::abs(ProductsMomentum.x()-theInitialState.GetMomentum().x()) > 1.0*keV || 
531        std::abs(ProductsMomentum.y()-theInitialState.GetMomentum().y()) > 1.0*keV ||
532        std::abs(ProductsMomentum.z()-theInitialState.GetMomentum().z()) > 1.0*keV) {
533        G4cout << "!!!!!!!!!! Momentum Conservation Violation !!!!!!!!!!" << G4endl;
534        G4cout << "G4CompetitiveFission.cc: Momentum Conservation test for fission fragments" 
535               << G4endl; 
536        G4cout << "Initial P = " << theInitialState.GetMomentum().vect() << " MeV"
537               << "   Fragments P = " << ProductsMomentum  << " MeV   Diference --> " 
538               << theInitialState.GetMomentum().vect() - ProductsMomentum << " MeV" << G4endl;
539    }
540    return;
541}
542#endif
543
544
545
546
Note: See TracBrowser for help on using the repository browser.