Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FissionProductYieldDist.cc
Go to the documentation of this file.
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 * File: G4FissionProductYieldDist.cc
28 * Author: B. Wendt ([email protected])
29 *
30 * Created on May 11, 2011, 12:04 PM
31 */
32
33/* * * * * * * * * * * * * * * * References * * * * * * * * * * * * * * * *
34 * *
35 * 1. "Systematics of fission fragment total kinetic energy release", *
36 * V. E. Viola, K. Kwiatkowski, and M. Walker, Physical Review C, 31.4, *
37 * April 1985 *
38 * 2. "Reactor Handbook", United States Atomic Energy Commission, *
39 * III.A:Physics, 1962 *
40 * 3. "Properties of the Alpha Particles Emitted in the Spontaneous Fission *
41 * of Cf252", Z. Fraenkel and S. G. Thompson, Physical Review Letters, *
42 * 13.14, October 1964 *
43 * *
44 * * * * * * * * * * * * * * * * References * * * * * * * * * * * * * * * */
45
46
47//#include <ios>
48//#include <iostream>
49
50#include "G4Alpha.hh"
51#include "G4Gamma.hh"
52#include "G4Ions.hh"
53#include "G4Neutron.hh"
54//#include "G4NeutronHPNames.hh"
55#include "G4NucleiProperties.hh"
56#include "G4ParticleTable.hh"
57#include "G4ThreeVector.hh"
58#include "Randomize.hh"
59#include "G4UImanager.hh"
60#include "globals.hh"
61#include "G4Exp.hh"
62#include "G4SystemOfUnits.hh"
63#include "G4HadFinalState.hh"
64#include "G4DynamicParticle.hh"
66#include "G4ReactionProduct.hh"
67
68#include "G4ArrayOps.hh"
69#include "G4ENDFTapeRead.hh"
72#include "G4FFGDefaultValues.hh"
73#include "G4FFGEnumerations.hh"
74#include "G4FPYNubarValues.hh"
75#include "G4FPYSamplingOps.hh"
78#include "G4Pow.hh"
79
80using CLHEP::pi;
81
82#ifdef G4MULTITHREADED
83#include "G4AutoLock.hh"
84G4Mutex G4FissionProductYieldDist::fissprodMutex = G4MUTEX_INITIALIZER;
85#endif
86
89 G4FFGEnumerations::MetaState WhichMetaState,
91 G4FFGEnumerations::YieldType WhichYieldType,
92 std::istringstream& dataStream )
93: Isotope_(WhichIsotope),
94 MetaState_(WhichMetaState),
95 Cause_(WhichCause),
96 YieldType_(WhichYieldType),
97 Verbosity_(G4FFGDefaultValues::Verbosity)
98{
100
101 try
102 {
103 // Initialize the class
104 Initialize(dataStream);
105 } catch (std::exception& e)
106 {
108 throw e;
109 }
110
112}
113
116 G4FFGEnumerations::MetaState WhichMetaState,
118 G4FFGEnumerations::YieldType WhichYieldType,
119 G4int Verbosity,
120 std::istringstream& dataStream)
121 : Isotope_(WhichIsotope), MetaState_(WhichMetaState), Cause_(WhichCause),
122 YieldType_(WhichYieldType), Verbosity_(Verbosity)
123{
125
126 try
127 {
128 // Initialize the class
129 Initialize(dataStream);
130 } catch (std::exception& e)
131 {
133 throw e;
134 }
135
137}
138
139
140void G4FissionProductYieldDist::Initialize(std::istringstream& dataStream)
141{
143
144 IncidentEnergy_ = 0.0;
147 SetNubar();
148
149 // Set miscellaneous variables
150 AlphaDefinition_ = reinterpret_cast<G4Ions*>(G4Alpha::Definition());
151 NeutronDefinition_ = reinterpret_cast<G4Ions*>(G4Neutron::Definition());
154
155 // Construct G4NeutronHPNames: provides access to the element names
157 // Get the pointer to G4ParticleTable: stores all G4Ions
159 // Construct the pointer to the random engine
160 // TODO Make G4FPSamplingOps a singleton so that only one instance is used across all classes
162
163 try
164 {
165 // Read in and sort the probability data
166 ENDFData_ = new G4ENDFTapeRead(dataStream,
168 Cause_,
169 Verbosity_);
170// ENDFData_ = new G4ENDFTapeRead(MakeDirectoryName(),
171// MakeFileName(Isotope_, MetaState_),
172// YieldType_,
173// Cause_);
179 MakeTrees();
181 } catch (std::exception& e)
182 {
183 delete ElementNames_;
184 delete RandomEngine_;
185
187 throw e;
188 }
189
191}
192
194{
196
197#ifdef G4MULTITHREADED
198 G4AutoLock lk(&G4FissionProductYieldDist::fissprodMutex);
199#endif
200
201 // Check to see if the user has set the alpha production to a somewhat
202 // reasonable level
204
205 // Generate the new G4DynamicParticle pointers to identify key locations in
206 // the G4DynamicParticle chain that will be passed to the G4FissionEvent
207 G4ReactionProduct* FirstDaughter = NULL;
208 G4ReactionProduct* SecondDaughter = NULL;
209 std::vector< G4ReactionProduct* >* Alphas = new std::vector< G4ReactionProduct* >;
210 std::vector< G4ReactionProduct* >* Neutrons = new std::vector< G4ReactionProduct* >;
211 std::vector< G4ReactionProduct* >* Gammas = new std::vector< G4ReactionProduct* >;
212
213 // Generate all the nucleonic fission products
214 // How many nucleons do we have to work with?
215 //TK modified 131108
216 const G4int ParentA = (Isotope_/10) % 1000;
217 const G4int ParentZ = ((Isotope_/10) - ParentA) / 1000;
218 RemainingA_ = ParentA;
219 RemainingZ_ = ParentZ;
220
221 // Don't forget the extra nucleons depending on the fission cause
222 switch(Cause_)
223 {
225 ++RemainingA_;
226 break;
227
229 ++RemainingZ_;
230 break;
231
234 default:
235 // Nothing to do here
236 break;
237 }
238
239 // Ternary fission can be set by the user. Thus, it is necessary to
240 // sample the alpha particle first and the first daughter product
241 // second. See the discussion in
242 // G4FissionProductYieldDist::G4GetFissionProduct() for more information
243 // as to why the fission events are sampled this way.
244 GenerateAlphas(Alphas);
245
246 // Generate the first daughter product
247 FirstDaughter = new G4ReactionProduct(GetFissionProduct());
248 RemainingA_ -= FirstDaughter->GetDefinition()->GetAtomicMass();
249 RemainingZ_ -= FirstDaughter->GetDefinition()->GetAtomicNumber();
253
254 G4cout << " -- First daughter product sampled" << G4endl;
256 G4cout << " Name: " << FirstDaughter->GetDefinition()->GetParticleName() << G4endl;
258 G4cout << " Z: " << FirstDaughter->GetDefinition()->GetAtomicNumber() << G4endl;
260 G4cout << " A: " << FirstDaughter->GetDefinition()->GetAtomicMass() << G4endl;
262 G4cout << " Meta State: " << (FirstDaughter->GetDefinition()->GetPDGEncoding() % 10) << G4endl;
263 }
264
265 GenerateNeutrons(Neutrons);
266
267 // Now that all the nucleonic particles have been generated, we can
268 // calculate the composition of the second daughter product.
269 G4int NewIsotope = RemainingZ_ * 1000 + RemainingA_;
274
275 G4cout << " -- Second daughter product sampled" << G4endl;
277 G4cout << " Name: " << SecondDaughter->GetDefinition()->GetParticleName() << G4endl;
279 G4cout << " Z: " << SecondDaughter->GetDefinition()->GetAtomicNumber() << G4endl;
281 G4cout << " A: " << SecondDaughter->GetDefinition()->GetAtomicMass() << G4endl;
283 G4cout << " Meta State: " << (SecondDaughter->GetDefinition()->GetPDGEncoding() % 10) << G4endl;
284 }
285
286 // Calculate how much kinetic energy will be available
287 // 195 to 205 MeV are available in a fission reaction, but about 20 MeV
288 // are from delayed sources. We are concerned only with prompt sources,
289 // so sample a Gaussian distribution about 20 MeV and subtract the
290 // result from the total available energy. Also, the energy of fission
291 // neutrinos is neglected. Fission neutrinos would add ~11 MeV
292 // additional energy to the fission. (Ref 2)
293 // Finally, add in the kinetic energy contribution of the fission
294 // inducing particle, if any.
295 const G4double TotalKE = RandomEngine_->G4SampleUniform(195.0, 205.0) * MeV
296 - RandomEngine_->G4SampleGaussian(20.0, 3.0) * MeV
298 RemainingEnergy_ = TotalKE;
299
300 // Calculate the energies of the alpha particles and neutrons
301 // Once again, since the alpha particles are user defined, we must
302 // sample their kinetic energy first. SampleAlphaEnergies() returns the
303 // amount of energy consumed by the alpha particles, so remove the total
304 // energy alloted to the alpha particles from the available energy
305 SampleAlphaEnergies(Alphas);
306
307 // Second, the neutrons are sampled from the Watt fission spectrum.
308 SampleNeutronEnergies(Neutrons);
309
310 // Calculate the total energy available to the daughter products
311 // A Gaussian distribution about the average calculated energy with
312 // a standard deviation of 1.5 MeV (Ref. 2) is used. Since the energy
313 // distribution is dependant on the alpha particle generation and the
314 // Watt fission sampling for neutrons, we only have the left-over energy
315 // to work with for the fission daughter products.
316 G4double FragmentsKE;
317 G4int icounter=0;
318 G4int icounter_max=1024;
319 do
320 {
321 icounter++;
322 if ( icounter > icounter_max ) {
323 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
324 break;
325 }
326 FragmentsKE = RandomEngine_->G4SampleGaussian(RemainingEnergy_, 1.5 *MeV);
327 } while(FragmentsKE > RemainingEnergy_); // Loop checking, 11.05.2015, T. Koi
328
329 // Make sure that we don't produce any sub-gamma photons
330 if((RemainingEnergy_ - FragmentsKE) / (100 * keV) < 1.0)
331 {
332 FragmentsKE = RemainingEnergy_;
333 }
334
335 // This energy has now been allotted to the fission fragments.
336 // Subtract FragmentsKE from the total available fission energy.
337 RemainingEnergy_ -= FragmentsKE;
338
339 // Sample the energies of the gamma rays
340 // Assign the remainder, if any, of the energy to the gamma rays
341 SampleGammaEnergies(Gammas);
342
343 // Prepare to balance the momenta of the system
344 // We will need these for sampling the angles and balancing the momenta
345 // of all the fission products
346 G4double NumeratorSqrt, NumeratorOther, Denominator;
347 G4double Theta, Phi, Magnitude;
348 G4ThreeVector Direction;
349 G4ParticleMomentum ResultantVector(0, 0, 0);
350
351 if(Alphas->size() > 0)
352 {
353 // Sample the angles of the alpha particles and neutrons, then calculate
354 // the total moment contribution to the system
355 // The average angle of the alpha particles with respect to the
356 // light fragment is dependent on the ratio of the kinetic energies.
357 // This equation was determined by performing a fit on data from
358 // Ref. 3 using the website:
359 // http://soft.arquimedex.com/linear_regression.php
360 //
361 // RatioOfKE Angle (rad) Angle (degrees)
362 // 1.05 1.257 72
363 // 1.155 1.361 78
364 // 1.28 1.414 81
365 // 1.5 1.518 87
366 // 1.75 1.606 92
367 // 1.9 1.623 93
368 // 2.2 1.728 99
369 // This equation generates the angle in radians. If the RatioOfKE is
370 // greater than 2.25 the angle defaults to 1.3963 rad (100 degrees)
371 G4double MassRatio = FirstDaughter->GetDefinition()->GetPDGMass() / SecondDaughter->GetDefinition()->GetPDGMass();
372
373 // Invert the mass ratio if the first daughter product is the lighter fragment
374 if(MassRatio < 1)
375 {
376 MassRatio = 1 / MassRatio;
377 }
378
379 // The empirical equation is valid for mass ratios up to 2.75
380 if(MassRatio > 2.75)
381 {
382 MassRatio = 2.75;
383 }
384 const G4double MeanAlphaAngle = 0.3644 * MassRatio * MassRatio * MassRatio
385 - 1.9766 * MassRatio * MassRatio
386 + 3.8207 * MassRatio
387 - 0.9917;
388
389 // Sample the directions of the alpha particles with respect to the
390 // light fragment. For the moment we will assume that the light
391 // fragment is traveling along the z-axis in the positive direction.
392 const G4double MeanAlphaAngleStdDev = 0.0523598776;
393 G4double PlusMinus;
394
395 for(unsigned int i = 0; i < Alphas->size(); ++i)
396 {
397 PlusMinus = std::acos(RandomEngine_->G4SampleGaussian(0, MeanAlphaAngleStdDev)) - (pi / 2);
398 Theta = MeanAlphaAngle + PlusMinus;
399 if(Theta < 0)
400 {
401 Theta = 0.0 - Theta;
402 } else if(Theta > pi)
403 {
404 Theta = (2 * pi - Theta);
405 }
406 Phi = RandomEngine_->G4SampleUniform(-pi, pi);
407
408 Direction.setRThetaPhi(1.0, Theta, Phi);
409 Magnitude = std::sqrt(2 * Alphas->at(i)->GetKineticEnergy() * Alphas->at(i)->GetDefinition()->GetPDGMass());
410 Alphas->at(i)->SetMomentum(Direction * Magnitude);
411 ResultantVector += Alphas->at(i)->GetMomentum();
412 }
413 }
414
415 // Sample the directions of the neutrons.
416 if(Neutrons->size() != 0)
417 {
418 for(unsigned int i = 0; i < Neutrons->size(); ++i)
419 {
420 Theta = std::acos(RandomEngine_->G4SampleUniform(-1, 1));
421 Phi = RandomEngine_->G4SampleUniform(-pi, pi);
422
423 Direction.setRThetaPhi(1.0, Theta, Phi);
424 Magnitude = std::sqrt(2 * Neutrons->at(i)->GetKineticEnergy() * Neutrons->at(i)->GetDefinition()->GetPDGMass());
425 Neutrons->at(i)->SetMomentum(Direction * Magnitude);
426 ResultantVector += Neutrons->at(i)->GetMomentum();
427 }
428 }
429
430 // Sample the directions of the gamma rays
431 if(Gammas->size() != 0)
432 {
433 for(unsigned int i = 0; i < Gammas->size(); ++i)
434 {
435 Theta = std::acos(RandomEngine_->G4SampleUniform(-1, 1));
436 Phi = RandomEngine_->G4SampleUniform(-pi, pi);
437
438 Direction.setRThetaPhi(1.0, Theta, Phi);
439 Magnitude = Gammas->at(i)->GetKineticEnergy() / CLHEP::c_light;
440 Gammas->at(i)->SetMomentum(Direction * Magnitude);
441 ResultantVector += Gammas->at(i)->GetMomentum();
442 }
443 }
444
445 // Calculate the momenta of the two daughter products
446 G4ReactionProduct* LightFragment;
447 G4ReactionProduct* HeavyFragment;
448 G4ThreeVector LightFragmentDirection;
449 G4ThreeVector HeavyFragmentDirection;
450 G4double ResultantX, ResultantY, ResultantZ;
451 ResultantX = ResultantVector.getX();
452 ResultantY = ResultantVector.getY();
453 ResultantZ = ResultantVector.getZ();
454
455 if(FirstDaughter->GetDefinition()->GetPDGMass() < SecondDaughter->GetDefinition()->GetPDGMass())
456 {
457 LightFragment = FirstDaughter;
458 HeavyFragment = SecondDaughter;
459 } else
460 {
461 LightFragment = SecondDaughter;
462 HeavyFragment = FirstDaughter;
463 }
464 const G4double LightFragmentMass = LightFragment->GetDefinition()->GetPDGMass();
465 const G4double HeavyFragmentMass = HeavyFragment->GetDefinition()->GetPDGMass();
466
467 LightFragmentDirection.setRThetaPhi(1.0, 0, 0);
468
469 // Fit the momenta of the daughter products to the resultant vector of
470 // the remaining fission products. This will be done in the Cartesian
471 // coordinate system, not spherical. This is done using the following
472 // table listing the system momenta and the corresponding equations:
473 // X Y Z
474 //
475 // A 0 0 P
476 //
477 // B -R_x -R_y -P - R_z
478 //
479 // R R_x R_y R_z
480 //
481 // v = sqrt(2*m*k) -> k = v^2/(2*m)
482 // tk = k_A + k_B
483 // k_L = P^2/(2*m_L)
484 // k_H = ((-R_x)^2 + (-R_y)^2 + (-P - R_z)^2)/(2*m_H)
485 // where:
486 // P: momentum of the light daughter product
487 // R: the remaining fission products' resultant vector
488 // v: momentum
489 // m: mass
490 // k: kinetic energy
491 // tk: total kinetic energy available to the daughter products
492 //
493 // Below is the solved form for P, with the solution generated using
494 // the WolframAlpha website:
495 // http://www.wolframalpha.com/input/?i=
496 // solve+((-x)^2+%2B+(-y)^2+%2B+(-P-z)^2)%2F(2*B)+%2B+L^2%2F(2*A)+%3D+k+
497 // for+P
498 //
499 //
500 // nsqrt = sqrt(m_L*(m_L*(2*m_H*tk - R_x^2 - R_y^2) +
501 // m_H*(2*m_H*tk - R_x^2 - R_y^2 - R_z^2))
502 NumeratorSqrt = std::sqrt(LightFragmentMass *
503 (LightFragmentMass * (2 * HeavyFragmentMass * FragmentsKE
504 - ResultantX * ResultantX
505 - ResultantY * ResultantY)
506 + HeavyFragmentMass * (2 * HeavyFragmentMass * FragmentsKE
507 - ResultantX * ResultantX
508 - ResultantY * ResultantY
509 - ResultantZ - ResultantZ)));
510
511 // nother = m_L*R_z
512 NumeratorOther = LightFragmentMass * ResultantZ;
513
514 // denom = m_L + m_H
515 Denominator = LightFragmentMass + HeavyFragmentMass;
516
517 // P = (nsqrt + nother) / denom
518 const G4double LightFragmentMomentum = (NumeratorSqrt + NumeratorOther) / Denominator;
519 const G4double HeavyFragmentMomentum = std::sqrt(ResultantX * ResultantX
520 + ResultantY * ResultantY
521 + G4Pow::GetInstance()->powN(LightFragmentMomentum + ResultantZ, 2));
522
523 // Finally! We now have everything we need for the daughter products
524 LightFragment->SetMomentum(LightFragmentDirection * LightFragmentMomentum);
525 HeavyFragmentDirection.setX(-ResultantX);
526 HeavyFragmentDirection.setY(-ResultantY);
527 HeavyFragmentDirection.setZ((-LightFragmentMomentum - ResultantZ));
528 // Don't forget to normalize the vector to the unit sphere
529 HeavyFragmentDirection.setR(1.0);
530 HeavyFragment->SetMomentum(HeavyFragmentDirection * HeavyFragmentMomentum);
531
533 {
536
537 G4cout << " -- Daugher product momenta finalized" << G4endl;
539 }
540
541 // Load all the particles into a contiguous set
542 //TK modifed 131108
543 //G4DynamicParticleVector* FissionProducts = new G4DynamicParticleVector(2 + Alphas->size() + Neutrons->size() + Gammas->size());
544 G4DynamicParticleVector* FissionProducts = new G4DynamicParticleVector();
545 // Load the fission fragments
546 FissionProducts->push_back(MakeG4DynamicParticle(LightFragment));
547 FissionProducts->push_back(MakeG4DynamicParticle(HeavyFragment));
548 // Load the neutrons
549 for(unsigned int i = 0; i < Neutrons->size(); i++)
550 {
551 FissionProducts->push_back(MakeG4DynamicParticle(Neutrons->at(i)));
552 }
553 // Load the gammas
554 for(unsigned int i = 0; i < Gammas->size(); i++)
555 {
556 FissionProducts->push_back(MakeG4DynamicParticle(Gammas->at(i)));
557 }
558 // Load the alphas
559 for(unsigned int i = 0; i < Alphas->size(); i++)
560 {
561 FissionProducts->push_back(MakeG4DynamicParticle(Alphas->at(i)));
562 }
563
564 // Rotate the system to a random location so that the light fission fragment
565 // is not always traveling along the positive z-axis
566 // Sample Theta and Phi.
567 G4ThreeVector RotationAxis;
568
569 Theta = std::acos(RandomEngine_->G4SampleUniform(-1, 1));
570 Phi = RandomEngine_->G4SampleUniform(-pi, pi);
571 RotationAxis.setRThetaPhi(1.0, Theta, Phi);
572
573 // We will also check the net momenta
574 ResultantVector.set(0.0, 0.0, 0.0);
575 for(unsigned int i = 0; i < FissionProducts->size(); i++)
576 {
577 Direction = FissionProducts->at(i)->GetMomentumDirection();
578 Direction.rotateUz(RotationAxis);
579 FissionProducts->at(i)->SetMomentumDirection(Direction);
580 ResultantVector += FissionProducts->at(i)->GetMomentum();
581 }
582
583 // Warn if the sum momenta of the system is not within a reasonable
584 // tolerance
585 G4double PossibleImbalance = ResultantVector.mag();
586 if(PossibleImbalance > 0.01)
587 {
588 std::ostringstream Temp;
589 Temp << "Momenta imbalance of ";
590 Temp << PossibleImbalance / (MeV / CLHEP::c_light);
591 Temp << " MeV/c in the system";
592 G4Exception("G4FissionProductYieldDist::G4GetFission()",
593 Temp.str().c_str(),
595 "Results may not be valid");
596 }
597
598 // Clean up
599 delete FirstDaughter;
600 delete SecondDaughter;
604
606 return FissionProducts;
607}
608
611{
613
615
617 return Product;
618}
619
621G4SetAlphaProduction(G4double WhatAlphaProduction)
622{
624
625 AlphaProduction_ = WhatAlphaProduction;
626
628}
629
631G4SetEnergy( G4double WhatIncidentEnergy )
632{
634
636 {
637 IncidentEnergy_ = WhatIncidentEnergy;
638 } else
639 {
640 IncidentEnergy_ = 0 * GeV;
641 }
642
644}
645
647G4SetTernaryProbability( G4double WhatTernaryProbability )
648{
650
651 TernaryProbability_ = WhatTernaryProbability;
652
654}
655
657G4SetVerbosity(G4int Verbosity)
658{
660
661 Verbosity_ = Verbosity;
662
665
667}
668
670CheckAlphaSanity( void )
671{
673
674 // This provides comfortable breathing room at 16 MeV per alpha
675 if(AlphaProduction_ > 10)
676 {
677 AlphaProduction_ = 10;
678 } else if(AlphaProduction_ < -7)
679 {
680 AlphaProduction_ = -7;
681 }
682
684}
685
687FindParticle( G4double RandomParticle )
688{
690
691 // Determine which energy group is currently in use
692 G4bool isExact = false;
693 G4bool lowerExists = false;
694 G4bool higherExists = false;
695 G4int energyGroup;
696 for(energyGroup = 0; energyGroup < YieldEnergyGroups_; energyGroup++)
697 {
698 if(IncidentEnergy_ == YieldEnergies_[energyGroup])
699 {
700 isExact = true;
701 break;
702 }
703
704 if(energyGroup == 0 && IncidentEnergy_ < YieldEnergies_[energyGroup])
705 {
706 // Break if the energy is less than the lowest energy
707 higherExists = true;
708 break;
709 } else if(energyGroup == YieldEnergyGroups_ - 1)
710 {
711 // The energy is greater than any values in the yield data.
712 lowerExists = true;
713 break;
714 } else
715 {
716 // Break if the energy is less than the lowest energy
717 if(IncidentEnergy_ > YieldEnergies_[energyGroup])
718 {
719 energyGroup--;
720 lowerExists = true;
721 higherExists = true;
722 break;
723 }
724 }
725 }
726
727 // Determine which particle it is
728 G4Ions* FoundParticle = NULL;
729 if(isExact || YieldEnergyGroups_ == 1)
730 {
731 // Determine which tree contains the random value
732 G4int tree;
733 for(tree = 0; tree < TreeCount_; tree++)
734 {
735 // Break if a tree is identified as containing the random particle
736 if(RandomParticle <= Trees_[tree].ProbabilityRangeEnd[energyGroup])
737 {
738 break;
739 }
740 }
741 ProbabilityBranch* Branch = Trees_[tree].Trunk;
742
743 // Iteratively traverse the tree until the particle addressed by the random
744 // variable is found
745 G4bool RangeIsSmaller;
746 G4bool RangeIsGreater;
747 while((RangeIsSmaller = (RandomParticle < Branch->ProbabilityRangeBottom[energyGroup]))
748 || (RangeIsGreater = (RandomParticle > Branch->ProbabilityRangeTop[energyGroup])))
749 // Loop checking, 11.05.2015, T. Koi
750 {
751 if(RangeIsSmaller)
752 {
753 Branch = Branch->Left;
754 } else {
755 Branch = Branch->Right;
756 }
757 }
758
759 FoundParticle = Branch->Particle;
760 } else if(lowerExists && higherExists)
761 {
762 // We need to do some interpolation
763 FoundParticle = FindParticleInterpolation(RandomParticle, energyGroup);
764 } else
765 {
766 // We need to do some extrapolation
767 FoundParticle = FindParticleExtrapolation(RandomParticle, lowerExists);
768 }
769
770 // Return the particle
772 return FoundParticle;
773}
774
776FindParticleExtrapolation( G4double RandomParticle,
777 G4bool LowerEnergyGroupExists )
778{
780
781 G4Ions* FoundParticle = NULL;
782 G4int NearestEnergy;
783 G4int NextNearestEnergy;
784
785 // Check to see if we are extrapolating above or below the data set
786 if(LowerEnergyGroupExists == true)
787 {
788 NearestEnergy = YieldEnergyGroups_ - 1;
789 NextNearestEnergy = NearestEnergy - 1;
790 } else
791 {
792 NearestEnergy = 0;
793 NextNearestEnergy = 1;
794 }
795
796 for(G4int Tree = 0; Tree < TreeCount_ && FoundParticle == NULL; Tree++)
797 {
798 FoundParticle = FindParticleBranchSearch(Trees_[Tree].Trunk,
799 RandomParticle,
800 NearestEnergy,
801 NextNearestEnergy);
802 }
803
805 return FoundParticle;
806}
807
809FindParticleInterpolation( G4double RandomParticle,
810 G4int LowerEnergyGroup )
811{
813
814 G4Ions* FoundParticle = NULL;
815 G4int HigherEnergyGroup = LowerEnergyGroup + 1;
816
817 for(G4int Tree = 0; Tree < TreeCount_ && FoundParticle == NULL; Tree++)
818 {
819 FoundParticle = FindParticleBranchSearch(Trees_[Tree].Trunk,
820 RandomParticle,
821 LowerEnergyGroup,
822 HigherEnergyGroup);
823 }
824
826 return FoundParticle;
827}
828
831 G4double RandomParticle,
832 G4int EnergyGroup1,
833 G4int EnergyGroup2 )
834{
836
837 G4Ions* Particle;
838
839 // Verify that the branch exists
840 if(Branch == NULL)
841 {
842 Particle = NULL;
843 } else if(EnergyGroup1 >= Branch->IncidentEnergiesCount
844 || EnergyGroup2 >= Branch->IncidentEnergiesCount
845 || EnergyGroup1 == EnergyGroup2
846 || Branch->IncidentEnergies[EnergyGroup1] == Branch->IncidentEnergies[EnergyGroup2])
847 {
848 // Set NULL if any invalid conditions exist
849 Particle = NULL;
850 } else
851 {
852 // Everything check out - proceed
853 G4Ions* FoundParticle = NULL;
854 G4double Intercept;
855 G4double Slope;
856 G4double RangeAtIncidentEnergy;
857 G4double Denominator = Branch->IncidentEnergies[EnergyGroup1] - Branch->IncidentEnergies[EnergyGroup2];
858
859 // Calculate the lower probability bounds
860 Slope = (Branch->ProbabilityRangeBottom[EnergyGroup1] - Branch->ProbabilityRangeBottom[EnergyGroup2]) / Denominator;
861 Intercept = Branch->ProbabilityRangeBottom[EnergyGroup1] - Slope * Branch->IncidentEnergies[EnergyGroup1];
862 RangeAtIncidentEnergy = Slope * IncidentEnergy_ + Intercept;
863
864 // Go right if the particle is below the probability bounds
865 if(RandomParticle < RangeAtIncidentEnergy)
866 {
867 FoundParticle = FindParticleBranchSearch(Branch->Left,
868 RandomParticle,
869 EnergyGroup1,
870 EnergyGroup2);
871 } else
872 {
873 // Calculate the upper probability bounds
874 Slope = (Branch->ProbabilityRangeTop[EnergyGroup1] - Branch->ProbabilityRangeTop[EnergyGroup2]) / Denominator;
875 Intercept = Branch->ProbabilityRangeTop[EnergyGroup1] - Slope * Branch->IncidentEnergies[EnergyGroup1];
876 RangeAtIncidentEnergy = Slope * IncidentEnergy_ + Intercept;
877
878 // Go left if the particle is above the probability bounds
879 if(RandomParticle > RangeAtIncidentEnergy)
880 {
881 FoundParticle = FindParticleBranchSearch(Branch->Right,
882 RandomParticle,
883 EnergyGroup1,
884 EnergyGroup2);
885 } else
886 {
887 // If the particle is bounded then we found it!
888 FoundParticle = Branch->Particle;
889 }
890 }
891
892 Particle = FoundParticle;
893 }
894
896 return Particle;
897}
898
900GenerateAlphas( std::vector< G4ReactionProduct* >* Alphas )
901{
903
904 // Throw the dice to determine if ternary fission occurs
906 if(MakeAlphas)
907 {
908 G4int NumberOfAlphasToProduce;
909
910 // Determine how many alpha particles to produce for the ternary fission
911 if(AlphaProduction_ < 0)
912 {
913 NumberOfAlphasToProduce = RandomEngine_->G4SampleIntegerGaussian(AlphaProduction_ * -1,
914 1,
916 } else
917 {
918 NumberOfAlphasToProduce = (G4int)AlphaProduction_;
919 }
920
921 //TK modifed 131108
922 //Alphas->resize(NumberOfAlphasToProduce);
923 for(int i = 0; i < NumberOfAlphasToProduce; i++)
924 {
925 // Set the G4Ions as an alpha particle
926 Alphas->push_back(new G4ReactionProduct(AlphaDefinition_));
927
928 // Remove 4 nucleons (2 protons and 2 neutrons) for each alpha added
929 RemainingZ_ -= 2;
930 RemainingA_ -= 4;
931 }
932 }
933
935}
936
938GenerateNeutrons( std::vector< G4ReactionProduct* >* Neutrons )
939{
941
942 G4int NeutronProduction;
944
945 //TK modifed 131108
946 //Neutrons->resize(NeutronProduction);
947 for(int i = 0; i < NeutronProduction; i++)
948 {
949 // Define the fragment as a neutron
950 Neutrons->push_back(new G4ReactionProduct(NeutronDefinition_));
951
952 // Remove 1 nucleon for each neutron added
953 RemainingA_--;
954 }
955
957}
958
961 //TK modified 131108
962 //G4FFGEnumerations::MetaState MetaState )
963 G4FFGEnumerations::MetaState /*MetaState*/ )
964{
966
967 G4Ions* Temp;
968
969 // Break Product down into its A and Z components
970 G4int A = Product % 1000; // Extract A
971 G4int Z = (Product - A) / 1000; // Extract Z
972
973 // Check to see if the particle is registered using the PDG code
974 // TODO Add metastable state when supported by G4IonTable::GetIon()
975 Temp = reinterpret_cast<G4Ions*>(IonTable_->GetIon(Z, A));
976
977 // Removed in favor of the G4IonTable::GetIon() method
978// // Register the particle if it does not exist
979// if(Temp == NULL)
980// {
981// // Define the particle properties
982// G4String Name = MakeIsotopeName(Product, MetaState);
983// // Calculate the rest mass using a function already in Geant4
984// G4double Mass = G4NucleiProperties::
985// GetNuclearMass((double)A, (double)Z );
986// G4double Charge = Z*eplus;
987// G4int BaryonNum = A;
988// G4bool Stable = TRUE;
989//
990// // I am unsure about the following properties:
991// // 2*Spin, Parity, C-conjugation, 2*Isospin, 2*Isospin3, G-parity.
992// // Perhaps is would be a good idea to have a physicist familiar with
993// // Geant4 nomenclature to review and correct these properties.
994// Temp = new G4Ions (
995// // Name Mass Width Charge
996// Name, Mass, 0.0, Charge,
997//
998// // 2*Spin Parity C-conjugation 2*Isospin
999// 0, 1, 0, 0,
1000//
1001// // 2*Isospin3 G-parity Type Lepton number
1002// 0, 0, "nucleus", 0,
1003//
1004// // Baryon number PDG encoding Stable Lifetime
1005// BaryonNum, PDGCode, Stable, -1,
1006//
1007// // Decay table Shortlived SubType Anti_encoding
1008// NULL, FALSE, "generic", 0,
1009//
1010// // Excitation
1011// 0.0);
1012// Temp->SetPDGMagneticMoment(0.0);
1013//
1014// // Declare that there is no anti-particle
1015// Temp->SetAntiPDGEncoding(0);
1016//
1017// // Define the processes to use in transporting the particles
1018// std::ostringstream osAdd;
1019// osAdd << "/run/particle/addProcManager " << Name;
1020// G4String cmdAdd = osAdd.str();
1021//
1022// // set /control/verbose 0
1023// G4int tempVerboseLevel = G4UImanager::GetUIpointer()->GetVerboseLevel();
1024// G4UImanager::GetUIpointer()->SetVerboseLevel(0);
1025//
1026// // issue /run/particle/addProcManage
1027// G4UImanager::GetUIpointer()->ApplyCommand(cmdAdd);
1028//
1029// // retrieve /control/verbose
1030// G4UImanager::GetUIpointer()->SetVerboseLevel(tempVerboseLevel);
1031// }
1032
1034 return Temp;
1035}
1036
1038MakeDirectoryName( void )
1039{
1041
1042 // Generate the file location starting in the Geant4 data directory
1043 std::ostringstream DirectoryName;
1044 DirectoryName << std::getenv("G4NEUTRONHPDATA") << G4FFGDefaultValues::ENDFFissionDataLocation;
1045
1046 // Return the directory structure
1048 return DirectoryName.str();
1049}
1050
1052MakeFileName( G4int Isotope,
1054{
1056
1057
1058 // Create the unique identifying name for the particle
1059 std::ostringstream FileName;
1060
1061 // Determine if a leading 0 is needed (ZZZAAA or 0ZZAAA)
1062 if(Isotope < 100000)
1063 {
1064 FileName << "0";
1065 }
1066
1067 // Add the name of the element and the extension
1068 FileName << MakeIsotopeName(Isotope, MetaState) << ".fpy";
1069
1071 return FileName.str();
1072}
1073
1075MakeG4DynamicParticle( G4ReactionProduct* ReactionProduct )
1076{
1078
1079 G4DynamicParticle* DynamicParticle = new G4DynamicParticle(ReactionProduct->GetDefinition(), ReactionProduct->GetMomentum());
1080
1082 return DynamicParticle;
1083}
1084
1086MakeIsotopeName( G4int Isotope,
1088{
1090
1091 // Break Product down into its A and Z components
1092 G4int A = Isotope % 1000;
1093 G4int Z = (Isotope - A) / 1000;
1094
1095 // Create the unique identifying name for the particle
1096 std::ostringstream IsotopeName;
1097
1098 IsotopeName << Z << "_" << A;
1099
1100 // If it is metastable then append "m" to the name
1101 if(MetaState != G4FFGEnumerations::GROUND_STATE)
1102 {
1103 IsotopeName << "m";
1104
1105 // If it is a second isomeric state then append "2" to the name
1106 if(MetaState == G4FFGEnumerations::META_2)
1107 {
1108 IsotopeName << "2";
1109 }
1110 }
1111 // Add the name of the element and the extension
1112 IsotopeName << "_" << ElementNames_->theString[Z - 1];
1113
1115 return IsotopeName.str();
1116}
1117
1119MakeTrees( void )
1120{
1122
1123 // Allocate the space
1124 // We will make each tree a binary search
1125 // The maximum number of iterations required to find a single fission product
1126 // based on it's probability is defined by the following:
1127 // x = number of fission products
1128 // Trees = T(x) = ceil( ln(x) )
1129 // Rows/Tree = R(x) = ceil(( sqrt( (8 * x / T(x)) + 1) - 1) / 2)
1130 // Maximum = M(x) = T(x) + R(x)
1131 // Results: x => M(x)
1132 // 10 5
1133 // 100 10
1134 // 1000 25
1135 // 10000 54
1136 // 100000 140
1139
1140 // Initialize the range of each node
1141 for(G4int i = 0; i < TreeCount_; i++)
1142 {
1144 Trees_[i].Trunk = NULL;
1145 Trees_[i].BranchCount = 0;
1146 Trees_[i].IsEnd = FALSE;
1147 }
1148 // Mark the last tree as the ending tree
1149 Trees_[TreeCount_ - 1].IsEnd = TRUE;
1150
1152}
1153
1155ReadProbabilities( void )
1156{
1158
1160 BranchCount_ = 0;
1162
1163 // Loop through all the products
1164 for(G4int i = 0; i < ProductCount; i++)
1165 {
1166 // Acquire the data and sort it
1168 }
1169
1170 // Generate the true normalization factor, since round-off errors may result
1171 // in non-singular normalization of the data files. Also, reset DataTotal_
1172 // since it is used by Renormalize() to set the probability segments.
1175
1176 // Go through all the trees one at a time
1177 for(G4int i = 0; i < TreeCount_; i++)
1178 {
1179 Renormalize(Trees_[i].Trunk);
1180 // Set the max range of the tree to DataTotal
1181 G4ArrayOps::Copy(YieldEnergyGroups_, Trees_[i].ProbabilityRangeEnd, DataTotal_);
1182 }
1183
1185}
1186
1189{
1191
1192 // Check to see if Branch exists. Branch will be a null pointer if it
1193 // doesn't exist
1194 if(Branch != NULL)
1195 {
1196 // Call the lower branch to set the probability segment first, since it
1197 // supposed to have a lower probability segment that this node
1198 Renormalize(Branch->Left);
1199
1200 // Set this node as the next sequential probability segment
1205
1206 // Now call the upper branch to set those probability segments
1207 Renormalize(Branch->Right);
1208 }
1209
1211}
1212
1214SampleAlphaEnergies( std::vector< G4ReactionProduct* >* Alphas )
1215{
1217
1218 // The condition of sampling more energy from the fission products than is
1219 // alloted is statistically unfavorable, but it could still happen. The
1220 // do-while loop prevents such an occurrence from happening
1221 G4double MeanAlphaEnergy = 16.0;
1222 G4double TotalAlphaEnergy;
1223
1224 do
1225 {
1226 G4double AlphaEnergy;
1227 TotalAlphaEnergy = 0;
1228
1229 // Walk through the alpha particles one at a time and sample each's
1230 // energy
1231 for(unsigned int i = 0; i < Alphas->size(); i++)
1232 {
1233 AlphaEnergy = RandomEngine_->G4SampleGaussian(MeanAlphaEnergy,
1234 2.35,
1236 // Assign the energy to the alpha particle
1237 Alphas->at(i)->SetKineticEnergy(AlphaEnergy);
1238
1239 // Add up the total amount of kinetic energy consumed.
1240 TotalAlphaEnergy += AlphaEnergy;
1241 }
1242
1243 // If true, decrement the mean alpha energy by 0.1 and try again.
1244 MeanAlphaEnergy -= 0.1;
1245 } while(TotalAlphaEnergy >= RemainingEnergy_); // Loop checking, 11.05.2015, T. Koi
1246
1247 // Subtract the total amount of energy that was assigned.
1248 RemainingEnergy_ -= TotalAlphaEnergy;
1249
1251}
1252
1254SampleGammaEnergies( std::vector< G4ReactionProduct* >* Gammas )
1255{
1257
1258 // Make sure that there is energy to assign to the gamma rays
1259 if(RemainingEnergy_ != 0)
1260 {
1261 G4double SampleEnergy;
1262
1263 // Sample from RemainingEnergy until it is all gone. Also,
1264 // RemainingEnergy should not be smaller than
1265 // G4FFGDefaultValues::MeanGammaEnergy. This will prevent the
1266 // sampling of a fractional portion of the Gaussian distribution
1267 // in an attempt to find a new gamma ray energy.
1268 G4int icounter=0;
1269 G4int icounter_max=1024;
1270 while(RemainingEnergy_ >= G4FFGDefaultValues::MeanGammaEnergy ) // Loop checking, 11.05.2015, T. Koi
1271 {
1272 icounter++;
1273 if ( icounter > icounter_max ) {
1274 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
1275 break;
1276 }
1277 SampleEnergy = RandomEngine_->
1278 G4SampleGaussian(G4FFGDefaultValues::MeanGammaEnergy, 1.0 * MeV, G4FFGEnumerations::POSITIVE);
1279 // Make sure that we didn't sample more energy than was available
1280 if(SampleEnergy <= RemainingEnergy_)
1281 {
1282 // If this energy assignment would leave less energy than the
1283 // 'intrinsic' minimal energy of a gamma ray then just assign
1284 // all of the remaining energy
1285 if(RemainingEnergy_ - SampleEnergy < 100 * keV)
1286 {
1287 SampleEnergy = RemainingEnergy_;
1288 }
1289
1290 // Create the new particle
1291 Gammas->push_back(new G4ReactionProduct());
1292
1293 // Set the properties
1294 Gammas->back()->SetDefinition(GammaDefinition_);
1295 Gammas->back()->SetTotalEnergy(SampleEnergy);
1296
1297 // Calculate how much is left
1298 RemainingEnergy_ -= SampleEnergy;
1299 }
1300 }
1301
1302 // If there is anything left over, the energy must be above 100 keV but
1303 // less than G4FFGDefaultValues::MeanGammaEnergy. Arbitrarily assign
1304 // RemainingEnergy to a new particle
1305 if(RemainingEnergy_ > 0)
1306 {
1307 SampleEnergy = RemainingEnergy_;
1308 Gammas->push_back(new G4ReactionProduct());
1309
1310 // Set the properties
1311 Gammas->back()->SetDefinition(GammaDefinition_);
1312 Gammas->back()->SetTotalEnergy(SampleEnergy);
1313
1314 // Calculate how much is left
1315 RemainingEnergy_ -= SampleEnergy;
1316 }
1317 }
1318
1320}
1321
1323SampleNeutronEnergies( std::vector< G4ReactionProduct* >* Neutrons )
1324{
1326
1327 // The condition of sampling more energy from the fission products than is
1328 // alloted is statistically unfavorable, but it could still happen. The
1329 // do-while loop prevents such an occurrence from happening
1330 G4double TotalNeutronEnergy;
1331 G4double NeutronEnergy;
1332
1333 // Make sure that we don't sample more energy than is available
1334 G4int icounter=0;
1335 G4int icounter_max=1024;
1336 do
1337 {
1338 icounter++;
1339 if ( icounter > icounter_max ) {
1340 G4cout << "Loop-counter exceeded the threshold value at " << __LINE__ << "th line of " << __FILE__ << "." << G4endl;
1341 break;
1342 }
1343 TotalNeutronEnergy = 0;
1344
1345 // Walk through the neutrons one at a time and sample the energies.
1346 // The gamma rays have not yet been sampled, so the last neutron will
1347 // have a NULL value for NextFragment
1348 for(unsigned int i = 0; i < Neutrons->size(); i++)
1349 {
1350 // Assign the energy to the neutron
1352 Neutrons->at(i)->SetKineticEnergy(NeutronEnergy);
1353
1354 // Add up the total amount of kinetic energy consumed.
1355 TotalNeutronEnergy +=NeutronEnergy;
1356 }
1357 } while (TotalNeutronEnergy > RemainingEnergy_); // Loop checking, 11.05.2015, T. Koi
1358
1359 // Subtract the total amount of energy that was assigned.
1360 RemainingEnergy_ -= TotalNeutronEnergy;
1361
1363}
1364
1366SetNubar( void )
1367{
1369
1370 G4int* WhichNubar;
1371 G4int* NubarWidth;
1372 G4double XFactor, BFactor;
1373
1375 {
1376 WhichNubar = const_cast<G4int*>(&SpontaneousNubar_[0][0]);
1377 NubarWidth = const_cast<G4int*>(&SpontaneousNubarWidth_[0][0]);
1378 } else
1379 {
1380 WhichNubar = const_cast<G4int*>(&NeutronInducedNubar_[0][0]);
1381 NubarWidth = const_cast<G4int*>(&NeutronInducedNubarWidth_[0][0]);
1382 }
1383
1384 XFactor = G4Pow::GetInstance()->powA(10.0, -13.0);
1385 BFactor = G4Pow::GetInstance()->powA(10.0, -4.0);
1386 Nubar_ = *(WhichNubar + 1) * IncidentEnergy_ * XFactor
1387 + *(WhichNubar + 2) * BFactor;
1388 while(*WhichNubar != -1) // Loop checking, 11.05.2015, T. Koi
1389 {
1390 if(*WhichNubar == Isotope_)
1391 {
1392 Nubar_ = *(WhichNubar + 1) * IncidentEnergy_ * XFactor
1393 + *(WhichNubar + 2) * BFactor;
1394
1395 break;
1396 }
1397 WhichNubar += 3;
1398 }
1399
1400 XFactor = G4Pow::GetInstance()->powN((G4double)10, -6);
1401 NubarWidth_ = *(NubarWidth + 1) * XFactor;
1402 while(*WhichNubar != -1) // Loop checking, 11.05.2015, T. Koi
1403 {
1404 if(*WhichNubar == Isotope_)
1405 {
1406 NubarWidth_ = *(NubarWidth + 1) * XFactor;
1407
1408 break;
1409 }
1410 WhichNubar += 2;
1411 }
1412
1414}
1415
1418{
1420
1421 // Initialize the new branch
1422 ProbabilityBranch* NewBranch = new ProbabilityBranch;
1424 NewBranch->Left = NULL;
1425 NewBranch->Right = NULL;
1426 NewBranch->Particle = GetParticleDefinition(YieldData->GetProduct(), YieldData->GetMetaState());
1427 NewBranch->IncidentEnergies = new G4double[YieldEnergyGroups_];
1428 NewBranch->ProbabilityRangeTop = new G4double[YieldEnergyGroups_];
1429 NewBranch->ProbabilityRangeBottom = new G4double[YieldEnergyGroups_];
1430 G4ArrayOps::Copy(YieldEnergyGroups_, NewBranch->ProbabilityRangeTop, YieldData->GetYieldProbability());
1431 G4ArrayOps::Copy(YieldEnergyGroups_, NewBranch->IncidentEnergies, YieldEnergies_);
1433
1434 // Check to see if the this is the smallest/largest particle. First, check
1435 // to see if this is the first particle in the system
1436 if(SmallestZ_ == NULL)
1437 {
1438 SmallestZ_ = SmallestA_ = LargestZ_ = LargestA_ = NewBranch->Particle;
1439 } else
1440 {
1441 G4bool IsSmallerZ = NewBranch->Particle->GetAtomicNumber() < SmallestZ_->GetAtomicNumber();
1442 G4bool IsSmallerA = NewBranch->Particle->GetAtomicMass() < SmallestA_->GetAtomicMass();
1443 G4bool IsLargerZ = NewBranch->Particle->GetAtomicNumber() > LargestZ_->GetAtomicNumber();
1444 G4bool IsLargerA = NewBranch->Particle->GetAtomicMass() > LargestA_->GetAtomicMass();
1445
1446 if(IsSmallerZ)
1447 {
1448 SmallestZ_ = NewBranch->Particle;
1449 }
1450
1451 if(IsLargerZ)
1452 {
1453 LargestA_ = NewBranch->Particle;
1454 }
1455
1456 if(IsSmallerA)
1457 {
1458 SmallestA_ = NewBranch->Particle;
1459 }
1460
1461 if(IsLargerA)
1462 {
1463 LargestA_ = NewBranch->Particle;
1464 }
1465 }
1466
1467 // Place the new branch
1468 // Determine which tree the new branch goes into
1469 G4int WhichTree = (G4int)floor((G4double)(BranchCount_ % TreeCount_));
1470 ProbabilityBranch** WhichBranch = &(Trees_[WhichTree].Trunk);
1471 Trees_[WhichTree].BranchCount++;
1472
1473 // Search for the position
1474 // Determine where the branch goes
1475 G4int BranchPosition = (G4int)floor((G4double)(BranchCount_ / TreeCount_)) + 1;
1476
1477 // Run through the tree until the end branch is reached
1478 while(BranchPosition > 1) // Loop checking, 11.05.2015, T. Koi
1479 {
1480 if(BranchPosition & 1)
1481 {
1482 // If the 1's bit is on then move to the next 'right' branch
1483 WhichBranch = &((*WhichBranch)->Right);
1484 } else
1485 {
1486 // If the 1's bit is off then move to the next 'down' branch
1487 WhichBranch = &((*WhichBranch)->Left);
1488 }
1489
1490 BranchPosition >>= 1;
1491 }
1492
1493 *WhichBranch = NewBranch;
1494 BranchCount_++;
1495
1497}
1498
1501{
1503
1504 // Burn each tree, one by one
1505 G4int WhichTree = 0;
1506 while(Trees_[WhichTree].IsEnd != TRUE) // Loop checking, 11.05.2015, T. Koi
1507 {
1508 BurnTree(Trees_[WhichTree].Trunk);
1509 delete Trees_[WhichTree].Trunk;
1510 delete[] Trees_[WhichTree].ProbabilityRangeEnd;
1511 WhichTree++;
1512 }
1513
1514 // Delete each dynamically allocated variable
1515 delete ENDFData_;
1516 delete[] Trees_;
1517 delete[] DataTotal_;
1518 delete[] MaintainNormalizedData_;
1519 delete ElementNames_;
1520 delete RandomEngine_;
1521
1523}
1524
1526BurnTree( ProbabilityBranch* Branch )
1527{
1529
1530 // Check to see it Branch exists. Branch will be a null pointer if it
1531 // doesn't exist
1532 if(Branch)
1533 {
1534 // Burn down before you burn up
1535 BurnTree(Branch->Left);
1536 delete Branch->Left;
1537 BurnTree(Branch->Right);
1538 delete Branch->Right;
1539
1540 delete[] Branch->IncidentEnergies;
1541 delete[] Branch->ProbabilityRangeTop;
1542 delete[] Branch->ProbabilityRangeBottom;
1543 }
1544
1546}
1547
double A(double temperature)
std::vector< G4DynamicParticle * > G4DynamicParticleVector
@ JustWarning
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4FFG_DATA_FUNCTIONENTER__
#define G4FFG_DATA_FUNCTIONLEAVE__
#define G4FFG_FUNCTIONLEAVE__
#define G4FFG_FUNCTIONENTER__
#define G4FFG_RECURSIVE_FUNCTIONENTER__
#define G4FFG_RECURSIVE_FUNCTIONLEAVE__
#define G4FFG_LOCATION__
#define G4FFG_SPACING__
#define G4MUTEX_INITIALIZER
Definition: G4Threading.hh:85
std::mutex G4Mutex
Definition: G4Threading.hh:81
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define TRUE
Definition: Globals.hh:27
#define FALSE
Definition: Globals.hh:23
double getZ() const
void setY(double)
void setRThetaPhi(double r, double theta, double phi)
void setR(double s)
void setZ(double)
double mag() const
void set(double x, double y, double z)
double getX() const
void setX(double)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
double getY() const
static G4Alpha * Definition()
Definition: G4Alpha.cc:48
G4int G4GetNumberOfFissionProducts(void)
void G4SetVerbosity(G4int WhatVerbosity)
G4double * G4GetEnergyGroupValues(void)
G4ENDFYieldDataContainer * G4GetYield(G4int WhichYield)
G4int G4GetNumberOfEnergyGroups(void)
G4FFGEnumerations::MetaState GetMetaState(void)
G4double G4SampleWatt(G4int WhatIsotope, G4FFGEnumerations::FissionCause WhatCause, G4double WhatEnergy)
G4double G4SampleGaussian(G4double Mean, G4double StdDev)
G4int G4SampleIntegerGaussian(G4double Mean, G4double StdDev)
void G4SetVerbosity(G4int WhatVerbosity)
G4double G4SampleUniform(void)
void Renormalize(ProbabilityBranch *Branch)
void BurnTree(ProbabilityBranch *Branch)
void G4SetEnergy(G4double WhatIncidentEnergy)
void SampleAlphaEnergies(std::vector< G4ReactionProduct * > *Alphas)
void SampleNeutronEnergies(std::vector< G4ReactionProduct * > *Neutrons)
G4Ions * FindParticleBranchSearch(ProbabilityBranch *Branch, G4double RandomParticle, G4int EnergyGroup1, G4int EnergyGroup2)
G4FissionProductYieldDist(G4int WhichIsotope, G4FFGEnumerations::MetaState WhichMetaState, G4FFGEnumerations::FissionCause WhichCause, G4FFGEnumerations::YieldType WhichYieldType, std::istringstream &dataStream)
G4DynamicParticleVector * G4GetFission(void)
virtual void SortProbability(G4ENDFYieldDataContainer *YieldData)
G4String MakeFileName(G4int Isotope, G4FFGEnumerations::MetaState MetaState)
G4String MakeIsotopeName(G4int Isotope, G4FFGEnumerations::MetaState MetaState)
const G4FFGEnumerations::FissionCause Cause_
G4Ions * FindParticle(G4double RandomParticle)
G4Ions * FindParticleInterpolation(G4double RandomParticle, G4int LowerEnergyGroup)
virtual void GenerateNeutrons(std::vector< G4ReactionProduct * > *Neutrons)
G4Ions * GetParticleDefinition(G4int Product, G4FFGEnumerations::MetaState MetaState)
void G4SetAlphaProduction(G4double WhatAlphaProduction)
virtual void GenerateAlphas(std::vector< G4ReactionProduct * > *Alphas)
G4DynamicParticle * MakeG4DynamicParticle(G4ReactionProduct *)
virtual G4Ions * GetFissionProduct(void)=0
void SampleGammaEnergies(std::vector< G4ReactionProduct * > *Gammas)
void G4SetVerbosity(G4int WhatVerbosity)
const G4FFGEnumerations::YieldType YieldType_
void G4SetTernaryProbability(G4double TernaryProbability)
G4Ions * FindParticleExtrapolation(G4double RandomParticle, G4bool LowerEnergyGroupExists)
static G4Gamma * Definition()
Definition: G4Gamma.cc:48
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4IonTable * GetIonTable()
Definition: G4IonTable.cc:170
Definition: G4Ions.hh:52
static G4Neutron * Definition()
Definition: G4Neutron.cc:53
G4int GetAtomicNumber() const
G4int GetAtomicMass() const
const G4String & GetParticleName() const
static const G4String theString[100]
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powN(G4double x, G4int n) const
Definition: G4Pow.cc:166
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
void SetMomentum(const G4double x, const G4double y, const G4double z)
const G4ParticleDefinition * GetDefinition() const
G4ThreeVector GetMomentum() const
void Add(G4int Elements, T *To, T *A1, T *A2=NULL)
Definition: G4ArrayOps.hh:77
void Multiply(G4int Elements, T *To, T *M1, T *M2=NULL)
Definition: G4ArrayOps.hh:138
void Divide(G4int Elements, T *To, T *Numerator, T *Denominator=NULL)
Definition: G4ArrayOps.hh:178
void DeleteVectorOfPointers(std::vector< T > &Vector)
Definition: G4ArrayOps.hh:216
void Copy(G4int Elements, T *To, T *From)
Definition: G4ArrayOps.hh:63
void Set(G4int Elements, T *To, T Value)
Definition: G4ArrayOps.hh:51
ProbabilityBranch * Left
G4double * ProbabilityRangeTop
ProbabilityBranch * Right
G4double * ProbabilityRangeBottom
G4double * ProbabilityRangeEnd
ProbabilityBranch * Trunk