G4PiMinusAbsorptionAtRest Class Reference

#include <G4PiMinusAbsorptionAtRest.hh>

Inheritance diagram for G4PiMinusAbsorptionAtRest:

G4VRestProcess G4VProcess

Public Member Functions

 G4PiMinusAbsorptionAtRest (const G4String &processName="PiMinusAbsorptionAtRest", G4ProcessType aType=fHadronic)
 ~G4PiMinusAbsorptionAtRest ()
G4bool IsApplicable (const G4ParticleDefinition &particle)
void PreparePhysicsTable (const G4ParticleDefinition &)
void BuildPhysicsTable (const G4ParticleDefinition &)
G4VParticleChangeAtRestDoIt (const G4Track &aTrack, const G4Step &aStep)
void SetDeexcitationAlgorithm (G4int index)

Protected Member Functions

G4double GetMeanLifeTime (const G4Track &aTrack, G4ForceCondition *)

Detailed Description

Definition at line 54 of file G4PiMinusAbsorptionAtRest.hh.


Constructor & Destructor Documentation

G4PiMinusAbsorptionAtRest::G4PiMinusAbsorptionAtRest ( const G4String processName = "PiMinusAbsorptionAtRest",
G4ProcessType  aType = fHadronic 
)

Definition at line 65 of file G4PiMinusAbsorptionAtRest.cc.

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

00066                                                                           :
00067   G4VRestProcess (processName, aType)
00068 {
00069   G4HadronicDeprecate("G4PiMinusAbsorptionAtRest");
00070 
00071   SetProcessSubType(fHadronAtRest);
00072 
00073   _indexDeexcitation = 0;
00074 
00075   if (verboseLevel>0) 
00076     { G4cout << GetProcessName() << " is created "<< G4endl; }
00077   G4HadronicProcessStore::Instance()->RegisterExtraProcess(this);
00078 }

G4PiMinusAbsorptionAtRest::~G4PiMinusAbsorptionAtRest (  ) 

Definition at line 83 of file G4PiMinusAbsorptionAtRest.cc.

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


Member Function Documentation

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

Reimplemented from G4VRestProcess.

Definition at line 98 of file G4PiMinusAbsorptionAtRest.cc.

