Geant4-11
Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4WilsonAblationModel Class Reference

#include <G4WilsonAblationModel.hh>

Inheritance diagram for G4WilsonAblationModel:
G4VEvaporation

Public Types

typedef std::vector< G4ParticleDefinition * > VectorOfFragmentTypes
 

Public Member Functions

virtual void BreakFragment (G4FragmentVector *, G4Fragment *theNucleus)
 
G4FragmentVectorBreakItUp (const G4Fragment &theNucleus)
 
 G4WilsonAblationModel ()
 
G4VFermiBreakUpGetFermiBreakUp () const
 
G4VEvaporationChannelGetFissionChannel ()
 
size_t GetNumberOfChannels () const
 
G4VEvaporationChannelGetPhotonEvaporation ()
 
G4bool GetProduceSecondaries ()
 
G4int GetVerboseLevel ()
 
virtual void InitialiseChannels ()
 
void SetFermiBreakUp (G4VFermiBreakUp *ptr)
 
void SetOPTxs (G4int opt)
 
virtual void SetPhotonEvaporation (G4VEvaporationChannel *ptr)
 
void SetProduceSecondaries (G4bool)
 
void SetVerboseLevel (G4int)
 
void UseSICB (G4bool use)
 
virtual ~G4WilsonAblationModel ()
 

Protected Member Functions

void CleanChannels ()
 

Protected Attributes

G4int OPTxs
 
G4VEvaporationFactorytheChannelFactory
 
std::vector< G4VEvaporationChannel * > * theChannels
 
G4VFermiBreakUptheFBU
 
G4VEvaporationChannelthePhotonEvaporation
 
G4bool useSICB
 

Private Member Functions

G4bool operator!= (const G4VEvaporation &right) const =delete
 
G4bool operator== (const G4VEvaporation &right) const =delete
 
void PrintWelcomeMessage ()
 
void SelectSecondariesByDefault (G4ThreeVector)
 
void SelectSecondariesByEvaporation (G4Fragment *)
 

Private Attributes

G4double B
 
VectorOfFragmentTypes evapType
 
G4FragmentVectorfragmentVector
 
G4ParticleDefinitionfragType [6]
 
G4double fSig [200]
 
G4int nFragTypes
 
G4bool produceSecondaries
 
G4int secID
 
G4int verboseLevel
 

Detailed Description

Definition at line 84 of file G4WilsonAblationModel.hh.

Member Typedef Documentation

◆ VectorOfFragmentTypes

Definition at line 90 of file G4WilsonAblationModel.hh.

Constructor & Destructor Documentation

◆ G4WilsonAblationModel()

G4WilsonAblationModel::G4WilsonAblationModel ( )

Definition at line 117 of file G4WilsonAblationModel.cc.

118{
119//
120//
121// Send message to stdout to advise that the G4Abrasion model is being used.
122//
124//
125//
126// Set the default verbose level to 0 - no output.
127//
128 verboseLevel = 0;
129//
130//
131// Set the binding energy per nucleon .... did I mention that this is a crude
132// model for nuclear de-excitation?
133//
134 B = 10.0 * MeV;
135//
136//
137// It is possuble to switch off secondary particle production (other than the
138// final nuclear fragment). The default is on.
139//
140 produceSecondaries = true;
141//
142//
143// Now we need to define the decay modes. We're using the G4Evaporation model
144// to help determine the kinematics of the decay.
145//
146 nFragTypes = 6;
148 fragType[1] = G4He3::He3();
153 for(G4int i=0; i<200; ++i) { fSig[i] = 0.0; }
154//
155//
156// Set verboseLevel default to no output.
157//
158 verboseLevel = 0;
161//
162//
163// Set defaults for evaporation classes. These can be overridden by user
164// "set" methods.
165//
166 OPTxs = 3;
167 useSICB = false;
168 fragmentVector = 0;
169
170 secID = G4PhysicsModelCatalog::GetModelID("model_G4WilsonAblationModel");
171}
static constexpr double MeV
Definition: G4SIunits.hh:200
int G4int
Definition: G4Types.hh:85
static G4Alpha * Alpha()
Definition: G4Alpha.cc:88
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:93
static G4He3 * He3()
Definition: G4He3.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:103
static G4int GetModelID(const G4int modelIndex)
static G4Proton * Proton()
Definition: G4Proton.cc:92
static G4Triton * Triton()
Definition: G4Triton.cc:93
virtual std::vector< G4VEvaporationChannel * > * GetChannel()=0
G4VEvaporationFactory * theChannelFactory
std::vector< G4VEvaporationChannel * > * theChannels
G4FragmentVector * fragmentVector
G4ParticleDefinition * fragType[6]

