Geant4 11.4.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4FermiBreakUpAN Class Reference

#include <G4FermiBreakUpAN.hh>

Inheritance diagram for G4FermiBreakUpAN:

Public Member Functions

 G4FermiBreakUpAN (G4int verbosity=0)
 ~G4FermiBreakUpAN () override=default
void Initialise () override
G4bool IsApplicable (G4int Z, G4int A, G4double eexc) const override
void BreakFragment (G4FragmentVector *results, G4Fragment *theNucleus) override
std::vector< G4FermiParticleBreakItUp (const G4FermiParticle &nucleus) const
Public Member Functions inherited from G4VFermiBreakUp
 G4VFermiBreakUp ()
virtual ~G4VFermiBreakUp ()=default
 G4VFermiBreakUp (const G4VFermiBreakUp &right)=delete
const G4VFermiBreakUpoperator= (const G4VFermiBreakUp &right)=delete
G4bool operator== (const G4VFermiBreakUp &right) const =delete
G4bool operator!= (const G4VFermiBreakUp &right) const =delete
void SetVerbose (G4int val)

Additional Inherited Members

Protected Attributes inherited from G4VFermiBreakUp
G4int verbose {0}

Detailed Description

Definition at line 48 of file G4FermiBreakUpAN.hh.

Constructor & Destructor Documentation

◆ G4FermiBreakUpAN()

G4FermiBreakUpAN::G4FermiBreakUpAN ( G4int verbosity = 0)
explicit

Definition at line 162 of file G4FermiBreakUpAN.cc.

163 : splits_(G4FermiAtomicMass(MAX_A)),
164 secID_(G4PhysicsModelCatalog::GetModelID("model_G4FermiBreakUpVI")),
165 verbosity_(verbosity)
166{}
static G4int GetModelID(const G4int modelIndex)

◆ ~G4FermiBreakUpAN()

G4FermiBreakUpAN::~G4FermiBreakUpAN ( )
overridedefault

Member Function Documentation

◆ BreakFragment()

void G4FermiBreakUpAN::BreakFragment ( G4FragmentVector * results,
G4Fragment * theNucleus )
overridevirtual

Reimplemented from G4VFermiBreakUp.

Definition at line 239 of file G4FermiBreakUpAN.cc.

240{
241 if (theNucleus == nullptr || results == nullptr) {
243 ed << "G4Fragment or result G4FragmentVector is not set in FermiBreakUp";
244 G4Exception("G4FermiBreakUpAN::BreakFragment()", "Fermi003", FatalErrorInArgument, ed);
245 return;
246 }
247
248 const auto particle =
249 G4FermiParticle(G4FermiAtomicMass(theNucleus->GetA_asInt()),
250 G4FermiChargeNumber(theNucleus->GetZ_asInt()), theNucleus->GetMomentum());
251 const auto fragments = BreakItUp(particle);
252
253 // decay impossible
254 if (fragments.size() <= 1) { return; }
255
256 const auto creationTime = theNucleus->GetCreationTime();
257 // primary should be deleted
258 delete theNucleus;
259
260 for (const auto& fragment : fragments) {
261 auto fr = new G4Fragment(static_cast<G4int>(fragment.GetAtomicMass()),
262 static_cast<G4int>(fragment.GetChargeNumber()),
263 fragment.GetMomentum());
264 results->push_back(fr);
265 fr->SetCreationTime(creationTime);
266 fr->SetCreatorModelID(secID_);
267 }
268}
@ FatalErrorInArgument
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::ostringstream G4ExceptionDescription
int G4int
Definition G4Types.hh:85
std::vector< G4FermiParticle > BreakItUp(const G4FermiParticle &nucleus) const
const G4LorentzVector & GetMomentum() const
G4double GetCreationTime() const
G4int GetZ_asInt() const
G4int GetA_asInt() const

◆ BreakItUp()

