CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_1/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/HadronPhysicsQGSPCMS_FTFP_BERT.h"
00007 #include "SimG4Core/PhysicsLists/interface/CMSEmStandardPhysics92.h"
00008 #include "SimG4Core/PhysicsLists/interface/CMSMonopolePhysics.h"
00009 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00010 
00011 #include "G4DecayPhysics.hh"
00012 #include "G4EmExtraPhysics.hh"
00013 #include "G4IonPhysics.hh"
00014 #include "G4QStoppingPhysics.hh"
00015 #include "G4HadronElasticPhysics.hh" 
00016 #include "G4NeutronTrackingCut.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 CMSEmStandardPhysics92("standard EM EML",ver,region));
00053 
00054     // Synchroton Radiation & GN Physics
00055     RegisterPhysics( new G4EmExtraPhysics("extra EM"));
00056   }
00057 
00058   // Decays
00059   RegisterPhysics( new G4DecayPhysics("decay",ver) );
00060 
00061   if (hadPhys) {
00062     // Hadron Elastic scattering
00063     RegisterPhysics( new G4HadronElasticPhysics("elastic",ver,false));
00064 
00065     // Hadron Physics
00066     G4bool quasiElastic=true;
00067     RegisterPhysics( new HadronPhysicsQGSPCMS_FTFP_BERT("hadron",quasiElastic));   
00068     // Stopping Physics
00069     RegisterPhysics( new G4QStoppingPhysics("stopping"));
00070 
00071     // Ion Physics
00072     RegisterPhysics( new G4IonPhysics("ion"));
00073 
00074     // Neutron tracking cut
00075     if (tracking) 
00076       RegisterPhysics( new G4NeutronTrackingCut("Neutron tracking cut", ver));
00077 
00078     /*
00079     if(hadronPhysics=="QGSP_FTFP_BERT") {
00080       RegisterPhysics( new HadronPhysicsQGSPCMS_FTFP_BERT_WP("hadron",quasiElastic)); 
00081     }
00082     else if(hadronPhysics=="QGSP_BERT") {
00083       RegisterPhysics( new HadronPhysicsQGSP_BERT_WP("hadron",quasiElastic));
00084     }
00085     else if (hadronPhysics=="QGSP") {
00086       RegisterPhysics( new HadronPhysicsQGSP_WP("hadron",quasiElastic));
00087     }
00088     else {
00089       edm::LogInfo("PhysicsList") << hadronPhysics << " is not available for GflashHadronPhysics!"
00090                                   << "... Using QGSP_FTFP_BERT\n";
00091       RegisterPhysics( new HadronPhysicsQGSPCMS_FTFP_BERT_WP("hadron",quasiElastic));
00092     }
00093     // Stopping Physics
00094     RegisterPhysics( new G4QStoppingPhysics("stopping"));
00095 
00096     // Ion Physics
00097     RegisterPhysics( new G4IonPhysics("ion"));
00098 
00099     // Neutron tracking cut
00100     if (tracking) 
00101       RegisterPhysics( new G4NeutronTrackingCut("Neutron tracking cut", ver));
00102     */
00103   }
00104 
00105   // Monopoles
00106   RegisterPhysics( new CMSMonopolePhysics(table_,fieldBuilder_,p));
00107 
00108   // singleton histogram object
00109   if(thePar.getParameter<bool>("GflashHistogram")) {
00110     theHisto = GflashHistogram::instance();
00111     theHisto->setStoreFlag(true);
00112     theHisto->bookHistogram(thePar.getParameter<std::string>("GflashHistogramName"));
00113   }
00114 
00115 }
00116 
00117 GFlash::~GFlash() {
00118 
00119   if(thePar.getParameter<bool>("GflashHistogram")) {
00120     if(theHisto) delete theHisto;
00121   }
00122 
00123 }
00124 
00125 //define this as a plug-in
00126 #include "FWCore/Framework/interface/MakerMacros.h"
00127 #include "SimG4Core/Physics/interface/PhysicsListFactory.h"
00128 
00129 
00130 DEFINE_PHYSICSLIST(GFlash);
00131