G4KaonMinusAbsorptionAtRest Class Reference

#include <G4KaonMinusAbsorptionAtRest.hh>

Inheritance diagram for G4KaonMinusAbsorptionAtRest:

G4VRestProcess G4VProcess

Public Member Functions

 G4KaonMinusAbsorptionAtRest (const G4String &processName="KaonMinusAbsorptionAtRest", G4ProcessType aType=fHadronic)
 ~G4KaonMinusAbsorptionAtRest ()
G4bool IsApplicable (const G4ParticleDefinition &particle)
void PreparePhysicsTable (const G4ParticleDefinition &)
void BuildPhysicsTable (const G4ParticleDefinition &)
G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)

Protected Member Functions

G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *)

Detailed Description

Definition at line 54 of file G4KaonMinusAbsorptionAtRest.hh.


Constructor & Destructor Documentation

G4KaonMinusAbsorptionAtRest::G4KaonMinusAbsorptionAtRest ( const G4String processName = "KaonMinusAbsorptionAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 48 of file G4KaonMinusAbsorptionAtRest.cc.

References fHadronAtRest, G4cout, G4endl, G4HadronicDeprecate, G4VProcess::GetProcessName(), G4HadronicProcessStore::Instance(), G4HadronicProcessStore::RegisterExtraProcess(), G4VProcess::SetProcessSubType(), and G4VProcess::verboseLevel.

00049                                                               :
00050   G4VRestProcess (processName, aType)
00051 {
00052   G4HadronicDeprecate("G4KaonMinusAbsorptionAtRest");
00053   if (verboseLevel>0) {
00054     G4cout << GetProcessName() << " is created "<< G4endl;
00055   }
00056   SetProcessSubType(fHadronAtRest);
00057 
00058   // see Cohn et al, PLB27(1968) 527;
00059   //     Davis et al, PLB1(1967) 434; 
00060   
00061   pionAbsorptionRate = 0.07;
00062   
00063   // see  VanderVelde-Wilquet et al, Nuov.Cim.39A(1978)538;
00064   // see  VanderVelde-Wilquet et al, Nuov.Cim.38A(1977)178;
00065   // see  VanderVelde-Wilquet et al, Nucl.Phys.A241(1975)511;
00066   // primary production rates ( for absorption on Carbon)
00067   // .. other elements are extrapolated by the halo factor.
00068   
00069   rateLambdaZeroPiZero = 0.052;
00070   rateSigmaMinusPiPlus = 0.199;
00071   rateSigmaPlusPiMinus = 0.446;
00072   rateSigmaZeroPiZero  = 0.303;
00073   rateLambdaZeroPiMinus = 0.568;
00074   rateSigmaZeroPiMinus  = 0.216;
00075   rateSigmaMinusPiZero  = 0.216;
00076   
00077   // for sigma- p -> lambda n
00078   //     sigma+ n -> lambda p
00079   //     sigma- n -> lambda 
00080   // all values compatible with 0.55 same literature as above.
00081   
00082   sigmaPlusLambdaConversionRate = 0.55; 
00083   sigmaMinusLambdaConversionRate = 0.55;
00084   sigmaZeroLambdaConversionRate = 0.55;
00085 
00086   G4HadronicProcessStore::Instance()->RegisterExtraProcess(this);
00087 }

G4KaonMinusAbsorptionAtRest::~G4KaonMinusAbsorptionAtRest (  ) 

Definition at line 90 of file G4KaonMinusAbsorptionAtRest.cc.

References G4HadronicProcessStore::DeRegisterExtraProcess(), and G4HadronicProcessStore::Instance().


Member Function Documentation

G4VParticleChange * G4KaonMinusAbsorptionAtRest::AtRestDoIt ( const G4Track aTrack,
const G4Step aStep 
) [virtual]

Reimplemented from G4VRestProcess.

Definition at line 106 of file G4KaonMinusAbsorptionAtRest.cc.

References G4StopDeexcitation::DoBreakUp(), FatalException, fStopAndKill, G4cerr, G4cout, G4endl, G4Exception(), G4NucleiProperties::GetBindingEnergy(), G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4Track::GetMaterial(), G4PionMinus::PionMinus(), G4PionPlus::PionPlus(), G4PionZero::PionZero(), G4SigmaMinus::SigmaMinus(), G4SigmaPlus::SigmaPlus(), and G4SigmaZero::SigmaZero().

00107 {
00108   stoppedHadron = track.GetDynamicParticle();
00109   
00110   // Check applicability
00111 
00112   if (!IsApplicable(*(stoppedHadron->GetDefinition()))) 
00113     {
00114       G4cerr  <<"G4KaonMinusAbsorptionAtRest:ERROR, particle must be a Kaon!" <<G4endl;
00115       return 0;
00116     }
00117   
00118   G4Material* material;
00119   material = track.GetMaterial();
00120   nucleus = 0;
00121   do
00122     {
00123       // Select the nucleus, get nucleon
00124       nucleus = new G4Nucleus(material);
00125       if (nucleus->GetA_asInt() < 1.5)
00126         {
00127           delete nucleus;
00128           nucleus = 0;
00129         }
00130     }  while(nucleus == 0);
00131     
00132   G4double Z = nucleus->GetZ_asInt();
00133   G4double A = nucleus->GetA_asInt();
00134 
00135   // Do the interaction with the nucleon
00136   G4DynamicParticleVector* absorptionProducts = KaonNucleonReaction();
00137 
00138   //A.R. 26-Jul-2012 Coverity fix
00139   if ( ! absorptionProducts ) {
00140     G4Exception("G4KaonMinusAbsorptionAtRest::AtRestDoIt()", "HAD_STOP_0001",
00141                 FatalException, "NULL absorptionProducts");
00142     return 0;
00143   }
00144   
00145   // Secondary interactions
00146   
00147   G4DynamicParticle* thePion;
00148   unsigned int i;
00149   for(i = 0; i < absorptionProducts->size(); i++)
00150     {
00151       thePion = (*absorptionProducts)[i];
00152       if (thePion->GetDefinition() == G4PionMinus::PionMinus()
00153           || thePion->GetDefinition() == G4PionPlus::PionPlus()
00154           || thePion->GetDefinition() == G4PionZero::PionZero()) 
00155         {
00156           if (AbsorbPionByNucleus(thePion))
00157             {
00158               absorptionProducts->erase(absorptionProducts->begin()+i);
00159               i--;
00160               delete thePion;
00161               if (verboseLevel > 1) 
00162                 G4cout << "G4KaonMinusAbsorption::AtRestDoIt: Pion absorbed in Nucleus" 
00163                        << G4endl;
00164             }                 
00165         }
00166     }
00167   
00168   G4DynamicParticle* theSigma;
00169   G4DynamicParticle* theLambda;
00170   for (i = 0; i < absorptionProducts->size(); i++)
00171     {
00172       theSigma = (*absorptionProducts)[i];
00173       if (theSigma->GetDefinition() == G4SigmaMinus::SigmaMinus()
00174           || theSigma->GetDefinition() == G4SigmaPlus::SigmaPlus()
00175           || theSigma->GetDefinition() == G4SigmaZero::SigmaZero()) 
00176         {
00177           theLambda = SigmaLambdaConversion(theSigma);
00178           if (theLambda  != 0){
00179             absorptionProducts->erase(absorptionProducts->begin()+i);
00180             i--;
00181             delete theSigma;
00182             absorptionProducts->push_back(theLambda);
00183 
00184             if (verboseLevel > 1) 
00185               G4cout << "G4KaonMinusAbsorption::AtRestDoIt: SigmaLambdaConversion Done" 
00186                      << G4endl;
00187           }                 
00188         }
00189     }
00190   
00191   // Nucleus deexcitation
00192   
00193   G4double productEnergy = 0.;
00194   G4ThreeVector pProducts(0.,0.,0.);
00195 
00196   unsigned int nAbsorptionProducts = 0;
00197   if (absorptionProducts != 0) nAbsorptionProducts = absorptionProducts->size();
00198   
00199   for ( i = 0; i<nAbsorptionProducts; i++)
00200     {
00201       pProducts += (*absorptionProducts)[i]->GetMomentum();
00202       productEnergy += (*absorptionProducts)[i]->GetKineticEnergy();
00203     }
00204 
00205   G4double newZ = nucleus->GetZ_asInt();
00206   G4double newA = nucleus->GetA_asInt();
00207 
00208   G4double bDiff = G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(A),static_cast<G4int>(Z)) - 
00209     G4NucleiProperties::GetBindingEnergy(static_cast<G4int>(newA), static_cast<G4int>(newZ));
00210   
00211   G4StopDeexcitationAlgorithm* nucleusAlgorithm = new G4StopTheoDeexcitation();
00212   G4StopDeexcitation stopDeexcitation(nucleusAlgorithm);
00213 
00214   nucleus->AddExcitationEnergy(bDiff);
00215    
00216   // returns excitation energy for the moment ..
00217   G4double energyDeposit = nucleus->GetEnergyDeposit(); 
00218   if (verboseLevel>0)  
00219     {
00220       G4cout << " -- KaonAtRest -- excitation = " 
00221              << energyDeposit 
00222              << ", pNucleus = "
00223              << pProducts 
00224              << ", A: "
00225              << A 
00226              << ", "
00227              << newA 
00228              << ", Z: "
00229              << Z
00230              << ", "
00231              << newZ
00232              << G4endl; 
00233     }
00234 
00235   if (energyDeposit < 0.)
00236     G4Exception("G4KaonMinusAbsorptionAtRest::AtRestDoIt()", "HAD_STOP_0002",
00237                 FatalException, "Excitation energy < 0");
00238   delete nucleus;    
00239 
00240   G4ReactionProductVector* fragmentationProducts = stopDeexcitation.DoBreakUp(newA,newZ,energyDeposit,pProducts);
00241   
00242   unsigned int nFragmentationProducts = 0;
00243   if (fragmentationProducts != 0) nFragmentationProducts = fragmentationProducts->size();
00244   
00245   //Initialize ParticleChange
00246    aParticleChange.Initialize(track);
00247   aParticleChange.SetNumberOfSecondaries(G4int(nAbsorptionProducts+nFragmentationProducts) ); 
00248   
00249   // update List of alive particles. put energy deposit at the right place ...
00250   for (i = 0; i < nAbsorptionProducts; i++)
00251     {aParticleChange.AddSecondary((*absorptionProducts)[i]); }
00252   if (absorptionProducts != 0) delete absorptionProducts;
00253   
00254 //  for (i = 0; i < nFragmentationProducts; i++)
00255 //    { aParticleChange.AddSecondary(fragmentationProducts->at(i)); }
00256   for(i=0; i<nFragmentationProducts; i++)
00257   {
00258     G4DynamicParticle * aNew = 
00259        new G4DynamicParticle((*fragmentationProducts)[i]->GetDefinition(),
00260                              (*fragmentationProducts)[i]->GetTotalEnergy(),
00261                              (*fragmentationProducts)[i]->GetMomentum());
00262     G4double newTime = aParticleChange.GetGlobalTime((*fragmentationProducts)[i]->GetFormationTime());
00263     aParticleChange.AddSecondary(aNew, newTime);
00264     delete (*fragmentationProducts)[i];
00265   }
00266   if (fragmentationProducts != 0) delete fragmentationProducts;
00267   
00268   // finally ...
00269   aParticleChange.ProposeTrackStatus(fStopAndKill); // Kill the incident Kaon
00270   return &aParticleChange;
00271 }