std::vector< G4FermiParticle > G4FermiBreakUpAN::BreakItUp ( const G4FermiParticle & nucleus) const

Definition at line 168 of file G4FermiBreakUpAN.cc.

169{
170 FERMI_LOG_DEBUG(verbosity_, "Breaking up particle: " << particle);
171
172 if (particle.GetExcitationEnergy() < 0.) {
173 FERMI_LOG_DEBUG(verbosity_, "G4FermiParticle is stable with excitation energy = "
174 << particle.GetExcitationEnergy());
175 return {particle};
176 }
177
178 const auto& splits = splits_.GetSplits(particle.GetAtomicMass(), particle.GetChargeNumber());
179 FERMI_LOG_DEBUG(verbosity_,
180 "Selecting Split for " << particle << " from " << splits.size() << " splits");
181 if (splits.empty()) {
182 FERMI_LOG_DEBUG(verbosity_, "No splits found");
183 return {particle};
184 }
185
186 // get phase space weights for every split
187 // we can't cache them, because calculations is probabilistic
188 weights_.resize(splits.size());
189 std::transform(splits.begin(), splits.end(), weights_.begin(),
190 [atomicMass = particle.GetAtomicMass(),
191 totalEnergy = particle.GetMomentum().m()](const auto& split) {
192 return G4FermiSplitter::DecayWeight(split, atomicMass, totalEnergy);
193 });
194
195 if (std::all_of(weights_.begin(), weights_.end(), [](auto weight) { return weight == 0.; })) {
196 FERMI_LOG_DEBUG(verbosity_, "Every split has zero weight");
197 return {particle};
198 }
199
200 const auto& chosenSplit = splits[SampleWeightDistribution(weights_)];
201 FERMI_LOG_DEBUG(verbosity_,
202 "From " << splits.size() << " splits chosen split: " << LogSplit(chosenSplit));
203
204 return SplitToParticles(particle, chosenSplit);
205}
#define FERMI_LOG_DEBUG(verbosity, msg)

Referenced by BreakFragment().

◆ Initialise()

void G4FermiBreakUpAN::Initialise ( )
overridevirtual

Reimplemented from G4VFermiBreakUp.

Definition at line 207 of file G4FermiBreakUpAN.cc.

208{
209 if (G4NucleiProperties::GetNuclearMass(2, 0) <= 0.) {
210 G4BaryonConstructor pCBar;
211 pCBar.ConstructParticle();
212 }
214
215 {
216 auto pool = G4FermiFragmentPoolAN::DefaultPoolANSource();
217 pool.Initialize();
219 }
220
221 // order is important here, we use G4FermiFragmentPool to create splits!
222 splits_ = PossibleSplits();
223 for (auto a = 1; a < MAX_A; ++a) {
224 for (auto z = 0; z <= a; ++z) {
225 const auto atomicMass = G4FermiAtomicMass(a);
226 const auto chargeNumber = G4FermiChargeNumber(z);
227
228 splits_.InsertSplits(atomicMass, chargeNumber,
229 G4FermiSplitter::GenerateSplits({atomicMass, chargeNumber}));
230 }
231 }
232}
static void ConstructParticle()
void Initialize(const DataSource &dataSource)
static G4FermiFragmentPoolAN & Instance()
static G4FermiNucleiProperties & Instance()
static void GenerateSplits(G4FermiNucleiData nucleiData, std::vector< G4FermiFragmentVector > &splits)
static G4double GetNuclearMass(const G4double A, const G4double Z)

◆ IsApplicable()

G4bool G4FermiBreakUpAN::IsApplicable ( G4int Z,
G4int A,
G4double eexc ) const
overridevirtual

Reimplemented from G4VFermiBreakUp.

Definition at line 234 of file G4FermiBreakUpAN.cc.

235{
236 return Z < MAX_Z && A < MAX_A;
237}
const G4double A[17]

The documentation for this class was generated from the following files: