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

#include <G4ExcitationHandler.hh>

Public Member Functions

 G4ExcitationHandler ()
 
 ~G4ExcitationHandler ()
 
G4ReactionProductVectorBreakItUp (const G4Fragment &theInitialState)
 
void ModelDescription (std::ostream &outFile) const
 
void Initialise ()
 
void SetEvaporation (G4VEvaporation *ptr, G4bool isLocal=false)
 
void SetMultiFragmentation (G4VMultiFragmentation *ptr)
 
void SetFermiModel (G4VFermiBreakUp *ptr)
 
void SetPhotonEvaporation (G4VEvaporationChannel *ptr)
 
void SetDeexChannelsType (G4DeexChannelType val)
 
void SetMaxZForFermiBreakUp (G4int aZ)
 
void SetMaxAForFermiBreakUp (G4int anA)
 
void SetMaxAandZForFermiBreakUp (G4int anA, G4int aZ)
 
void SetMinEForMultiFrag (G4double anE)
 
G4VEvaporationGetEvaporation ()
 
G4VMultiFragmentationGetMultiFragmentation ()
 
G4VFermiBreakUpGetFermiModel ()
 
G4VEvaporationChannelGetPhotonEvaporation ()
 
void SetOPTxs (G4int opt)
 
void UseSICB ()
 

Detailed Description

Definition at line 63 of file G4ExcitationHandler.hh.

Constructor & Destructor Documentation

◆ G4ExcitationHandler()

G4ExcitationHandler::G4ExcitationHandler ( )
explicit

Definition at line 88 of file G4ExcitationHandler.cc.

89 : icID(0),maxZForFermiBreakUp(9),maxAForFermiBreakUp(17),
90 fVerbose(1),fWarnings(0),minEForMultiFrag(1.*CLHEP::TeV),
91 minExcitation(1.*CLHEP::eV),maxExcitation(100.*CLHEP::MeV),
92 isInitialised(false),isEvapLocal(true),isActive(true)
93{
94 thePartTable = G4ParticleTable::GetParticleTable();
95 theTableOfIons = thePartTable->GetIonTable();
97
98 theMultiFragmentation = nullptr;
99 theFermiModel = nullptr;
100 theEvaporation = nullptr;
101 thePhotonEvaporation = nullptr;
102 theResults.reserve(60);
103 results.reserve(30);
104 theEvapList.reserve(30);
105
107 theElectron = G4Electron::Electron();
108 theNeutron = G4Neutron::NeutronDefinition();
109 theProton = G4Proton::ProtonDefinition();
110 theDeuteron = G4Deuteron::DeuteronDefinition();
111 theTriton = G4Triton::TritonDefinition();
112 theHe3 = G4He3::He3Definition();
113 theAlpha = G4Alpha::AlphaDefinition();
114 theLambda = G4Lambda::Lambda();
115
116 fLambdaMass = theLambda->GetPDGMass();
117
118 if(fVerbose > 1) { G4cout << "### New handler " << this << G4endl; }
119}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Alpha * AlphaDefinition()
Definition: G4Alpha.cc:83
static G4Deuteron * DeuteronDefinition()
Definition: G4Deuteron.cc:88
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4He3 * He3Definition()
Definition: G4He3.cc:88
static G4Lambda * Lambda()
Definition: G4Lambda.cc:107
static G4Neutron * NeutronDefinition()
Definition: G4Neutron.cc:98
static G4NistManager * Instance()
G4IonTable * GetIonTable() const
static G4ParticleTable * GetParticleTable()
static G4Pow * GetInstance()
Definition: G4Pow.cc:41
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87
static G4Triton * TritonDefinition()
Definition: G4Triton.cc:88

◆ ~G4ExcitationHandler()

G4ExcitationHandler::~G4ExcitationHandler ( )

Definition at line 121 of file G4ExcitationHandler.cc.

122{
123 delete theMultiFragmentation;
124 delete theFermiModel;
125 if(isEvapLocal) { delete theEvaporation; }
126}

