G4EmStandardPhysics_option4 Class Reference

#include <G4EmStandardPhysics_option4.hh>

Inheritance diagram for G4EmStandardPhysics_option4:

G4VPhysicsConstructor

Public Member Functions

 G4EmStandardPhysics_option4 (G4int ver=1)
 G4EmStandardPhysics_option4 (G4int ver, const G4String &name)
virtual ~G4EmStandardPhysics_option4 ()
virtual void ConstructParticle ()
virtual void ConstructProcess ()

Detailed Description

Definition at line 53 of file G4EmStandardPhysics_option4.hh.


Constructor & Destructor Documentation

G4EmStandardPhysics_option4::G4EmStandardPhysics_option4 ( G4int  ver = 1  ) 

Definition at line 123 of file G4EmStandardPhysics_option4.cc.

References bElectromagnetic, G4LossTableManager::Instance(), and G4VPhysicsConstructor::SetPhysicsType().

00124   : G4VPhysicsConstructor("G4EmStandard_opt4"), verbose(ver)
00125 {
00126   G4LossTableManager::Instance();
00127   SetPhysicsType(bElectromagnetic);
00128 }

G4EmStandardPhysics_option4::G4EmStandardPhysics_option4 ( G4int  ver,
const G4String name 
)

Definition at line 132 of file G4EmStandardPhysics_option4.cc.

References bElectromagnetic, G4LossTableManager::Instance(), and G4VPhysicsConstructor::SetPhysicsType().

00133   : G4VPhysicsConstructor("G4EmStandard_opt3"), verbose(ver)
00134 {
00135   G4LossTableManager::Instance();
00136   SetPhysicsType(bElectromagnetic);
00137 }

G4EmStandardPhysics_option4::~G4EmStandardPhysics_option4 (  )  [virtual]

Definition at line 141 of file G4EmStandardPhysics_option4.cc.

00142 {}


Member Function Documentation

void G4EmStandardPhysics_option4::ConstructParticle (  )  [virtual]

Implements G4VPhysicsConstructor.

Definition at line 146 of file G4EmStandardPhysics_option4.cc.

References G4Alpha::Alpha(), G4AntiProton::AntiProton(), G4Deuteron::Deuteron(), G4Electron::Electron(), G4Gamma::Gamma(), G4GenericIon::GenericIonDefinition(), G4He3::He3(), G4KaonMinus::KaonMinusDefinition(), G4KaonPlus::KaonPlusDefinition(), G4MuonMinus::MuonMinus(), G4MuonPlus::MuonPlus(), G4PionMinus::PionMinusDefinition(), G4PionPlus::PionPlusDefinition(), G4Positron::Positron(), G4Proton::Proton(), and G4Triton::Triton().

00147 {
00148 // gamma
00149   G4Gamma::Gamma();
00150 
00151 // leptons
00152   G4Electron::Electron();
00153   G4Positron::Positron();
00154   G4MuonPlus::MuonPlus();
00155   G4MuonMinus::MuonMinus();
00156 
00157 // mesons
00158   G4PionPlus::PionPlusDefinition();
00159   G4PionMinus::PionMinusDefinition();
00160   G4KaonPlus::KaonPlusDefinition();
00161   G4KaonMinus::KaonMinusDefinition();
00162 
00163 // barions
00164   G4Proton::Proton();
00165   G4AntiProton::AntiProton();
00166 
00167 // ions
00168   G4Deuteron::Deuteron();
00169   G4Triton::Triton();
00170   G4He3::He3();
00171   G4Alpha::Alpha();
00172   G4GenericIon::GenericIonDefinition();
00173 }

void G4EmStandardPhysics_option4::ConstructProcess (  )  [virtual]

Implements G4VPhysicsConstructor.

Definition at line 177 of file G4EmStandardPhysics_option4.cc.