References G4Alpha::Alpha(), B, G4Deuteron::Deuteron(), fragmentVector, fragType, fSig, G4VEvaporationFactory::GetChannel(), G4PhysicsModelCatalog::GetModelID(), G4He3::He3(), MeV, G4Neutron::Neutron(), nFragTypes, G4VEvaporation::OPTxs, PrintWelcomeMessage(), produceSecondaries, G4Proton::Proton(), secID, G4VEvaporation::theChannelFactory, G4VEvaporation::theChannels, G4Triton::Triton(), G4VEvaporation::useSICB, and verboseLevel.

◆ ~G4WilsonAblationModel()

G4WilsonAblationModel::~G4WilsonAblationModel ( )
virtual

Definition at line 174 of file G4WilsonAblationModel.cc.

175{}

Member Function Documentation

◆ BreakFragment()

void G4VEvaporation::BreakFragment ( G4FragmentVector ,
G4Fragment theNucleus 
)
virtualinherited

Reimplemented in G4Evaporation.

Definition at line 75 of file G4VEvaporation.cc.

76{}

Referenced by G4ExcitationHandler::BreakItUp().

◆ BreakItUp()

G4FragmentVector * G4WilsonAblationModel::BreakItUp ( const G4Fragment theNucleus)

Definition at line 179 of file G4WilsonAblationModel.cc.