Member Function Documentation

◆ BreakItUp()

G4ReactionProductVector * G4ExcitationHandler::BreakItUp ( const G4Fragment theInitialState)

Definition at line 285 of file G4ExcitationHandler.cc.

286{
287 // Variables existing until end of method
288 G4Fragment * theInitialStatePtr = new G4Fragment(theInitialState);
289 if(fVerbose > 1) {
290 G4cout << "@@@@@@@@@@ Start G4Excitation Handler @@@@@@@@@@@@@ " << G4endl;
291 G4cout << theInitialState << G4endl;
292 }
293 if(!isInitialised) { Initialise(); }
294
295 // pointer to fragment vector which receives temporal results
296 G4FragmentVector * theTempResult = nullptr;
297
298 theResults.clear();
299 theEvapList.clear();
300
301 // Variables to describe the excited configuration
302 G4double exEnergy = theInitialState.GetExcitationEnergy();
303 G4int A = theInitialState.GetA_asInt();
304 G4int Z = theInitialState.GetZ_asInt();
305 G4int nL = theInitialState.GetNumberOfLambdas();
306
307 // too much excitation
308 if(exEnergy > A*maxExcitation && A > 0) {
309 ++fWarnings;
310 if(fWarnings < 0) {
312 ed << "High excitation Fragment Z= " << Z << " A= " << A
313 << " Eex/A(MeV)= " << exEnergy/A;
314 G4Exception("G4ExcitationHandler::BreakItUp()","had0034",JustWarning,ed,"");
315 }
316 }
317
318 // for hyper-nuclei subtract lambdas from the projectile fragment
319 G4double lambdaF = 0.0;
320 G4LorentzVector lambdaLV = theInitialStatePtr->GetMomentum();
321 if(0 < nL) {
322
323 // is it a stable hyper-nuclei?
324 if(A >= 3 && A <= 5 && nL <= 2) {
325 G4int pdg = 0;
326 if(3 == A && 1 == nL) {
327 pdg = 1010010030;
328 } else if(5 == A && 2 == Z && 1 == nL) {
329 pdg = 1010020050;
330 } else if(4 == A) {
331 if(1 == Z && 1 == nL) {
332 pdg = 1010010040;
333 } else if(2 == Z && 1 == nL) {
334 pdg = 1010020040;
335 } else if(0 == Z && 2 == nL) {
336 pdg = 1020000040;
337 } else if(1 == Z && 2 == nL) {
338 pdg = 1020010040;
339 }
340 }
341 // initial state is one of hyper-nuclei
342 if(0 < pdg) {
343 const G4ParticleDefinition* part = thePartTable->FindParticle(pdg);
344 if(nullptr != part) {
345 G4ReactionProduct* theNew = new G4ReactionProduct(part);
346 G4ThreeVector dir = G4ThreeVector( 0.0, 0.0, 0.0 );
347 if ( lambdaLV.vect().mag() > CLHEP::eV ) {
348 dir = lambdaLV.vect().unit();
349 }
350 G4double mass = part->GetPDGMass();
351 G4double etot = std::max(lambdaLV.e(), mass);
352 dir *= std::sqrt((etot - mass)*(etot + mass));
353 theNew->SetMomentum(dir);
354 theNew->SetTotalEnergy(etot);
355 theNew->SetFormationTime(theInitialState.GetCreationTime());
356 theNew->SetCreatorModelID(theInitialState.GetCreatorModelID());
358 v->push_back(theNew);
359 return v;
360 }
361 }
362 }
363 G4double mass = theInitialStatePtr->GetGroundStateMass();
364 lambdaF = nL*(fLambdaMass - CLHEP::neutron_mass_c2)/mass;
365
366 // de-excitation with neutrons instead of lambda inside the fragment
367 theInitialStatePtr->SetZAandMomentum(lambdaLV*(1. - lambdaF), Z, A, 0);
368
369 // 4-momentum not used in de-excitation
370 lambdaLV *= lambdaF;
371 } else if(0 > nL) {
372 ++fWarnings;
373 if(fWarnings < 0) {
375 ed << "Fragment with negative L: Z=" << Z << " A=" << A << " L=" << nL
376 << " Eex/A(MeV)= " << exEnergy/A;
377 G4Exception("G4ExcitationHandler::BreakItUp()","had0034",JustWarning,ed,"");
378 }
379 }
380
381 // In case A <= 1 the fragment will not perform any nucleon emission
382 if (A <= 1 || !isActive) {
383 theResults.push_back( theInitialStatePtr );
384
385 // check if a fragment is stable
386 } else if(exEnergy < minExcitation &&
387 nist->GetIsotopeAbundance(Z, A) > 0.0) {
388 theResults.push_back( theInitialStatePtr );
389
390 // JMQ 150909: first step in de-excitation is treated separately
391 // Fragments after the first step are stored in theEvapList
392 } else {
393 if((A<maxAForFermiBreakUp && Z<maxZForFermiBreakUp)
394 || exEnergy <= minEForMultiFrag*A) {
395 theEvapList.push_back(theInitialStatePtr);
396
397 // Statistical Multifragmentation will take place only once
398 } else {
399 theTempResult = theMultiFragmentation->BreakItUp(theInitialState);
400 if(!theTempResult) {
401 theEvapList.push_back(theInitialStatePtr);
402 } else {
403 size_t nsec = theTempResult->size();
404
405 // no fragmentation
406 if(0 == nsec) {
407 theEvapList.push_back(theInitialStatePtr);
408
409 // secondary are produced - sort out secondary fragments
410 } else {
411 G4bool deletePrimary = true;
412 for (auto ptr : *theTempResult) {
413 if(ptr == theInitialStatePtr) { deletePrimary = false; }
414 SortSecondaryFragment(ptr);
415 }
416 if( deletePrimary ) { delete theInitialStatePtr; }
417 }
418 delete theTempResult; // end multifragmentation
419 }
420 }
421 }
422 if(fVerbose > 2) {
423 G4cout << "## After first step of handler " << theEvapList.size()
424 << " for evap; "
425 << theResults.size() << " results. " << G4endl;
426 }
427 // -----------------------------------
428 // FermiBreakUp and De-excitation loop
429 // -----------------------------------
430
431 static const G4int countmax = 1000;
432 size_t kk;
433 for (kk=0; kk<theEvapList.size(); ++kk) {
434 G4Fragment* frag = theEvapList[kk];
435 if(fVerbose > 3) {
436 G4cout << "Next evaporate: " << G4endl;
437 G4cout << *frag << G4endl;
438 }
439 if(kk >= countmax) {
441 ed << "Infinite loop in the de-excitation module: " << kk
442 << " iterations \n"
443 << " Initial fragment: \n" << theInitialState
444 << "\n Current fragment: \n" << *frag;
445 G4Exception("G4ExcitationHandler::BreakItUp","had0333",FatalException,
446 ed,"Stop execution");
447
448 }
449 A = frag->GetA_asInt();
450 Z = frag->GetZ_asInt();
451 results.clear();
452 if(fVerbose > 2) {
453 G4cout << "G4ExcitationHandler# " << kk << " Z= " << Z << " A= " << A
454 << " Eex(MeV)= " << frag->GetExcitationEnergy() << G4endl;
455 }
456 // Fermi Break-Up
457 if(theFermiModel->IsApplicable(Z, A, frag->GetExcitationEnergy())) {
458 theFermiModel->BreakFragment(&results, frag);
459 size_t nsec = results.size();
460 if(fVerbose > 2) { G4cout << "FermiBreakUp Nsec= " << nsec << G4endl; }
461
462 // FBU takes care to delete input fragment or add it to the results
463 // The secondary may be excited - photo-evaporation should be applied
464 if(1 < nsec) {
465 for(auto & res : results) {
466 SortSecondaryFragment(res);
467 }
468 continue;
469 }
470 // evaporation will be applied
471 }
472 // apply Evaporation, residual nucleus is always added to the results
473 // photon evaporation is possible
474 theEvaporation->BreakFragment(&results, frag);
475 if(fVerbose > 3) {
476 G4cout << "Evaporation Nsec= " << results.size() << G4endl;
477 }
478 if(0 == results.size()) {
479 theResults.push_back(frag);
480 } else {
481 SortSecondaryFragment(frag);
482 }
483
484 // Sort out secondary fragments
485 for (auto & res : results) {
486 if(fVerbose > 4) {
487 G4cout << "Evaporated product #" << *res << G4endl;
488 }
489 SortSecondaryFragment(res);
490 } // end of loop on secondary
491 } // end of the loop over theEvapList
492 if(fVerbose > 2) {
493 G4cout << "## After 2nd step of handler " << theEvapList.size()
494 << " was evap; "
495 << theResults.size() << " results. " << G4endl;
496 }
497 G4ReactionProductVector * theReactionProductVector =
499
500 // MAC (24/07/08)
501 // To optimise the storing speed, we reserve space
502 // in memory for the vector
503 theReactionProductVector->reserve( theResults.size() );
504
505 if(fVerbose > 2) {
506 G4cout << "### ExcitationHandler provides " << theResults.size()
507 << " evaporated products:" << G4endl;
508 }
509 G4LorentzVector partOfLambdaLV;
510 if ( nL > 0 ) partOfLambdaLV = lambdaLV/(G4double)nL;
511 for (auto & frag : theResults) {
512 G4LorentzVector lv0 = frag->GetMomentum();
513 G4double etot = lv0.e();
514
515 // in the case of dummy de-excitation, excitation energy is transfered
516 // into kinetic energy of output ion
517 if(!isActive) {
518 G4double mass = frag->GetGroundStateMass();
519 G4double ptot = lv0.vect().mag();
520 G4double fac = (etot <= mass || 0.0 == ptot) ? 0.0
521 : std::sqrt((etot - mass)*(etot + mass))/ptot;
522 G4LorentzVector lv((frag->GetMomentum()).px()*fac,
523 (frag->GetMomentum()).py()*fac,
524 (frag->GetMomentum()).pz()*fac, etot);
525 frag->SetMomentum(lv);
526 }
527 if(fVerbose > 3) {
528 G4cout << *frag;
529 if(frag->NuclearPolarization()) {
530 G4cout << " " << frag->NuclearPolarization();
531 }
532 G4cout << G4endl;
533 }
534
535 G4int fragmentA = frag->GetA_asInt();
536 G4int fragmentZ = frag->GetZ_asInt();
537 G4double eexc = 0.0;
538 const G4ParticleDefinition* theKindOfFragment = nullptr;
539 G4bool isHyperN = false;
540 if (fragmentA == 0) { // photon or e-
541 theKindOfFragment = frag->GetParticleDefinition();
542 } else if (fragmentA == 1 && fragmentZ == 0) { // neutron
543 theKindOfFragment = theNeutron;
544 } else if (fragmentA == 1 && fragmentZ == 1) { // proton
545 theKindOfFragment = theProton;
546 } else if (fragmentA == 2 && fragmentZ == 1) { // deuteron
547 theKindOfFragment = theDeuteron;
548 } else if (fragmentA == 3 && fragmentZ == 1) { // triton
549 theKindOfFragment = theTriton;
550 if(0 < nL) {
551 const G4ParticleDefinition* p = thePartTable->FindParticle(1010010030);
552 if(nullptr != p) {
553 theKindOfFragment = p;
554 isHyperN = true;
555 --nL;
556 }
557 }
558 } else if (fragmentA == 3 && fragmentZ == 2) { // helium3
559 theKindOfFragment = theHe3;
560 } else if (fragmentA == 4 && fragmentZ == 2) { // alpha
561 theKindOfFragment = theAlpha;
562 if(0 < nL) {
563 const G4ParticleDefinition* p = thePartTable->FindParticle(1010020040);
564 if(nullptr != p) {
565 theKindOfFragment = p;
566 isHyperN = true;
567 --nL;
568 }
569 }
570 } else {
571
572 // fragment
573 eexc = frag->GetExcitationEnergy();
574 G4int idxf = frag->GetFloatingLevelNumber();
575 if(eexc < minExcitation) {
576 eexc = 0.0;
577 idxf = 0;
578 }
579
580 theKindOfFragment = theTableOfIons->GetIon(fragmentZ, fragmentA, eexc,
582 if(fVerbose > 3) {
583 G4cout << "### EXCH: Find ion Z= " << fragmentZ
584 << " A= " << fragmentA
585 << " Eexc(MeV)= " << eexc/MeV << " idx= " << idxf
586 << G4endl;
587 }
588 }
589 // fragment identified
590 if(nullptr != theKindOfFragment) {
591 G4ReactionProduct * theNew = new G4ReactionProduct(theKindOfFragment);
592 if(isHyperN) {
593 G4LorentzVector lv = lv0 + partOfLambdaLV;
594 G4ThreeVector dir = lv.vect().unit();
595 G4double mass = theKindOfFragment->GetPDGMass();
596 etot = std::max(lv.e(), mass);
597 G4double ptot = std::sqrt((etot - mass)*(etot + mass));
598 dir *= ptot;
599 theNew->SetMomentum(dir);
600 // remaining not compensated 4-momentum
601 lambdaLV += (lv0 - G4LorentzVector(dir, etot));
602 } else {
603 theNew->SetMomentum(lv0.vect());
604 }
605 theNew->SetTotalEnergy(etot);
606 theNew->SetFormationTime(frag->GetCreationTime());
607 if(theKindOfFragment == theElectron) {
608 theNew->SetCreatorModelID(icID);
609 } else {
610 theNew->SetCreatorModelID(frag->GetCreatorModelID());
611 }
612 theReactionProductVector->push_back(theNew);
613
614 // fragment not found out ground state is created
615 } else {
616 theKindOfFragment =
617 theTableOfIons->GetIon(fragmentZ,fragmentA,0.0,noFloat,0);
618 if(theKindOfFragment) {
619 G4ThreeVector mom(0.0,0.0,0.0);
620 G4double ionmass = theKindOfFragment->GetPDGMass();
621 if(etot <= ionmass) {
622 etot = ionmass;
623 } else {
624 G4double ptot = std::sqrt((etot - ionmass)*(etot + ionmass));
625 mom = (frag->GetMomentum().vect().unit())*ptot;
626 }
627 G4ReactionProduct * theNew = new G4ReactionProduct(theKindOfFragment);
628 theNew->SetMomentum(mom);
629 theNew->SetTotalEnergy(etot);
630 theNew->SetFormationTime(frag->GetCreationTime());
631 theNew->SetCreatorModelID(frag->GetCreatorModelID());
632 theReactionProductVector->push_back(theNew);
633 if(fVerbose > 3) {
634 G4cout << " ground state, energy corrected E(MeV)= "
635 << etot << G4endl;
636 }
637 }
638 }
639 delete frag;
640 }
641 // remaining lambdas are free; conserve quantum numbers but
642 // not 4-momentum
643 if(0 < nL) {
644 G4ThreeVector dir = G4ThreeVector( 0.0, 0.0, 0.0 );
645 if ( lambdaLV.vect().mag() > CLHEP::eV ) {
646 dir = lambdaLV.vect().unit();
647 }
648 G4double etot = std::max(lambdaLV.e()/(G4double)nL, fLambdaMass);
649 dir *= std::sqrt((etot - fLambdaMass)*(etot + fLambdaMass));
650 for(G4int i=0; i<nL; ++i) {
651 G4ReactionProduct* theNew = new G4ReactionProduct(theLambda);
652 theNew->SetMomentum(dir);
653 theNew->SetTotalEnergy(etot);
654 theNew->SetFormationTime(theInitialState.GetCreationTime());
655 theNew->SetCreatorModelID(theInitialState.GetCreatorModelID());
656 theReactionProductVector->push_back(theNew);
657 }
658 }
659 if(fVerbose > 3) {
660 G4cout << "@@@@@@@@@@ End G4Excitation Handler "<< G4endl;
661 }
662 return theReactionProductVector;
663}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
std::ostringstream G4ExceptionDescription
Definition: G4Exception.hh:40
std::vector< G4Fragment * > G4FragmentVector
Definition: G4Fragment.hh:65
#define noFloat
Definition: G4Ions.hh:112
CLHEP::HepLorentzVector G4LorentzVector
std::vector< G4ReactionProduct * > G4ReactionProductVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
const G4double A[17]
Hep3Vector unit() const
double mag() const
Hep3Vector vect() const
G4double GetGroundStateMass() const
Definition: G4Fragment.hh:317
G4int GetCreatorModelID() const
Definition: G4Fragment.hh:428
G4double GetExcitationEnergy() const
Definition: G4Fragment.hh:312
const G4LorentzVector & GetMomentum() const
Definition: G4Fragment.hh:322
G4double GetCreationTime() const
Definition: G4Fragment.hh:479
G4int GetZ_asInt() const
Definition: G4Fragment.hh:289
G4int GetNumberOfLambdas() const
Definition: G4Fragment.hh:307
void SetMomentum(const G4LorentzVector &value)
Definition: G4Fragment.hh:327
G4int GetA_asInt() const
Definition: G4Fragment.hh:284
void SetZAandMomentum(const G4LorentzVector &, G4int Z, G4int A, G4int nLambdas=0)
Definition: G4Fragment.hh:334
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4Ions::G4FloatLevelBase FloatLevelBase(char flbChar)
Definition: G4Ions.cc:110
G4double GetIsotopeAbundance(G4int Z, G4int N) const
G4ParticleDefinition * FindParticle(G4int PDGEncoding)
void SetMomentum(const G4double x, const G4double y, const G4double z)
void SetTotalEnergy(const G4double en)
void SetCreatorModelID(const G4int mod)
void SetFormationTime(G4double aTime)
virtual void BreakFragment(G4FragmentVector *, G4Fragment *theNucleus)
virtual G4bool IsApplicable(G4int Z, G4int A, G4double mass) const =0
virtual void BreakFragment(G4FragmentVector *results, G4Fragment *theNucleus)=0
virtual G4FragmentVector * BreakItUp(const G4Fragment &theNucleus)=0

Referenced by G4WilsonAbrasionModel::ApplyYourself(), G4EMDissociation::ApplyYourself(), G4QMDReaction::ApplyYourself(), G4LowEIonFragmentation::ApplyYourself(), and G4PreCompoundDeexcitation::deExcite().

◆ GetEvaporation()

G4VEvaporation * G4ExcitationHandler::GetEvaporation ( )

Definition at line 260 of file G4ExcitationHandler.cc.

261{
262 if(!theEvaporation) { SetParameters(); }
263 return theEvaporation;
264}

Referenced by G4INCLXXInterface::G4INCLXXInterface().

◆ GetFermiModel()

G4VFermiBreakUp * G4ExcitationHandler::GetFermiModel ( )

Definition at line 272 of file G4ExcitationHandler.cc.

273{
274 if(!theFermiModel) { SetParameters(); }
275 return theFermiModel;
276}

◆ GetMultiFragmentation()

G4VMultiFragmentation * G4ExcitationHandler::GetMultiFragmentation ( )

Definition at line 266 of file G4ExcitationHandler.cc.

267{
268 if(!theMultiFragmentation) { SetParameters(); }
269 return theMultiFragmentation;
270}

◆ GetPhotonEvaporation()

G4VEvaporationChannel * G4ExcitationHandler::GetPhotonEvaporation ( )

Definition at line 278 of file G4ExcitationHandler.cc.

279{
280 if(!thePhotonEvaporation) { SetParameters(); }
281 return thePhotonEvaporation;
282}

◆ Initialise()

void G4ExcitationHandler::Initialise ( )

Definition at line 164 of file G4ExcitationHandler.cc.

165{
166 if(isInitialised) { return; }
167 if(fVerbose > 1) {
168 G4cout << "G4ExcitationHandler::Initialise() started " << this << G4endl;
169 }
170 G4DeexPrecoParameters* param =
172 isInitialised = true;
173 SetParameters();
174 if(isActive) {
175 theFermiModel->Initialise();
176 theEvaporation->InitialiseChannels();
177 }
178 // dump level is controlled by parameter class
179 param->Dump();
180}
G4DeexPrecoParameters * GetParameters()
static G4NuclearLevelData * GetInstance()
virtual void InitialiseChannels()
virtual void Initialise()=0

Referenced by BreakItUp(), G4AblaInterface::InitialiseModel(), and G4PreCompoundModel::InitialiseModel().

◆ ModelDescription()

void G4ExcitationHandler::ModelDescription ( std::ostream &  outFile) const

Definition at line 665 of file G4ExcitationHandler.cc.

666{
667 outFile << "G4ExcitationHandler description\n"
668 << "This class samples de-excitation of excited nucleus using\n"
669 << "Fermi Break-up model for light fragments (Z < 9, A < 17), "
670 << "evaporation, fission, and photo-evaporation models. Evaporated\n"
671 << "particle may be proton, neutron, and other light fragment \n"
672 << "(Z < 13, A < 29). During photon evaporation produced gamma \n"
673 << "or electrons due to internal conversion \n";
674}

Referenced by G4BinaryCascade::ModelDescription(), and G4BinaryCascade::PropagateModelDescription().

◆ SetDeexChannelsType()

void G4ExcitationHandler::SetDeexChannelsType ( G4DeexChannelType  val)

Definition at line 228 of file G4ExcitationHandler.cc.

229{
230 G4Evaporation* evap = static_cast<G4Evaporation*>(theEvaporation);
231 if(fVerbose > 1) {
232 G4cout << "G4ExcitationHandler::SetDeexChannelsType " << val
233 << " for " << this << G4endl;
234 }
235 if(val == fDummy) {
236 isActive = false;
237 return;
238 }
239 if(!evap) { return; }
240 if(val == fEvaporation) {
241 evap->SetDefaultChannel();
242 } else if(val == fCombined) {
243 evap->SetCombinedChannel();
244 } else if(val == fGEM) {
245 evap->SetGEMChannel();
246 } else if(val == fGEMVI) {
247 evap->SetGEMVIChannel();
248 }
249 evap->InitialiseChannels();
250 if(fVerbose > 1) {
252 G4cout << "Number of de-excitation channels is changed to: "
253 << theEvaporation->GetNumberOfChannels();
254 G4cout << " " << this;
255 }
256 G4cout << G4endl;
257 }
258}
@ fEvaporation
virtual void InitialiseChannels() final
void SetGEMChannel()
void SetDefaultChannel()
void SetCombinedChannel()
void SetGEMVIChannel()
size_t GetNumberOfChannels() const
G4bool IsMasterThread()
Definition: G4Threading.cc:124

