G4EmLivermorePolarizedPhysics.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id$
00027 
00028 #include "G4EmLivermorePolarizedPhysics.hh"
00029 #include "G4ParticleDefinition.hh"
00030 #include "G4SystemOfUnits.hh"
00031 
00032 // *** Processes and models
00033 
00034 // gamma
00035 #include "G4PhotoElectricEffect.hh"
00036 #include "G4LivermorePolarizedPhotoElectricModel.hh"
00037 
00038 #include "G4ComptonScattering.hh"
00039 #include "G4LivermorePolarizedComptonModel.hh"
00040 
00041 #include "G4GammaConversion.hh"
00042 #include "G4LivermorePolarizedGammaConversionModel.hh"
00043 
00044 #include "G4RayleighScattering.hh" 
00045 #include "G4LivermorePolarizedRayleighModel.hh"
00046 
00047 // e+-
00048 #include "G4eMultipleScattering.hh"
00049 #include "G4UniversalFluctuation.hh"
00050 
00051 #include "G4eIonisation.hh"
00052 #include "G4LivermoreIonisationModel.hh"
00053 
00054 #include "G4eBremsstrahlung.hh"
00055 #include "G4LivermoreBremsstrahlungModel.hh"
00056 
00057 // e+
00058 #include "G4eplusAnnihilation.hh"
00059 
00060 // mu+-
00061 #include "G4MuMultipleScattering.hh"
00062 #include "G4MuIonisation.hh"
00063 #include "G4MuBremsstrahlung.hh"
00064 #include "G4MuPairProduction.hh"
00065 
00066 #include "G4MuBremsstrahlungModel.hh"
00067 #include "G4MuPairProductionModel.hh"
00068 #include "G4hBremsstrahlungModel.hh"
00069 #include "G4hPairProductionModel.hh"
00070 
00071 // hadrons
00072 #include "G4hMultipleScattering.hh"
00073 #include "G4MscStepLimitType.hh"
00074 
00075 #include "G4hBremsstrahlung.hh"
00076 #include "G4hPairProduction.hh"
00077 
00078 #include "G4hIonisation.hh"
00079 #include "G4ionIonisation.hh"
00080 #include "G4alphaIonisation.hh"
00081 #include "G4IonParametrisedLossModel.hh"
00082 #include "G4NuclearStopping.hh"
00083 
00084 // msc models
00085 #include "G4UrbanMscModel95.hh"
00086 #include "G4WentzelVIModel.hh"
00087 #include "G4GoudsmitSaundersonMscModel.hh"
00088 #include "G4CoulombScattering.hh"
00089 #include "G4eCoulombScatteringModel.hh"
00090 
00091 // interfaces
00092 #include "G4LossTableManager.hh"
00093 #include "G4EmProcessOptions.hh"
00094 #include "G4UAtomicDeexcitation.hh"
00095 
00096 // particles
00097 #include "G4Gamma.hh"
00098 #include "G4Electron.hh"
00099 #include "G4Positron.hh"
00100 #include "G4MuonPlus.hh"
00101 #include "G4MuonMinus.hh"
00102 #include "G4PionPlus.hh"
00103 #include "G4PionMinus.hh"
00104 #include "G4KaonPlus.hh"
00105 #include "G4KaonMinus.hh"
00106 #include "G4Proton.hh"
00107 #include "G4AntiProton.hh"
00108 #include "G4Deuteron.hh"
00109 #include "G4Triton.hh"
00110 #include "G4He3.hh"
00111 #include "G4Alpha.hh"
00112 #include "G4GenericIon.hh"
00113 
00114 //
00115 #include "G4PhysicsListHelper.hh"
00116 #include "G4BuilderType.hh"
00117 
00118 // factory
00119 #include "G4PhysicsConstructorFactory.hh"
00120 //
00121 G4_DECLARE_PHYSCONSTR_FACTORY(G4EmLivermorePolarizedPhysics);
00122 
00123 
00124 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00125 
00126 G4EmLivermorePolarizedPhysics::G4EmLivermorePolarizedPhysics(G4int ver)
00127   : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
00128 {
00129   G4LossTableManager::Instance();
00130   SetPhysicsType(bElectromagnetic);
00131 }
00132 
00133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00134 
00135 G4EmLivermorePolarizedPhysics::G4EmLivermorePolarizedPhysics(G4int ver, const G4String&)
00136   : G4VPhysicsConstructor("G4EmLivermorePolarizedPhysics"), verbose(ver)
00137 {
00138   G4LossTableManager::Instance();
00139   SetPhysicsType(bElectromagnetic);
00140 }
00141 
00142 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00143 
00144 G4EmLivermorePolarizedPhysics::~G4EmLivermorePolarizedPhysics()
00145 {}
00146 
00147 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00148 
00149 void G4EmLivermorePolarizedPhysics::ConstructParticle()
00150 {
00151 // gamma
00152   G4Gamma::Gamma();
00153 
00154 // leptons
00155   G4Electron::Electron();
00156   G4Positron::Positron();
00157   G4MuonPlus::MuonPlus();
00158   G4MuonMinus::MuonMinus();
00159 
00160 // mesons
00161   G4PionPlus::PionPlusDefinition();
00162   G4PionMinus::PionMinusDefinition();
00163   G4KaonPlus::KaonPlusDefinition();
00164   G4KaonMinus::KaonMinusDefinition();
00165 
00166 // baryons
00167   G4Proton::Proton();
00168   G4AntiProton::AntiProton();
00169 
00170 // ions
00171   G4Deuteron::Deuteron();
00172   G4Triton::Triton();
00173   G4He3::He3();
00174   G4Alpha::Alpha();
00175   G4GenericIon::GenericIonDefinition();
00176 }
00177 
00178 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00179 
00180 void G4EmLivermorePolarizedPhysics::ConstructProcess()
00181 {
00182   G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
00183 
00184   // muon & hadron bremsstrahlung and pair production
00185   G4MuBremsstrahlung* mub = new G4MuBremsstrahlung();
00186   G4MuPairProduction* mup = new G4MuPairProduction();
00187   G4hBremsstrahlung* pib = new G4hBremsstrahlung();
00188   G4hPairProduction* pip = new G4hPairProduction();
00189   G4hBremsstrahlung* kb = new G4hBremsstrahlung();
00190   G4hPairProduction* kp = new G4hPairProduction();
00191   G4hBremsstrahlung* pb = new G4hBremsstrahlung();
00192   G4hPairProduction* pp = new G4hPairProduction();
00193 
00194   // muon & hadron multiple scattering
00195   G4MuMultipleScattering* mumsc = new G4MuMultipleScattering();
00196   mumsc->AddEmModel(0, new G4WentzelVIModel());
00197   G4MuMultipleScattering* pimsc = new G4MuMultipleScattering();
00198   pimsc->AddEmModel(0, new G4WentzelVIModel());
00199   G4MuMultipleScattering* kmsc = new G4MuMultipleScattering();
00200   kmsc->AddEmModel(0, new G4WentzelVIModel());
00201   G4MuMultipleScattering* pmsc = new G4MuMultipleScattering();
00202   pmsc->AddEmModel(0, new G4WentzelVIModel());
00203   G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
00204 
00205   // high energy limit for e+- scattering models
00206   G4double highEnergyLimit = 100*MeV;
00207 
00208   // nuclear stopping
00209   G4NuclearStopping* ionnuc = new G4NuclearStopping();
00210   G4NuclearStopping* pnuc = new G4NuclearStopping();
00211 
00212   // Add Livermore EM Processes
00213   theParticleIterator->reset();
00214 
00215   while( (*theParticleIterator)() ){
00216   
00217     G4ParticleDefinition* particle = theParticleIterator->value();
00218     G4String particleName = particle->GetParticleName();
00219     
00220     if(verbose > 1)
00221       G4cout << "### " << GetPhysicsName() << " instantiates for " 
00222              << particleName << G4endl;
00223 
00224     //Applicability range for Livermore models
00225     //for higher energies, the Standard models are used   
00226     G4double LivermoreHighEnergyLimit = GeV;
00227 
00228     if (particleName == "gamma") {
00229 
00230       G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
00231       G4LivermorePolarizedPhotoElectricModel* theLivermorePhotoElectricModel = new G4LivermorePolarizedPhotoElectricModel();
00232       theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00233       thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
00234       ph->RegisterProcess(thePhotoElectricEffect, particle);
00235 
00236       G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
00237       G4LivermorePolarizedComptonModel* theLivermoreComptonModel = new G4LivermorePolarizedComptonModel();
00238       theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00239       theComptonScattering->AddEmModel(0, theLivermoreComptonModel);
00240       ph->RegisterProcess(theComptonScattering, particle);
00241 
00242       G4GammaConversion* theGammaConversion = new G4GammaConversion();
00243       G4LivermorePolarizedGammaConversionModel* theLivermoreGammaConversionModel = new G4LivermorePolarizedGammaConversionModel();
00244       theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00245       theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
00246       ph->RegisterProcess(theGammaConversion, particle);
00247 
00248       G4RayleighScattering* theRayleigh = new G4RayleighScattering();
00249       G4LivermorePolarizedRayleighModel* theRayleighModel = new G4LivermorePolarizedRayleighModel();
00250       theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
00251       theRayleigh->AddEmModel(0, theRayleighModel);
00252       ph->RegisterProcess(theRayleigh, particle);
00253 
00254     } else if (particleName == "e-") {
00255 
00256       // multiple scattering
00257       G4eMultipleScattering* msc = new G4eMultipleScattering;
00258       msc->SetStepLimitType(fUseDistanceToBoundary);
00259       G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
00260       G4WentzelVIModel* msc2 = new G4WentzelVIModel();
00261       msc1->SetHighEnergyLimit(highEnergyLimit);
00262       msc2->SetLowEnergyLimit(highEnergyLimit);
00263       msc->AddEmModel(0, msc1);
00264       msc->AddEmModel(0, msc2);
00265 
00266       G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 
00267       G4CoulombScattering* ss = new G4CoulombScattering();
00268       ss->SetEmModel(ssm, 1); 
00269       ss->SetMinKinEnergy(highEnergyLimit);
00270       ssm->SetLowEnergyLimit(highEnergyLimit);
00271       ssm->SetActivationLowEnergyLimit(highEnergyLimit);
00272       ph->RegisterProcess(msc, particle);
00273       ph->RegisterProcess(ss, particle);
00274       
00275       // Ionisation
00276       G4eIonisation* eIoni = new G4eIonisation();
00277       G4LivermoreIonisationModel* theIoniLivermore = new
00278         G4LivermoreIonisationModel();
00279       theIoniLivermore->SetHighEnergyLimit(0.1*MeV); 
00280       eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
00281       eIoni->SetStepFunction(0.2, 100*um); //     
00282       ph->RegisterProcess(eIoni, particle);
00283       
00284       // Bremsstrahlung from standard
00285       G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
00286       ph->RegisterProcess(eBrem, particle);
00287 
00288     } else if (particleName == "e+") {
00289 
00290       // multiple scattering
00291       G4eMultipleScattering* msc = new G4eMultipleScattering;
00292       msc->SetStepLimitType(fUseDistanceToBoundary);
00293       G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
00294       G4WentzelVIModel* msc2 = new G4WentzelVIModel();
00295       msc1->SetHighEnergyLimit(highEnergyLimit);
00296       msc2->SetLowEnergyLimit(highEnergyLimit);
00297       msc->AddEmModel(0, msc1);
00298       msc->AddEmModel(0, msc2);
00299 
00300       G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 
00301       G4CoulombScattering* ss = new G4CoulombScattering();
00302       ss->SetEmModel(ssm, 1); 
00303       ss->SetMinKinEnergy(highEnergyLimit);
00304       ssm->SetLowEnergyLimit(highEnergyLimit);
00305       ssm->SetActivationLowEnergyLimit(highEnergyLimit);
00306 
00307       // Ionisation
00308       G4eIonisation* eIoni = new G4eIonisation();
00309       eIoni->SetStepFunction(0.2, 100*um);      
00310 
00311       ph->RegisterProcess(msc, particle);
00312       ph->RegisterProcess(eIoni, particle);
00313       ph->RegisterProcess(new G4eBremsstrahlung(), particle);
00314       ph->RegisterProcess(new G4eplusAnnihilation(), particle);
00315       ph->RegisterProcess(ss, particle);
00316 
00317     } else if (particleName == "mu+" ||
00318                particleName == "mu-"    ) {
00319 
00320       G4MuIonisation* muIoni = new G4MuIonisation();
00321       muIoni->SetStepFunction(0.2, 50*um);          
00322 
00323       ph->RegisterProcess(mumsc, particle);
00324       ph->RegisterProcess(muIoni, particle);
00325       ph->RegisterProcess(mub, particle);
00326       ph->RegisterProcess(mup, particle);
00327       ph->RegisterProcess(new G4CoulombScattering(), particle);
00328 
00329     } else if (particleName == "alpha" ||
00330                particleName == "He3" ) {
00331 
00332       // Identical to G4EmStandardPhysics_option3
00333       
00334       G4hMultipleScattering* msc = new G4hMultipleScattering();
00335       G4ionIonisation* ionIoni = new G4ionIonisation();
00336       ionIoni->SetStepFunction(0.1, 10*um);
00337 
00338       ph->RegisterProcess(msc, particle);
00339       ph->RegisterProcess(ionIoni, particle);
00340       ph->RegisterProcess(ionnuc, particle);
00341 
00342     } else if (particleName == "GenericIon") {
00343 
00344       // Identical to G4EmStandardPhysics_option3
00345       
00346       G4ionIonisation* ionIoni = new G4ionIonisation();
00347       ionIoni->SetEmModel(new G4IonParametrisedLossModel());
00348       ionIoni->SetStepFunction(0.1, 1*um);
00349 
00350       ph->RegisterProcess(hmsc, particle);
00351       ph->RegisterProcess(ionIoni, particle);
00352       ph->RegisterProcess(ionnuc, particle);
00353 
00354     } else if (particleName == "pi+" ||
00355                particleName == "pi-" ) {
00356 
00357       //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
00358       G4hIonisation* hIoni = new G4hIonisation();
00359       hIoni->SetStepFunction(0.2, 50*um);
00360 
00361       ph->RegisterProcess(pimsc, particle);
00362       ph->RegisterProcess(hIoni, particle);
00363       ph->RegisterProcess(pib, particle);
00364       ph->RegisterProcess(pip, particle);
00365 
00366     } else if (particleName == "kaon+" ||
00367                particleName == "kaon-" ) {
00368 
00369       //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
00370       G4hIonisation* hIoni = new G4hIonisation();
00371       hIoni->SetStepFunction(0.2, 50*um);
00372 
00373       ph->RegisterProcess(kmsc, particle);
00374       ph->RegisterProcess(hIoni, particle);
00375       ph->RegisterProcess(kb, particle);
00376       ph->RegisterProcess(kp, particle);
00377 
00378     } else if (particleName == "proton" ||
00379                particleName == "anti_proton") {
00380 
00381       //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
00382       G4hIonisation* hIoni = new G4hIonisation();
00383       hIoni->SetStepFunction(0.2, 50*um);
00384 
00385       ph->RegisterProcess(pmsc, particle);
00386       ph->RegisterProcess(hIoni, particle);
00387       ph->RegisterProcess(pb, particle);
00388       ph->RegisterProcess(pp, particle);
00389       ph->RegisterProcess(pnuc, particle);
00390 
00391     } else if (particleName == "B+" ||
00392                particleName == "B-" ||
00393                particleName == "D+" ||
00394                particleName == "D-" ||
00395                particleName == "Ds+" ||
00396                particleName == "Ds-" ||
00397                particleName == "anti_He3" ||
00398                particleName == "anti_alpha" ||
00399                particleName == "anti_deuteron" ||
00400                particleName == "anti_lambda_c+" ||
00401                particleName == "anti_omega-" ||
00402                particleName == "anti_sigma_c+" ||
00403                particleName == "anti_sigma_c++" ||
00404                particleName == "anti_sigma+" ||
00405                particleName == "anti_sigma-" ||
00406                particleName == "anti_triton" ||
00407                particleName == "anti_xi_c+" ||
00408                particleName == "anti_xi-" ||
00409                particleName == "deuteron" ||
00410                particleName == "lambda_c+" ||
00411                particleName == "omega-" ||
00412                particleName == "sigma_c+" ||
00413                particleName == "sigma_c++" ||
00414                particleName == "sigma+" ||
00415                particleName == "sigma-" ||
00416                particleName == "tau+" ||
00417                particleName == "tau-" ||
00418                particleName == "triton" ||
00419                particleName == "xi_c+" ||
00420                particleName == "xi-" ) {
00421 
00422       // Identical to G4EmStandardPhysics_option3
00423       
00424       ph->RegisterProcess(hmsc, particle);
00425       ph->RegisterProcess(new G4hIonisation(), particle);
00426       ph->RegisterProcess(pnuc, particle);
00427     }
00428   }
00429     
00430   // Em options
00431   //      
00432   G4EmProcessOptions opt;
00433   opt.SetVerbose(verbose);
00434   
00435   // Multiple Coulomb scattering
00436   //
00437   opt.SetPolarAngleLimit(CLHEP::pi);
00438     
00439   // Physics tables
00440   //
00441 
00442   opt.SetMinEnergy(100*eV);
00443   opt.SetMaxEnergy(10*TeV);
00444   opt.SetDEDXBinning(220);
00445   opt.SetLambdaBinning(220);
00446 
00447   // Nuclear stopping
00448   pnuc->SetMaxKinEnergy(MeV);
00449 
00450   // Ionization
00451   //
00452   //opt.SetSubCutoff(true);    
00453 
00454   // Deexcitation
00455   //
00456   G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
00457   G4LossTableManager::Instance()->SetAtomDeexcitation(de);
00458   de->SetFluo(true);
00459 }
00460 
00461 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

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