File indexing completed on 2023-03-17 11:25:01
0001 #include "QGSPCMS_FTFP_BERT.h"
0002 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0003 #include "SimG4Core/PhysicsLists/interface/HadronPhysicsQGSPCMS_FTFP_BERT.h"
0004
0005 #include "G4EmStandardPhysics.hh"
0006 #include "G4DecayPhysics.hh"
0007 #include "G4EmExtraPhysics.hh"
0008 #include "G4IonPhysics.hh"
0009 #include "G4StoppingPhysics.hh"
0010 #include "G4HadronElasticPhysics.hh"
0011 #include "G4NeutronTrackingCut.hh"
0012 #include "G4HadronicProcessStore.hh"
0013
0014 #include "G4HadronPhysicsQGSP_FTFP_BERT.hh"
0015
0016 QGSPCMS_FTFP_BERT::QGSPCMS_FTFP_BERT(const edm::ParameterSet& p) : PhysicsList(p) {
0017 int ver = p.getUntrackedParameter<int>("Verbosity", 0);
0018 bool emPhys = p.getUntrackedParameter<bool>("EMPhysics", true);
0019 bool hadPhys = p.getUntrackedParameter<bool>("HadPhysics", true);
0020 bool tracking = p.getParameter<bool>("TrackingCut");
0021 double timeLimit = p.getParameter<double>("MaxTrackTime") * CLHEP::ns;
0022 double minFTFP = p.getParameter<double>("EminFTFP") * CLHEP::GeV;
0023 double maxBERT = p.getParameter<double>("EmaxBERT") * CLHEP::GeV;
0024 double minQGSP = p.getParameter<double>("EminQGSP") * CLHEP::GeV;
0025 double maxFTFP = p.getParameter<double>("EmaxFTFP") * CLHEP::GeV;
0026 double maxBERTpi = p.getParameter<double>("EmaxBERTpi") * CLHEP::GeV;
0027 edm::LogInfo("PhysicsList") << "You are using the simulation engine: "
0028 << "QGSP_FTFP_BERT \n Flags for EM Physics " << emPhys << ", for Hadronic Physics "
0029 << hadPhys << " and tracking cut " << tracking << " t(ns)= " << timeLimit / CLHEP::ns
0030 << "\n transition energy Bertini/FTFP from " << minFTFP / CLHEP::GeV << " to "
0031 << maxBERTpi / CLHEP::GeV << ":" << maxBERT / CLHEP::GeV << " GeV"
0032 << "\n transition energy FTFP/QGSP from " << minQGSP / CLHEP::GeV << " to "
0033 << maxFTFP / CLHEP::GeV << " GeV";
0034
0035 if (emPhys) {
0036
0037 RegisterPhysics(new G4EmStandardPhysics(ver));
0038
0039
0040 G4EmExtraPhysics* gn = new G4EmExtraPhysics(ver);
0041 RegisterPhysics(gn);
0042 }
0043
0044
0045 this->RegisterPhysics(new G4DecayPhysics(ver));
0046
0047 if (hadPhys) {
0048 G4HadronicProcessStore::Instance()->SetVerbose(ver);
0049
0050
0051 RegisterPhysics(new G4HadronElasticPhysics(ver));
0052
0053
0054 RegisterPhysics(new HadronPhysicsQGSPCMS_FTFP_BERT(minFTFP, maxBERT, minQGSP, maxFTFP, maxBERTpi));
0055
0056
0057 RegisterPhysics(new G4StoppingPhysics(ver));
0058
0059
0060 RegisterPhysics(new G4IonPhysics(ver));
0061
0062
0063 if (tracking) {
0064 G4NeutronTrackingCut* ncut = new G4NeutronTrackingCut(ver);
0065 ncut->SetTimeLimit(timeLimit);
0066 RegisterPhysics(ncut);
0067 }
0068 }
0069 }