Referenced by G4QMDReaction::G4QMDReaction().

◆ SetEvaporation()

void G4ExcitationHandler::SetEvaporation ( G4VEvaporation ptr,
G4bool  isLocal = false 
)

Definition at line 182 of file G4ExcitationHandler.cc.

183{
184 if(ptr && ptr != theEvaporation) {
185 delete theEvaporation;
186 theEvaporation = ptr;
188 theEvaporation->SetFermiBreakUp(theFermiModel);
189 isEvapLocal = flag;
190 if(fVerbose > 1) {
191 G4cout << "G4ExcitationHandler::SetEvaporation() for " << this << G4endl;
192 }
193 }
194}
void SetPhotonEvaporation(G4VEvaporationChannel *ptr)
G4VEvaporationChannel * GetPhotonEvaporation()
void SetFermiBreakUp(G4VFermiBreakUp *ptr)

Referenced by G4QMDReaction::G4QMDReaction(), G4WilsonAbrasionModel::G4WilsonAbrasionModel(), and G4WilsonAbrasionModel::SetUseAblation().

◆ SetFermiModel()

void G4ExcitationHandler::SetFermiModel ( G4VFermiBreakUp ptr)

Definition at line 205 of file G4ExcitationHandler.cc.

206{
207 if(ptr && ptr != theFermiModel) {
208 delete theFermiModel;
209 theFermiModel = ptr;
210 if(theEvaporation) { theEvaporation->SetFermiBreakUp(theFermiModel); }
211 }
212}