void G4KaonMinusAbsorptionAtRest::BuildPhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 100 of file G4KaonMinusAbsorptionAtRest.cc.

References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::PrintInfo().

00101 {
00102   G4HadronicProcessStore::Instance()->PrintInfo(&p);
00103 }

G4double G4KaonMinusAbsorptionAtRest::GetMeanLifeTime ( const G4Track aTrack,
G4ForceCondition  
) [inline, protected, virtual]

Implements G4VRestProcess.

Definition at line 81 of file G4KaonMinusAbsorptionAtRest.hh.

References DBL_MAX, G4Track::GetMaterial(), G4Material::GetNumberOfElements(), and G4Material::GetZ().

00083      {
00084      G4double result = 0;
00085      if(aTrack.GetMaterial()->GetNumberOfElements() == 1)
00086         if(aTrack.GetMaterial()->GetZ()<1.5) result = DBL_MAX;
00087      return result;
00088      }

G4bool G4KaonMinusAbsorptionAtRest::IsApplicable ( const G4ParticleDefinition particle  )  [inline, virtual]

Reimplemented from G4VProcess.

Definition at line 68 of file G4KaonMinusAbsorptionAtRest.hh.

References G4KaonMinus::KaonMinus().

00068                                                             {
00069                return( particle == *(G4KaonMinus::KaonMinus()) );
00070   }

void G4KaonMinusAbsorptionAtRest::PreparePhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VProcess.

Definition at line 95 of file G4KaonMinusAbsorptionAtRest.cc.

References G4HadronicProcessStore::Instance(), and G4HadronicProcessStore::RegisterParticleForExtraProcess().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:20 2013 for Geant4 by  doxygen 1.4.7