181{
182//
183//
184// Initilise the pointer to the G4FragmentVector used to return the information
185// about the breakup.
186//
188 fragmentVector->clear();
189//
190//
191// Get the A, Z and excitation of the nucleus.
192//
193 G4int A = theNucleus.GetA_asInt();
194 G4int Z = theNucleus.GetZ_asInt();
195 G4double ex = theNucleus.GetExcitationEnergy();
196 if (verboseLevel >= 2)
197 {
198 G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
199 <<"oooooooooooooooooooooooooooooooooooooooo"
200 <<G4endl;
201 G4cout.precision(6);
202 G4cout <<"IN G4WilsonAblationModel" <<G4endl;
203 G4cout <<"Initial prefragment A=" <<A
204 <<", Z=" <<Z
205 <<", excitation energy = " <<ex/MeV <<" MeV"
206 <<G4endl;
207 }
208//
209//
210// Check that there is a nucleus to speak of. It's possible there isn't one
211// or its just a proton or neutron. In either case, the excitation energy
212// (from the Lorentz vector) is not used.
213//
214 if (A == 0)
215 {
216 if (verboseLevel >= 2)
217 {
218 G4cout <<"No nucleus to decay" <<G4endl;
219 G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
220 <<"oooooooooooooooooooooooooooooooooooooooo"
221 <<G4endl;
222 }
223 return fragmentVector;
224 }
225 else if (A == 1)
226 {
227 G4LorentzVector lorentzVector = theNucleus.GetMomentum();
228 lorentzVector.setE(lorentzVector.e()-ex+10.0*eV);
229 if (Z == 0)
230 {
231 G4Fragment *fragment = new G4Fragment(lorentzVector,G4Neutron::Neutron());
232 if (fragment != nullptr) { fragment->SetCreatorModelID(secID); }
233 fragmentVector->push_back(fragment);
234 }
235 else
236 {
237 G4Fragment *fragment = new G4Fragment(lorentzVector,G4Proton::Proton());
238 if (fragment != nullptr) { fragment->SetCreatorModelID(secID); }
239 fragmentVector->push_back(fragment);
240 }
241 if (verboseLevel >= 2)
242 {
243 G4cout <<"Final fragment is in fact only a nucleon) :" <<G4endl;
244 G4cout <<(*fragmentVector)[0] <<G4endl;
245 G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
246 <<"oooooooooooooooooooooooooooooooooooooooo"
247 <<G4endl;
248 }
249 return fragmentVector;
250 }
251//
252//
253// Then the number of nucleons ablated (either as nucleons or light nuclear
254// fragments) is based on a simple argument for the binding energy per nucleon.
255//
256 G4int DAabl = (G4int) (ex / B);
257 if (DAabl > A) DAabl = A;
258// The following lines are no longer accurate given we now treat the final fragment
259// if (verboseLevel >= 2)
260// G4cout <<"Number of nucleons ejected = " <<DAabl <<G4endl;
261
262//
263//
264// Determine the nuclear fragment from the ablation process by sampling the
265// Rudstam equation.
266//
267 G4int AF = A - DAabl;
268 G4int ZF = 0;
269
270 if (AF > 0)
271 {
272 G4Pow* g4calc = G4Pow::GetInstance();
273 G4double AFd = (G4double) AF;
274 G4double R = 11.8 / g4calc->powZ(AF, 0.45);
275 G4int minZ = std::max(1, Z - DAabl);
276//
277//
278// Here we define an integral probability distribution based on the Rudstam
279// equation assuming a constant AF.
280//
281 G4int zmax = std::min(199, Z);
282 G4double sum = 0.0;
283 for (ZF=minZ; ZF<=zmax; ++ZF)
284 {
285 sum += G4Exp(-R*g4calc->powA(std::abs(ZF - 0.486*AFd + 3.8E-04*AFd*AFd),1.5));
286 fSig[ZF] = sum;
287 }
288//
289//
290// Now sample that distribution to determine a value for ZF.
291//
292 sum *= G4UniformRand();
293 for (ZF=minZ; ZF<=zmax; ++ZF) {
294 if(sum <= fSig[ZF]) { break; }
295 }
296 }
297 G4int DZabl = Z - ZF;
298//
299//
300// Now determine the nucleons or nuclei which have bee ablated. The preference
301// is for the production of alphas, then other nuclei in order of decreasing
302// binding energy. The energies assigned to the products of the decay are
303// provisional for the moment (the 10eV is just to avoid errors with negative
304// excitation energies due to rounding).
305//
306 G4double totalEpost = 0.0;
307 evapType.clear();
308 for (G4int ift=0; ift<nFragTypes; ift++)
309 {
310 G4ParticleDefinition *type = fragType[ift];
311 G4double n = std::floor((G4double) DAabl / type->GetBaryonNumber() + 1.0E-10);
312 G4double n1 = 1.0E+10;
313 if (fragType[ift]->GetPDGCharge() > 0.0)
314 n1 = std::floor((G4double) DZabl / type->GetPDGCharge() + 1.0E-10);
315 if (n > n1) n = n1;
316 if (n > 0.0)
317 {
318 G4double mass = type->GetPDGMass();
319 for (G4int j=0; j<(G4int) n; j++)
320 {
321 totalEpost += mass;
322 evapType.push_back(type);
323 }
324 DAabl -= (G4int) (n * type->GetBaryonNumber() + 1.0E-10);
325 DZabl -= (G4int) (n * type->GetPDGCharge() + 1.0E-10);
326 }
327 }
328//
329//
330// Determine the properties of the final nuclear fragment. Note that if
331// the final fragment is predicted to have a nucleon number of zero, then
332// really it's the particle last in the vector evapType which becomes the
333// final fragment. Therefore delete this from the vector if this is the
334// case.
335//
336 G4double massFinalFrag = 0.0;
337 if (AF > 0)
339 GetIonMass(ZF,AF);
340 else
341 {
342 G4ParticleDefinition *type = evapType[evapType.size()-1];
343 AF = type->GetBaryonNumber();
344 ZF = (G4int) (type->GetPDGCharge() + 1.0E-10);
345 evapType.erase(evapType.end()-1);
346 }
347 totalEpost += massFinalFrag;
348//
349//
350// Provide verbose output on the nuclear fragment if requested.
351//
352 if (verboseLevel >= 2)
353 {
354 G4cout <<"Final fragment A=" <<AF
355 <<", Z=" <<ZF
356 <<G4endl;
357 for (G4int ift=0; ift<nFragTypes; ift++)
358 {
359 G4ParticleDefinition *type = fragType[ift];
360 G4int n = std::count(evapType.begin(),evapType.end(),type);
361 if (n > 0)
362 G4cout <<"Particle type: " <<std::setw(10) <<type->GetParticleName()
363 <<", number of particles emitted = " <<n <<G4endl;
364 }
365 }
366//
367// Add the total energy from the fragment. Note that the fragment is assumed
368// to be de-excited and does not undergo photo-evaporation .... I did mention
369// this is a bit of a crude model?
370//
371 G4double massPreFrag = theNucleus.GetGroundStateMass();
372 G4double totalEpre = massPreFrag + ex;
373 G4double excess = totalEpre - totalEpost;
374// G4Fragment *resultNucleus(theNucleus);
375 G4Fragment *resultNucleus = new G4Fragment(A, Z, theNucleus.GetMomentum());
376 G4ThreeVector boost(0.0,0.0,0.0);
377 G4int nEvap = 0;
378 if (produceSecondaries && evapType.size()>0)
379 {
380 if (excess > 0.0)
381 {
382 SelectSecondariesByEvaporation (resultNucleus);
383 nEvap = fragmentVector->size();
384 boost = resultNucleus->GetMomentum().findBoostToCM();
385 if (evapType.size() > 0)
387 }
388 else
390 }
391
392 if (AF > 0)
393 {
395 GetIonMass(ZF,AF);
396 G4double e = mass + 10.0*eV;
397 G4double p = std::sqrt(e*e-mass*mass);
398 G4ThreeVector direction(0.0,0.0,1.0);
399 G4LorentzVector lorentzVector = G4LorentzVector(direction*p, e);
400 lorentzVector.boost(-boost);
401 G4Fragment* frag = new G4Fragment(AF, ZF, lorentzVector);
402 if (frag != nullptr) { frag->SetCreatorModelID(secID); }
403 fragmentVector->push_back(frag);
404 }
405 delete resultNucleus;
406//
407//
408// Provide verbose output on the ablation products if requested.
409//
410 if (verboseLevel >= 2)
411 {
412 if (nEvap > 0)
413 {
414 G4cout <<"----------------------" <<G4endl;
415 G4cout <<"Evaporated particles :" <<G4endl;
416 G4cout <<"----------------------" <<G4endl;
417 }
418 G4int ie = 0;
419 G4FragmentVector::iterator iter;
420 for (iter = fragmentVector->begin(); iter != fragmentVector->end(); iter++)
421 {
422 if (ie == nEvap)
423 {
424// G4cout <<*iter <<G4endl;
425 G4cout <<"---------------------------------" <<G4endl;
426 G4cout <<"Particles from default emission :" <<G4endl;
427 G4cout <<"---------------------------------" <<G4endl;
428 }
429 G4cout <<*iter <<G4endl;
430 }
431 G4cout <<"oooooooooooooooooooooooooooooooooooooooo"
432 <<"oooooooooooooooooooooooooooooooooooooooo"
433 <<G4endl;
434 }
435
436 return fragmentVector;
437}
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:64
CLHEP::HepLorentzVector G4LorentzVector
static constexpr double eV
Definition: G4SIunits.hh:201
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
const G4int Z[17]
const G4double A[17]
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:52
HepLorentzVector & boost(double, double, double)
Hep3Vector findBoostToCM() const
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:304
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:299
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:323
void SetCreatorModelID(G4int value)
Definition: G4Fragment.hh:428
G4int GetZ_asInt() const
Definition: G4Fragment.hh:276
G4int GetA_asInt() const
Definition: G4Fragment.hh:271
G4double GetPDGCharge() const
const G4String & GetParticleName() const
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
Definition: G4Pow.hh:49
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
G4double powZ(G4int Z, G4double y) const
Definition: G4Pow.hh:225
G4double powA(G4double A, G4double y) const
Definition: G4Pow.hh:230
void SelectSecondariesByEvaporation(G4Fragment *)
VectorOfFragmentTypes evapType
void SelectSecondariesByDefault(G4ThreeVector)
T max(const T t1, const T t2)
brief Return the largest of the two arguments
T min(const T t1, const T t2)
brief Return the smallest of the two arguments

