CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/SimG4Core/GFlash/plugins/GFlash.cc

Go to the documentation of this file.
00001 #include "SimG4Core/GFlash/interface/GFlash.h"
00002 #include "SimG4Core/GFlash/interface/ParametrisedPhysics.h"
00003 #include "SimG4Core/GFlash/interface/HadronPhysicsQGSP_WP.h"
00004 #include "SimG4Core/GFlash/interface/HadronPhysicsQGSP_BERT_WP.h"
00005 #include "SimG4Core/GFlash/interface/HadronPhysicsQGSPCMS_FTFP_BERT_WP.h"
00006 #include "SimG4Core/PhysicsLists/interface/CMSEmStandardPhysics95msc93.h"
00007 #include "SimG4Core/PhysicsLists/interface/CMSMonopolePhysics.h"
00008 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00009 
00010 #include "G4DecayPhysics.hh"
00011 #include "G4EmExtraPhysics.hh"
00012 #include "G4IonPhysics.hh"
00013 #include "G4StoppingPhysics.hh"
00014 #include "G4HadronElasticPhysics.hh" 
00015 #include "G4NeutronTrackingCut.hh"
00016 #include "HadronPhysicsQGSP_FTFP_BERT.hh"
00017 
00018 #include "G4DataQuestionaire.hh"
00019 #include "SimGeneral/GFlash/interface/GflashHistogram.h"
00020 
00021 #include <string>
00022 
00023 GFlash::GFlash(G4LogicalVolumeToDDLogicalPartMap& map, 
00024                const HepPDT::ParticleDataTable * table_, 
00025                sim::FieldBuilder *fieldBuilder_, 
00026                const edm::ParameterSet & p) : PhysicsList(map, table_, fieldBuilder_, p), 
00027                                               thePar(p.getParameter<edm::ParameterSet>("GFlash")) {
00028 
00029   G4DataQuestionaire it(photon);
00030 
00031   //std::string hadronPhysics = thePar.getParameter<std::string>("GflashHadronPhysics");
00032 
00033   int  ver     = p.getUntrackedParameter<int>("Verbosity",0);
00034   bool emPhys  = p.getUntrackedParameter<bool>("EMPhysics",true);
00035   bool hadPhys = p.getUntrackedParameter<bool>("HadPhysics",true);
00036   bool tracking= p.getParameter<bool>("TrackingCut");
00037   std::string region = p.getParameter<std::string>("Region");
00038   bool gem  = thePar.getParameter<bool>("GflashEMShowerModel");
00039   bool ghad = thePar.getParameter<bool>("GflashHadronShowerModel");
00040 
00041   edm::LogInfo("PhysicsList") << "You are using the simulation engine: "
00042                               << " + CMS GFLASH with Flags for EM Physics "
00043                               << gem << ", for Hadronic Physics "
00044                               << ghad 
00045                               << " and tracking cut " << tracking
00046                               << " with special region " << region;
00047 
00048   RegisterPhysics(new ParametrisedPhysics("parametrised",thePar)); 
00049 
00050   if (emPhys) {
00051     // EM Physics
00052     RegisterPhysics( new CMSEmStandardPhysics95msc93("EM standard msc93",ver,region));
00053 
00054     // Synchroton Radiation & GN Physics
00055     RegisterPhysics( new G4EmExtraPhysics(ver));
00056   }
00057 
00058   // Decays
00059   RegisterPhysics( new G4DecayPhysics(ver) );
00060 
00061   if (hadPhys) {
00062     // Hadron Elastic scattering
00063     RegisterPhysics( new G4HadronElasticPhysics(ver));
00064 
00065     // Hadron Physics
00066     RegisterPhysics( new HadronPhysicsQGSP_FTFP_BERT(ver));   
00067     // Stopping Physics
00068     RegisterPhysics( new G4StoppingPhysics(ver));
00069 
00070     // Ion Physics
00071     RegisterPhysics( new G4IonPhysics(ver));
00072 
00073     // Neutron tracking cut
00074     if (tracking) {
00075       RegisterPhysics( new G4NeutronTrackingCut(ver));
00076     }
00077     /*
00078     if(hadronPhysics=="QGSP_FTFP_BERT") {
00079       RegisterPhysics( new HadronPhysicsQGSPCMS_FTFP_BERT_WP("hadron",quasiElastic)); 
00080     }
00081     else if(hadronPhysics=="QGSP_BERT") {
00082       RegisterPhysics( new HadronPhysicsQGSP_BERT_WP("hadron",quasiElastic));
00083     }
00084     else if (hadronPhysics=="QGSP") {
00085       RegisterPhysics( new HadronPhysicsQGSP_WP("hadron",quasiElastic));
00086     }
00087     else {
00088       edm::LogInfo("PhysicsList") << hadronPhysics << " is not available for GflashHadronPhysics!"
00089                                   << "... Using QGSP_FTFP_BERT\n";
00090       RegisterPhysics( new HadronPhysicsQGSPCMS_FTFP_BERT_WP("hadron",quasiElastic));
00091     }
00092     // Stopping Physics
00093     RegisterPhysics( new G4QStoppingPhysics("stopping"));
00094 
00095     // Ion Physics
00096     RegisterPhysics( new G4IonPhysics("ion"));
00097 
00098     // Neutron tracking cut
00099     if (tracking) 
00100       RegisterPhysics( new G4NeutronTrackingCut("Neutron tracking cut", ver));
00101     */
00102   }
00103 
00104   // Monopoles
00105   RegisterPhysics( new CMSMonopolePhysics(table_,fieldBuilder_,p));
00106 
00107   // singleton histogram object
00108   if(thePar.getParameter<bool>("GflashHistogram")) {
00109     theHisto = GflashHistogram::instance();
00110     theHisto->setStoreFlag(true);
00111     theHisto->bookHistogram(thePar.getParameter<std::string>("GflashHistogramName"));
00112   }
00113 
00114 }
00115 
00116 GFlash::~GFlash() {
00117 
00118   if(thePar.getParameter<bool>("GflashHistogram")) {
00119     if(theHisto) delete theHisto;
00120   }
00121 
00122 }
00123 
00124 //define this as a plug-in
00125 #include "FWCore/Framework/interface/MakerMacros.h"
00126 #include "SimG4Core/Physics/interface/PhysicsListFactory.h"
00127 
00128 
00129 DEFINE_PHYSICSLIST(GFlash);
00130