References G4VEnergyLossProcess::AddEmModel(), G4VEmProcess::AddEmModel(), G4VMultipleScattering::AddEmModel(), fUseDistanceToBoundary, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), G4PhysicsListHelper::GetPhysicsListHelper(), G4VPhysicsConstructor::GetPhysicsName(), G4LossTableManager::Instance(), G4INCL::Math::pi, G4InuclParticleNames::pip, G4InuclParticleNames::pp, G4PhysicsListHelper::RegisterProcess(), G4ParticleTableIterator< K, V >::reset(), G4VEmModel::SetActivationLowEnergyLimit(), G4LossTableManager::SetAtomDeexcitation(), G4EmProcessOptions::SetDEDXBinning(), G4VEnergyLossProcess::SetEmModel(), G4VEmProcess::SetEmModel(), G4VAtomDeexcitation::SetFluo(), G4VEmModel::SetHighEnergyLimit(), G4EmProcessOptions::SetLambdaBinning(), G4VEmModel::SetLowEnergyLimit(), G4EmProcessOptions::SetMaxEnergy(), G4VEmProcess::SetMaxKinEnergy(), G4EmProcessOptions::SetMinEnergy(), G4VEmProcess::SetMinKinEnergy(), G4EmProcessOptions::SetPolarAngleLimit(), G4VEnergyLossProcess::SetStepFunction(), G4VMultipleScattering::SetStepLimitType(), G4EmProcessOptions::SetVerbose(), G4VPhysicsConstructor::theParticleIterator, and G4ParticleTableIterator< K, V >::value().

