G4EmPenelopePhysics Class Reference

#include <G4EmPenelopePhysics.hh>

Inheritance diagram for G4EmPenelopePhysics:

G4VPhysicsConstructor

Public Member Functions

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

Detailed Description

Definition at line 36 of file G4EmPenelopePhysics.hh.


Constructor & Destructor Documentation

G4EmPenelopePhysics::G4EmPenelopePhysics ( G4int  ver = 1  ) 

Definition at line 134 of file G4EmPenelopePhysics.cc.

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

00135   : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
00136 {
00137   G4LossTableManager::Instance();
00138   SetPhysicsType(bElectromagnetic);
00139 }

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

Definition at line 143 of file G4EmPenelopePhysics.cc.

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

00144   : G4VPhysicsConstructor("G4EmPenelopePhysics"), verbose(ver)
00145 {
00146   G4LossTableManager::Instance();
00147   SetPhysicsType(bElectromagnetic);
00148 }

G4EmPenelopePhysics::~G4EmPenelopePhysics (  )  [virtual]

Definition at line 152 of file G4EmPenelopePhysics.cc.

00153 {}


Member Function Documentation

void G4EmPenelopePhysics::ConstructParticle (  )  [virtual]

Implements G4VPhysicsConstructor.

Definition at line 157 of file G4EmPenelopePhysics.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().

00158 {
00159 // gamma
00160   G4Gamma::Gamma();
00161 
00162 // leptons
00163   G4Electron::Electron();
00164   G4Positron::Positron();
00165   G4MuonPlus::MuonPlus();
00166   G4MuonMinus::MuonMinus();
00167 
00168 // mesons
00169   G4PionPlus::PionPlusDefinition();
00170   G4PionMinus::PionMinusDefinition();
00171   G4KaonPlus::KaonPlusDefinition();
00172   G4KaonMinus::KaonMinusDefinition();
00173 
00174 // baryons
00175   G4Proton::Proton();
00176   G4AntiProton::AntiProton();
00177 
00178 // ions
00179   G4Deuteron::Deuteron();
00180   G4Triton::Triton();
00181   G4He3::He3();
00182   G4Alpha::Alpha();
00183   G4GenericIon::GenericIonDefinition();
00184 }

void G4EmPenelopePhysics::ConstructProcess (  )  [virtual]

Implements G4VPhysicsConstructor.

Definition at line 188 of file G4EmPenelopePhysics.cc.

References G4VEmProcess::AddEmModel(), G4VEnergyLossProcess::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().