◆ SetMaxAandZForFermiBreakUp()

void G4ExcitationHandler::SetMaxAandZForFermiBreakUp ( G4int  anA,
G4int  aZ 
)
inline

Definition at line 172 of file G4ExcitationHandler.hh.

173{
176}
void SetMaxZForFermiBreakUp(G4int aZ)
void SetMaxAForFermiBreakUp(G4int anA)

◆ SetMaxAForFermiBreakUp()

void G4ExcitationHandler::SetMaxAForFermiBreakUp ( G4int  anA)
inline

Definition at line 167 of file G4ExcitationHandler.hh.

168{
169 maxAForFermiBreakUp = anA;
170}

Referenced by SetMaxAandZForFermiBreakUp().

◆ SetMaxZForFermiBreakUp()

void G4ExcitationHandler::SetMaxZForFermiBreakUp ( G4int  aZ)
inline

Definition at line 162 of file G4ExcitationHandler.hh.

163{
164 maxZForFermiBreakUp = aZ;
165}

Referenced by SetMaxAandZForFermiBreakUp().

◆ SetMinEForMultiFrag()

void G4ExcitationHandler::SetMinEForMultiFrag ( G4double  anE)
inline

Definition at line 178 of file G4ExcitationHandler.hh.

179{
180 minEForMultiFrag = anE;
181}

