Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Types | Public Member Functions
G4WilsonAblationModel Class Reference

#include <G4WilsonAblationModel.hh>

Inheritance diagram for G4WilsonAblationModel:
G4VEvaporation

Public Types

typedef std::vector
< G4ParticleDefinition * > 
VectorOfFragmentTypes
 

Public Member Functions

 G4WilsonAblationModel ()
 
 ~G4WilsonAblationModel ()
 
G4FragmentVectorBreakItUp (const G4Fragment &theNucleus)
 
void SetProduceSecondaries (G4bool)
 
G4bool GetProduceSecondaries ()
 
void SetVerboseLevel (G4int)
 
G4int GetVerboseLevel ()
 
- Public Member Functions inherited from G4VEvaporation
 G4VEvaporation ()
 
virtual ~G4VEvaporation ()
 
virtual void Initialise ()
 
virtual void SetPhotonEvaporation (G4VEvaporationChannel *ptr)
 
G4VEvaporationChannelGetPhotonEvaporation ()
 
void SetOPTxs (G4int opt)
 
void UseSICB (G4bool use)
 

Additional Inherited Members

- Protected Attributes inherited from G4VEvaporation
G4VEvaporationChannelthePhotonEvaporation
 
G4int OPTxs
 
G4bool useSICB
 

Detailed Description

Definition at line 83 of file G4WilsonAblationModel.hh.

Member Typedef Documentation

Definition at line 89 of file G4WilsonAblationModel.hh.

Constructor & Destructor Documentation

G4WilsonAblationModel::G4WilsonAblationModel ( )

Definition at line 105 of file G4WilsonAblationModel.cc.

References G4Alpha::Alpha(), G4Deuteron::Deuteron(), G4VEvaporationFactory::GetChannel(), G4He3::He3(), python.hepunit::MeV, G4Neutron::Neutron(), G4VEvaporation::OPTxs, G4Proton::Proton(), G4Triton::Triton(), and G4VEvaporation::useSICB.

106 {
107 //
108 //
109 // Send message to stdout to advise that the G4Abrasion model is being used.
110 //
111  PrintWelcomeMessage();
112 //
113 //
114 // Set the default verbose level to 0 - no output.
115 //
116  verboseLevel = 0;
117 //
118 //
119 // Set the binding energy per nucleon .... did I mention that this is a crude
120 // model for nuclear de-excitation?
121 //
122  B = 10.0 * MeV;
123 //
124 //
125 // It is possuble to switch off secondary particle production (other than the
126 // final nuclear fragment). The default is on.
127 //
128  produceSecondaries = true;
129 //
130 //
131 // Now we need to define the decay modes. We're using the G4Evaporation model
132 // to help determine the kinematics of the decay.
133 //
134  nFragTypes = 6;
135  fragType[0] = G4Alpha::Alpha();
136  fragType[1] = G4He3::He3();
137  fragType[2] = G4Triton::Triton();
138  fragType[3] = G4Deuteron::Deuteron();
139  fragType[4] = G4Proton::Proton();
140  fragType[5] = G4Neutron::Neutron();
141 //
142 //
143 // Set verboseLevel default to no output.
144 //
145  verboseLevel = 0;
146  theChannelFactory = new G4EvaporationFactory();
147  theChannels = theChannelFactory->GetChannel();
148 //
149 //
150 // Set defaults for evaporation classes. These can be overridden by user
151 // "set" methods.
152 //
153  OPTxs = 3;
154  useSICB = false;
155  fragmentVector = 0;
156 }
virtual std::vector< G4VEvaporationChannel * > * GetChannel()=0
static G4Triton * Triton()
Definition: G4Triton.cc:95
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Neutron * Neutron()
Definition: G4Neutron.cc:104
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4He3 * He3()
Definition: G4He3.cc:94
G4WilsonAblationModel::~G4WilsonAblationModel ( )

Definition at line 159 of file G4WilsonAblationModel.cc.

160 {
161  if (theChannels != 0) theChannels = 0;
162  if (theChannelFactory != 0) delete theChannelFactory;
163 }

Member Function Documentation

G4FragmentVector * G4WilsonAblationModel::BreakItUp ( const G4Fragment theNucleus)
virtual

Implements G4VEvaporation.

Definition at line 167 of file G4WilsonAblationModel.cc.

References CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::e(), python.hepunit::eV, CLHEP::HepLorentzVector::findBoostToCM(), G4cout, G4endl, G4UniformRand, G4Fragment::GetA_asInt(), G4ParticleDefinition::GetBaryonNumber(), G4Fragment::GetExcitationEnergy(), G4Fragment::GetGroundStateMass(), G4ParticleTable::GetIonTable(), G4Fragment::GetMomentum(), G4ParticleDefinition::GetParticleName(), G4ParticleTable::GetParticleTable(), G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4Fragment::GetZ_asInt(), iz, python.hepunit::MeV, n, G4Neutron::Neutron(), G4Proton::Proton(), and CLHEP::HepLorentzVector::setE().

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

Definition at line 122 of file G4WilsonAblationModel.hh.

123  {return produceSecondaries;}
G4int G4WilsonAblationModel::GetVerboseLevel ( )
inline

Definition at line 130 of file G4WilsonAblationModel.hh.

131  {return verboseLevel;}
void G4WilsonAblationModel::SetProduceSecondaries ( G4bool  produceSecondaries1)
inline

Definition at line 118 of file G4WilsonAblationModel.hh.

119  {produceSecondaries = produceSecondaries1;}
void G4WilsonAblationModel::SetVerboseLevel ( G4int  verboseLevel1)
inline

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