References A, B, CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::e(), eV, evapType, CLHEP::HepLorentzVector::findBoostToCM(), fragmentVector, fragType, fSig, G4cout, G4endl, G4Exp(), G4UniformRand, G4Fragment::GetA_asInt(), G4ParticleDefinition::GetBaryonNumber(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetGroundStateMass(), G4Pow::GetInstance(), G4ParticleTable::GetIonTable(), G4Fragment::GetMomentum(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4Fragment::GetZ_asInt(), G4INCL::Math::max(), MeV, G4INCL::Math::min(), CLHEP::detail::n, G4Neutron::Neutron(), nFragTypes, G4Pow::powA(), G4Pow::powZ(), produceSecondaries, G4Proton::Proton(), secID, SelectSecondariesByDefault(), SelectSecondariesByEvaporation(), G4Fragment::SetCreatorModelID(), CLHEP::HepLorentzVector::setE(), verboseLevel, and Z.

◆ CleanChannels()

void G4VEvaporation::CleanChannels ( )
protectedinherited

Definition at line 49 of file G4VEvaporation.cc.

50{
51 // clean all except photon evaporation
52 if(theChannels) {
53 for (size_t i=1; i<theChannels->size(); ++i) {
54 delete (*theChannels)[i];
55 }
56 delete theChannels;
57 theChannels = nullptr;
58 }
59}

References G4VEvaporation::theChannels.

Referenced by G4Evaporation::SetCombinedChannel(), G4Evaporation::SetDefaultChannel(), G4Evaporation::SetGEMChannel(), G4Evaporation::SetGEMVIChannel(), and G4VEvaporation::~G4VEvaporation().

◆ GetFermiBreakUp()

G4VFermiBreakUp * G4VEvaporation::GetFermiBreakUp ( ) const
inlineinherited

Definition at line 108 of file G4VEvaporation.hh.

109{
110 return theFBU;
111}
G4VFermiBreakUp * theFBU

References G4VEvaporation::theFBU.

◆ GetFissionChannel()

G4VEvaporationChannel * G4VEvaporation::GetFissionChannel ( )
inlineinherited

Definition at line 118 of file G4VEvaporation.hh.

119{
120 return (theChannels && theChannels->size() > 1) ? (*theChannels)[1] : nullptr;
121}

References G4VEvaporation::theChannels.

Referenced by G4INCLXXInterface::G4INCLXXInterface().

◆ GetNumberOfChannels()

size_t G4VEvaporation::GetNumberOfChannels ( ) const
inlineinherited

Definition at line 133 of file G4VEvaporation.hh.

134{
135 return theChannels ? theChannels->size() : 0;
136}

References G4VEvaporation::theChannels.

Referenced by G4ExcitationHandler::SetDeexChannelsType().

◆ GetPhotonEvaporation()

G4VEvaporationChannel * G4VEvaporation::GetPhotonEvaporation ( )
inlineinherited

Definition at line 113 of file G4VEvaporation.hh.

114{
116}
G4VEvaporationChannel * thePhotonEvaporation

References G4VEvaporation::thePhotonEvaporation.

Referenced by G4ExcitationHandler::SetEvaporation().

◆ GetProduceSecondaries()

G4bool G4WilsonAblationModel::GetProduceSecondaries ( )
inline

Definition at line 121 of file G4WilsonAblationModel.hh.

122 {return produceSecondaries;}

References produceSecondaries.

◆ GetVerboseLevel()

G4int G4WilsonAblationModel::GetVerboseLevel ( )
inline

Definition at line 129 of file G4WilsonAblationModel.hh.

130 {return verboseLevel;}

References verboseLevel.

◆ InitialiseChannels()

void G4VEvaporation::InitialiseChannels ( )
virtualinherited

Reimplemented in G4Evaporation.

Definition at line 61 of file G4VEvaporation.cc.

62{}

Referenced by G4ExcitationHandler::Initialise().

◆ operator!=()

G4bool G4VEvaporation::operator!= ( const G4VEvaporation right) const
privatedeleteinherited

◆ operator==()

G4bool G4VEvaporation::operator== ( const G4VEvaporation right) const
privatedeleteinherited

◆ PrintWelcomeMessage()

void G4WilsonAblationModel::PrintWelcomeMessage ( )
private

Definition at line 594 of file G4WilsonAblationModel.cc.

595{
596 G4cout <<G4endl;
597 G4cout <<" *****************************************************************"
598 <<G4endl;
599 G4cout <<" Nuclear ablation model for nuclear-nuclear interactions activated"
600 <<G4endl;
601 G4cout <<" (Written by QinetiQ Ltd for the European Space Agency)"
602 <<G4endl;
603 G4cout <<" !!! WARNING: This model is not well validation and should not be used for accurate simulation !!!"
604 <<G4endl;
605 G4cout <<" *****************************************************************"
606 <<G4endl;
607 G4cout << G4endl;
608
609 return;
610}

References G4cout, and G4endl.

Referenced by G4WilsonAblationModel().

◆ SelectSecondariesByDefault()

void G4WilsonAblationModel::SelectSecondariesByDefault ( G4ThreeVector  boost)
private

Definition at line 566 of file G4WilsonAblationModel.cc.

567{
568 for (unsigned i=0; i<evapType.size(); i++)
569 {
571 G4double mass = type->GetPDGMass();
572 G4double e = mass + 10.0*eV;
573 G4double p = std::sqrt(e*e-mass*mass);
574 G4double costheta = 2.0*G4UniformRand() - 1.0;
575 G4double sintheta = std::sqrt((1.0 - costheta)*(1.0 + costheta));
576 G4double phi = twopi * G4UniformRand() * rad;
577 G4ThreeVector direction(sintheta*std::cos(phi),sintheta*std::sin(phi),costheta);
578 G4LorentzVector lorentzVector = G4LorentzVector(direction*p, e);
579 lorentzVector.boost(-boost);
580// Possibility that the following line is not correctly carrying over A and Z
581// from particle definition. Force values. PRT 03/12/2009.
582// G4Fragment *fragment =
583// new G4Fragment(lorentzVector, type);
584 G4int A = type->GetBaryonNumber();
585 G4int Z = (G4int) (type->GetPDGCharge() + 1.0E-10);
586 G4Fragment *fragment =
587 new G4Fragment(A, Z, lorentzVector);
588 if (fragment != nullptr) { fragment->SetCreatorModelID(secID); }
589 fragmentVector->push_back(fragment);
590 }
591}
static constexpr double twopi
Definition: G4SIunits.hh:56
static constexpr double rad
Definition: G4SIunits.hh:129

References A, CLHEP::HepLorentzVector::boost(), eV, evapType, fragmentVector, G4UniformRand, G4ParticleDefinition::GetBaryonNumber(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), rad, secID, G4Fragment::SetCreatorModelID(), twopi, and Z.

Referenced by BreakItUp().

◆ SelectSecondariesByEvaporation()

void G4WilsonAblationModel::SelectSecondariesByEvaporation ( G4Fragment intermediateNucleus)
private

Definition at line 440 of file G4WilsonAblationModel.cc.

442{
443 G4Fragment theResidualNucleus = *intermediateNucleus;
444 G4bool evaporate = true;
445 // Loop checking, 05-Aug-2015, Vladimir Ivanchenko
446 while (evaporate && evapType.size() != 0)
447 {
448//
449//
450// Here's the cheaky bit. We're hijacking the G4Evaporation model, in order to
451// more accurately sample to kinematics, but the species of the nuclear
452// fragments will be the ones of our choosing as above.
453//
454 std::vector <G4VEvaporationChannel*> theChannels1;
455 theChannels1.clear();
456 std::vector <G4VEvaporationChannel*>::iterator i;
457 VectorOfFragmentTypes::iterator iter;
458 std::vector <VectorOfFragmentTypes::iterator> iters;
459 iters.clear();
460 iter = std::find(evapType.begin(), evapType.end(), G4Alpha::Alpha());
461 if (iter != evapType.end())
462 {
463 theChannels1.push_back(new G4AlphaEvaporationChannel);
464 i = theChannels1.end() - 1;
465 (*i)->SetOPTxs(OPTxs);
466 (*i)->UseSICB(useSICB);
467// (*i)->Initialize(theResidualNucleus);
468 iters.push_back(iter);
469 }
470 iter = std::find(evapType.begin(), evapType.end(), G4He3::He3());
471 if (iter != evapType.end())
472 {
473 theChannels1.push_back(new G4He3EvaporationChannel);
474 i = theChannels1.end() - 1;
475 (*i)->SetOPTxs(OPTxs);
476 (*i)->UseSICB(useSICB);
477// (*i)->Initialize(theResidualNucleus);
478 iters.push_back(iter);
479 }
480 iter = std::find(evapType.begin(), evapType.end(), G4Triton::Triton());
481 if (iter != evapType.end())
482 {
483 theChannels1.push_back(new G4TritonEvaporationChannel);
484 i = theChannels1.end() - 1;
485 (*i)->SetOPTxs(OPTxs);
486 (*i)->UseSICB(useSICB);
487// (*i)->Initialize(theResidualNucleus);
488 iters.push_back(iter);
489 }
490 iter = std::find(evapType.begin(), evapType.end(), G4Deuteron::Deuteron());
491 if (iter != evapType.end())
492 {
493 theChannels1.push_back(new G4DeuteronEvaporationChannel);
494 i = theChannels1.end() - 1;
495 (*i)->SetOPTxs(OPTxs);
496 (*i)->UseSICB(useSICB);
497// (*i)->Initialize(theResidualNucleus);
498 iters.push_back(iter);
499 }
500 iter = std::find(evapType.begin(), evapType.end(), G4Proton::Proton());
501 if (iter != evapType.end())
502 {
503 theChannels1.push_back(new G4ProtonEvaporationChannel);
504 i = theChannels1.end() - 1;
505 (*i)->SetOPTxs(OPTxs);
506 (*i)->UseSICB(useSICB);
507// (*i)->Initialize(theResidualNucleus);
508 iters.push_back(iter);
509 }
510 iter = std::find(evapType.begin(), evapType.end(), G4Neutron::Neutron());
511 if (iter != evapType.end())
512 {
513 theChannels1.push_back(new G4NeutronEvaporationChannel);
514 i = theChannels1.end() - 1;
515 (*i)->SetOPTxs(OPTxs);
516 (*i)->UseSICB(useSICB);
517// (*i)->Initialize(theResidualNucleus);
518 iters.push_back(iter);
519 }
520 G4int nChannels = theChannels1.size();
521
522 G4double totalProb = 0.0;
523 G4int ich = 0;
524 G4double probEvapType[6] = {0.0};
525 std::vector<G4VEvaporationChannel*>::iterator iterEv;
526 for (iterEv=theChannels1.begin(); iterEv!=theChannels1.end(); iterEv++) {
527 totalProb += (*iterEv)->GetEmissionProbability(intermediateNucleus);
528 probEvapType[ich] = totalProb;
529 ++ich;
530 }
531 if (totalProb > 0.0) {
532//
533//
534// The emission probability for at least one of the evaporation channels is
535// positive, therefore work out which one should be selected and decay
536// the nucleus.
537//
538 G4double xi = totalProb*G4UniformRand();
539 G4int ii = 0;
540 for (ii=0; ii<nChannels; ii++) {
541 if (xi < probEvapType[ii]) { break; }
542 }
543 if (ii >= nChannels) { ii = nChannels - 1; }
544 G4FragmentVector *evaporationResult = theChannels1[ii]->
545 BreakUpFragment(intermediateNucleus);
546 if ((*evaporationResult)[0] != nullptr) { (*evaporationResult)[0]->SetCreatorModelID(secID); }
547 fragmentVector->push_back((*evaporationResult)[0]);
548 intermediateNucleus = (*evaporationResult)[1];
549 delete evaporationResult;
550 }
551 else
552 {
553//
554//
555// Probability for further evaporation is nil so have to escape from this
556// routine and set the energies of the secondaries to 10eV.
557//
558 evaporate = false;
559 }
560 }
561
562 return;
563}
bool G4bool
Definition: G4Types.hh:86

References G4Alpha::Alpha(), G4Deuteron::Deuteron(), evapType, fragmentVector, G4UniformRand, G4He3::He3(), G4Neutron::Neutron(), G4VEvaporation::OPTxs, G4Proton::Proton(), secID, G4Triton::Triton(), and G4VEvaporation::useSICB.

Referenced by BreakItUp().

◆ SetFermiBreakUp()

void G4VEvaporation::SetFermiBreakUp ( G4VFermiBreakUp ptr)
inlineinherited

Definition at line 103 of file G4VEvaporation.hh.

104{
105 theFBU = ptr;
106}

References G4VEvaporation::theFBU.

Referenced by G4ExcitationHandler::SetEvaporation(), and G4ExcitationHandler::SetFermiModel().

◆ SetOPTxs()

void G4VEvaporation::SetOPTxs ( G4int  opt)
inlineinherited

Definition at line 123 of file G4VEvaporation.hh.

124{
125 OPTxs = opt;
126}

References G4VEvaporation::OPTxs.

◆ SetPhotonEvaporation()

void G4VEvaporation::SetPhotonEvaporation ( G4VEvaporationChannel ptr)
virtualinherited

Definition at line 64 of file G4VEvaporation.cc.

65{
66 // photon evaporation channel is the first
67 // G4VEvaporation is responsible for its deletion
68 if(thePhotonEvaporation != ptr) {
71 if(theChannels && 0 < theChannels->size()) { (*theChannels)[0] = ptr; }
72 }
73}

References G4VEvaporation::theChannels, and G4VEvaporation::thePhotonEvaporation.

Referenced by G4Evaporation::G4Evaporation(), and G4ExcitationHandler::SetPhotonEvaporation().

◆ SetProduceSecondaries()

void G4WilsonAblationModel::SetProduceSecondaries ( G4bool  produceSecondaries1)
inline

Definition at line 116 of file G4WilsonAblationModel.hh.

118 {produceSecondaries = produceSecondaries1;}

References produceSecondaries.

◆ SetVerboseLevel()

void G4WilsonAblationModel::SetVerboseLevel ( G4int  verboseLevel1)
inline

◆ UseSICB()

void G4VEvaporation::UseSICB ( G4bool  use)
inlineinherited

Definition at line 128 of file G4VEvaporation.hh.

129{
130 useSICB = use;
131}

References G4VEvaporation::useSICB.

Field Documentation

◆ B

G4double G4WilsonAblationModel::B
private

Definition at line 106 of file G4WilsonAblationModel.hh.

Referenced by BreakItUp(), and G4WilsonAblationModel().

◆ evapType

VectorOfFragmentTypes G4WilsonAblationModel::evapType
private

◆ fragmentVector

G4FragmentVector* G4WilsonAblationModel::fragmentVector
private

◆ fragType

G4ParticleDefinition* G4WilsonAblationModel::fragType[6]
private

Definition at line 108 of file G4WilsonAblationModel.hh.

Referenced by BreakItUp(), and G4WilsonAblationModel().

◆ fSig

G4double G4WilsonAblationModel::fSig[200]
private

Definition at line 111 of file G4WilsonAblationModel.hh.

Referenced by BreakItUp(), and G4WilsonAblationModel().

◆ nFragTypes

G4int G4WilsonAblationModel::nFragTypes
private

Definition at line 107 of file G4WilsonAblationModel.hh.

Referenced by BreakItUp(), and G4WilsonAblationModel().

◆ OPTxs

G4int G4VEvaporation::OPTxs
protectedinherited

◆ produceSecondaries

G4bool G4WilsonAblationModel::produceSecondaries
private

◆ secID

G4int G4WilsonAblationModel::secID
private

◆ theChannelFactory

G4VEvaporationFactory* G4VEvaporation::theChannelFactory
protectedinherited

◆ theChannels

std::vector<G4VEvaporationChannel*>* G4VEvaporation::theChannels
protectedinherited

◆ theFBU

G4VFermiBreakUp* G4VEvaporation::theFBU
protectedinherited

◆ thePhotonEvaporation

G4VEvaporationChannel* G4VEvaporation::thePhotonEvaporation
protectedinherited

◆ useSICB

G4bool G4VEvaporation::useSICB
protectedinherited

◆ verboseLevel

G4int G4WilsonAblationModel::verboseLevel
private

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