LBE.icc

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 //
00027 // --------------------------------------------------------------
00028 //
00029 //      For information related to this code contact: Alex Howard
00030 //      e-mail: alexander.howard@cern.ch
00031 // --------------------------------------------------------------
00032 // Comments
00033 //
00034 //                  Underground Advanced
00035 //
00036 // This physics list is taken from the underground_physics example with small
00037 // modifications.  It is an example of a "flat" physics list with no dependence
00038 // on builders.  The physics covered would be suitable for a low background
00039 // experiment including the neutron_hp package
00040 //
00041 //
00042 //
00043 // PhysicsList program
00044 //
00045 // Modified:
00046 //
00047 // 14-02-03 Fix bugs in msc and hIon instanciation + cut per region
00048 // 16-08-10 Remove inclusion of obsolete class of G4ParticleWithCuts 
00049 // 20-10-10 Migrate LowEnergy process to Livermore models, LP
00050 // --------------------------------------------------------------
00051 
00052 #include <iomanip>  
00053 
00054 #include "globals.hh"
00055 #include "G4ios.hh"
00056 #include "G4ProcessManager.hh"
00057 #include "G4ProcessVector.hh"
00058 
00059 #include "G4ParticleTypes.hh"
00060 #include "G4ParticleTable.hh"
00061 #include "G4ProductionCutsTable.hh"
00062 
00063 #include "G4UserLimits.hh"
00064 #include "G4DataQuestionaire.hh"
00065 #include "G4WarnPLStatus.hh"
00066 
00067 
00068 // Constructor /////////////////////////////////////////////////////////////
00069 template<class T> TLBE<T>::TLBE() : G4VUserPhysicsList() 
00070 {
00071 
00072   G4DataQuestionaire it(photon, lowenergy, neutron, radioactive);
00073   G4cout << "You are using the simulation engine: LBE 5.3"<<G4endl;
00074   G4cout <<G4endl<<G4endl;
00075   this->defaultCutValue     = 1.0*CLHEP::micrometer; //
00076   cutForGamma         = this->defaultCutValue;
00077   cutForElectron      = 1.0*CLHEP::nanometer;
00078   cutForPositron      = this->defaultCutValue;
00079   cutForProton        = this->defaultCutValue;
00080   cutForAlpha         = 1.0*CLHEP::nanometer;
00081   cutForGenericIon    = 1.0*CLHEP::nanometer;
00082 
00083   VerboseLevel = 1;
00084   OpVerbLevel = 0;
00085 
00086   this->SetVerboseLevel(VerboseLevel);
00087 }
00088 
00089 
00090 // Destructor //////////////////////////////////////////////////////////////
00091 template<class T> TLBE<T>::~TLBE() 
00092 {;}
00093 
00094 
00095 // Construct Particles /////////////////////////////////////////////////////
00096  template<class T> void TLBE<T>::ConstructParticle() 
00097 {
00098 
00099   // In this method, static member functions should be called
00100   // for all particles which you want to use.
00101   // This ensures that objects of these particle types will be
00102   // created in the program. 
00103 
00104   ConstructMyBosons();
00105   ConstructMyLeptons();
00106   ConstructMyMesons();
00107   ConstructMyBaryons();
00108   ConstructMyIons();
00109   ConstructMyShortLiveds();
00110 
00111 }
00112 
00113 
00114 // construct Bosons://///////////////////////////////////////////////////
00115  template<class T> void TLBE<T>::ConstructMyBosons()
00116 {
00117   // pseudo-particles
00118   G4Geantino::GeantinoDefinition();
00119   G4ChargedGeantino::ChargedGeantinoDefinition();
00120   
00121   // gamma
00122   G4Gamma::GammaDefinition();
00123 
00124   //OpticalPhotons
00125   G4OpticalPhoton::OpticalPhotonDefinition();
00126 
00127 }
00128 
00129 
00130 // construct Leptons://///////////////////////////////////////////////////
00131  template<class T> void TLBE<T>::ConstructMyLeptons()
00132 {
00133   // leptons
00134   G4Electron::ElectronDefinition();
00135   G4Positron::PositronDefinition();
00136   G4MuonPlus::MuonPlusDefinition();
00137   G4MuonMinus::MuonMinusDefinition();
00138 
00139   G4NeutrinoE::NeutrinoEDefinition();
00140   G4AntiNeutrinoE::AntiNeutrinoEDefinition();
00141   G4NeutrinoMu::NeutrinoMuDefinition();
00142   G4AntiNeutrinoMu::AntiNeutrinoMuDefinition();
00143 }
00144 
00145 #include "G4MesonConstructor.hh"
00146 #include "G4BaryonConstructor.hh"
00147 #include "G4IonConstructor.hh"
00148 
00149 
00150 // construct Mesons://///////////////////////////////////////////////////
00151  template<class T> void TLBE<T>::ConstructMyMesons()
00152 {
00153  //  mesons
00154   G4MesonConstructor mConstructor;
00155   mConstructor.ConstructParticle();
00156 
00157 }
00158 
00159 
00160 // construct Baryons://///////////////////////////////////////////////////
00161  template<class T> void TLBE<T>::ConstructMyBaryons()
00162 {
00163  //  baryons
00164   G4BaryonConstructor bConstructor;
00165   bConstructor.ConstructParticle();
00166 
00167 }
00168 
00169 
00170 // construct Ions://///////////////////////////////////////////////////
00171  template<class T> void TLBE<T>::ConstructMyIons()
00172 {
00173  //  ions
00174   G4IonConstructor iConstructor;
00175   iConstructor.ConstructParticle();
00176 
00177 }
00178 
00179 // construct Shortliveds://///////////////////////////////////////////////////
00180  template<class T> void TLBE<T>::ConstructMyShortLiveds()
00181 {
00182   // ShortLiveds
00183   ;
00184 }
00185 
00186 
00187 
00188 
00189 // Construct Processes //////////////////////////////////////////////////////
00190  template<class T> void TLBE<T>::ConstructProcess() 
00191 {
00192 
00193   AddTransportation();
00194 
00195   ConstructEM();
00196 
00197   ConstructOp();
00198 
00199   ConstructHad();
00200 
00201   ConstructGeneral();
00202 
00203 }
00204 
00205 
00206 // Transportation ///////////////////////////////////////////////////////////
00207 #include "MaxTimeCuts.hh"
00208 #include "MinEkineCuts.hh"
00209 
00210  template<class T> void TLBE<T>::AddTransportation() {
00211   
00212   G4VUserPhysicsList::AddTransportation();
00213   
00214   this->theParticleIterator->reset();
00215   while( (*(this->theParticleIterator))() ){
00216     G4ParticleDefinition* particle = this->theParticleIterator->value();
00217     G4ProcessManager* pmanager = particle->GetProcessManager();
00218     G4String particleName = particle->GetParticleName();
00219     // time cuts for ONLY neutrons:
00220     if(particleName == "neutron") 
00221     pmanager->AddDiscreteProcess(new MaxTimeCuts());
00222     // Energy cuts to kill charged (embedded in method) particles:
00223     pmanager->AddDiscreteProcess(new MinEkineCuts());
00224   }                   
00225 }
00226 
00227 
00228 // Electromagnetic Processes ////////////////////////////////////////////////
00229 // all charged particles
00230 
00231 #include "G4eMultipleScattering.hh"
00232 #include "G4MuMultipleScattering.hh"
00233 #include "G4hMultipleScattering.hh"
00234 
00235 // gamma. Use Livermore models
00236 #include "G4PhotoElectricEffect.hh"
00237 #include "G4LivermorePhotoElectricModel.hh"
00238 #include "G4ComptonScattering.hh"
00239 #include "G4LivermoreComptonModel.hh"
00240 #include "G4GammaConversion.hh"
00241 #include "G4LivermoreGammaConversionModel.hh"
00242 #include "G4RayleighScattering.hh" 
00243 #include "G4LivermoreRayleighModel.hh"
00244 
00245 
00246 // e-
00247 #include "G4eMultipleScattering.hh"
00248 #include "G4UniversalFluctuation.hh"
00249 #include "G4UrbanMscModel93.hh"
00250 
00251 #include "G4eIonisation.hh"
00252 #include "G4LivermoreIonisationModel.hh"
00253 
00254 #include "G4eBremsstrahlung.hh"
00255 #include "G4LivermoreBremsstrahlungModel.hh"
00256 
00257 // e+
00258 #include "G4eplusAnnihilation.hh"
00259 
00260 
00261 // alpha and GenericIon and deuterons, triton, He3:
00262 #include "G4ionIonisation.hh"
00263 #include "G4hIonisation.hh"
00264 #include "G4hBremsstrahlung.hh"
00265 //
00266 #include "G4IonParametrisedLossModel.hh"
00267 #include "G4NuclearStopping.hh"
00268 #include "G4EnergyLossTables.hh"
00269 
00270 
00271 // msc models
00272 #include "G4UrbanMscModel93.hh"
00273 
00274 
00275 //muon:
00276 #include "G4MuIonisation.hh"
00277 #include "G4MuBremsstrahlung.hh"
00278 #include "G4MuPairProduction.hh"
00279 #include "G4MuonMinusCaptureAtRest.hh"
00280 
00281 //OTHERS:
00282 //#include "G4hIonisation.hh" // standard hadron ionisation
00283 
00284  template<class T> void TLBE<T>::ConstructEM() {
00285 
00286  // models & processes:
00287  // Use Livermore models up to 20 MeV, and standard 
00288  // models for higher energy
00289  G4double LivermoreHighEnergyLimit = 20*CLHEP::MeV;
00290  //
00291   this->theParticleIterator->reset();
00292   while( (*(this->theParticleIterator))() ){
00293     G4ParticleDefinition* particle = this->theParticleIterator->value();
00294     G4ProcessManager* pmanager = particle->GetProcessManager();
00295     G4String particleName = particle->GetParticleName();
00296     G4String particleType = particle->GetParticleType();
00297     G4double charge = particle->GetPDGCharge();
00298     
00299     if (particleName == "gamma") 
00300       {
00301       G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
00302       G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel = 
00303         new G4LivermorePhotoElectricModel();
00304       theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00305       thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
00306       pmanager->AddDiscreteProcess(thePhotoElectricEffect);
00307 
00308       G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
00309       G4LivermoreComptonModel* theLivermoreComptonModel = 
00310         new G4LivermoreComptonModel();
00311       theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00312       theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
00313       pmanager->AddDiscreteProcess(theComptonScattering);
00314 
00315       G4GammaConversion* theGammaConversion = new G4GammaConversion();
00316       G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel = 
00317         new G4LivermoreGammaConversionModel();
00318       theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00319       theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
00320       pmanager->AddDiscreteProcess(theGammaConversion);
00321 
00322       G4RayleighScattering* theRayleigh = new G4RayleighScattering();
00323       G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
00324       theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00325       theRayleigh->AddEmModel(0, theRayleighModel);
00326       pmanager->AddDiscreteProcess(theRayleigh);
00327 
00328       } 
00329     else if (particleName == "e-") 
00330       {
00331        //electron
00332        // process ordering: AddProcess(name, at rest, along step, post step)
00333        // -1 = not implemented, then ordering
00334         G4eMultipleScattering* msc = new G4eMultipleScattering();     
00335         msc->AddEmModel(0, new G4UrbanMscModel93());
00336         msc->SetStepLimitType(fUseDistanceToBoundary);
00337         pmanager->AddProcess(msc,                   -1, 1, 1);
00338       
00339        // Ionisation
00340        G4eIonisation* eIoni = new G4eIonisation();
00341        G4LivermoreIonisationModel* theIoniLivermore = new
00342         G4LivermoreIonisationModel();
00343        theIoniLivermore->SetHighEnergyLimit(1*CLHEP::MeV); 
00344        eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
00345        eIoni->SetStepFunction(0.2, 100*CLHEP::um); //     
00346        pmanager->AddProcess(eIoni,                 -1, 2, 2);
00347       
00348        // Bremsstrahlung
00349        G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
00350        G4LivermoreBremsstrahlungModel* theBremLivermore = new
00351          G4LivermoreBremsstrahlungModel();
00352        theBremLivermore->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00353        eBrem->AddEmModel(0, theBremLivermore);
00354        pmanager->AddProcess(eBrem, -1,-3, 3);   
00355       } 
00356     else if (particleName == "e+") 
00357       {
00358         //positron
00359       G4eMultipleScattering* msc = new G4eMultipleScattering();
00360       msc->AddEmModel(0, new G4UrbanMscModel93());      
00361       msc->SetStepLimitType(fUseDistanceToBoundary);
00362       pmanager->AddProcess(msc,                   -1, 1, 1);
00363       G4eIonisation* eIoni = new G4eIonisation();
00364       eIoni->SetStepFunction(0.2, 100*CLHEP::um);      
00365       pmanager->AddProcess(eIoni,                 -1, 2, 2);
00366       pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3);      
00367       pmanager->AddProcess(new G4eplusAnnihilation,0,-1, 4);
00368       } 
00369     else if( particleName == "mu+" || 
00370              particleName == "mu-"    ) 
00371       {
00372         //muon  
00373         G4MuMultipleScattering* aMultipleScattering = new G4MuMultipleScattering();
00374         pmanager->AddProcess(aMultipleScattering,           -1, 1, 1);
00375         pmanager->AddProcess(new G4MuIonisation(),          -1, 2, 2);
00376         pmanager->AddProcess(new G4MuBremsstrahlung(),      -1,-1, 3);
00377         pmanager->AddProcess(new G4MuPairProduction(),      -1,-1, 4);
00378         if( particleName == "mu-" )
00379           pmanager->AddProcess(new G4MuonMinusCaptureAtRest(), 0,-1,-1);
00380       } 
00381     else if (particleName == "GenericIon")
00382     {
00383       pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
00384       G4ionIonisation* ionIoni = new G4ionIonisation();
00385       ionIoni->SetEmModel(new G4IonParametrisedLossModel());
00386       ionIoni->SetStepFunction(0.1, 10*CLHEP::um);
00387       pmanager->AddProcess(ionIoni,                   -1, 2, 2);
00388       pmanager->AddProcess(new G4NuclearStopping(),   -1, 3,-1);        
00389     }
00390     else if (particleName == "alpha" || particleName == "He3")
00391     {
00392       //MSC, ion-Ionisation, Nuclear Stopping
00393       pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
00394 
00395       G4ionIonisation* ionIoni = new G4ionIonisation();
00396       ionIoni->SetStepFunction(0.1, 20*CLHEP::um);
00397       pmanager->AddProcess(ionIoni,                   -1, 2, 2);
00398       pmanager->AddProcess(new G4NuclearStopping(),   -1, 3,-1);
00399     }
00400     else if (particleName == "proton"     ||      
00401              particleName == "deuteron"   ||
00402              particleName == "triton"     ||
00403              particleName == "pi+" ||
00404              particleName == "pi-" ||
00405              particleName == "kaon+" ||
00406              particleName == "kaon-") 
00407       {
00408        //MSC, h-ionisation, bremsstrahlung
00409        pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);      
00410        G4hIonisation* hIoni = new G4hIonisation();
00411        hIoni->SetStepFunction(0.2, 50*CLHEP::um);
00412        pmanager->AddProcess(hIoni,                     -1, 2, 2);      
00413        pmanager->AddProcess(new G4hBremsstrahlung,     -1,-3, 3);    
00414       } 
00415     else if ((!particle->IsShortLived()) &&
00416              (charge != 0.0) && 
00417              (particle->GetParticleName() != "chargedgeantino")) 
00418       {
00419         //all others charged particles except geantino
00420         pmanager->AddProcess(new G4hMultipleScattering, -1, 1, 1);
00421         pmanager->AddProcess(new G4hIonisation,         -1, 2, 2);
00422       }
00423     
00424   }
00425 }
00426 
00427 
00428 // Optical Processes ////////////////////////////////////////////////////////
00429 #include "G4Scintillation.hh"
00430 #include "G4OpAbsorption.hh"
00431 //#include "G4OpRayleigh.hh"
00432 #include "G4OpBoundaryProcess.hh"
00433 
00434  template<class T> void TLBE<T>::ConstructOp() 
00435 {
00436   // default scintillation process
00437   //Coverity report: check that the process is actually used, if not must delete
00438   G4bool theScintProcessDefNeverUsed = true;
00439   G4Scintillation* theScintProcessDef = new G4Scintillation("Scintillation");
00440   // theScintProcessDef->DumpPhysicsTable();
00441   theScintProcessDef->SetTrackSecondariesFirst(true);
00442   theScintProcessDef->SetScintillationYieldFactor(1.0); //
00443   theScintProcessDef->SetScintillationExcitationRatio(0.0); //
00444   theScintProcessDef->SetVerboseLevel(OpVerbLevel);
00445 
00446   // scintillation process for alpha:
00447   G4bool theScintProcessAlphaNeverUsed = true;
00448   G4Scintillation* theScintProcessAlpha = new G4Scintillation("Scintillation");
00449   // theScintProcessNuc->DumpPhysicsTable();
00450   theScintProcessAlpha->SetTrackSecondariesFirst(true);
00451   theScintProcessAlpha->SetScintillationYieldFactor(1.1);
00452   theScintProcessAlpha->SetScintillationExcitationRatio(1.0);
00453   theScintProcessAlpha->SetVerboseLevel(OpVerbLevel);
00454 
00455   // scintillation process for heavy nuclei
00456   G4bool theScintProcessNucNeverUsed = true;  
00457   G4Scintillation* theScintProcessNuc = new G4Scintillation("Scintillation");
00458   // theScintProcessNuc->DumpPhysicsTable();
00459   theScintProcessNuc->SetTrackSecondariesFirst(true);
00460   theScintProcessNuc->SetScintillationYieldFactor(0.2);
00461   theScintProcessNuc->SetScintillationExcitationRatio(1.0);
00462   theScintProcessNuc->SetVerboseLevel(OpVerbLevel);
00463 
00464   // optical processes
00465   G4bool theAbsorptionProcessNeverUsed = true;
00466   G4OpAbsorption* theAbsorptionProcess = new G4OpAbsorption();
00467   //  G4OpRayleigh* theRayleighScatteringProcess = new G4OpRayleigh();
00468   G4bool theBoundaryProcessNeverUsed = true;
00469   G4OpBoundaryProcess* theBoundaryProcess = new G4OpBoundaryProcess();
00470   //  theAbsorptionProcess->DumpPhysicsTable();
00471   //  theRayleighScatteringProcess->DumpPhysicsTable();
00472   theAbsorptionProcess->SetVerboseLevel(OpVerbLevel);
00473   // theRayleighScatteringProcess->SetVerboseLevel(OpVerbLevel);
00474   theBoundaryProcess->SetVerboseLevel(OpVerbLevel);
00475 
00476   this->theParticleIterator->reset();
00477   while( (*(this->theParticleIterator))() )
00478     {
00479       G4ParticleDefinition* particle = this->theParticleIterator->value();
00480       G4ProcessManager* pmanager = particle->GetProcessManager();
00481       G4String particleName = particle->GetParticleName();
00482       if (theScintProcessDef->IsApplicable(*particle)) {
00483         //      if(particle->GetPDGMass() > 5.0*CLHEP::GeV) 
00484         if(particle->GetParticleName() == "GenericIon") {
00485           pmanager->AddProcess(theScintProcessNuc); // AtRestDiscrete
00486           pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxAtRest);
00487           pmanager->SetProcessOrderingToLast(theScintProcessNuc,idxPostStep);
00488           theScintProcessNucNeverUsed = false;
00489         }         
00490         else if(particle->GetParticleName() == "alpha") {
00491           pmanager->AddProcess(theScintProcessAlpha);
00492           pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxAtRest);
00493           pmanager->SetProcessOrderingToLast(theScintProcessAlpha,idxPostStep);
00494           theScintProcessAlphaNeverUsed = false;
00495         }
00496         else {
00497           pmanager->AddProcess(theScintProcessDef);
00498           pmanager->SetProcessOrderingToLast(theScintProcessDef,idxAtRest);
00499           pmanager->SetProcessOrderingToLast(theScintProcessDef,idxPostStep);
00500           theScintProcessDefNeverUsed = false;
00501         }         
00502       }
00503       
00504       if (particleName == "opticalphoton") {
00505         pmanager->AddDiscreteProcess(theAbsorptionProcess);
00506         theAbsorptionProcessNeverUsed = false;
00507         //      pmanager->AddDiscreteProcess(theRayleighScatteringProcess);
00508         theBoundaryProcessNeverUsed = false;
00509         pmanager->AddDiscreteProcess(theBoundaryProcess);
00510       }
00511     }
00512     if ( theScintProcessDefNeverUsed ) delete theScintProcessDef;
00513     if ( theScintProcessAlphaNeverUsed ) delete theScintProcessAlpha;
00514     if ( theScintProcessNucNeverUsed ) delete theScintProcessNuc;
00515     if ( theBoundaryProcessNeverUsed ) delete theBoundaryProcess;
00516     if ( theAbsorptionProcessNeverUsed ) delete theAbsorptionProcess;
00517     if ( theBoundaryProcessNeverUsed ) delete theBoundaryProcess;
00518 }
00519 
00520 
00521 // Hadronic processes ////////////////////////////////////////////////////////
00522 
00523 // Elastic processes:
00524 #include "G4HadronElasticProcess.hh"
00525 
00526 // Inelastic processes:
00527 #include "G4PionPlusInelasticProcess.hh"
00528 #include "G4PionMinusInelasticProcess.hh"
00529 #include "G4KaonPlusInelasticProcess.hh"
00530 #include "G4KaonZeroSInelasticProcess.hh"
00531 #include "G4KaonZeroLInelasticProcess.hh"
00532 #include "G4KaonMinusInelasticProcess.hh"
00533 #include "G4ProtonInelasticProcess.hh"
00534 #include "G4AntiProtonInelasticProcess.hh"
00535 #include "G4NeutronInelasticProcess.hh"
00536 #include "G4AntiNeutronInelasticProcess.hh"
00537 #include "G4DeuteronInelasticProcess.hh"
00538 #include "G4TritonInelasticProcess.hh"
00539 #include "G4AlphaInelasticProcess.hh"
00540 
00541 // Low-energy Models: < 20GeV
00542 #include "G4LElastic.hh"
00543 #include "G4LEPionPlusInelastic.hh"
00544 #include "G4LEPionMinusInelastic.hh"
00545 #include "G4LEKaonPlusInelastic.hh"
00546 #include "G4LEKaonZeroSInelastic.hh"
00547 #include "G4LEKaonZeroLInelastic.hh"
00548 #include "G4LEKaonMinusInelastic.hh"
00549 #include "G4LEProtonInelastic.hh"
00550 #include "G4LEAntiProtonInelastic.hh"
00551 #include "G4LENeutronInelastic.hh"
00552 #include "G4LEAntiNeutronInelastic.hh"
00553 #include "G4LEDeuteronInelastic.hh"
00554 #include "G4LETritonInelastic.hh"
00555 #include "G4LEAlphaInelastic.hh"
00556 
00557 // High-energy Models: >20 GeV
00558 #include "G4HEPionPlusInelastic.hh"
00559 #include "G4HEPionMinusInelastic.hh"
00560 #include "G4HEKaonPlusInelastic.hh"
00561 #include "G4HEKaonZeroInelastic.hh"
00562 #include "G4HEKaonZeroInelastic.hh"
00563 #include "G4HEKaonMinusInelastic.hh"
00564 #include "G4HEProtonInelastic.hh"
00565 #include "G4HEAntiProtonInelastic.hh"
00566 #include "G4HENeutronInelastic.hh"
00567 #include "G4HEAntiNeutronInelastic.hh"
00568 
00569 // Neutron high-precision models: <20 MeV
00570 #include "G4NeutronHPElastic.hh"
00571 #include "G4NeutronHPElasticData.hh"
00572 #include "G4NeutronHPCapture.hh"
00573 #include "G4NeutronHPCaptureData.hh"
00574 #include "G4NeutronHPInelastic.hh"
00575 #include "G4NeutronHPInelasticData.hh"
00576 #include "G4LCapture.hh"
00577 
00578 // Stopping processes
00579 #include "G4PiMinusAbsorptionAtRest.hh"
00580 #include "G4KaonMinusAbsorptionAtRest.hh"
00581 #include "G4AntiProtonAnnihilationAtRest.hh"
00582 #include "G4AntiNeutronAnnihilationAtRest.hh"
00583 
00584 
00585 // ConstructHad()
00586 // Makes discrete physics processes for the hadrons, at present limited
00587 // to those particles with GHEISHA interactions (INTRC > 0).
00588 // The processes are: Elastic scattering and Inelastic scattering.
00589 // F.W.Jones  09-JUL-1998
00590  template<class T> void TLBE<T>::ConstructHad() 
00591 {
00592   G4HadronElasticProcess* theElasticProcess = new G4HadronElasticProcess;
00593   G4LElastic* theElasticModel = new G4LElastic;
00594   theElasticProcess->RegisterMe(theElasticModel);
00595   
00596   this->theParticleIterator->reset();
00597   while ((*(this->theParticleIterator))()) 
00598     {
00599       G4ParticleDefinition* particle = this->theParticleIterator->value();
00600       G4ProcessManager* pmanager = particle->GetProcessManager();
00601       G4String particleName = particle->GetParticleName();
00602 
00603       if (particleName == "pi+") 
00604         {
00605           pmanager->AddDiscreteProcess(theElasticProcess);
00606           G4PionPlusInelasticProcess* theInelasticProcess = 
00607             new G4PionPlusInelasticProcess("inelastic");
00608           G4LEPionPlusInelastic* theLEInelasticModel = 
00609             new G4LEPionPlusInelastic;
00610           theInelasticProcess->RegisterMe(theLEInelasticModel);
00611           G4HEPionPlusInelastic* theHEInelasticModel = 
00612             new G4HEPionPlusInelastic;
00613           theInelasticProcess->RegisterMe(theHEInelasticModel);
00614           pmanager->AddDiscreteProcess(theInelasticProcess);
00615         } 
00616 
00617       else if (particleName == "pi-") 
00618         {
00619           pmanager->AddDiscreteProcess(theElasticProcess);
00620           G4PionMinusInelasticProcess* theInelasticProcess = 
00621             new G4PionMinusInelasticProcess("inelastic");
00622           G4LEPionMinusInelastic* theLEInelasticModel = 
00623             new G4LEPionMinusInelastic;
00624           theInelasticProcess->RegisterMe(theLEInelasticModel);
00625           G4HEPionMinusInelastic* theHEInelasticModel = 
00626             new G4HEPionMinusInelastic;
00627           theInelasticProcess->RegisterMe(theHEInelasticModel);
00628           pmanager->AddDiscreteProcess(theInelasticProcess);
00629           G4String prcNam;
00630           pmanager->AddRestProcess(new G4PiMinusAbsorptionAtRest, ordDefault);
00631         }
00632       
00633       else if (particleName == "kaon+") 
00634         {
00635           pmanager->AddDiscreteProcess(theElasticProcess);
00636           G4KaonPlusInelasticProcess* theInelasticProcess = 
00637             new G4KaonPlusInelasticProcess("inelastic");
00638           G4LEKaonPlusInelastic* theLEInelasticModel = 
00639             new G4LEKaonPlusInelastic;
00640           theInelasticProcess->RegisterMe(theLEInelasticModel);
00641           G4HEKaonPlusInelastic* theHEInelasticModel = 
00642             new G4HEKaonPlusInelastic;
00643           theInelasticProcess->RegisterMe(theHEInelasticModel);
00644           pmanager->AddDiscreteProcess(theInelasticProcess);
00645         }
00646       
00647       else if (particleName == "kaon0S") 
00648         {
00649           pmanager->AddDiscreteProcess(theElasticProcess);
00650           G4KaonZeroSInelasticProcess* theInelasticProcess = 
00651             new G4KaonZeroSInelasticProcess("inelastic");
00652           G4LEKaonZeroSInelastic* theLEInelasticModel = 
00653             new G4LEKaonZeroSInelastic;
00654           theInelasticProcess->RegisterMe(theLEInelasticModel);
00655           G4HEKaonZeroInelastic* theHEInelasticModel = 
00656             new G4HEKaonZeroInelastic;
00657           theInelasticProcess->RegisterMe(theHEInelasticModel);
00658           pmanager->AddDiscreteProcess(theInelasticProcess);
00659         }
00660 
00661       else if (particleName == "kaon0L") 
00662         {
00663           pmanager->AddDiscreteProcess(theElasticProcess);
00664           G4KaonZeroLInelasticProcess* theInelasticProcess = 
00665             new G4KaonZeroLInelasticProcess("inelastic");
00666           G4LEKaonZeroLInelastic* theLEInelasticModel = 
00667             new G4LEKaonZeroLInelastic;
00668           theInelasticProcess->RegisterMe(theLEInelasticModel);
00669           G4HEKaonZeroInelastic* theHEInelasticModel = 
00670             new G4HEKaonZeroInelastic;
00671           theInelasticProcess->RegisterMe(theHEInelasticModel);
00672           pmanager->AddDiscreteProcess(theInelasticProcess);
00673         }
00674 
00675       else if (particleName == "kaon-") 
00676         {
00677           pmanager->AddDiscreteProcess(theElasticProcess);
00678           G4KaonMinusInelasticProcess* theInelasticProcess = 
00679             new G4KaonMinusInelasticProcess("inelastic");
00680           G4LEKaonMinusInelastic* theLEInelasticModel = 
00681             new G4LEKaonMinusInelastic;
00682           theInelasticProcess->RegisterMe(theLEInelasticModel);
00683           G4HEKaonMinusInelastic* theHEInelasticModel = 
00684             new G4HEKaonMinusInelastic;
00685           theInelasticProcess->RegisterMe(theHEInelasticModel);
00686           pmanager->AddDiscreteProcess(theInelasticProcess);
00687           pmanager->AddRestProcess(new G4KaonMinusAbsorptionAtRest, ordDefault);
00688         }
00689 
00690       else if (particleName == "proton") 
00691         {
00692           pmanager->AddDiscreteProcess(theElasticProcess);
00693           G4ProtonInelasticProcess* theInelasticProcess = 
00694             new G4ProtonInelasticProcess("inelastic");
00695           G4LEProtonInelastic* theLEInelasticModel = new G4LEProtonInelastic;
00696           theInelasticProcess->RegisterMe(theLEInelasticModel);
00697           G4HEProtonInelastic* theHEInelasticModel = new G4HEProtonInelastic;
00698           theInelasticProcess->RegisterMe(theHEInelasticModel);
00699           pmanager->AddDiscreteProcess(theInelasticProcess);
00700         }
00701 
00702       else if (particleName == "anti_proton") 
00703         {
00704           pmanager->AddDiscreteProcess(theElasticProcess);
00705           G4AntiProtonInelasticProcess* theInelasticProcess = 
00706             new G4AntiProtonInelasticProcess("inelastic");
00707           G4LEAntiProtonInelastic* theLEInelasticModel = 
00708             new G4LEAntiProtonInelastic;
00709           theInelasticProcess->RegisterMe(theLEInelasticModel);
00710           G4HEAntiProtonInelastic* theHEInelasticModel = 
00711             new G4HEAntiProtonInelastic;
00712           theInelasticProcess->RegisterMe(theHEInelasticModel);
00713           pmanager->AddDiscreteProcess(theInelasticProcess);
00714         }
00715 
00716       else if (particleName == "neutron") {
00717         // elastic scattering
00718         G4HadronElasticProcess* theNeutronElasticProcess = 
00719           new G4HadronElasticProcess;
00720         G4LElastic* theElasticModel1 = new G4LElastic;
00721         G4NeutronHPElastic * theElasticNeutron = new G4NeutronHPElastic;
00722         theNeutronElasticProcess->RegisterMe(theElasticModel1);
00723         theElasticModel1->SetMinEnergy(19*CLHEP::MeV);
00724         theNeutronElasticProcess->RegisterMe(theElasticNeutron);
00725         G4NeutronHPElasticData * theNeutronData = new G4NeutronHPElasticData;
00726         theNeutronElasticProcess->AddDataSet(theNeutronData);
00727         pmanager->AddDiscreteProcess(theNeutronElasticProcess);
00728         // inelastic scattering
00729         G4NeutronInelasticProcess* theInelasticProcess =
00730           new G4NeutronInelasticProcess("inelastic");
00731         G4LENeutronInelastic* theInelasticModel = new G4LENeutronInelastic;
00732         theInelasticModel->SetMinEnergy(19*CLHEP::MeV);
00733         theInelasticProcess->RegisterMe(theInelasticModel);
00734         G4NeutronHPInelastic * theLENeutronInelasticModel =
00735           new G4NeutronHPInelastic;
00736         theInelasticProcess->RegisterMe(theLENeutronInelasticModel);
00737         G4NeutronHPInelasticData * theNeutronData1 = 
00738           new G4NeutronHPInelasticData;
00739         theInelasticProcess->AddDataSet(theNeutronData1);
00740         pmanager->AddDiscreteProcess(theInelasticProcess);
00741         // capture
00742         G4HadronCaptureProcess* theCaptureProcess =
00743           new G4HadronCaptureProcess;
00744         G4LCapture* theCaptureModel = new G4LCapture;
00745         theCaptureModel->SetMinEnergy(19*CLHEP::MeV);
00746         theCaptureProcess->RegisterMe(theCaptureModel);
00747         G4NeutronHPCapture * theLENeutronCaptureModel = new G4NeutronHPCapture;
00748         theCaptureProcess->RegisterMe(theLENeutronCaptureModel);
00749         G4NeutronHPCaptureData * theNeutronData3 = new G4NeutronHPCaptureData;
00750         theCaptureProcess->AddDataSet(theNeutronData3);
00751         pmanager->AddDiscreteProcess(theCaptureProcess);
00752         //  G4ProcessManager* pmanager = G4Neutron::Neutron->GetProcessManager();
00753         //  pmanager->AddProcess(new G4UserSpecialCuts(),-1,-1,1);
00754       }
00755       else if (particleName == "anti_neutron") 
00756         {
00757           pmanager->AddDiscreteProcess(theElasticProcess);
00758           G4AntiNeutronInelasticProcess* theInelasticProcess = 
00759             new G4AntiNeutronInelasticProcess("inelastic");
00760           G4LEAntiNeutronInelastic* theLEInelasticModel = 
00761             new G4LEAntiNeutronInelastic;
00762           theInelasticProcess->RegisterMe(theLEInelasticModel);
00763           G4HEAntiNeutronInelastic* theHEInelasticModel = 
00764             new G4HEAntiNeutronInelastic;
00765           theInelasticProcess->RegisterMe(theHEInelasticModel);
00766           pmanager->AddDiscreteProcess(theInelasticProcess);
00767         }
00768 
00769       else if (particleName == "deuteron") 
00770         {
00771           pmanager->AddDiscreteProcess(theElasticProcess);
00772           G4DeuteronInelasticProcess* theInelasticProcess = 
00773             new G4DeuteronInelasticProcess("inelastic");
00774           G4LEDeuteronInelastic* theLEInelasticModel = 
00775             new G4LEDeuteronInelastic;
00776           theInelasticProcess->RegisterMe(theLEInelasticModel);
00777           pmanager->AddDiscreteProcess(theInelasticProcess);
00778         }
00779       
00780       else if (particleName == "triton") 
00781         {
00782           pmanager->AddDiscreteProcess(theElasticProcess);
00783           G4TritonInelasticProcess* theInelasticProcess = 
00784             new G4TritonInelasticProcess("inelastic");
00785           G4LETritonInelastic* theLEInelasticModel = 
00786             new G4LETritonInelastic;
00787           theInelasticProcess->RegisterMe(theLEInelasticModel);
00788           pmanager->AddDiscreteProcess(theInelasticProcess);
00789         }
00790 
00791       else if (particleName == "alpha") 
00792         {
00793           pmanager->AddDiscreteProcess(theElasticProcess);
00794           G4AlphaInelasticProcess* theInelasticProcess = 
00795             new G4AlphaInelasticProcess("inelastic");
00796           G4LEAlphaInelastic* theLEInelasticModel = 
00797             new G4LEAlphaInelastic;
00798           theInelasticProcess->RegisterMe(theLEInelasticModel);
00799           pmanager->AddDiscreteProcess(theInelasticProcess);
00800         }
00801 
00802     }
00803 }
00804 
00805 
00806 // Decays ///////////////////////////////////////////////////////////////////
00807 #include "G4Decay.hh"
00808 #include "G4RadioactiveDecay.hh"
00809 #include "G4IonTable.hh"
00810 #include "G4Ions.hh"
00811 
00812  template<class T> void TLBE<T>::ConstructGeneral() {
00813 
00814   // Add Decay Process
00815   G4Decay* theDecayProcess = new G4Decay();
00816   G4bool theDecayProcessNeverUsed = true; //Check if theDecayProcess will be used
00817   this->theParticleIterator->reset();
00818   while( (*(this->theParticleIterator))() )
00819     {
00820       G4ParticleDefinition* particle = this->theParticleIterator->value();
00821       G4ProcessManager* pmanager = particle->GetProcessManager();
00822       
00823       if (theDecayProcess->IsApplicable(*particle) && !particle->IsShortLived()) 
00824         { 
00825           theDecayProcessNeverUsed = false;
00826           pmanager ->AddProcess(theDecayProcess);
00827           // set ordering for PostStepDoIt and AtRestDoIt
00828           pmanager ->SetProcessOrdering(theDecayProcess, idxPostStep);
00829           pmanager ->SetProcessOrdering(theDecayProcess, idxAtRest);
00830         }
00831     }
00832 
00833   // Declare radioactive decay to the GenericIon in the IonTable.
00834   const G4IonTable *theIonTable = 
00835     G4ParticleTable::GetParticleTable()->GetIonTable();
00836   G4RadioactiveDecay *theRadioactiveDecay = new G4RadioactiveDecay();
00837 
00838   for (G4int i=0; i<theIonTable->Entries(); i++) 
00839     {
00840       G4String particleName = theIonTable->GetParticle(i)->GetParticleName();
00841       G4String particleType = theIonTable->GetParticle(i)->GetParticleType();
00842       
00843       if (particleName == "GenericIon") 
00844         {
00845           G4ProcessManager* pmanager = 
00846             theIonTable->GetParticle(i)->GetProcessManager();
00847           pmanager->SetVerboseLevel(VerboseLevel);
00848           pmanager ->AddProcess(theRadioactiveDecay);
00849           pmanager ->SetProcessOrdering(theRadioactiveDecay, idxPostStep);
00850           pmanager ->SetProcessOrdering(theRadioactiveDecay, idxAtRest);
00851         } 
00852     }
00853     //If we actually never used the process, delete it
00854     //From Coverity report
00855     if ( theDecayProcessNeverUsed ) delete theDecayProcess;
00856 }
00857 
00858 // Cuts /////////////////////////////////////////////////////////////////////
00859 template<class T> void TLBE<T>::SetCuts() 
00860 {
00861   
00862   if (this->verboseLevel >1)
00863     G4cout << "LBE::SetCuts:";
00864   
00865   if (this->verboseLevel>0){
00866     G4cout << "LBE::SetCuts:";
00867     G4cout << "CutLength : " 
00868            << G4BestUnit(this->defaultCutValue,"Length") << G4endl;
00869   }
00870 
00871   //special for low energy physics
00872   G4double lowlimit=250*CLHEP::eV;  
00873   G4ProductionCutsTable * aPCTable = G4ProductionCutsTable::GetProductionCutsTable();
00874   aPCTable->SetEnergyRange(lowlimit,100*CLHEP::GeV);
00875    
00876   // set cut values for gamma at first and for e- second and next for e+,
00877   // because some processes for e+/e- need cut values for gamma 
00878   this->SetCutValue(cutForGamma, "gamma");
00879   this->SetCutValue(cutForElectron, "e-");
00880   this->SetCutValue(cutForPositron, "e+");
00881   
00882   //  this->SetCutValue(cutForProton, "proton");
00883   //  this->SetCutValue(cutForProton, "anti_proton");
00884   //  this->SetCutValue(cutForAlpha,  "alpha");
00885   //  this->SetCutValue(cutForGenericIon,  "GenericIon");
00886   
00887   //  this->SetCutValueForOthers(this->defaultCutValue);
00888   
00889   if (this->verboseLevel>0) this->DumpCutValuesTable();
00890 }
00891 

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