G4HadronStoppingProcess.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 //
00028 //---------------------------------------------------------------------
00029 //
00030 // GEANT4 Class 
00031 //
00032 // File name:     G4HadronStoppingProcess
00033 //
00034 // Author V.Ivanchenko 21 April 2012 
00035 //
00036 //
00037 // Class Description:
00038 //
00039 // Base process class for nuclear capture of negatively charged particles
00040 //
00041 // Modifications:
00042 //
00043 //  20120522  M. Kelsey -- Set enableAtRestDoIt flag for G4ProcessManager
00044 //  20120914  M. Kelsey -- Pass subType in base ctor, remove enable flags
00045 //  20121004  K. Genser -- use G4HadronicProcessType in the constructor
00046 //  20121016  K. Genser -- Reverting to use one argument c'tor
00047 //
00048 //------------------------------------------------------------------------
00049 
00050 #include "G4HadronStoppingProcess.hh"
00051 #include "G4HadronicProcessStore.hh"
00052 #include "G4HadronicProcessType.hh"
00053 #include "G4EmCaptureCascade.hh"
00054 #include "G4Nucleus.hh"
00055 #include "G4HadFinalState.hh"
00056 #include "G4HadProjectile.hh"
00057 #include "G4HadSecondary.hh"
00058 #include "G4Material.hh"
00059 #include "G4Element.hh"
00060 
00061 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00062 
00063 G4HadronStoppingProcess::G4HadronStoppingProcess(const G4String& name)
00064   : G4HadronicProcess(name, fHadronAtRest)
00065 {
00066   // Modify G4VProcess flags to emulate G4VRest instead of G4VDiscrete
00067   enableAtRestDoIt = true;
00068   enablePostStepDoIt = false;
00069 
00070   fElementSelector = new G4ElementSelector();
00071   fEmCascade = new G4EmCaptureCascade();        // Owned by InteractionRegistry
00072   fBoundDecay = 0;
00073   G4HadronicProcessStore::Instance()->RegisterExtraProcess(this);
00074 }
00075 
00076 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00077 
00078 G4HadronStoppingProcess::~G4HadronStoppingProcess()
00079 {
00080   G4HadronicProcessStore::Instance()->DeRegisterExtraProcess(this);
00081   delete fElementSelector;
00082   // NOTE: fEmCascade and fEmBoundDecay owned by registry, not locally
00083 }
00084 
00085 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00086 
00087 G4bool G4HadronStoppingProcess::IsApplicable(const G4ParticleDefinition& p)
00088 {
00089   return (p.GetPDGCharge() < 0.0);
00090 }
00091 
00092 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00093 
00094 void G4HadronStoppingProcess::PreparePhysicsTable(const G4ParticleDefinition& p)
00095 {
00096   G4HadronicProcessStore::Instance()->RegisterParticleForExtraProcess(this, &p);
00097 }
00098 
00099 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00100 
00101 void G4HadronStoppingProcess::BuildPhysicsTable(const G4ParticleDefinition& p) 
00102 {
00103   G4HadronicProcessStore::Instance()->PrintInfo(&p);
00104 }
00105 
00106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00107 
00108 G4double G4HadronStoppingProcess::AtRestGetPhysicalInteractionLength(
00109     const G4Track&, G4ForceCondition* condition)
00110 {
00111   *condition = NotForced;
00112   return 0.0;
00113 }
00114 
00115 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00116 
00117 G4double G4HadronStoppingProcess::PostStepGetPhysicalInteractionLength(
00118     const G4Track&, G4double, G4ForceCondition* condition)
00119 {
00120   *condition = NotForced;
00121   return DBL_MAX;
00122 }
00123 
00124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00125 
00126 G4VParticleChange* G4HadronStoppingProcess::AtRestDoIt(const G4Track& track, 
00127                                                        const G4Step&)
00128 {
00129   // if primary is not Alive then do nothing
00130   theTotalResult->Initialize(track);
00131 
00132   G4Nucleus* nucleus = GetTargetNucleusPointer();
00133   G4Element* elm = fElementSelector->SelectZandA(track, nucleus);
00134 
00135   G4HadFinalState* result = 0;
00136   thePro.Initialise(track);
00137   G4double time0 = track.GetGlobalTime();
00138   G4bool nuclearCapture = true;
00139 
00140   // Do the electromagnetic cascade in the nuclear field.
00141   // EM cascade should keep G4HadFinalState object,
00142   // because it will not be deleted at the end of this method
00143   //
00144   result = fEmCascade->ApplyYourself(thePro, *nucleus);
00145   G4double ebound = result->GetLocalEnergyDeposit(); 
00146   G4double edep = 0.0; 
00147   G4int nSecondaries = result->GetNumberOfSecondaries();
00148 
00149   // Try decay from bound level 
00150   // For mu- the time of projectile should be changed.
00151   // Decay should keep G4HadFinalState object,
00152   // because it will not be deleted at the end of this method
00153   //
00154   thePro.SetBoundEnergy(ebound);
00155   if(fBoundDecay) {
00156     G4HadFinalState* resultDecay = 
00157       fBoundDecay->ApplyYourself(thePro, *nucleus);
00158     G4int n = resultDecay->GetNumberOfSecondaries();
00159     if(0 < n) {
00160       nSecondaries += n;
00161       result->AddSecondaries(resultDecay);
00162     } 
00163     if(resultDecay->GetStatusChange() == stopAndKill) {
00164       nuclearCapture = false; 
00165     }
00166     resultDecay->Clear();
00167   }
00168 
00169   if(nuclearCapture) {
00170 
00171     // select model
00172     G4HadronicInteraction* model = 0;
00173     try {
00174       model = ChooseHadronicInteraction(0.0, track.GetMaterial(), elm);
00175     }
00176     catch(G4HadronicException & aE) {
00177       G4ExceptionDescription ed;
00178       ed << "Target element "<<elm->GetName()<<"  Z= " 
00179          << nucleus->GetZ_asInt() << "  A= " 
00180          << nucleus->GetA_asInt() << G4endl;
00181       DumpState(track,"ChooseHadronicInteraction",ed);
00182       ed << " No HadronicInteraction found out" << G4endl;
00183       G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had005", 
00184                   FatalException, ed);
00185     }
00186 
00187     G4HadFinalState* resultNuc = 0;
00188     G4int reentryCount = 0;
00189     do {
00190       // sample final state
00191       // nuclear interaction should keep G4HadFinalState object
00192       // model should define time of each secondary particle
00193       try {
00194         resultNuc = model->ApplyYourself(thePro, *nucleus);
00195         ++reentryCount;
00196       }
00197       catch(G4HadronicException aR) {
00198         G4ExceptionDescription ed;
00199         ed << "Call for " << model->GetModelName() << G4endl;
00200         ed << "Target element "<<elm->GetName()<<"  Z= " 
00201            << nucleus->GetZ_asInt() 
00202            << "  A= " << nucleus->GetA_asInt() << G4endl;
00203         DumpState(track,"ApplyYourself",ed);
00204         ed << " ApplyYourself failed" << G4endl;
00205         G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006", 
00206                     FatalException, ed);
00207       }
00208 
00209       // Check the result for catastrophic energy non-conservation
00210       resultNuc = CheckResult(thePro, *nucleus, resultNuc);
00211 
00212       if(reentryCount>100) {
00213         G4ExceptionDescription ed;
00214         ed << "Call for " << model->GetModelName() << G4endl;
00215         ed << "Target element "<<elm->GetName()<<"  Z= " 
00216            << nucleus->GetZ_asInt() 
00217            << "  A= " << nucleus->GetA_asInt() << G4endl;
00218         DumpState(track,"ApplyYourself",ed);
00219         ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
00220         G4Exception("G4HadronStoppingProcess::AtRestDoIt", "had006", 
00221                     FatalException, ed);  
00222       }
00223     }
00224     while(!resultNuc);
00225 
00226     edep = resultNuc->GetLocalEnergyDeposit();
00227     nSecondaries += resultNuc->GetNumberOfSecondaries();
00228     result->AddSecondaries(resultNuc); 
00229     resultNuc->Clear();
00230   }
00231 
00232   // Fill results
00233   //
00234   theTotalResult->ProposeTrackStatus(fStopAndKill);
00235   theTotalResult->ProposeLocalEnergyDeposit(edep);  
00236   theTotalResult->SetNumberOfSecondaries(nSecondaries);
00237   G4double w  = track.GetWeight();
00238   theTotalResult->ProposeWeight(w);
00239   for(G4int i=0; i<nSecondaries; ++i) {
00240     G4HadSecondary* sec = result->GetSecondary(i);
00241 
00242     // add track global time to the reaction time
00243     G4double time = sec->GetTime();
00244     if(time < 0.0) { time = 0.0; }
00245     time += time0;
00246 
00247     // create secondary track
00248     G4Track* t = new G4Track(sec->GetParticle(),
00249                              time, 
00250                              track.GetPosition());
00251     t->SetWeight(w*sec->GetWeight());
00252     t->SetTouchableHandle(track.GetTouchableHandle());
00253     theTotalResult->AddSecondary(t);
00254   }
00255   result->Clear();
00256 
00257   if (epReportLevel != 0) { 
00258     CheckEnergyMomentumConservation(track, *nucleus);
00259   }
00260   return theTotalResult;
00261 }
00262 
00263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00264 
00265 void G4HadronStoppingProcess::ProcessDescription(std::ostream& outFile) const
00266 {
00267   outFile << "Base process for negatively charged particle capture at rest.\n";
00268 }
00269 
00270 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....

Generated on Mon May 27 17:48:27 2013 for Geant4 by  doxygen 1.4.7