src/BDSPhysicsList.cc

00001       
00002 //   BDSIM, (C) 2001-2007 
00003 //    
00004 //   version 0.4 
00005 //   last modified : 10 Sept 2007 by malton@pp.rhul.ac.uk
00006 //  
00007 
00008 //
00009 //    Physics lists
00010 //
00011 
00012 
00013 #include "BDSExecOptions.hh"
00014 #include "BDSGlobalConstants.hh" 
00015 #include "BDSPhysicsList.hh"
00016 
00017 #include "globals.hh"
00018 #include "G4ParticleDefinition.hh"
00019 #include "G4ParticleWithCuts.hh"
00020 #include "G4ProcessManager.hh"
00021 #include "G4ProcessVector.hh"
00022 #include "G4ParticleTypes.hh"
00023 #include "G4ParticleTable.hh"
00024 #include "G4FastSimulationManagerProcess.hh"
00025 
00026 #include "G4Version.hh"
00027 
00028 #include "G4Material.hh"
00029 #include "G4MaterialTable.hh"
00030 #include "G4ios.hh"
00031 #include <iomanip>   
00032 #include "QGSP_BERT.hh"
00033 #include "QGSP_BERT_HP.hh"
00034 #include "HadronPhysicsQGSP_BERT.hh"
00035 #include "HadronPhysicsQGSP_BERT_HP.hh"
00036 #include "HadronPhysicsFTFP_BERT.hh"
00037 #include "G4Decay.hh"
00038 #include "G4eeToHadrons.hh"
00039 
00040 #include "HadronPhysicsLHEP.hh"
00041 #include "G4EmStandardPhysics.hh"
00042 #include "G4EmLivermorePhysics.hh"
00043 
00044 //#include "IonPhysics.hh"
00045 
00046 
00047 // physics processes
00048 
00049 //
00050 // EM
00051 //
00052 
00053 // gamma
00054 #include "G4ComptonScattering.hh"
00055 #include "G4GammaConversion.hh"
00056 #include "G4GammaConversionToMuons.hh"
00057 #include "G4PhotoElectricEffect.hh"
00058 #include "BDSXSBias.hh" //added by L.D. 16/11/09
00059 #include "BDSGammaConversion_LPB.hh" //added by M.D. Salt, R.B. Appleby, 15/10/09
00060 
00061 // charged particles
00062 #if G4VERSION_NUMBER > 819
00063 #include "G4eMultipleScattering.hh"
00064 #include "G4MuMultipleScattering.hh"
00065 #include "G4hMultipleScattering.hh"
00066 #else
00067 #include "G4MultipleScattering.hh"
00068 #endif
00069 #include "G4Cerenkov.hh"
00070 
00071 #include "G4CoulombScattering.hh"
00072 
00073 // e
00074 #include "G4eIonisation.hh"
00075 #include "G4eBremsstrahlung.hh"
00076 #include "G4eplusAnnihilation.hh"
00077 #include "BDSeBremsstrahlung_LPB.hh" //added by M.D. Salt, R.B.Appleby,15/10/09
00078 
00079 // mu
00080 #include "G4MuIonisation.hh"
00081 #include "G4MuBremsstrahlung.hh"
00082 #include "G4MuPairProduction.hh"
00083 #if G4VERSION_NUMBER < 950
00084 #include "G4MuonNucleusProcess.hh"
00085 #elif G4VERSION_NUMBER < 953
00086 #include "G4MuNuclearInteraction.hh"
00087 #else 
00088 #include "G4MuonVDNuclearModel.hh"
00089 #endif
00090 
00091 //ions
00092 #include "G4hIonisation.hh"
00093 #include "G4ionIonisation.hh"
00094 
00095 
00096 //Decay of unstable particles
00097 #include "G4Decay.hh"
00098 
00099 //
00100 // Low EM
00101 //
00102 
00103 //gamma
00104 #if G4VERSION_NUMBER < 950
00105 #include "G4LowEnergyRayleigh.hh"
00106 #include "G4LowEnergyPhotoElectric.hh"
00107 #include "G4LowEnergyCompton.hh"
00108 #include "G4LowEnergyGammaConversion.hh"
00109 #else
00110 #include "G4RayleighScattering.hh"
00111 #include "G4LivermoreRayleighModel.hh"
00112 #include "G4PhotoElectricEffect.hh"
00113 #include "G4LivermorePhotoElectricModel.hh"
00114 #include "G4ComptonScattering.hh"
00115 #include "G4LivermoreComptonModel.hh"
00116 #include "G4GammaConversion.hh"
00117 #include "G4LivermoreGammaConversionModel.hh"
00118 #endif
00119 
00120 //e
00121 #if G4VERSION_NUMBER < 950
00122 #include "G4LowEnergyIonisation.hh"
00123 #include "G4LowEnergyBremsstrahlung.hh"
00124 #else
00125 #include "G4eIonisation.hh"
00126 #include "G4LivermoreIonisationModel.hh"
00127 #include "G4UniversalFluctuation.hh"
00128 
00129 #include "G4eBremsstrahlung.hh"
00130 #include "G4LivermoreBremsstrahlungModel.hh"
00131 #endif
00132 
00133 #include "G4AnnihiToMuPair.hh"
00134 
00135 //ions
00136 #if G4VERSION_NUMBER < 950
00137 #include "G4hLowEnergyIonisation.hh"
00138 #endif
00139 
00140 #include "BDSLaserCompton.hh"
00141 #include "BDSSynchrotronRadiation.hh"
00142 #include "BDSContinuousSR.hh"
00143 #include "G4StepLimiter.hh"
00144 #include "G4UserSpecialCuts.hh"
00145 
00146 
00147 //
00148 // Hadronic
00149 //
00150 #include "G4TheoFSGenerator.hh"
00151 #include "G4GeneratorPrecompoundInterface.hh"
00152 #include "G4QGSModel.hh"
00153 #include "G4GammaParticipants.hh"
00154 #include "G4QGSMFragmentation.hh"
00155 #include "G4ExcitedStringDecay.hh"
00156 
00157 #include "G4GammaNuclearReaction.hh"
00158 #include "G4ElectroNuclearReaction.hh"
00159 #include "G4PhotoNuclearProcess.hh"
00160 #include "G4ElectronNuclearProcess.hh"
00161 #include "G4PositronNuclearProcess.hh"
00162 
00163 //Planck scattering
00164 #include "BDSPlanckScatterBuilder.hh"
00165 
00166 //
00167 // particle definition
00168 //
00169 
00170 // Bosons
00171 #include "G4ChargedGeantino.hh"
00172 #include "G4Geantino.hh"
00173 #include "G4Gamma.hh"
00174 #include "G4OpticalPhoton.hh"
00175 
00176 // leptons
00177 #include "G4MuonPlus.hh"
00178 #include "G4MuonMinus.hh"
00179 #include "G4NeutrinoMu.hh"
00180 #include "G4AntiNeutrinoMu.hh"
00181 
00182 #include "G4Electron.hh"
00183 #include "G4Positron.hh"
00184 #include "G4NeutrinoE.hh"
00185 #include "G4AntiNeutrinoE.hh"
00186 
00187 // Hadrons
00188 #include "G4MesonConstructor.hh"
00189 #include "G4BaryonConstructor.hh"
00190 #include "G4IonConstructor.hh"
00191 
00192 //ShortLived
00193 #include "G4ShortLivedConstructor.hh"
00194 
00195 
00196 
00197 
00198 BDSPhysicsList::BDSPhysicsList():  G4VUserPhysicsList()
00199 {
00200   verbose = BDSExecOptions::Instance()->GetVerbose();
00201 
00202   // construct particles
00203 
00204   //defaultCutValue = 0.7*mm;  
00205   defaultCutValue = BDSGlobalConstants::Instance()->GetDefaultRangeCut()*m;  
00206   SetDefaultCutValue(BDSGlobalConstants::Instance()->GetDefaultRangeCut()*m);
00207 
00208   G4cout  << __METHOD_NAME__ << "Charged Thresholdcut=" 
00209           << BDSGlobalConstants::Instance()->GetThresholdCutCharged()/GeV<<" GeV"<<G4endl;
00210   G4cout  << __METHOD_NAME__ << "Photon Thresholdcut=" 
00211           << BDSGlobalConstants::Instance()->GetThresholdCutPhotons()/GeV<<" GeV"<<G4endl;
00212   G4cout  << __METHOD_NAME__ << "Default range cut=" 
00213           << BDSGlobalConstants::Instance()->GetDefaultRangeCut()/m<<" m"<<G4endl;
00214 
00215   SetVerboseLevel(1);   
00216 }
00217 
00218 BDSPhysicsList::~BDSPhysicsList()
00219 {
00220 }
00221 
00222 void BDSPhysicsList::ConstructProcess()
00223 { 
00224 #if DEBUG
00225   G4cout << __METHOD_NAME__ << G4endl;
00226 #endif 
00227 
00228   bool plistFound=false;
00229   //standard physics lists
00230   if(BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_HP"){
00231     QGSP_BERT_HP* physList = new QGSP_BERT_HP;
00232     physList->ConstructProcess();
00233     plistFound=true;
00234   } else if(BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT"){
00235     QGSP_BERT* physList = new QGSP_BERT;
00236     physList->ConstructProcess();
00237     plistFound=true;
00238   } else if (BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_HP_muon"){ //Modified standard physics lists
00239     QGSP_BERT_HP* physList = new QGSP_BERT_HP;
00240     physList->ConstructProcess();
00241     ConstructMuon();
00242     plistFound=true;
00243   } else if (BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_muon"){
00244     QGSP_BERT* physList = new QGSP_BERT;
00245     physList->ConstructProcess();
00246     ConstructMuon();
00247     plistFound=true;
00248   } else if(BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_HP_muon_em_low"){
00249     QGSP_BERT_HP* physList = new QGSP_BERT_HP;
00250     physList->ConstructProcess();
00251     ConstructMuon();
00252     ConstructEM_Low_Energy();
00253     plistFound=true;
00254   } else if(BDSGlobalConstants::Instance()->GetPhysListName() == "livermore"){
00255     G4EmLivermorePhysics* physList = new G4EmLivermorePhysics;
00256     physList->ConstructProcess();
00257     plistFound=true;
00258   }
00259 
00260 
00261   if(!plistFound){
00262     //Need to add transportation if non-standard physics list
00263     AddTransportation();
00264   }
00265   
00266   //Apply the following in all cases - step limiter etc.
00267   theParticleIterator->reset();
00268   while( (*theParticleIterator)() ){
00269     G4ParticleDefinition* particle = theParticleIterator->value();
00270     if((particle->GetParticleName()=="gamma")||(particle->GetParticleName()=="e-")||(particle->GetParticleName()=="e+")){
00271       particle->SetApplyCutsFlag(true);
00272     }
00273     G4ProcessManager *pmanager = particle->GetProcessManager();
00274     pmanager->AddProcess(new G4StepLimiter,-1,-1,1);
00275 #ifndef NOUSERSPECIALCUTS
00276     pmanager->AddDiscreteProcess(new G4UserSpecialCuts);
00277 #endif
00278   }
00279   
00280   //===========================================
00281   //Some options
00282   //-------------------------------------------
00283   //Build planck scattering if option is set
00284   if(BDSGlobalConstants::Instance()->GetDoPlanckScattering()){
00285     BDSPlanckScatterBuilder* psbuild = new BDSPlanckScatterBuilder();
00286     psbuild->Build();
00287     // psbuild is leaked. This can't be deleted as its BDSPlanckScatterProcess is registered
00288   }
00289   //A flag to switch on hadronic lead particle biasing
00290   if (BDSGlobalConstants::Instance()->GetUseHadLPB() ){
00291     setenv("SwitchLeadBiasOn","1",1); 
00292   }
00293   //Synchrotron radiation
00294   if(BDSGlobalConstants::Instance()->GetSynchRadOn()) {
00295 #ifdef DEBUG
00296     G4cout << __METHOD_NAME__ << "synch. rad. is turned on" << G4endl;
00297 #endif
00298     ConstructSR();
00299   } else {
00300 #ifdef DEBUG
00301     G4cout << __METHOD_NAME__ << "synch. rad. is turned OFF!" << G4endl;
00302 #endif
00303   }
00304   //Particle decay
00305   if(BDSGlobalConstants::Instance()->GetDecayOn()) ConstructDecay();
00306   //============================================
00307   
00308   if(!plistFound){ //Search BDSIM physics lists
00309     if (BDSGlobalConstants::Instance()->GetPhysListName() != "standard"){ // register physics processes here
00310       
00311       //Always add parameterisation
00312       AddParameterisation();
00313 
00314       // standard e+/e-/gamma electromagnetic interactions
00315       if(BDSGlobalConstants::Instance()->GetPhysListName() == "em_standard") 
00316         {
00317           ConstructEM();
00318           return;
00319         }
00320 
00321       if(BDSGlobalConstants::Instance()->GetPhysListName() == "em_single_scatter") 
00322         {
00323           ConstructEMSingleScatter();
00324           return;
00325         }
00326       
00327       if(BDSGlobalConstants::Instance()->GetPhysListName() == "merlin") 
00328         {
00329           ConstructMerlin();
00330           return;
00331         }
00332       
00333       // low energy em processes
00334       if(BDSGlobalConstants::Instance()->GetPhysListName() == "em_low") 
00335         {
00336           ConstructEM_Low_Energy();
00337           return;
00338         }
00339       
00340       // standard electromagnetic + muon
00341       if(BDSGlobalConstants::Instance()->GetPhysListName() == "em_muon") 
00342         {
00343           ConstructEM();
00344           ConstructMuon();
00345           return;
00346         }
00347       // standard hadronic - photo-nuclear etc.
00348       if(BDSGlobalConstants::Instance()->GetPhysListName() == "hadronic_standard") 
00349         {
00350           ConstructEM();
00351           ConstructHadronic();
00352           return;
00353         }
00354       
00355       // standard electromagnetic + muon + hadronic
00356       if(BDSGlobalConstants::Instance()->GetPhysListName() == "hadronic_muon") 
00357         {
00358           ConstructEM();
00359           ConstructMuon();
00360           ConstructHadronic();
00361           return;
00362         }
00363       
00364       if(BDSGlobalConstants::Instance()->GetPhysListName() == "hadronic_QGSP_BERT") {
00365         ConstructEM();
00366         G4VPhysicsConstructor* hadPhysList = new HadronPhysicsQGSP_BERT("hadron");
00367         hadPhysList -> ConstructProcess();
00368         return;
00369       }
00370       
00371       if(BDSGlobalConstants::Instance()->GetPhysListName() == "hadronic_QGSP_BERT_muon") {
00372         ConstructEM();
00373         ConstructMuon();
00374         G4VPhysicsConstructor* hadPhysList = new HadronPhysicsQGSP_BERT("hadron");
00375         hadPhysList -> ConstructProcess();
00376         return;
00377       }
00378       
00379       if(BDSGlobalConstants::Instance()->GetPhysListName() == "hadronic_FTFP_BERT"){
00380         ConstructEM();
00381         HadronPhysicsFTFP_BERT *myHadPhysList = new HadronPhysicsFTFP_BERT;
00382         myHadPhysList->ConstructProcess();
00383         return;
00384       }
00385       
00386       if(BDSGlobalConstants::Instance()->GetPhysListName() == "hadronic_QGSP_BERT_HP_muon"){
00387         ConstructEM();
00388         ConstructMuon();
00389         ConstructHadronic();
00390         //      ConstructPhotolepton_Hadron();
00391         HadronPhysicsQGSP_BERT_HP *myHadPhysList = new HadronPhysicsQGSP_BERT_HP;
00392         myHadPhysList->ConstructProcess();
00393         return;
00394       }
00395       
00396       if(BDSGlobalConstants::Instance()->GetPhysListName() == "hadronic_FTFP_BERT_muon"){
00397         G4cout << __METHOD_NAME__ << "Using hadronic_FTFP_BERT_muon" << G4endl;
00398         ConstructEM();
00399         ConstructMuon();
00400         HadronPhysicsFTFP_BERT *myHadPhysList = new HadronPhysicsFTFP_BERT;
00401         myHadPhysList->ConstructProcess();
00402         
00403         return;
00404       }
00405       // physics list for laser wire - standard em stuff +
00406       // weighted compton scattering from laser wire
00407       if(BDSGlobalConstants::Instance()->GetPhysListName() == "lw") {
00408         ConstructEM();
00409         ConstructLaserWire();
00410         return;
00411       } 
00412       //default - standard (only transportation)
00413       G4cerr<<"WARNING : Unknown physics list "<<BDSGlobalConstants::Instance()->GetPhysListName()<<
00414         "  using transportation only (standard) "<<G4endl;
00415       return;
00416     }
00417   }
00418 
00419 }
00420 
00421 void BDSPhysicsList::ConstructParticle()
00422 {
00423   //standard physics lists
00424   if(BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_HP"){
00425     QGSP_BERT_HP* physList = new QGSP_BERT_HP;
00426     physList->ConstructParticle();
00427   } else if(BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT"){
00428     QGSP_BERT* physList = new QGSP_BERT;
00429     physList->ConstructParticle();
00430   } else if (BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_HP_muon"){
00431     QGSP_BERT_HP* physList = new QGSP_BERT_HP;
00432     physList->ConstructParticle();
00433   } else if (BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_muon"){
00434     QGSP_BERT* physList = new QGSP_BERT;
00435     physList->ConstructParticle();
00436   } else if(BDSGlobalConstants::Instance()->GetPhysListName() == "QGSP_BERT_HP_muon_em_low"){
00437     QGSP_BERT_HP* physList = new QGSP_BERT_HP;
00438     physList->ConstructParticle();
00439   } else {
00440     // pseudo-particles
00441     G4Geantino::GeantinoDefinition();
00442     G4ChargedGeantino::ChargedGeantinoDefinition();
00443     
00444     // gamma
00445     G4Gamma::GammaDefinition();
00446     
00447     // optical photon
00448     G4OpticalPhoton::OpticalPhotonDefinition();
00449     
00450     // leptons
00451     G4Electron::ElectronDefinition();
00452     G4Positron::PositronDefinition();
00453     G4MuonPlus::MuonPlusDefinition();
00454     G4MuonMinus::MuonMinusDefinition();
00455     
00456     G4NeutrinoE::NeutrinoEDefinition();
00457     G4AntiNeutrinoE::AntiNeutrinoEDefinition();
00458     G4NeutrinoMu::NeutrinoMuDefinition();
00459     G4AntiNeutrinoMu::AntiNeutrinoMuDefinition();  
00460     
00461     // mesons
00462     G4MesonConstructor mConstructor;
00463     mConstructor.ConstructParticle();
00464     
00465     // barions
00466     G4BaryonConstructor bConstructor;
00467     bConstructor.ConstructParticle();
00468     
00469     // ions
00470     G4IonConstructor iConstructor;
00471     iConstructor.ConstructParticle();
00472     
00473     //  Construct  resonaces and quarks
00474     G4ShortLivedConstructor pShortLivedConstructor;
00475     pShortLivedConstructor.ConstructParticle();
00476   }
00477   
00478   // set primary particle definition and kinetic beam parameters other than total energy
00479   G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
00480   BDSGlobalConstants::Instance()->SetParticleDefinition(particleTable->
00481                                     FindParticle(BDSGlobalConstants::Instance()->GetParticleName()) );  
00482   
00483   if(!BDSGlobalConstants::Instance()->GetParticleDefinition()) 
00484     {
00485       G4Exception("Particle not found, quitting!", "-1", FatalException, "");
00486       exit(1);
00487     }
00488   
00489   // set kinetic beam parameters other than total energy
00490   BDSGlobalConstants::Instance()->SetBeamMomentum( sqrt(pow(BDSGlobalConstants::Instance()->GetBeamTotalEnergy(),2)-
00491                                     pow(BDSGlobalConstants::Instance()->GetParticleDefinition()->GetPDGMass(),2)) );
00492   
00493   BDSGlobalConstants::Instance()->SetBeamKineticEnergy(BDSGlobalConstants::Instance()->GetBeamTotalEnergy() - 
00494                                    BDSGlobalConstants::Instance()->GetParticleDefinition()->GetPDGMass() );
00495   
00496   G4cout << __METHOD_NAME__ << "Beam properties:"<<G4endl;
00497   G4cout << __METHOD_NAME__ << "Particle : " 
00498          << BDSGlobalConstants::Instance()->GetParticleDefinition()->GetParticleName()<<G4endl;
00499   G4cout << __METHOD_NAME__ << "Mass : " 
00500          << BDSGlobalConstants::Instance()->GetParticleDefinition()->GetPDGMass()/GeV<< " GeV"<<G4endl;
00501   G4cout << __METHOD_NAME__ << "Charge : " 
00502          << BDSGlobalConstants::Instance()->GetParticleDefinition()->GetPDGCharge()<< " e"<<G4endl;
00503   G4cout << __METHOD_NAME__ << "Total Energy : "
00504          << BDSGlobalConstants::Instance()->GetBeamTotalEnergy()/GeV<<" GeV"<<G4endl;
00505   G4cout << __METHOD_NAME__ << "Kinetic Energy : "
00506          << BDSGlobalConstants::Instance()->GetBeamKineticEnergy()/GeV<<" GeV"<<G4endl;
00507   G4cout << __METHOD_NAME__ << "Momentum : "
00508          << BDSGlobalConstants::Instance()->GetBeamMomentum()/GeV<<" GeV"<<G4endl;
00509 }
00510 
00511 #include "G4Region.hh"
00512 #include "G4ProductionCuts.hh"
00513 void BDSPhysicsList::SetCuts()
00514 {
00515   if (verbose){
00516     G4cout << __METHOD_NAME__ << " setting cuts\n";
00517     
00518   }
00519   
00520   SetCutsWithDefault();   
00521 
00522 
00523   
00524     if(BDSGlobalConstants::Instance()->GetProdCutPhotons()>0)
00525       SetCutValue(BDSGlobalConstants::Instance()->GetProdCutPhotons(),G4ProductionCuts::GetIndex("gamma"));
00526   
00527    if(BDSGlobalConstants::Instance()->GetProdCutElectrons()>0)
00528      SetCutValue(BDSGlobalConstants::Instance()->GetProdCutElectrons(),G4ProductionCuts::GetIndex("e-"));
00529   
00530   if(BDSGlobalConstants::Instance()->GetProdCutPositrons()>0)
00531     SetCutValue(BDSGlobalConstants::Instance()->GetProdCutPositrons(),G4ProductionCuts::GetIndex("e+"));
00532   
00533 
00534     
00535   if(1)
00536     DumpCutValuesTable(); 
00537 
00538 }
00539 
00540 
00541 // particular physics process constructors
00542 
00543 void BDSPhysicsList::ConstructEM(){
00544   ConstructEMMisc();
00545   ConstructMultipleScattering();
00546 }
00547 
00548 void BDSPhysicsList::ConstructEMSingleScatter(){
00549   ConstructEMMisc();
00550   ConstructCoulombScattering();
00551 }
00552 
00553 void BDSPhysicsList::ConstructEMMisc()
00554 {
00555   theParticleIterator->reset();
00556   while( (*theParticleIterator)() ){
00557     G4ParticleDefinition* particle = theParticleIterator->value();
00558     G4ProcessManager* pmanager = particle->GetProcessManager();
00559     G4String particleName = particle->GetParticleName();
00560     if (particleName == "gamma") {
00561       // gamma         
00562       pmanager->AddDiscreteProcess(new G4PhotoElectricEffect);
00563       pmanager->AddDiscreteProcess(new G4ComptonScattering);
00564 
00565       if(0){
00566         G4GammaConversion* gammaconversion = new G4GammaConversion();
00567         gammaconversion->SetLambdaFactor(1/1.0e-20);
00568         BDSXSBias* gammaconversion_xsbias = new BDSXSBias();
00569         gammaconversion_xsbias->RegisterProcess(gammaconversion);
00570         gammaconversion_xsbias->SetEnhanceFactor(1e-20);
00571         pmanager->AddDiscreteProcess(gammaconversion_xsbias);
00572         
00573       } else if (BDSGlobalConstants::Instance()->GetUseEMLPB()){ //added by M.D. Salt, R.B. Appleby, 15/10/09
00574           G4GammaConversion* gammaconversion = new G4GammaConversion();
00575           GammaConversion_LPB* gammaconversion_lpb = new GammaConversion_LPB();
00576           gammaconversion_lpb->RegisterProcess(gammaconversion);
00577           pmanager->AddDiscreteProcess(gammaconversion_lpb);
00578       } else {
00579         pmanager->AddDiscreteProcess(new G4GammaConversion);
00580       }
00581     
00582       
00583     } else if (particleName == "e-") {
00584       pmanager->AddProcess(new G4eIonisation,       -1, 2,2);
00585       if(0){
00586         G4eBremsstrahlung* ebremsstrahlung = new G4eBremsstrahlung();
00587         ebremsstrahlung->SetLambdaFactor(1/1.0e-20);
00588         BDSXSBias* ebremsstrahlung_xsbias = new BDSXSBias();
00589         ebremsstrahlung_xsbias->RegisterProcess(ebremsstrahlung);
00590         ebremsstrahlung_xsbias->SetEnhanceFactor(1e-20);
00591         pmanager->AddDiscreteProcess(ebremsstrahlung_xsbias);     
00592       } else if(BDSGlobalConstants::Instance()->GetUseEMLPB()){ //added by M.D. Salt, R.B. Appleby, 15/10/09
00593           
00594         G4eBremsstrahlung* ebremsstrahlung = new G4eBremsstrahlung();
00595         eBremsstrahlung_LPB* ebremsstrahlung_lpb = new eBremsstrahlung_LPB();
00596         ebremsstrahlung_lpb->RegisterProcess(ebremsstrahlung);
00597         pmanager->AddProcess(ebremsstrahlung_lpb,     -1,-1,3);
00598       } else {
00599         G4eBremsstrahlung* ebremsstrahlung = new G4eBremsstrahlung();
00600         pmanager->AddProcess(ebremsstrahlung,   -1, 3,3);     
00601       }
00602             
00603       if(BDSGlobalConstants::Instance()->GetTurnOnCerenkov()){
00604 #if G4VERSION_NUMBER > 909
00605         G4Cerenkov* theCerenkovProcess = new G4Cerenkov;
00606         pmanager->AddProcess(theCerenkovProcess);
00607         pmanager->SetProcessOrdering(theCerenkovProcess,idxPostStep);
00608 #else
00609         pmanager->AddProcess(new G4Cerenkov,          -1, 5,-1);
00610 #endif
00611       }
00612       
00613     } else if (particleName == "e+") {
00614       //positron
00615       pmanager->AddProcess(new G4eIonisation,       -1, 2,2);
00616       if(0){
00617         G4eBremsstrahlung* ebremsstrahlung = new G4eBremsstrahlung();
00618         ebremsstrahlung->SetLambdaFactor(1/1.0e-20);
00619         BDSXSBias* ebremsstrahlung_xsbias = new BDSXSBias();
00620         ebremsstrahlung_xsbias->RegisterProcess(ebremsstrahlung);
00621         ebremsstrahlung_xsbias->SetEnhanceFactor(1e-20);
00622         pmanager->AddDiscreteProcess(ebremsstrahlung_xsbias);      
00623       } else if (BDSGlobalConstants::Instance()->GetUseEMLPB()){
00624         G4eBremsstrahlung* ebremsstrahlung = new G4eBremsstrahlung();
00625         eBremsstrahlung_LPB* ebremsstrahlung_lpb = new eBremsstrahlung_LPB();
00626         ebremsstrahlung_lpb->RegisterProcess(ebremsstrahlung);
00627         pmanager->AddProcess(ebremsstrahlung_lpb,     -1,-1,3);
00628       } else {
00629         pmanager->AddProcess(new G4eBremsstrahlung,   -1, 3,3);
00630       }
00631       pmanager->AddProcess(new G4eplusAnnihilation,  0,-1,4);
00632       if(BDSGlobalConstants::Instance()->GetTurnOnCerenkov()){      
00633 #if G4VERSION_NUMBER > 909
00634         G4Cerenkov* theCerenkovProcess = new G4Cerenkov;
00635         pmanager->AddProcess(theCerenkovProcess);
00636         pmanager->SetProcessOrdering(theCerenkovProcess,idxPostStep);
00637 #else
00638         pmanager->AddProcess(new G4Cerenkov,          -1, 5,-1);
00639 #endif 
00640       }
00641     } else if ((!particle->IsShortLived()) &&
00642                (particle->GetPDGCharge() != 0.0) && 
00643                (particle->GetParticleName() != "chargedgeantino")) {
00644       //all others charged particles except geantino
00645       pmanager->AddProcess(new G4hIonisation,       -1, 2,2);
00646            if(BDSGlobalConstants::Instance()->GetTurnOnCerenkov()){
00647 #if  G4VERSION_NUMBER > 909
00648         G4Cerenkov* theCerenkovProcess = new G4Cerenkov;
00649         pmanager->AddProcess(theCerenkovProcess);
00650         pmanager->SetProcessOrdering(theCerenkovProcess,idxPostStep);
00651 #else
00652         pmanager->AddProcess(new G4Cerenkov,          -1, 3,-1);
00653 #endif
00654       }
00655     }
00656   }
00657 }
00658 
00659 void BDSPhysicsList::ConstructMultipleScattering(){
00660   theParticleIterator->reset();
00661   while( (*theParticleIterator)() ){
00662     G4ParticleDefinition* particle = theParticleIterator->value();
00663     G4ProcessManager* pmanager = particle->GetProcessManager();
00664     G4String particleName = particle->GetParticleName();
00665     if (particleName == "e-") {
00666       //electron
00667 #if G4VERSION_NUMBER>919
00668       pmanager->AddProcess(new G4eMultipleScattering,-1, 1,1);
00669 #else
00670       pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00671 #endif
00672       
00673     } else if (particleName == "e+") {
00674       
00675 #if G4VERSION_NUMBER>919
00676       pmanager->AddProcess(new G4eMultipleScattering,-1, 1,1);
00677 #else
00678       pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00679 #endif
00680     } else if ((!particle->IsShortLived()) &&
00681                (particle->GetPDGCharge() != 0.0) && 
00682                (particle->GetParticleName() != "chargedgeantino")) {
00683       
00684       //all others charged particles except geantino
00685 #if G4VERSION_NUMBER>919
00686       pmanager->AddProcess(new G4hMultipleScattering,-1, 1,1);
00687 #else
00688       pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00689 #endif 
00690     }
00691   }
00692 }
00693 
00694 void BDSPhysicsList::ConstructCoulombScattering(){
00695   theParticleIterator->reset();
00696   while( (*theParticleIterator)() ){
00697     G4ParticleDefinition* particle = theParticleIterator->value();
00698     G4ProcessManager* pmanager = particle->GetProcessManager();
00699     G4String particleName = particle->GetParticleName(); 
00700     if (particleName == "e-") {
00701       pmanager->AddDiscreteProcess(new G4CoulombScattering);
00702     } else if (particleName == "e+") {
00703       pmanager->AddDiscreteProcess(new G4CoulombScattering);
00704     }else if ((!particle->IsShortLived()) &&
00705               (particle->GetPDGCharge() != 0.0) && 
00706               (particle->GetParticleName() != "chargedgeantino")) {
00707       //all others charged particles except geantino
00708       pmanager->AddDiscreteProcess(new G4CoulombScattering);
00709     } 
00710   }
00711 }
00712 
00713 // particular physics process constructors
00714 void BDSPhysicsList::ConstructMuon()
00715 {
00716   theParticleIterator->reset();
00717   while( (*theParticleIterator)() ){
00718     G4ParticleDefinition* particle = theParticleIterator->value();
00719     G4ProcessManager* pmanager = particle->GetProcessManager();
00720     G4String particleName = particle->GetParticleName();
00721 
00722     if (particleName == "gamma") {
00723       // gamma         
00724       G4GammaConversionToMuons* gammaconversiontomuons = new G4GammaConversionToMuons();
00725       BDSXSBias* gammaconversiontomuon_xsbias = new BDSXSBias();
00726       gammaconversiontomuons->SetCrossSecFactor(BDSGlobalConstants::Instance()->GetGammaToMuFe());
00727       gammaconversiontomuon_xsbias->RegisterProcess(gammaconversiontomuons);
00728       gammaconversiontomuon_xsbias->SetEnhanceFactor(BDSGlobalConstants::Instance()->GetGammaToMuFe());
00729       pmanager->AddDiscreteProcess(gammaconversiontomuon_xsbias);
00730 #ifdef DEBUG
00731       G4cout << __METHOD_NAME__ << "GammaToMuFe = " << BDSGlobalConstants::Instance()->GetGammaToMuFe() << G4endl;
00732 #endif
00733     } else if (particleName == "e+") {
00734       //positron
00735       //===========ee to hadrons in development================
00736       G4eeToHadrons* eetohadrons = new G4eeToHadrons();
00737       // BDSXSBias* eetohadrons_xsbias = new BDSXSBias();
00738       // G4cout << "eeToHadronsXSBias = " << BDSGlobalConstants::Instance()->GetEeToHadronsFe() << G4endl;
00739       eetohadrons->SetCrossSecFactor(BDSGlobalConstants::Instance()->GetEeToHadronsFe());
00740       //eetohadrons_xsbias->RegisterProcess(eetohadrons);
00741       //eetohadrons_xsbias->SetEnhanceFactor(BDSGlobalConstants::Instance()->GetEeToHadronsFe());
00742       //pmanager->AddDiscreteProcess(eetohadrons_xsbias);
00743       pmanager->AddDiscreteProcess(eetohadrons);
00744       //-------------------------------------------------------
00745       G4AnnihiToMuPair* annihitomupair = new G4AnnihiToMuPair();
00746       BDSXSBias* annihitomupair_xsbias = new BDSXSBias();
00747       annihitomupair->SetCrossSecFactor(BDSGlobalConstants::Instance()->GetAnnihiToMuFe());
00748       annihitomupair_xsbias->RegisterProcess(annihitomupair);
00749       annihitomupair_xsbias->SetEnhanceFactor(BDSGlobalConstants::Instance()->GetAnnihiToMuFe());
00750       pmanager->AddDiscreteProcess(annihitomupair_xsbias); 
00751 #ifdef DEBUG
00752       G4cout << __METHOD_NAME__ << "AnnihiToMuFe = " << BDSGlobalConstants::Instance()->GetAnnihiToMuFe() << G4endl;
00753 #endif    
00754     } else if( particleName == "mu+" || 
00755                particleName == "mu-"    ) {
00756       //muon  
00757 #if  G4VERSION_NUMBER>919
00758       pmanager->AddProcess(new G4MuMultipleScattering,-1, 1,1);
00759 #else 
00760       pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00761 #endif
00762       pmanager->AddProcess(new G4MuIonisation,      -1, 2,2);
00763       pmanager->AddProcess(new G4MuBremsstrahlung,  -1, 3,3);
00764       pmanager->AddProcess(new G4MuPairProduction,  -1, 4,4);
00765       if(BDSGlobalConstants::Instance()->GetTurnOnCerenkov()){
00766 #if  G4VERSION_NUMBER > 909
00767         G4Cerenkov* theCerenkovProcess = new G4Cerenkov;
00768         pmanager->AddProcess(theCerenkovProcess);
00769         pmanager->SetProcessOrdering(theCerenkovProcess,idxPostStep);
00770 #else
00771         pmanager->AddProcess(new G4Cerenkov,          -1, 5,-1);
00772 #endif
00773 #if G4VERSION_NUMBER < 950
00774         pmanager->AddDiscreteProcess(new G4MuonNucleusProcess);     
00775 #elif G4VERSION_NUMBER < 953
00776         pmanager->AddDiscreteProcess(new G4MuNuclearInteraction);     
00777 #else 
00778         /*      pmanager->AddDiscreteProcess(new G4MuonVDNuclearModel); */
00779 #endif
00780       }
00781     }
00782   }
00783 }
00784    
00785 
00786 void BDSPhysicsList::ConstructDecay()
00787 {
00788   theParticleIterator->reset();
00789   G4Decay* theDecayProcess = new G4Decay();
00790   while( (*theParticleIterator)() ){
00791     G4ParticleDefinition* particle = theParticleIterator->value();
00792     G4ProcessManager* pmanager = particle->GetProcessManager();
00793     G4String particleName = particle->GetParticleName();
00794     
00795     if (theDecayProcess->IsApplicable(*particle)) { 
00796       pmanager -> AddProcess(theDecayProcess);
00797       pmanager -> SetProcessOrdering(theDecayProcess, idxPostStep);
00798       pmanager -> SetProcessOrdering(theDecayProcess, idxAtRest);
00799     }
00800   }
00801 }
00802 
00803 
00804 void BDSPhysicsList::ConstructMerlin()
00805 {
00806   theParticleIterator->reset();
00807   while( (*theParticleIterator)() ){
00808     G4ParticleDefinition* particle = theParticleIterator->value();
00809     G4ProcessManager* pmanager = particle->GetProcessManager();
00810     G4String particleName = particle->GetParticleName();
00811     
00812     if (particleName == "e-") {
00813       //electron
00814 #if G4VERSION_NUMBER>919
00815       pmanager->AddProcess(new G4eMultipleScattering,-1, 1,1);
00816 #else 
00817       pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00818 #endif
00819       pmanager->AddProcess(new G4eIonisation,       -1, 2,2);
00820       pmanager->AddProcess(new G4eBremsstrahlung,   -1, 3,3);      
00821     } 
00822   }
00823 }
00824 
00825 
00826 void BDSPhysicsList::ConstructEM_Low_Energy()
00827 {
00828 #if G4VERSION_NUMBER > 949
00829   //Applicability range for Livermore models                                                                                                                                
00830   //for higher energies, the Standard models are used                                                                                                                       
00831   G4double highEnergyLimit = 1*GeV;
00832 #endif
00833 
00834   theParticleIterator->reset();
00835   while( (*theParticleIterator)() ){
00836     G4ParticleDefinition* particle = theParticleIterator->value();
00837     G4ProcessManager* pmanager = particle->GetProcessManager();
00838     G4String particleName = particle->GetParticleName();
00839      
00840     if (particleName == "gamma") {
00841 #if G4VERSION_NUMBER < 950
00842       pmanager->AddDiscreteProcess(new G4LowEnergyRayleigh());
00843       pmanager->AddDiscreteProcess(new G4LowEnergyPhotoElectric);
00844       pmanager->AddDiscreteProcess(new G4LowEnergyCompton);
00845       pmanager->AddDiscreteProcess(new G4LowEnergyGammaConversion);
00846 #else
00847       G4RayleighScattering* rayl = new G4RayleighScattering();
00848       G4LivermoreRayleighModel*
00849         raylModel = new G4LivermoreRayleighModel();
00850       raylModel->SetHighEnergyLimit(highEnergyLimit);
00851       rayl->AddEmModel(0, raylModel);
00852       pmanager->AddDiscreteProcess(rayl);
00853 
00854       G4PhotoElectricEffect* phot = new G4PhotoElectricEffect();
00855       G4LivermorePhotoElectricModel*
00856         photModel = new G4LivermorePhotoElectricModel();
00857       photModel->SetHighEnergyLimit(highEnergyLimit);
00858       phot->AddEmModel(0, photModel);
00859       pmanager->AddDiscreteProcess(phot);
00860 
00861       G4ComptonScattering* compt = new G4ComptonScattering();
00862       G4LivermoreComptonModel*
00863         comptModel = new G4LivermoreComptonModel();
00864       comptModel->SetHighEnergyLimit(highEnergyLimit);
00865       compt->AddEmModel(0, comptModel);
00866       pmanager->AddDiscreteProcess(compt);
00867       
00868       G4GammaConversion* conv = new G4GammaConversion();
00869       G4LivermoreGammaConversionModel*
00870         convModel = new G4LivermoreGammaConversionModel();
00871       convModel->SetHighEnergyLimit(highEnergyLimit);
00872       conv->AddEmModel(0, convModel);
00873       pmanager->AddDiscreteProcess(conv);
00874 #endif
00875       
00876     } else if (particleName == "e-") {
00877       #if G4VERSION_NUMBER>919
00878         pmanager->AddProcess(new G4eMultipleScattering,-1, 1,1);
00879       #else
00880         pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00881       #endif
00882 #if G4VERSION_NUMBER < 950
00883         pmanager->AddProcess(new G4LowEnergyIonisation,        -1, 2,2);
00884         pmanager->AddProcess(new G4LowEnergyBremsstrahlung,    -1, 3,3);
00885 #else
00886         G4eIonisation* eIoni = new G4eIonisation();
00887         G4LivermoreIonisationModel*
00888           eIoniModel = new G4LivermoreIonisationModel();
00889         eIoniModel->SetHighEnergyLimit(highEnergyLimit);
00890         eIoni->AddEmModel(0, eIoniModel, new G4UniversalFluctuation() );
00891         pmanager->AddProcess(eIoni,                   -1,-1, 1);
00892         
00893         G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
00894         G4LivermoreBremsstrahlungModel*
00895           eBremModel = new G4LivermoreBremsstrahlungModel();
00896         eBremModel->SetHighEnergyLimit(highEnergyLimit);
00897         eBrem->AddEmModel(0, eBremModel);
00898         pmanager->AddProcess(eBrem,                   -1,-1, 2);
00899 #endif
00900             
00901     } else if (particleName == "e+") {
00902       #if G4VERSION_NUMBER>919
00903         pmanager->AddProcess(new G4eMultipleScattering,-1, 1,1);
00904       #else 
00905         pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00906       #endif
00907         pmanager->AddProcess(new G4eIonisation,        -1, 2,2);
00908         pmanager->AddProcess(new G4eBremsstrahlung,    -1, 3,3);
00909         pmanager->AddProcess(new G4eplusAnnihilation,   0,-1,4);
00910     } else if( particleName == "mu+" || 
00911                particleName == "mu-"    ) {
00912       #if G4VERSION_NUMBER>919
00913         pmanager->AddProcess(new G4MuMultipleScattering,-1, 1,1);
00914       #else 
00915         pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00916       #endif
00917       pmanager->AddProcess(new G4MuIonisation,      -1, 2,2);
00918       pmanager->AddProcess(new G4MuBremsstrahlung,  -1, 3,3);
00919       pmanager->AddProcess(new G4MuPairProduction,  -1, 4,4);       
00920 
00921     } else if (particleName == "GenericIon") {
00922       #if G4VERSION_NUMBER>919
00923         pmanager->AddProcess(new G4hMultipleScattering,-1, 1,1);
00924       #else 
00925         pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00926 #endif
00927 #if G4VERSION_NUMBER < 950
00928         pmanager->AddProcess(new G4hLowEnergyIonisation,       -1,2,2);
00929       //      pmanager->AddProcess(new G4ionIonisation,      -1, 2,2);
00930       // it dose not work here
00931 #else
00932         pmanager->AddProcess(new G4ionIonisation,     -1,-1, 1);
00933 #endif
00934 
00935     } else if ((!particle->IsShortLived()) &&
00936                (particle->GetPDGCharge() != 0.0) && 
00937                (particle->GetParticleName() != "chargedgeantino")) {
00938 
00939       #if G4VERSION_NUMBER>919
00940         pmanager->AddProcess(new G4hMultipleScattering,-1, 1,1);
00941       #else 
00942         pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00943       #endif
00944 #if G4VERSION_NUMBER < 950
00945       pmanager->AddProcess(new G4hLowEnergyIonisation,       -1,2,2);
00946 #else
00947       pmanager->AddProcess(new G4hIonisation,       -1,-1, 1);
00948 #endif
00949     }
00950   }
00951 }
00952 
00953 void BDSPhysicsList::ConstructLaserWire()
00954 {
00955   G4cout << "Constructing laser-wire" << G4endl;
00956 
00957   theParticleIterator->reset();
00958 
00959   BDSLaserCompton* lwProcess = new BDSLaserCompton;
00960 
00961   while( (*theParticleIterator)() ){
00962     G4ParticleDefinition* particle = theParticleIterator->value();
00963     G4ProcessManager* pmanager = particle->GetProcessManager();
00964     G4String particleName = particle->GetParticleName();
00965     
00966     if (particleName == "e-") {
00967       pmanager->AddProcess(lwProcess);
00968       pmanager->SetProcessOrderingToLast(lwProcess,idxPostStep);
00969     }
00970 
00971     if (particleName == "e+") {
00972       pmanager->AddProcess(lwProcess);
00973       pmanager->SetProcessOrderingToLast(lwProcess,idxPostStep);
00974     }
00975 
00976   }
00977 
00978 }
00979 
00980 
00981 // Hadron Processes
00982 
00983 #include "G4HadronElasticProcess.hh"
00984 #include "G4HadronFissionProcess.hh"
00985 #include "G4HadronCaptureProcess.hh"
00986 
00987 #include "G4PionPlusInelasticProcess.hh"
00988 #include "G4PionMinusInelasticProcess.hh"
00989 #include "G4KaonPlusInelasticProcess.hh"
00990 #include "G4KaonZeroSInelasticProcess.hh"
00991 #include "G4KaonZeroLInelasticProcess.hh"
00992 #include "G4KaonMinusInelasticProcess.hh"
00993 #include "G4ProtonInelasticProcess.hh"
00994 #include "G4AntiProtonInelasticProcess.hh"
00995 #include "G4NeutronInelasticProcess.hh"
00996 #include "G4AntiNeutronInelasticProcess.hh"
00997 #include "G4LambdaInelasticProcess.hh"
00998 #include "G4AntiLambdaInelasticProcess.hh"
00999 #include "G4SigmaPlusInelasticProcess.hh"
01000 #include "G4SigmaMinusInelasticProcess.hh"
01001 #include "G4AntiSigmaPlusInelasticProcess.hh"
01002 #include "G4AntiSigmaMinusInelasticProcess.hh"
01003 #include "G4XiZeroInelasticProcess.hh"
01004 #include "G4XiMinusInelasticProcess.hh"
01005 #include "G4AntiXiZeroInelasticProcess.hh"
01006 #include "G4AntiXiMinusInelasticProcess.hh"
01007 #include "G4DeuteronInelasticProcess.hh"
01008 #include "G4TritonInelasticProcess.hh"
01009 #include "G4AlphaInelasticProcess.hh"
01010 #include "G4OmegaMinusInelasticProcess.hh"
01011 #include "G4AntiOmegaMinusInelasticProcess.hh"
01012 
01013 // Low-energy Models
01014 
01015 #include "G4LElastic.hh"
01016 #include "G4LFission.hh"
01017 #include "G4LCapture.hh"
01018 
01019 #include "G4LEPionPlusInelastic.hh"
01020 #include "G4LEPionMinusInelastic.hh"
01021 #include "G4LEKaonPlusInelastic.hh"
01022 #include "G4LEKaonZeroSInelastic.hh"
01023 #include "G4LEKaonZeroLInelastic.hh"
01024 #include "G4LEKaonMinusInelastic.hh"
01025 #include "G4LEProtonInelastic.hh"
01026 #include "G4LEAntiProtonInelastic.hh"
01027 #include "G4LENeutronInelastic.hh"
01028 #include "G4LEAntiNeutronInelastic.hh"
01029 #include "G4LELambdaInelastic.hh"
01030 #include "G4LEAntiLambdaInelastic.hh"
01031 #include "G4LESigmaPlusInelastic.hh"
01032 #include "G4LESigmaMinusInelastic.hh"
01033 #include "G4LEAntiSigmaPlusInelastic.hh"
01034 #include "G4LEAntiSigmaMinusInelastic.hh"
01035 #include "G4LEXiZeroInelastic.hh"
01036 #include "G4LEXiMinusInelastic.hh"
01037 #include "G4LEAntiXiZeroInelastic.hh"
01038 #include "G4LEAntiXiMinusInelastic.hh"
01039 #include "G4LEDeuteronInelastic.hh"
01040 #include "G4LETritonInelastic.hh"
01041 #include "G4LEAlphaInelastic.hh"
01042 #include "G4LEOmegaMinusInelastic.hh"
01043 #include "G4LEAntiOmegaMinusInelastic.hh"
01044 
01045 // -- generator models
01046 #include "G4TheoFSGenerator.hh"
01047 #include "G4ExcitationHandler.hh"
01048 #include "G4Evaporation.hh"
01049 #include "G4CompetitiveFission.hh"
01050 #include "G4FermiBreakUp.hh"
01051 #include "G4StatMF.hh"
01052 #include "G4GeneratorPrecompoundInterface.hh"
01053 #include "G4Fancy3DNucleus.hh"
01054 #include "G4LEProtonInelastic.hh"
01055 #include "G4StringModel.hh"
01056 #include "G4PreCompoundModel.hh"
01057 #include "G4FTFModel.hh"
01058 #include "G4QGSMFragmentation.hh"
01059 #include "G4ExcitedStringDecay.hh"
01060 
01061 //
01062 // ConstructHad()
01063 //
01064 // Makes discrete physics processes for the hadrons, at present limited
01065 // to those particles with GHEISHA interactions (INTRC > 0).
01066 // The processes are: Elastic scattering, Inelastic scattering,
01067 // Fission (for neutron only), and Capture (neutron).
01068 //
01069 // F.W.Jones  06-JUL-1998
01070 //
01071 
01072 void BDSPhysicsList::ConstructHad()
01073 {
01074     // this will be the model class for high energies
01075     G4TheoFSGenerator * theTheoModel = new G4TheoFSGenerator;
01076        
01077     // all models for treatment of thermal nucleus 
01078     G4Evaporation * theEvaporation = new G4Evaporation;
01079     G4FermiBreakUp * theFermiBreakUp = new G4FermiBreakUp;
01080     G4StatMF * theMF = new G4StatMF;
01081 
01082     // Evaporation logic
01083     G4ExcitationHandler * theHandler = new G4ExcitationHandler;
01084         theHandler->SetEvaporation(theEvaporation);
01085         theHandler->SetFermiModel(theFermiBreakUp);
01086         theHandler->SetMultiFragmentation(theMF);
01087         theHandler->SetMaxAandZForFermiBreakUp(12, 6);
01088         theHandler->SetMinEForMultiFrag(3*MeV);
01089         
01090     // Pre equilibrium stage 
01091     G4PreCompoundModel * thePreEquilib = new G4PreCompoundModel(theHandler);
01092 
01093     
01094     // a no-cascade generator-precompound interaface
01095     G4GeneratorPrecompoundInterface * theCascade = new G4GeneratorPrecompoundInterface;
01096             theCascade->SetDeExcitation(thePreEquilib);  
01097         
01098     // here come the high energy parts
01099     // the string model; still not quite according to design - Explicite use of the forseen interfaces 
01100     // will be tested and documented in this program by beta-02 at latest.
01101     G4VPartonStringModel * theStringModel;
01102     theStringModel = new G4FTFModel;
01103     theTheoModel->SetTransport(theCascade);
01104     theTheoModel->SetHighEnergyGenerator(theStringModel);
01105     theTheoModel->SetMinEnergy(19*GeV);
01106     theTheoModel->SetMaxEnergy(100*TeV);
01107 
01108       G4VLongitudinalStringDecay * theFragmentation = new G4QGSMFragmentation;
01109       G4ExcitedStringDecay * theStringDecay = new G4ExcitedStringDecay(theFragmentation);
01110       theStringModel->SetFragmentationModel(theStringDecay);
01111 
01112 // done with the generator model (most of the above is also available as default)
01113    G4HadronElasticProcess* theElasticProcess = 
01114                                     new G4HadronElasticProcess;
01115    G4LElastic* theElasticModel = new G4LElastic;
01116    theElasticProcess->RegisterMe(theElasticModel);
01117    G4HadronElasticProcess* theElasticProcess1 = 
01118                                     new G4HadronElasticProcess;
01119    theParticleIterator->reset();
01120    while ((*theParticleIterator)()) {
01121       G4ParticleDefinition* particle = theParticleIterator->value();
01122       G4ProcessManager* pmanager = particle->GetProcessManager();
01123       G4String particleName = particle->GetParticleName();
01124      
01125       if (particleName == "pi+") {
01126          pmanager->AddDiscreteProcess(theElasticProcess);
01127          G4PionPlusInelasticProcess* theInelasticProcess = 
01128                                 new G4PionPlusInelasticProcess("inelastic");
01129          G4LEPionPlusInelastic* theInelasticModel = 
01130                                 new G4LEPionPlusInelastic;
01131          theInelasticProcess->RegisterMe(theInelasticModel);
01132          theInelasticProcess->RegisterMe(theTheoModel);
01133          pmanager->AddDiscreteProcess(theInelasticProcess);
01134       }
01135       else if (particleName == "pi-") {
01136          pmanager->AddDiscreteProcess(theElasticProcess);
01137          G4PionMinusInelasticProcess* theInelasticProcess = 
01138                                 new G4PionMinusInelasticProcess("inelastic");
01139          G4LEPionMinusInelastic* theInelasticModel = 
01140                                 new G4LEPionMinusInelastic;
01141          theInelasticProcess->RegisterMe(theInelasticModel);
01142          theInelasticProcess->RegisterMe(theTheoModel);
01143          pmanager->AddDiscreteProcess(theInelasticProcess);
01144       }
01145       else if (particleName == "kaon+") {
01146          pmanager->AddDiscreteProcess(theElasticProcess);
01147          G4KaonPlusInelasticProcess* theInelasticProcess = 
01148                                   new G4KaonPlusInelasticProcess("inelastic");
01149          G4LEKaonPlusInelastic* theInelasticModel = new G4LEKaonPlusInelastic;
01150          theInelasticProcess->RegisterMe(theInelasticModel);
01151          theInelasticProcess->RegisterMe(theTheoModel);
01152          pmanager->AddDiscreteProcess(theInelasticProcess);
01153       }
01154       else if (particleName == "kaon0S") {
01155          pmanager->AddDiscreteProcess(theElasticProcess);
01156          G4KaonZeroSInelasticProcess* theInelasticProcess = 
01157                              new G4KaonZeroSInelasticProcess("inelastic");
01158          G4LEKaonZeroSInelastic* theInelasticModel = 
01159                              new G4LEKaonZeroSInelastic;
01160          theInelasticProcess->RegisterMe(theInelasticModel);
01161          theInelasticProcess->RegisterMe(theTheoModel);
01162          pmanager->AddDiscreteProcess(theInelasticProcess);
01163       }
01164       else if (particleName == "kaon0L") {
01165          pmanager->AddDiscreteProcess(theElasticProcess);
01166          G4KaonZeroLInelasticProcess* theInelasticProcess = 
01167                              new G4KaonZeroLInelasticProcess("inelastic");
01168          G4LEKaonZeroLInelastic* theInelasticModel = 
01169                              new G4LEKaonZeroLInelastic;
01170          theInelasticProcess->RegisterMe(theInelasticModel);
01171          theInelasticProcess->RegisterMe(theTheoModel);
01172          pmanager->AddDiscreteProcess(theInelasticProcess);
01173       }
01174       else if (particleName == "kaon-") {
01175          pmanager->AddDiscreteProcess(theElasticProcess);
01176          G4KaonMinusInelasticProcess* theInelasticProcess = 
01177                                  new G4KaonMinusInelasticProcess("inelastic");
01178          G4LEKaonMinusInelastic* theInelasticModel = 
01179                                  new G4LEKaonMinusInelastic;
01180          theInelasticProcess->RegisterMe(theInelasticModel);
01181          theInelasticProcess->RegisterMe(theTheoModel);
01182          pmanager->AddDiscreteProcess(theInelasticProcess);
01183       }
01184       else if (particleName == "proton") {
01185          pmanager->AddDiscreteProcess(theElasticProcess);
01186          G4ProtonInelasticProcess* theInelasticProcess = 
01187                                     new G4ProtonInelasticProcess("inelastic");
01188          G4LEProtonInelastic* theInelasticModel = new G4LEProtonInelastic;
01189          theInelasticProcess->RegisterMe(theInelasticModel);
01190          theInelasticProcess->RegisterMe(theTheoModel);
01191          pmanager->AddDiscreteProcess(theInelasticProcess);
01192       }
01193       else if (particleName == "anti_proton") {
01194          pmanager->AddDiscreteProcess(theElasticProcess);
01195          G4AntiProtonInelasticProcess* theInelasticProcess = 
01196                                 new G4AntiProtonInelasticProcess("inelastic");
01197          G4LEAntiProtonInelastic* theInelasticModel = 
01198                                 new G4LEAntiProtonInelastic;
01199          theInelasticProcess->RegisterMe(theInelasticModel);
01200          theInelasticProcess->RegisterMe(theTheoModel);
01201          pmanager->AddDiscreteProcess(theInelasticProcess);
01202       }
01203       else if (particleName == "neutron") {
01204          
01205           // elastic scattering
01206          G4LElastic* theElasticModel1 = new G4LElastic;
01207          theElasticProcess1->RegisterMe(theElasticModel1);
01208          pmanager->AddDiscreteProcess(theElasticProcess1);
01209           // inelastic scattering
01210          G4NeutronInelasticProcess* theInelasticProcess = 
01211                                     new G4NeutronInelasticProcess("inelastic");
01212          G4LENeutronInelastic* theInelasticModel = new G4LENeutronInelastic;
01213          theInelasticProcess->RegisterMe(theInelasticModel);
01214          theInelasticProcess->RegisterMe(theTheoModel);
01215          pmanager->AddDiscreteProcess(theInelasticProcess);
01216           // fission
01217          G4HadronFissionProcess* theFissionProcess =
01218                                     new G4HadronFissionProcess;
01219          G4LFission* theFissionModel = new G4LFission;
01220          theFissionProcess->RegisterMe(theFissionModel);
01221          pmanager->AddDiscreteProcess(theFissionProcess);
01222          // capture
01223          G4HadronCaptureProcess* theCaptureProcess =
01224                                     new G4HadronCaptureProcess;
01225          G4LCapture* theCaptureModel = new G4LCapture;
01226          theCaptureProcess->RegisterMe(theCaptureModel);
01227          pmanager->AddDiscreteProcess(theCaptureProcess);
01228       }  
01229       else if (particleName == "anti_neutron") {
01230          pmanager->AddDiscreteProcess(theElasticProcess);
01231          G4AntiNeutronInelasticProcess* theInelasticProcess = 
01232                                new G4AntiNeutronInelasticProcess("inelastic");
01233          G4LEAntiNeutronInelastic* theInelasticModel = 
01234                                new G4LEAntiNeutronInelastic;
01235          theInelasticProcess->RegisterMe(theInelasticModel);
01236          theInelasticProcess->RegisterMe(theTheoModel);
01237          pmanager->AddDiscreteProcess(theInelasticProcess);
01238       }
01239       else if (particleName == "lambda") {
01240          pmanager->AddDiscreteProcess(theElasticProcess);
01241          G4LambdaInelasticProcess* theInelasticProcess = 
01242                                     new G4LambdaInelasticProcess("inelastic");
01243          G4LELambdaInelastic* theInelasticModel = new G4LELambdaInelastic;
01244          theInelasticProcess->RegisterMe(theInelasticModel);
01245          theInelasticProcess->RegisterMe(theTheoModel);
01246          pmanager->AddDiscreteProcess(theInelasticProcess);
01247       }
01248       else if (particleName == "anti_lambda") {
01249          pmanager->AddDiscreteProcess(theElasticProcess);
01250          G4AntiLambdaInelasticProcess* theInelasticProcess = 
01251                                 new G4AntiLambdaInelasticProcess("inelastic");
01252          G4LEAntiLambdaInelastic* theInelasticModel = 
01253                                 new G4LEAntiLambdaInelastic;
01254          theInelasticProcess->RegisterMe(theInelasticModel);
01255          theInelasticProcess->RegisterMe(theTheoModel);
01256          pmanager->AddDiscreteProcess(theInelasticProcess);
01257       }
01258       else if (particleName == "sigma+") {
01259          pmanager->AddDiscreteProcess(theElasticProcess);
01260          G4SigmaPlusInelasticProcess* theInelasticProcess = 
01261                                  new G4SigmaPlusInelasticProcess("inelastic");
01262          G4LESigmaPlusInelastic* theInelasticModel = 
01263                                  new G4LESigmaPlusInelastic;
01264          theInelasticProcess->RegisterMe(theInelasticModel);
01265          theInelasticProcess->RegisterMe(theTheoModel);
01266          pmanager->AddDiscreteProcess(theInelasticProcess);
01267       }
01268       else if (particleName == "sigma-") {
01269          pmanager->AddDiscreteProcess(theElasticProcess);
01270          G4SigmaMinusInelasticProcess* theInelasticProcess = 
01271                                  new G4SigmaMinusInelasticProcess("inelastic");
01272          G4LESigmaMinusInelastic* theInelasticModel = 
01273                                  new G4LESigmaMinusInelastic;
01274          theInelasticProcess->RegisterMe(theInelasticModel);
01275          theInelasticProcess->RegisterMe(theTheoModel);
01276          pmanager->AddDiscreteProcess(theInelasticProcess);
01277       }
01278       else if (particleName == "anti_sigma+") {
01279          pmanager->AddDiscreteProcess(theElasticProcess);
01280          G4AntiSigmaPlusInelasticProcess* theInelasticProcess = 
01281                              new G4AntiSigmaPlusInelasticProcess("inelastic");
01282          G4LEAntiSigmaPlusInelastic* theInelasticModel = 
01283                                  new G4LEAntiSigmaPlusInelastic;
01284          theInelasticProcess->RegisterMe(theInelasticModel);
01285          theInelasticProcess->RegisterMe(theTheoModel);
01286          pmanager->AddDiscreteProcess(theInelasticProcess);
01287       }
01288       else if (particleName == "anti_sigma-") {
01289          pmanager->AddDiscreteProcess(theElasticProcess);
01290          G4AntiSigmaMinusInelasticProcess* theInelasticProcess = 
01291                             new G4AntiSigmaMinusInelasticProcess("inelastic");
01292          G4LEAntiSigmaMinusInelastic* theInelasticModel = 
01293                                  new G4LEAntiSigmaMinusInelastic;
01294          theInelasticProcess->RegisterMe(theInelasticModel);
01295          theInelasticProcess->RegisterMe(theTheoModel);
01296          pmanager->AddDiscreteProcess(theInelasticProcess);
01297       }
01298       else if (particleName == "xi0") {
01299          pmanager->AddDiscreteProcess(theElasticProcess);
01300          G4XiZeroInelasticProcess* theInelasticProcess = 
01301                             new G4XiZeroInelasticProcess("inelastic");
01302          G4LEXiZeroInelastic* theInelasticModel = 
01303                                  new G4LEXiZeroInelastic;
01304          theInelasticProcess->RegisterMe(theInelasticModel);
01305          theInelasticProcess->RegisterMe(theTheoModel);
01306          pmanager->AddDiscreteProcess(theInelasticProcess);
01307       }
01308       else if (particleName == "xi-") {
01309          pmanager->AddDiscreteProcess(theElasticProcess);
01310          G4XiMinusInelasticProcess* theInelasticProcess = 
01311                             new G4XiMinusInelasticProcess("inelastic");
01312          G4LEXiMinusInelastic* theInelasticModel = 
01313                                  new G4LEXiMinusInelastic;
01314          theInelasticProcess->RegisterMe(theInelasticModel);
01315          theInelasticProcess->RegisterMe(theTheoModel);
01316          pmanager->AddDiscreteProcess(theInelasticProcess);
01317       }
01318       else if (particleName == "anti_xi0") {
01319          pmanager->AddDiscreteProcess(theElasticProcess);
01320          G4AntiXiZeroInelasticProcess* theInelasticProcess = 
01321                             new G4AntiXiZeroInelasticProcess("inelastic");
01322          G4LEAntiXiZeroInelastic* theInelasticModel = 
01323                                  new G4LEAntiXiZeroInelastic;
01324          theInelasticProcess->RegisterMe(theInelasticModel);
01325          theInelasticProcess->RegisterMe(theTheoModel);
01326          pmanager->AddDiscreteProcess(theInelasticProcess);
01327       }
01328       else if (particleName == "anti_xi-") {
01329          pmanager->AddDiscreteProcess(theElasticProcess);
01330          G4AntiXiMinusInelasticProcess* theInelasticProcess = 
01331                             new G4AntiXiMinusInelasticProcess("inelastic");
01332          G4LEAntiXiMinusInelastic* theInelasticModel = 
01333                                  new G4LEAntiXiMinusInelastic;
01334          theInelasticProcess->RegisterMe(theInelasticModel);
01335          theInelasticProcess->RegisterMe(theTheoModel);
01336          pmanager->AddDiscreteProcess(theInelasticProcess);
01337       }
01338       else if (particleName == "deuteron") {
01339          pmanager->AddDiscreteProcess(theElasticProcess);
01340          G4DeuteronInelasticProcess* theInelasticProcess = 
01341                             new G4DeuteronInelasticProcess("inelastic");
01342          G4LEDeuteronInelastic* theInelasticModel = 
01343                                  new G4LEDeuteronInelastic;
01344          theInelasticProcess->RegisterMe(theInelasticModel);
01345          theInelasticProcess->RegisterMe(theTheoModel);
01346          pmanager->AddDiscreteProcess(theInelasticProcess);
01347       }
01348       else if (particleName == "triton") {
01349          pmanager->AddDiscreteProcess(theElasticProcess);
01350          G4TritonInelasticProcess* theInelasticProcess = 
01351                             new G4TritonInelasticProcess("inelastic");
01352          G4LETritonInelastic* theInelasticModel = 
01353                                  new G4LETritonInelastic;
01354          theInelasticProcess->RegisterMe(theInelasticModel);
01355          theInelasticProcess->RegisterMe(theTheoModel);
01356          pmanager->AddDiscreteProcess(theInelasticProcess);
01357       }
01358       else if (particleName == "alpha") {
01359          pmanager->AddDiscreteProcess(theElasticProcess);
01360          G4AlphaInelasticProcess* theInelasticProcess = 
01361                             new G4AlphaInelasticProcess("inelastic");
01362          G4LEAlphaInelastic* theInelasticModel = 
01363                                  new G4LEAlphaInelastic;
01364          theInelasticProcess->RegisterMe(theInelasticModel);
01365          theInelasticProcess->RegisterMe(theTheoModel);
01366          pmanager->AddDiscreteProcess(theInelasticProcess);
01367       }
01368       else if (particleName == "omega-") {
01369          pmanager->AddDiscreteProcess(theElasticProcess);
01370          G4OmegaMinusInelasticProcess* theInelasticProcess = 
01371                             new G4OmegaMinusInelasticProcess("inelastic");
01372          G4LEOmegaMinusInelastic* theInelasticModel = 
01373                                  new G4LEOmegaMinusInelastic;
01374          theInelasticProcess->RegisterMe(theInelasticModel);
01375          theInelasticProcess->RegisterMe(theTheoModel);
01376          pmanager->AddDiscreteProcess(theInelasticProcess);
01377       }
01378       else if (particleName == "anti_omega-") {
01379          pmanager->AddDiscreteProcess(theElasticProcess);
01380          G4AntiOmegaMinusInelasticProcess* theInelasticProcess = 
01381                             new G4AntiOmegaMinusInelasticProcess("inelastic");
01382          G4LEAntiOmegaMinusInelastic* theInelasticModel = 
01383                                  new G4LEAntiOmegaMinusInelastic;
01384          theInelasticProcess->RegisterMe(theInelasticModel);
01385          theInelasticProcess->RegisterMe(theTheoModel);
01386          pmanager->AddDiscreteProcess(theInelasticProcess);
01387       }
01388    }
01389 }
01390 
01391 /*
01392 void BDSPhysicsList::ConstructPhotolepton_Hadron(){
01393   G4TheoFSGenerator * theModel;
01394   G4GeneratorPrecompoundInterface * theCascade;
01395   G4QGSModel< G4GammaParticipants > * theStringModel;
01396   G4QGSMFragmentation * theFragmentation;
01397   G4ExcitedStringDecay * theStringDecay;
01398 
01399   G4PhotoNuclearProcess * thePhotoNuclearProcess;
01400   G4ElectronNuclearProcess * theElectronNuclearProcess;
01401   G4PositronNuclearProcess * thePositronNuclearProcess;
01402   G4ElectroNuclearReaction * theElectroReaction;
01403   G4GammaNuclearReaction * theGammaReaction;  
01404 
01405   theModel = new G4TheoFSGenerator;
01406   
01407   theStringModel = new G4QGSModel< G4GammaParticipants >;
01408   theStringDecay = new G4ExcitedStringDecay(theFragmentation=new G4QGSMFragmentation);
01409   theStringModel->SetFragmentationModel(theStringDecay);
01410   
01411   theCascade = new G4GeneratorPrecompoundInterface;
01412   
01413   theModel->SetTransport(theCascade);
01414   theModel->SetHighEnergyGenerator(theStringModel);
01415 
01416   G4ProcessManager * aProcMan = 0;
01417   
01418   aProcMan = G4Gamma::Gamma()->GetProcessManager();
01419   theGammaReaction->SetMaxEnergy(3.5*GeV);
01420   thePhotoNuclearProcess->RegisterMe(theGammaReaction);
01421   theModel->SetMinEnergy(3.*GeV);
01422   theModel->SetMaxEnergy(100*TeV);
01423   thePhotoNuclearProcess->RegisterMe(theModel);
01424   aProcMan->AddDiscreteProcess(thePhotoNuclearProcess);
01425 
01426   aProcMan = G4Electron::Electron()->GetProcessManager();
01427   theElectronNuclearProcess->RegisterMe(theElectroReaction);
01428   aProcMan->AddDiscreteProcess(theElectronNuclearProcess);
01429   
01430   aProcMan = G4Positron::Positron()->GetProcessManager();
01431   thePositronNuclearProcess->RegisterMe(theElectroReaction);
01432   aProcMan->AddDiscreteProcess(thePositronNuclearProcess);
01433 }
01434 */
01435 
01436 void BDSPhysicsList::ConstructHadronic()
01437 {
01438 
01439   G4NeutronBuilder* theNeutrons=new G4NeutronBuilder;
01440   G4LHEPNeutronBuilder * theLHEPNeutron;
01441   theNeutrons->RegisterMe(theLHEPNeutron=new G4LHEPNeutronBuilder);
01442 
01443   G4ProtonBuilder * thePro;
01444   G4LHEPProtonBuilder * theLHEPPro;
01445 
01446   thePro=new G4ProtonBuilder;
01447   thePro->RegisterMe(theLHEPPro=new G4LHEPProtonBuilder);
01448 
01449   G4PiKBuilder * thePiK;
01450   G4LHEPPiKBuilder * theLHEPPiK;
01451 
01452   thePiK=new G4PiKBuilder;
01453   thePiK->RegisterMe(theLHEPPiK=new G4LHEPPiKBuilder);
01454 
01455   theNeutrons->Build();
01456   thePro->Build();
01457   thePiK->Build();
01458 
01459   // Photonuclear processes
01460 
01461   G4PhotoNuclearProcess * thePhotoNuclearProcess;
01462   G4ElectronNuclearProcess * theElectronNuclearProcess;
01463   G4PositronNuclearProcess * thePositronNuclearProcess;
01464   G4ElectroNuclearReaction * theElectroReaction;
01465   G4GammaNuclearReaction * theGammaReaction;  
01466   
01467   G4TheoFSGenerator * theModel;
01468   G4GeneratorPrecompoundInterface * theCascade;
01469   G4QGSModel< G4GammaParticipants > * theStringModel;
01470   G4QGSMFragmentation * theFragmentation;
01471   G4ExcitedStringDecay * theStringDecay;
01472 
01473   thePhotoNuclearProcess = new G4PhotoNuclearProcess;
01474   theGammaReaction = new G4GammaNuclearReaction;
01475   theElectronNuclearProcess = new G4ElectronNuclearProcess;
01476   thePositronNuclearProcess = new G4PositronNuclearProcess;
01477   theElectroReaction = new G4ElectroNuclearReaction;
01478 
01479   theModel = new G4TheoFSGenerator;
01480   
01481   theStringModel = new G4QGSModel< G4GammaParticipants >;
01482   theStringDecay = new G4ExcitedStringDecay(theFragmentation=new G4QGSMFragmentation);
01483   theStringModel->SetFragmentationModel(theStringDecay);
01484   
01485   theCascade = new G4GeneratorPrecompoundInterface;
01486   
01487   theModel->SetTransport(theCascade);
01488   theModel->SetHighEnergyGenerator(theStringModel);
01489 
01490   G4ProcessManager * aProcMan = 0;
01491   
01492   aProcMan = G4Gamma::Gamma()->GetProcessManager();
01493   theGammaReaction->SetMaxEnergy(3.5*GeV);
01494   thePhotoNuclearProcess->RegisterMe(theGammaReaction);
01495   theModel->SetMinEnergy(3.*GeV);
01496   theModel->SetMaxEnergy(100*TeV);
01497   thePhotoNuclearProcess->RegisterMe(theModel);
01498   aProcMan->AddDiscreteProcess(thePhotoNuclearProcess);
01499 
01500   aProcMan = G4Electron::Electron()->GetProcessManager();
01501   theElectronNuclearProcess->RegisterMe(theElectroReaction);
01502   aProcMan->AddDiscreteProcess(theElectronNuclearProcess);
01503   
01504   aProcMan = G4Positron::Positron()->GetProcessManager();
01505   thePositronNuclearProcess->RegisterMe(theElectroReaction);
01506   aProcMan->AddDiscreteProcess(thePositronNuclearProcess);
01507 }
01508 
01509 void BDSPhysicsList::ConstructSR()
01510 {
01511   // BDSIM's version of Synchrotron Radiation
01512   BDSSynchrotronRadiation* srProcess = new BDSSynchrotronRadiation;
01513   
01514   BDSContinuousSR *contSR = new BDSContinuousSR(); // contin. energy loss process
01515 
01516   // G4's version of Synchrotron Radiation - not used because does not have
01517   // Multiplicity or MeanFreeFactor capability
01518   // G4SynchrotronRadiation* srProcess = new G4SynchrotronRadiation;
01519 
01520   theParticleIterator->reset();
01521 
01522   while( (*theParticleIterator)() ){
01523     G4ParticleDefinition* particle = theParticleIterator->value();
01524     G4ProcessManager* pmanager = particle->GetProcessManager();
01525     G4String particleName = particle->GetParticleName();
01526     
01527     if (particleName == "e-") {
01528       pmanager->AddProcess(srProcess);
01529       pmanager->SetProcessOrderingToLast(srProcess,idxPostStep);
01530 
01531       G4int idx = pmanager->AddProcess(contSR);
01532       pmanager->SetProcessOrderingToLast(contSR,idxPostStep);
01533       pmanager->SetProcessActivation(idx, false);
01534     }
01535     
01536     if (particleName == "e+") {
01537       pmanager->AddProcess(srProcess);
01538       pmanager->SetProcessOrderingToLast(srProcess,idxPostStep);
01539 
01540       //G4int idx = pmanager->AddProcess(contSR);
01541       //      pmanager->SetProcessOrderingToLast(contSR,idxPostStep);
01542       //      pmanager->SetProcessActivation(idx, false);
01543     }
01544     
01545   }
01546   return; 
01547 }
01548 
01549 void BDSPhysicsList::AddParameterisation()
01550 {
01551   G4FastSimulationManagerProcess*
01552     theFastSimulationManagerProcess =
01553     new G4FastSimulationManagerProcess();
01554   G4cout << "FastSimulationManagerProcess" <<G4endl;
01555   theParticleIterator->reset();
01556   //G4cout<<"---"<<G4endl;                                                                                                                                              
01557   while( (*theParticleIterator)() ){
01558     //G4cout<<"+++"<<G4endl;                                                                                                                                            
01559 
01560     G4ParticleDefinition* particle = theParticleIterator->value();
01561     // G4cout<<"--- particle "<<particle->GetParticleName()<<G4endl;                                                                                                    
01562     G4ProcessManager* pmanager = particle->GetProcessManager();
01563     // The fast simulation process becomes a discrete process only since 9.0:                                                                                                 
01564     pmanager->AddDiscreteProcess(theFastSimulationManagerProcess);
01565   }
01566 }

Generated on 27 Aug 2013 for BDSIM by  doxygen 1.4.7