Referenced by G4EMDissociation::G4EMDissociation().

◆ SetMultiFragmentation()

void G4ExcitationHandler::SetMultiFragmentation ( G4VMultiFragmentation ptr)

Definition at line 197 of file G4ExcitationHandler.cc.

198{
199 if(ptr && ptr != theMultiFragmentation) {
200 delete theMultiFragmentation;
201 theMultiFragmentation = ptr;
202 }
203}

◆ SetOPTxs()

void G4ExcitationHandler::SetOPTxs ( G4int  opt)
inline

Definition at line 183 of file G4ExcitationHandler.hh.

184{}

◆ SetPhotonEvaporation()

void G4ExcitationHandler::SetPhotonEvaporation ( G4VEvaporationChannel ptr)

Definition at line 215 of file G4ExcitationHandler.cc.

216{
217 if(ptr && ptr != thePhotonEvaporation) {
218 delete thePhotonEvaporation;
219 thePhotonEvaporation = ptr;
220 if(theEvaporation) { theEvaporation->SetPhotonEvaporation(ptr); }
221 if(fVerbose > 1) {
222 G4cout << "G4ExcitationHandler::SetPhotonEvaporation() " << ptr
223 << " for handler " << this << G4endl;
224 }
225 }
226}
virtual void SetPhotonEvaporation(G4VEvaporationChannel *ptr)

Referenced by SetEvaporation().

◆ UseSICB()

void G4ExcitationHandler::UseSICB ( )
inline

Definition at line 186 of file G4ExcitationHandler.hh.

187{}

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