00189 {
00190   G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
00191 
00192   // muon & hadron bremsstrahlung and pair production
00193   G4MuBremsstrahlung* mub = new G4MuBremsstrahlung();
00194   G4MuPairProduction* mup = new G4MuPairProduction();
00195   G4hBremsstrahlung* pib = new G4hBremsstrahlung();
00196   G4hPairProduction* pip = new G4hPairProduction();
00197   G4hBremsstrahlung* kb = new G4hBremsstrahlung();
00198   G4hPairProduction* kp = new G4hPairProduction();
00199   G4hBremsstrahlung* pb = new G4hBremsstrahlung();
00200   G4hPairProduction* pp = new G4hPairProduction();
00201 
00202   // muon & hadron multiple scattering
00203   G4MuMultipleScattering* mumsc = new G4MuMultipleScattering();
00204   mumsc->AddEmModel(0, new G4WentzelVIModel());
00205   G4MuMultipleScattering* pimsc = new G4MuMultipleScattering();
00206   pimsc->AddEmModel(0, new G4WentzelVIModel());
00207   G4MuMultipleScattering* kmsc = new G4MuMultipleScattering();
00208   kmsc->AddEmModel(0, new G4WentzelVIModel());
00209   G4MuMultipleScattering* pmsc = new G4MuMultipleScattering();
00210   pmsc->AddEmModel(0, new G4WentzelVIModel());
00211   G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
00212 
00213   // high energy limit for e+- scattering models
00214   G4double highEnergyLimit = 100*MeV;
00215 
00216   // nuclear stopping
00217   G4NuclearStopping* ionnuc = new G4NuclearStopping();
00218   G4NuclearStopping* pnuc = new G4NuclearStopping();
00219 
00220   // Add Penelope EM Processes
00221   theParticleIterator->reset();
00222 
00223   while( (*theParticleIterator)() ){
00224   
00225     G4ParticleDefinition* particle = theParticleIterator->value();
00226     G4String particleName = particle->GetParticleName();
00227     
00228     if(verbose > 1)
00229       G4cout << "### " << GetPhysicsName() << " instantiates for " 
00230              << particleName << G4endl;
00231 
00232     //Applicability range for Penelope models
00233     //for higher energies, the Standard models are used   
00234     G4double PenelopeHighEnergyLimit = 1.0*GeV;
00235 
00236     if (particleName == "gamma") {
00237 
00238       //Photo-electric effect
00239       G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
00240       G4PenelopePhotoElectricModel* thePEPenelopeModel = new 
00241         G4PenelopePhotoElectricModel();   
00242       thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
00243       thePhotoElectricEffect->SetEmModel(thePEPenelopeModel, 1);
00244       ph->RegisterProcess(thePhotoElectricEffect, particle);
00245 
00246       //Compton scattering
00247       G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
00248       G4PenelopeComptonModel* theComptonPenelopeModel = 
00249         new G4PenelopeComptonModel();
00250       theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
00251       theComptonScattering->SetEmModel(theComptonPenelopeModel, 1);
00252       ph->RegisterProcess(theComptonScattering, particle);
00253 
00254       //Gamma conversion
00255       G4GammaConversion* theGammaConversion = new G4GammaConversion();
00256       G4PenelopeGammaConversionModel* theGCPenelopeModel = 
00257         new G4PenelopeGammaConversionModel();
00258       theGammaConversion->SetEmModel(theGCPenelopeModel,1);
00259       ph->RegisterProcess(theGammaConversion, particle);
00260 
00261       //Rayleigh scattering
00262       G4RayleighScattering* theRayleigh = new G4RayleighScattering();
00263       G4PenelopeRayleighModel* theRayleighPenelopeModel = 
00264         new G4PenelopeRayleighModel();
00265       //theRayleighPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
00266       theRayleigh->SetEmModel(theRayleighPenelopeModel, 1);
00267       ph->RegisterProcess(theRayleigh, particle);
00268 
00269     } else if (particleName == "e-") {
00270 
00271       // multiple scattering
00272       G4eMultipleScattering* msc = new G4eMultipleScattering;
00273       msc->SetStepLimitType(fUseDistanceToBoundary);
00274       G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
00275       G4WentzelVIModel* msc2 = new G4WentzelVIModel();
00276       msc1->SetHighEnergyLimit(highEnergyLimit);
00277       msc2->SetLowEnergyLimit(highEnergyLimit);
00278       msc->AddEmModel(0, msc1);
00279       msc->AddEmModel(0, msc2);
00280 
00281       G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 
00282       G4CoulombScattering* ss = new G4CoulombScattering();
00283       ss->SetEmModel(ssm, 1); 
00284       ss->SetMinKinEnergy(highEnergyLimit);
00285       ssm->SetLowEnergyLimit(highEnergyLimit);
00286       ssm->SetActivationLowEnergyLimit(highEnergyLimit);
00287       ph->RegisterProcess(msc, particle);
00288       ph->RegisterProcess(ss, particle);
00289       
00290       //Ionisation
00291       G4eIonisation* eIoni = new G4eIonisation();
00292       G4PenelopeIonisationModel* theIoniPenelope = 
00293         new G4PenelopeIonisationModel();
00294       theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);     
00295       eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
00296       eIoni->SetStepFunction(0.2, 100*um); //     
00297       ph->RegisterProcess(eIoni, particle);
00298       
00299       //Bremsstrahlung
00300       G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
00301       G4PenelopeBremsstrahlungModel* theBremPenelope = new 
00302         G4PenelopeBremsstrahlungModel();
00303       theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
00304       eBrem->AddEmModel(0,theBremPenelope);
00305       ph->RegisterProcess(eBrem, particle);
00306 
00307     } else if (particleName == "e+") {
00308     
00309       // multiple scattering
00310       G4eMultipleScattering* msc = new G4eMultipleScattering;
00311       msc->SetStepLimitType(fUseDistanceToBoundary);
00312       G4UrbanMscModel95* msc1 = new G4UrbanMscModel95();
00313       G4WentzelVIModel* msc2 = new G4WentzelVIModel();
00314       msc1->SetHighEnergyLimit(highEnergyLimit);
00315       msc2->SetLowEnergyLimit(highEnergyLimit);
00316       msc->AddEmModel(0, msc1);
00317       msc->AddEmModel(0, msc2);
00318 
00319       G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel(); 
00320       G4CoulombScattering* ss = new G4CoulombScattering();
00321       ss->SetEmModel(ssm, 1); 
00322       ss->SetMinKinEnergy(highEnergyLimit);
00323       ssm->SetLowEnergyLimit(highEnergyLimit);
00324       ssm->SetActivationLowEnergyLimit(highEnergyLimit);
00325       ph->RegisterProcess(msc, particle);
00326       ph->RegisterProcess(ss, particle);
00327       
00328       //Ionisation
00329       G4eIonisation* eIoni = new G4eIonisation();
00330       G4PenelopeIonisationModel* theIoniPenelope = 
00331         new G4PenelopeIonisationModel();
00332       theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
00333       eIoni->AddEmModel(0,theIoniPenelope,new G4UniversalFluctuation());
00334       eIoni->SetStepFunction(0.2, 100*um); //     
00335       ph->RegisterProcess(eIoni, particle);
00336 
00337        //Bremsstrahlung
00338       G4eBremsstrahlung* eBrem = new G4eBremsstrahlung();
00339       G4PenelopeBremsstrahlungModel* theBremPenelope = new 
00340         G4PenelopeBremsstrahlungModel();
00341       theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
00342       eBrem->AddEmModel(0,theBremPenelope);
00343       ph->RegisterProcess(eBrem, particle);
00344       
00345       //Annihilation
00346       G4eplusAnnihilation* eAnni = new G4eplusAnnihilation();
00347       G4PenelopeAnnihilationModel* theAnnPenelope = new 
00348         G4PenelopeAnnihilationModel();
00349       theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
00350       eAnni->AddEmModel(0,theAnnPenelope);
00351       ph->RegisterProcess(eAnni, particle);
00352 
00353     } else if (particleName == "mu+" ||
00354                particleName == "mu-"    ) {
00355 
00356       // Identical to G4EmStandardPhysics_option3
00357 
00358       G4MuIonisation* muIoni = new G4MuIonisation();
00359       muIoni->SetStepFunction(0.2, 50*um);          
00360 
00361       ph->RegisterProcess(mumsc, particle);
00362       ph->RegisterProcess(muIoni, particle);
00363       ph->RegisterProcess(mub, particle);
00364       ph->RegisterProcess(mup, particle);
00365       ph->RegisterProcess(new G4CoulombScattering(), particle);
00366 
00367     } else if (particleName == "alpha" ||
00368                particleName == "He3" ) {
00369 
00370       // Identical to G4EmStandardPhysics_option3
00371       
00372       G4hMultipleScattering* msc = new G4hMultipleScattering();
00373       G4ionIonisation* ionIoni = new G4ionIonisation();
00374       ionIoni->SetStepFunction(0.1, 10*um);
00375 
00376       ph->RegisterProcess(msc, particle);
00377       ph->RegisterProcess(ionIoni, particle);
00378       ph->RegisterProcess(ionnuc, particle);
00379 
00380     } else if (particleName == "GenericIon") {
00381 
00382       // Identical to G4EmStandardPhysics_option3
00383       
00384       G4ionIonisation* ionIoni = new G4ionIonisation();
00385       ionIoni->SetEmModel(new G4IonParametrisedLossModel());
00386       ionIoni->SetStepFunction(0.1, 1*um);
00387 
00388       ph->RegisterProcess(hmsc, particle);
00389       ph->RegisterProcess(ionIoni, particle);
00390       ph->RegisterProcess(ionnuc, particle);
00391 
00392     } else if (particleName == "pi+" ||
00393                particleName == "pi-" ) {
00394 
00395       G4hIonisation* hIoni = new G4hIonisation();
00396       hIoni->SetStepFunction(0.2, 50*um);
00397 
00398       ph->RegisterProcess(pimsc, particle);
00399       ph->RegisterProcess(hIoni, particle);
00400       ph->RegisterProcess(pib, particle);
00401       ph->RegisterProcess(pip, particle);
00402 
00403     } else if (particleName == "kaon+" ||
00404                particleName == "kaon-" ) {
00405 
00406       G4hIonisation* hIoni = new G4hIonisation();
00407       hIoni->SetStepFunction(0.2, 50*um);
00408 
00409       ph->RegisterProcess(kmsc, particle);
00410       ph->RegisterProcess(hIoni, particle);
00411       ph->RegisterProcess(kb, particle);
00412       ph->RegisterProcess(kp, particle);
00413 
00414     } else if (particleName == "proton" ||
00415                particleName == "anti_proton") {
00416 
00417       G4hIonisation* hIoni = new G4hIonisation();
00418       hIoni->SetStepFunction(0.2, 50*um);
00419 
00420       ph->RegisterProcess(pmsc, particle);
00421       ph->RegisterProcess(hIoni, particle);
00422       ph->RegisterProcess(pb, particle);
00423       ph->RegisterProcess(pp, particle);
00424       ph->RegisterProcess(pnuc, particle);
00425 
00426     } else if (particleName == "B+" ||
00427                particleName == "B-" ||
00428                particleName == "D+" ||
00429                particleName == "D-" ||
00430                particleName == "Ds+" ||
00431                particleName == "Ds-" ||
00432                particleName == "anti_He3" ||
00433                particleName == "anti_alpha" ||
00434                particleName == "anti_deuteron" ||
00435                particleName == "anti_lambda_c+" ||
00436                particleName == "anti_omega-" ||
00437                particleName == "anti_sigma_c+" ||
00438                particleName == "anti_sigma_c++" ||
00439                particleName == "anti_sigma+" ||
00440                particleName == "anti_sigma-" ||
00441                particleName == "anti_triton" ||
00442                particleName == "anti_xi_c+" ||
00443                particleName == "anti_xi-" ||
00444                particleName == "deuteron" ||
00445                particleName == "lambda_c+" ||
00446                particleName == "omega-" ||
00447                particleName == "sigma_c+" ||
00448                particleName == "sigma_c++" ||
00449                particleName == "sigma+" ||
00450                particleName == "sigma-" ||
00451                particleName == "tau+" ||
00452                particleName == "tau-" ||
00453                particleName == "triton" ||
00454                particleName == "xi_c+" ||
00455                particleName == "xi-" ) {
00456 
00457       // Identical to G4EmStandardPhysics_option3
00458       
00459       ph->RegisterProcess(hmsc, particle);
00460       ph->RegisterProcess(new G4hIonisation(), particle);
00461       ph->RegisterProcess(pnuc, particle);
00462     }
00463   }
00464     
00465   // Em options
00466   //      
00467   G4EmProcessOptions opt;
00468   opt.SetVerbose(verbose);
00469   
00470   // Multiple Coulomb scattering
00471   //
00472   //opt.SetMscStepLimitation(fUseDistanceToBoundary);
00473   //opt.SetMscRangeFactor(0.02);
00474     
00475   // Physics tables
00476   //
00477 
00478   opt.SetMinEnergy(100*eV);
00479   opt.SetMaxEnergy(10*TeV);
00480   opt.SetDEDXBinning(220);
00481   opt.SetLambdaBinning(220);
00482 
00483   // Nuclear stopping
00484   pnuc->SetMaxKinEnergy(MeV);
00485 
00486   //opt.SetSplineFlag(true);
00487   opt.SetPolarAngleLimit(CLHEP::pi);
00488     
00489   // Ionization
00490   //
00491   //opt.SetSubCutoff(true);    
00492 
00493   
00494   // Deexcitation
00495   //
00496   G4VAtomDeexcitation* deexcitation = new G4UAtomicDeexcitation();
00497   G4LossTableManager::Instance()->SetAtomDeexcitation(deexcitation);
00498   deexcitation->SetFluo(true); 
00499 }


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