00178 {
00179   G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
00180 
00181   // muon & hadron bremsstrahlung and pair production
00182   G4MuBremsstrahlung* mub = new G4MuBremsstrahlung();
00183   G4MuPairProduction* mup = new G4MuPairProduction();
00184   G4hBremsstrahlung* pib = new G4hBremsstrahlung();
00185   G4hPairProduction* pip = new G4hPairProduction();
00186   G4hBremsstrahlung* kb = new G4hBremsstrahlung();
00187   G4hPairProduction* kp = new G4hPairProduction();
00188   G4hBremsstrahlung* pb = new G4hBremsstrahlung();
00189   G4hPairProduction* pp = new G4hPairProduction();
00190 
00191   // muon & hadron multiple scattering
00192   G4MuMultipleScattering* mumsc = new G4MuMultipleScattering();
00193   mumsc->AddEmModel(0, new G4WentzelVIModel());
00194   G4hMultipleScattering* pimsc = new G4hMultipleScattering();
00195   pimsc->AddEmModel(0, new G4WentzelVIModel());
00196   G4hMultipleScattering* kmsc = new G4hMultipleScattering();
00197   kmsc->AddEmModel(0, new G4WentzelVIModel());
00198   G4hMultipleScattering* pmsc = new G4hMultipleScattering();
00199   pmsc->AddEmModel(0, new G4WentzelVIModel());
00200   G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
00201 
00202   // high energy limit for e+- scattering models
00203   G4double highEnergyLimit = 100*MeV;
00204 
00205   // nuclear stopping
00206   G4NuclearStopping* ionnuc = new G4NuclearStopping();
00207   G4NuclearStopping* pnuc = new G4NuclearStopping();
00208 
00209   // Add standard EM Processes
00210   theParticleIterator->reset();
00211   while( (*theParticleIterator)() ){
00212     G4ParticleDefinition* particle = theParticleIterator->value();
00213     G4String particleName = particle->GetParticleName();
00214     if(verbose > 1)
00215       G4cout << "### " << GetPhysicsName() << " instantiates for " 
00216              << particleName << G4endl;
00217 
00218     if (particleName == "gamma") {
00219 
00220       // Compton scattering
00221       G4ComptonScattering* cs = new G4ComptonScattering;
00222       cs->SetEmModel(new G4KleinNishinaModel(),1);
00223       G4VEmModel* theLowEPComptonModel = new G4LivermoreComptonModel();
00224       //G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
00225       theLowEPComptonModel->SetHighEnergyLimit(2*MeV);
00226       cs->AddEmModel(0, theLowEPComptonModel);
00227       ph->RegisterProcess(cs, particle);
00228 
00229       // Photoelectric
00230       G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
00231       G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
00232       theLivermorePEModel->SetHighEnergyLimit(10*GeV);
00233       pe->SetEmModel(theLivermorePEModel,1);
00234       ph->RegisterProcess(pe, particle);
00235 
00236       // Gamma conversion
00237       G4GammaConversion* gc = new G4GammaConversion();
00238       G4VEmModel* thePenelopeGCModel = new G4PenelopeGammaConversionModel();
00239       thePenelopeGCModel->SetHighEnergyLimit(1*GeV);
00240       gc->SetEmModel(thePenelopeGCModel,1);
00241       ph->RegisterProcess(gc, particle);
00242 
00243       // Rayleigh scattering
00244       ph->RegisterProcess(new G4RayleighScattering(), particle);
00245  
00246     } else if (particleName == "e-") {
00247 
00248       // multiple scattering
00249       G4eMultipleScattering* msc = new G4eMultipleScattering;
00250       msc->SetStepLimitType(fUseDistanceToBoundary);
00251       G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
00252       G4WentzelVIModel* msc2 = new G4WentzelVIModel();
00253       msc1->SetHighEnergyLimit(highEnergyLimit);
00254       msc2->SetLowEnergyLimit(highEnergyLimit);
00255       msc->AddEmModel(0, msc1);
00256       msc->AddEmModel(0, msc2);
00257 
00258       G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 
00259       G4CoulombScattering* ss = new G4CoulombScattering();
00260       ss->SetEmModel(ssm, 1); 
00261       ss->SetMinKinEnergy(highEnergyLimit);
00262       ssm->SetLowEnergyLimit(highEnergyLimit);
00263       ssm->SetActivationLowEnergyLimit(highEnergyLimit);
00264 
00265       // ionisation
00266       G4eIonisation* eIoni = new G4eIonisation();
00267       eIoni->SetStepFunction(0.2, 100*um);
00268       G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
00269       theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
00270       eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
00271 
00272       // bremsstrahlung
00273       G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
00274 
00275       ph->RegisterProcess(msc, particle);
00276       ph->RegisterProcess(eIoni, particle);
00277       ph->RegisterProcess(eBrem, particle);
00278       ph->RegisterProcess(ss, particle);
00279 
00280     } else if (particleName == "e+") {
00281 
00282       // multiple scattering
00283       G4eMultipleScattering* msc = new G4eMultipleScattering;
00284       msc->SetStepLimitType(fUseDistanceToBoundary);
00285       G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
00286       G4WentzelVIModel* msc2 = new G4WentzelVIModel();
00287       msc1->SetHighEnergyLimit(highEnergyLimit);
00288       msc2->SetLowEnergyLimit(highEnergyLimit);
00289       msc->AddEmModel(0, msc1);
00290       msc->AddEmModel(0, msc2);
00291 
00292       G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 
00293       G4CoulombScattering* ss = new G4CoulombScattering();
00294       ss->SetEmModel(ssm, 1); 
00295       ss->SetMinKinEnergy(highEnergyLimit);
00296       ssm->SetLowEnergyLimit(highEnergyLimit);
00297       ssm->SetActivationLowEnergyLimit(highEnergyLimit);
00298 
00299       // ionisation
00300       G4eIonisation* eIoni = new G4eIonisation();
00301       eIoni->SetStepFunction(0.2, 100*um);
00302       G4VEmModel* theIoniPenelope = new G4PenelopeIonisationModel();
00303       theIoniPenelope->SetHighEnergyLimit(0.1*MeV);
00304       eIoni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
00305 
00306       // bremsstrahlung
00307       G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
00308 
00309       ph->RegisterProcess(msc, particle);
00310       ph->RegisterProcess(eIoni, particle);
00311       ph->RegisterProcess(eBrem, particle);
00312       ph->RegisterProcess(ss, particle);
00313 
00314       // annihilation at rest and in flight
00315       ph->RegisterProcess(new G4eplusAnnihilation(), 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       G4hMultipleScattering* msc = new G4hMultipleScattering();
00333       G4ionIonisation* ionIoni = new G4ionIonisation();
00334       ionIoni->SetStepFunction(0.1, 10*um);
00335 
00336       ph->RegisterProcess(msc, particle);
00337       ph->RegisterProcess(ionIoni, particle);
00338       ph->RegisterProcess(ionnuc, particle);
00339 
00340     } else if (particleName == "GenericIon") {
00341 
00342       G4ionIonisation* ionIoni = new G4ionIonisation();
00343       ionIoni->SetEmModel(new G4IonParametrisedLossModel());
00344       ionIoni->SetStepFunction(0.1, 1*um);
00345 
00346       ph->RegisterProcess(hmsc, particle);
00347       ph->RegisterProcess(ionIoni, particle);
00348       ph->RegisterProcess(ionnuc, particle);
00349 
00350     } else if (particleName == "pi+" ||
00351                particleName == "pi-" ) {
00352 
00353       //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
00354       G4hIonisation* hIoni = new G4hIonisation();
00355       hIoni->SetStepFunction(0.2, 50*um);
00356 
00357       ph->RegisterProcess(pimsc, particle);
00358       ph->RegisterProcess(hIoni, particle);
00359       ph->RegisterProcess(pib, particle);
00360       ph->RegisterProcess(pip, particle);
00361 
00362     } else if (particleName == "kaon+" ||
00363                particleName == "kaon-" ) {
00364 
00365       //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
00366       G4hIonisation* hIoni = new G4hIonisation();
00367       hIoni->SetStepFunction(0.2, 50*um);
00368 
00369       ph->RegisterProcess(kmsc, particle);
00370       ph->RegisterProcess(hIoni, particle);
00371       ph->RegisterProcess(kb, particle);
00372       ph->RegisterProcess(kp, particle);
00373 
00374     } else if (particleName == "proton" ||
00375                particleName == "anti_proton") {
00376 
00377       //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
00378       G4hIonisation* hIoni = new G4hIonisation();
00379       hIoni->SetStepFunction(0.2, 50*um);
00380 
00381       ph->RegisterProcess(pmsc, particle);
00382       ph->RegisterProcess(hIoni, particle);
00383       ph->RegisterProcess(pb, particle);
00384       ph->RegisterProcess(pp, particle);
00385       ph->RegisterProcess(pnuc, particle);
00386 
00387     } else if (particleName == "B+" ||
00388                particleName == "B-" ||
00389                particleName == "D+" ||
00390                particleName == "D-" ||
00391                particleName == "Ds+" ||
00392                particleName == "Ds-" ||
00393                particleName == "anti_He3" ||
00394                particleName == "anti_alpha" ||
00395                particleName == "anti_deuteron" ||
00396                particleName == "anti_lambda_c+" ||
00397                particleName == "anti_omega-" ||
00398                particleName == "anti_sigma_c+" ||
00399                particleName == "anti_sigma_c++" ||
00400                particleName == "anti_sigma+" ||
00401                particleName == "anti_sigma-" ||
00402                particleName == "anti_triton" ||
00403                particleName == "anti_xi_c+" ||
00404                particleName == "anti_xi-" ||
00405                particleName == "deuteron" ||
00406                particleName == "lambda_c+" ||
00407                particleName == "omega-" ||
00408                particleName == "sigma_c+" ||
00409                particleName == "sigma_c++" ||
00410                particleName == "sigma+" ||
00411                particleName == "sigma-" ||
00412                particleName == "tau+" ||
00413                particleName == "tau-" ||
00414                particleName == "triton" ||
00415                particleName == "xi_c+" ||
00416                particleName == "xi-" ) {
00417 
00418       ph->RegisterProcess(hmsc, particle);
00419       ph->RegisterProcess(new G4hIonisation(), particle);
00420       ph->RegisterProcess(pnuc, particle);
00421     }
00422   }
00423     
00424   // Em options
00425   //      
00426   G4EmProcessOptions opt;
00427   opt.SetVerbose(verbose);
00428   
00429   // Multiple Coulomb scattering
00430   //
00431   opt.SetPolarAngleLimit(CLHEP::pi);
00432     
00433   // Physics tables
00434   //
00435   opt.SetMinEnergy(100*eV);
00436   opt.SetMaxEnergy(10*TeV);
00437   opt.SetDEDXBinning(220);
00438   opt.SetLambdaBinning(220);
00439 
00440   // Nuclear stopping
00441   pnuc->SetMaxKinEnergy(MeV);
00442     
00443   // Ionization
00444   //
00445   //opt.SetSubCutoff(true);    
00446 
00447   // Deexcitation
00448   G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
00449   G4LossTableManager::Instance()->SetAtomDeexcitation(de);
00450   de->SetFluo(true);
00451 }


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