CMS 3D CMS Logo

CustomPhysicsList.cc

Go to the documentation of this file.
00001 #include "SimG4Core/CustomPhysics/interface/CustomPhysicsList.h"
00002 #include "SimG4Core/CustomPhysics/interface/CustomParticleFactory.h"
00003 #include "SimG4Core/CustomPhysics/interface/DummyChargeFlipProcess.h"
00004 
00005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00006 #include "FWCore/ParameterSet/interface/FileInPath.h"
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 #include "G4Decay.hh"
00010 #include "G4MultipleScattering.hh"
00011 #include "G4hIonisation.hh"
00012 #include "G4ProcessManager.hh"
00013 
00014 
00015 #include "G4LeptonConstructor.hh"
00016 #include "G4MesonConstructor.hh"
00017 #include "G4BaryonConstructor.hh"
00018 #include "G4ShortLivedConstructor.hh"
00019 #include "G4IonConstructor.hh"
00020 
00021 #include "SimG4Core/CustomPhysics/interface/FullModelHadronicProcess.hh"
00022 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
00023  
00024 
00025 CustomPhysicsList::CustomPhysicsList(std::string name, const edm::ParameterSet & p)  :  G4VPhysicsConstructor(name)
00026 {
00027   
00028   myConfig = p;
00029   edm::FileInPath fp = p.getParameter<edm::FileInPath>("particlesDef");
00030   particleDefFilePath = fp.fullPath();
00031   edm::LogInfo("CustomPhysics")<<"Path for custom particle definition file: "<<particleDefFilePath<<std::endl;
00032   myHelper = 0;
00033   
00034  }
00035 
00036 CustomPhysicsList::~CustomPhysicsList() {
00037   delete myHelper;
00038 }
00039  
00040 void CustomPhysicsList::ConstructParticle(){
00041   CustomParticleFactory::loadCustomParticles(particleDefFilePath);     
00042 }
00043  
00044 void CustomPhysicsList::ConstructProcess() {
00045   addCustomPhysics();
00046 }
00047  
00048 void CustomPhysicsList::addCustomPhysics(){
00049     LogDebug("CustomPhysics") << " CustomPhysics: adding CustomPhysics processes  " <<std::endl;
00050     theParticleIterator->reset();
00051 
00052   while((*theParticleIterator)())    {
00053       int i = 0;
00054       G4ParticleDefinition* particle = theParticleIterator->value();
00055       CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
00056       if(CustomParticleFactory::isCustomParticle(particle))
00057         {
00058           LogDebug("CustomPhysics") << particle->GetParticleName()<<", "<<particle->GetPDGEncoding()
00059                        << " is Custom. Mass is "
00060                        <<particle->GetPDGMass()/GeV
00061                        <<" GeV."<<G4endl;
00062           if(cp->GetCloud()!=0)
00063             {
00064               LogDebug("CustomPhysics")<<"Cloud mass is "
00065                           <<cp->GetCloud()->GetPDGMass()/GeV
00066                           <<" GeV. Spectator mass is "
00067                           <<static_cast<CustomParticle*>(particle)->GetSpectator()->GetPDGMass()/GeV
00068                           <<" GeV."       << G4endl;
00069             }
00070           G4ProcessManager* pmanager = particle->GetProcessManager();
00071           if(pmanager)
00072             {
00073               if(cp!=0) {
00074                 if(particle->GetParticleType()=="rhadron" || particle->GetParticleType()=="mesonino" || particle->GetParticleType() == "sbaryon"){
00075                   if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
00076                   pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper)); //GHEISHA
00077                 }
00078               }
00079               if(particle->GetPDGCharge()/eplus != 0)
00080                 {
00081                   pmanager->AddProcess(new G4MultipleScattering,-1, 1,i+1);
00082                   pmanager->AddProcess(new G4hIonisation,       -1, 2,i+2);
00083                 }
00084             }
00085           else      LogDebug("CustomPhysics") << "   No pmanager" << G4endl;
00086         }
00087 
00088 
00089 
00090 
00091 
00092         
00093     }
00094 }
00095 
00096 
00097 void CustomPhysicsList::setupRHadronPhycis(G4ParticleDefinition* particle){
00098 
00099   //    LogDebug("CustomPhysics")<<"Configuring rHadron: "
00100   //    <<cp->
00101 
00102   CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
00103   if(cp->GetCloud()!=0) 
00104     LogDebug("CustomPhysics")<<"Cloud mass is "
00105                 <<cp->GetCloud()->GetPDGMass()/GeV
00106                 <<" GeV. Spectator mass is "
00107                 <<static_cast<CustomParticle*>(particle)->GetSpectator()->GetPDGMass()/GeV
00108                 <<" GeV."       << G4endl;
00109   
00110   G4ProcessManager* pmanager = particle->GetProcessManager();
00111   if(pmanager){
00112     if(!myHelper) myHelper = new G4ProcessHelper(myConfig);
00113     pmanager->AddDiscreteProcess(new FullModelHadronicProcess(myHelper)); //GHEISHA
00114     if(particle->GetPDGCharge()/eplus != 0){
00115       pmanager->AddProcess(new G4MultipleScattering,-1, 1,1);
00116       pmanager->AddProcess(new G4hIonisation,       -1, 2,2);
00117     }
00118   }
00119   else      LogDebug("CustomPhysics") << "   No pmanager" << G4endl;
00120 }
00121                                                
00122 
00123 void CustomPhysicsList::setupSUSYPhycis(G4ParticleDefinition* particle){
00124 
00125   CustomParticle* cp = dynamic_cast<CustomParticle*>(particle);
00126   G4ProcessManager* pmanager = particle->GetProcessManager();
00127   if(pmanager){
00128     pmanager->AddProcess(new G4Decay,1, 1,1);
00129     if(particle->GetPDGCharge()/eplus != 0){
00130       pmanager->AddProcess(new G4MultipleScattering,-1, 2,2);
00131       pmanager->AddProcess(new G4hIonisation,       -1, 3,3);
00132     }
00133   }
00134   else      LogDebug("CustomPhysics") << "   No pmanager" << G4endl;
00135 }

Generated on Tue Jun 9 17:47:00 2009 for CMSSW by  doxygen 1.5.4