References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, G4PiMinusStopAbsorption::DoAbsorption(), G4StopDeexcitation::DoBreakUp(), G4PiMinusStopAbsorption::Energy(), FatalException, fStopAndKill, G4cerr, G4cout, G4endl, G4Exception(), G4UniformRand, G4DynamicParticle::GetDefinition(), G4Track::GetDynamicParticle(), G4Material::GetElementVector(), G4Material::GetFractionVector(), G4ParticleChange::GetGlobalTime(), G4Track::GetMaterial(), G4Material::GetNumberOfElements(), G4DynamicParticle::GetTotalEnergy(), G4ParticleChange::Initialize(), IsApplicable(), G4PiMinusStopAbsorption::NNeutrons(), G4PiMinusStopAbsorption::NProtons(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4VParticleChange::ProposeTrackStatus(), G4PiMinusStopAbsorption::RecoilMomentum(), G4VParticleChange::SetNumberOfSecondaries(), G4PiMinusStopAbsorption::SetVerboseLevel(), and G4VProcess::verboseLevel.

00099 {
00100   const G4DynamicParticle* stoppedHadron = track.GetDynamicParticle();
00101       
00102   // Check applicability
00103   if (! IsApplicable(*(stoppedHadron->GetDefinition())))
00104     {
00105       G4cerr  
00106       << "G4PiMinusAbsorptionAtRest: ERROR, particle must be a pion minus!" 
00107       << G4endl;
00108       return NULL;
00109     }
00110 
00111   // Get the current material
00112   const G4Material* material = track.GetMaterial();
00113 
00114   G4double A=-1;
00115   G4double Z=-1;
00116   G4double random = G4UniformRand();
00117   const G4ElementVector* theElementVector = material->GetElementVector();
00118   unsigned int i;
00119   G4double sum = 0;
00120   G4double totalsum=0;
00121   for(i=0; i<material->GetNumberOfElements(); ++i)
00122     {
00123       if((*theElementVector)[i]->GetZ()!=1) totalsum+=material->GetFractionVector()[i];
00124     }
00125   for (i = 0; i<material->GetNumberOfElements(); ++i)
00126     {
00127       if((*theElementVector)[i]->GetZ()!=1) sum += material->GetFractionVector()[i];
00128       if ( sum/totalsum > random )  
00129         { 
00130           A = (*theElementVector)[i]->GetA()*mole/g;
00131           Z = (*theElementVector)[i]->GetZ();
00132           break;
00133         }
00134     }
00135 
00136   // Do the interaction with the nucleon cluster
00137   
00138   G4PiMinusStopMaterial* algorithm = LoadAlgorithm(static_cast<G4int>(Z));
00139   G4PiMinusStopAbsorption stopAbsorption(algorithm,Z,A);
00140   stopAbsorption.SetVerboseLevel(verboseLevel);
00141 
00142   G4DynamicParticleVector* absorptionProducts = stopAbsorption.DoAbsorption();
00143 
00144   // Deal with the leftover nucleus
00145 
00146   G4double pionEnergy = stoppedHadron->GetTotalEnergy();
00147   G4double excitation = pionEnergy - stopAbsorption.Energy();
00148   if (excitation < 0.) 
00149   {
00150     G4Exception("G4PiMinusAbsorptionAtRest::AtRestDoIt()", "HAD_STOP_0000",
00151                 FatalException, "Excitation energy < 0");
00152   }
00153   if (verboseLevel>0) { G4cout << " excitation " << excitation << G4endl; }
00154 
00155   G4StopDeexcitationAlgorithm* nucleusAlgorithm = LoadNucleusAlgorithm();
00156   G4StopDeexcitation stopDeexcitation(nucleusAlgorithm);
00157 
00158   G4double newZ = Z - stopAbsorption.NProtons();
00159   G4double newN = A - Z - stopAbsorption.NNeutrons();
00160   G4double newA = newZ + newN;
00161   G4ReactionProductVector* fragmentationProducts = stopDeexcitation.DoBreakUp(newA,newZ,excitation,stopAbsorption.RecoilMomentum());
00162 
00163   unsigned int nAbsorptionProducts = 0;
00164   if (absorptionProducts != 0)     
00165     { nAbsorptionProducts  =  absorptionProducts->size(); }
00166 
00167   unsigned int nFragmentationProducts = 0;
00168   if (fragmentationProducts != 0) 
00169     { nFragmentationProducts = fragmentationProducts->size(); }
00170   
00171   if (verboseLevel>0) 
00172     {
00173       G4cout << "nAbsorptionProducts = " << nAbsorptionProducts
00174              << "  nFragmentationProducts = " << nFragmentationProducts
00175              << G4endl;
00176     }
00177 
00178   // Deal with ParticleChange final state
00179 
00180   aParticleChange.Initialize(track);
00181   aParticleChange.SetNumberOfSecondaries(G4int(nAbsorptionProducts + nFragmentationProducts)); 
00182      
00183   for (i = 0; i<nAbsorptionProducts; i++)
00184     { aParticleChange.AddSecondary((*absorptionProducts)[i]); }
00185 
00186 //  for (i = 0; i<nFragmentationProducts; i++)
00187 //    { aParticleChange.AddSecondary(fragmentationProducts->at(i)); }
00188   for(i=0; i<nFragmentationProducts; i++)
00189   {
00190     G4DynamicParticle * aNew = 
00191        new G4DynamicParticle((*fragmentationProducts)[i]->GetDefinition(),
00192                              (*fragmentationProducts)[i]->GetMomentum());
00193     G4double newTime = aParticleChange.GetGlobalTime((*fragmentationProducts)[i]->GetFormationTime());
00194     aParticleChange.AddSecondary(aNew, newTime);
00195     delete (*fragmentationProducts)[i];
00196   }
00197 
00198   if (fragmentationProducts != 0)   delete fragmentationProducts;   
00199 
00200   if (_indexDeexcitation == 1) aParticleChange.ProposeLocalEnergyDeposit(excitation);
00201 
00202   // Kill the absorbed pion
00203   aParticleChange.ProposeTrackStatus(fStopAndKill); 
00204 
00205   return &aParticleChange;
00206 
00207 }

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

Reimplemented from G4VProcess.

Definition at line 93 of file G4PiMinusAbsorptionAtRest.cc.

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

00094 {
00095   G4HadronicProcessStore::Instance()->PrintInfo(&p);
00096 }

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

Implements G4VRestProcess.

Definition at line 88 of file G4PiMinusAbsorptionAtRest.hh.

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

00090   {
00091      G4double result = 0;
00092      if(aTrack.GetMaterial()->GetNumberOfElements() == 1)
00093         if(aTrack.GetMaterial()->GetZ()<1.5) result = DBL_MAX;
00094      return result;
00095    }

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

Reimplemented from G4VProcess.

Definition at line 74 of file G4PiMinusAbsorptionAtRest.hh.

References G4PionMinus::PionMinus().

Referenced by AtRestDoIt().

00075     { return ( particle == *(G4PionMinus::PionMinus()) ); }

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

Reimplemented from G4VProcess.

Definition at line 88 of file G4PiMinusAbsorptionAtRest.cc.

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

void G4PiMinusAbsorptionAtRest::SetDeexcitationAlgorithm ( G4int  index  ) 

Definition at line 287 of file G4PiMinusAbsorptionAtRest.cc.

00288 {  
00289   _indexDeexcitation = index;
00290 }


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