CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Attributes
ParametrisedEMPhysics Class Reference

#include <ParametrisedEMPhysics.h>

Inheritance diagram for ParametrisedEMPhysics:

Public Member Functions

virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
 ParametrisedEMPhysics (std::string name, const edm::ParameterSet &p)
 
virtual ~ParametrisedEMPhysics ()
 

Private Attributes

GFlashEMShowerModeltheEcalEMShowerModel
 
GFlashHadronShowerModeltheEcalHadShowerModel
 
ElectronLimitertheElectronLimiter
 
GFlashEMShowerModeltheHcalEMShowerModel
 
GFlashHadronShowerModeltheHcalHadShowerModel
 
edm::ParameterSet theParSet
 
ElectronLimiterthePositronLimiter
 

Detailed Description

Definition at line 18 of file ParametrisedEMPhysics.h.

Constructor & Destructor Documentation

ParametrisedEMPhysics::ParametrisedEMPhysics ( std::string  name,
const edm::ParameterSet p 
)

Definition at line 43 of file ParametrisedEMPhysics.cc.

References theEcalEMShowerModel, theEcalHadShowerModel, theHcalEMShowerModel, and theHcalHadShowerModel.

45  : G4VPhysicsConstructor(name), theParSet(p)
46 {
51 }
edm::ParameterSet theParSet
GFlashHadronShowerModel * theEcalHadShowerModel
GFlashHadronShowerModel * theHcalHadShowerModel
GFlashEMShowerModel * theHcalEMShowerModel
GFlashEMShowerModel * theEcalEMShowerModel
ParametrisedEMPhysics::~ParametrisedEMPhysics ( )
virtual

Definition at line 53 of file ParametrisedEMPhysics.cc.

References theEcalEMShowerModel, theEcalHadShowerModel, theHcalEMShowerModel, and theHcalHadShowerModel.

53  {
54  delete theEcalEMShowerModel;
55  delete theEcalHadShowerModel;
56  delete theHcalEMShowerModel;
57  delete theHcalHadShowerModel;
58 }
GFlashHadronShowerModel * theEcalHadShowerModel
GFlashHadronShowerModel * theHcalHadShowerModel
GFlashEMShowerModel * theHcalEMShowerModel
GFlashEMShowerModel * theEcalEMShowerModel

Member Function Documentation

void ParametrisedEMPhysics::ConstructParticle ( )
virtual

Definition at line 60 of file ParametrisedEMPhysics.cc.

61 {
62  G4LeptonConstructor pLeptonConstructor;
63  pLeptonConstructor.ConstructParticle();
64 
65  G4MesonConstructor pMesonConstructor;
66  pMesonConstructor.ConstructParticle();
67 
68  G4BaryonConstructor pBaryonConstructor;
69  pBaryonConstructor.ConstructParticle();
70 
71  G4ShortLivedConstructor pShortLivedConstructor;
72  pShortLivedConstructor.ConstructParticle();
73 
74  G4IonConstructor pConstructor;
75  pConstructor.ConstructParticle();
76 }
void ParametrisedEMPhysics::ConstructProcess ( )
virtual

Definition at line 78 of file ParametrisedEMPhysics.cc.

References configurableAnalysis::Electron, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), GeV, i, MeV, ElectronLimiter::SetFieldCheckFlag(), ElectronLimiter::SetRangeCheckFlag(), theEcalEMShowerModel, theEcalHadShowerModel, theElectronLimiter, theHcalEMShowerModel, theHcalHadShowerModel, theParSet, and thePositronLimiter.

78  {
79 
80  // GFlash part
81  bool gem = theParSet.getParameter<bool>("GflashEcal");
82  bool ghad = theParSet.getParameter<bool>("GflashHcal");
83  bool gemHad = theParSet.getParameter<bool>("GflashEcalHad");
84  bool ghadHad = theParSet.getParameter<bool>("GflashHcalHad");
85 
86  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
87  if(gem || ghad || gemHad || ghadHad) {
88  edm::LogInfo("SimG4CoreApplication")
89  << "ParametrisedEMPhysics: GFlash Construct for e+-: "
90  << gem << " " << ghad << " for hadrons: " << gemHad << " " << ghadHad;
91 
92  G4FastSimulationManagerProcess * theFastSimulationManagerProcess =
93  new G4FastSimulationManagerProcess();
94 
95  if(gem || ghad) {
96  ph->RegisterProcess(theFastSimulationManagerProcess, G4Electron::Electron());
97  ph->RegisterProcess(theFastSimulationManagerProcess, G4Positron::Positron());
98  }
99  if(gemHad || ghadHad) {
100  ph->RegisterProcess(theFastSimulationManagerProcess, G4Proton::Proton());
101  ph->RegisterProcess(theFastSimulationManagerProcess, G4AntiProton::AntiProton());
102  ph->RegisterProcess(theFastSimulationManagerProcess, G4PionPlus::PionPlus());
103  ph->RegisterProcess(theFastSimulationManagerProcess, G4PionMinus::PionMinus());
104  ph->RegisterProcess(theFastSimulationManagerProcess, G4KaonPlus::KaonPlus());
105  ph->RegisterProcess(theFastSimulationManagerProcess, G4KaonMinus::KaonMinus());
106  }
107 
108  if(gem || gemHad) {
109  G4Region* aRegion =
110  G4RegionStore::GetInstance()->GetRegion("EcalRegion");
111 
112  if(!aRegion){
113  edm::LogInfo("SimG4CoreApplication")
114  << "ParametrisedEMPhysics::ConstructProcess: "
115  << "EcalRegion is not defined, GFlash will not be enabled for ECAL!";
116 
117  } else {
118  if(gem) {
119 
120  //Electromagnetic Shower Model for ECAL
122  new GFlashEMShowerModel("GflashEcalEMShowerModel",aRegion,theParSet);
123  //std::cout << "GFlash is defined for EcalRegion" << std::endl;
124  }
125  if(gemHad) {
126 
127  //Electromagnetic Shower Model for ECAL
129  new GFlashHadronShowerModel("GflashEcalHadShowerModel",aRegion,theParSet);
130  //std::cout << "GFlash is defined for EcalRegion" << std::endl;
131  }
132  }
133  }
134  if(ghad || ghadHad) {
135  G4Region* aRegion =
136  G4RegionStore::GetInstance()->GetRegion("HcalRegion");
137  if(!aRegion) {
138  edm::LogInfo("SimG4CoreApplication")
139  << "ParametrisedEMPhysics::ConstructProcess: "
140  << "HcalRegion is not defined, GFlash will not be enabled for HCAL!";
141 
142  } else {
143  if(ghad) {
144 
145  //Electromagnetic Shower Model for HCAL
147  new GFlashEMShowerModel("GflashHcalEMShowerModel",aRegion,theParSet);
148  //std::cout << "GFlash is defined for HcalRegion" << std::endl;
149  }
150  if(ghadHad) {
151 
152  //Electromagnetic Shower Model for ECAL
154  new GFlashHadronShowerModel("GflashHcalHadShowerModel",aRegion,theParSet);
155  //std::cout << "GFlash is defined for EcalRegion" << std::endl;
156  }
157  }
158  }
159  }
160  // bremsstrahlung threshold and EM verbosity
161  G4EmProcessOptions opt;
162  G4int verb = theParSet.getUntrackedParameter<int>("Verbosity",0);
163  opt.SetVerbose(verb - 1);
164 
165  G4double bremth = theParSet.getParameter<double>("G4BremsstrahlungThreshold")*GeV;
166  edm::LogInfo("SimG4CoreApplication")
167  << "ParametrisedEMPhysics::ConstructProcess: bremsstrahlung threshold Eth= "
168  << bremth/GeV << " GeV";
169  opt.SetBremsstrahlungTh(bremth);
170 
171  // Russian roulette and tracking cut for e+-
172  const G4int NREG = 6;
173  const G4String rname[NREG] = {"EcalRegion", "HcalRegion", "MuonIron",
174  "PreshowerRegion","CastorRegion",
175  "DefaultRegionForTheWorld"};
176  G4double rrfact[NREG] = { 1.0 };
177 
178  double energyLim =
179  theParSet.getParameter<double>("RusRoElectronEnergyLimit")*MeV;
180  if(energyLim > 0.0) {
181  rrfact[0] = theParSet.getParameter<double>("RusRoEcalElectron");
182  rrfact[1] = theParSet.getParameter<double>("RusRoHcalElectron");
183  rrfact[2] = theParSet.getParameter<double>("RusRoMuonIronElectron");
184  rrfact[3] = theParSet.getParameter<double>("RusRoPreShowerElectron");
185  rrfact[4] = theParSet.getParameter<double>("RusRoCastorElectron");
186  rrfact[5] = theParSet.getParameter<double>("RusRoWorldElectron");
187  for(int i=0; i<NREG; ++i) {
188  if(rrfact[i] < 1.0) {
189  opt.ActivateSecondaryBiasing("eIoni",rname[i],rrfact[i],energyLim);
190  opt.ActivateSecondaryBiasing("hIoni",rname[i],rrfact[i],energyLim);
191  edm::LogInfo("SimG4CoreApplication")
192  << "ParametrisedEMPhysics: Russian Roulette"
193  << " for e- Prob= " << rrfact[i]
194  << " Elimit(MeV)= " << energyLim/CLHEP::MeV
195  << " inside " << rname[i];
196  }
197  }
198  }
199 
200  // Step limiters for e+-
201  bool eLimiter = theParSet.getParameter<bool>("ElectronStepLimit");
202  bool rLimiter = theParSet.getParameter<bool>("ElectronRangeTest");
203  bool pLimiter = theParSet.getParameter<bool>("PositronStepLimit");
204 
205  if(eLimiter || rLimiter) {
209  ph->RegisterProcess(theElectronLimiter, G4Electron::Electron());
210  }
211 
212  if(pLimiter){
215  ph->RegisterProcess(theElectronLimiter, G4Positron::Positron());
216  }
217  // enable fluorescence
218  bool fluo = theParSet.getParameter<bool>("FlagFluo");
219  if(fluo) {
220  G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
221  G4LossTableManager::Instance()->SetAtomDeexcitation(de);
222  de->SetFluo(true);
223  }
224  // enable muon nuclear (valid option for Geant4 10.0pX only)
225  bool munuc = theParSet.getParameter<bool>("FlagMuNucl");
226  if(munuc) {
227  G4MuonNuclearProcess* muNucProcess = new G4MuonNuclearProcess();
228  muNucProcess->RegisterMe(new G4MuonVDNuclearModel());
229  ph->RegisterProcess(muNucProcess, G4MuonPlus::MuonPlus());
230  ph->RegisterProcess(muNucProcess, G4MuonMinus::MuonMinus());
231  }
232 
233 }
ElectronLimiter * theElectronLimiter
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
int i
Definition: DBlmapReader.cc:9
const double GeV
Definition: MathUtil.h:16
void SetRangeCheckFlag(G4bool)
edm::ParameterSet theParSet
ElectronLimiter * thePositronLimiter
GFlashHadronShowerModel * theEcalHadShowerModel
const double MeV
void SetFieldCheckFlag(G4bool)
GFlashHadronShowerModel * theHcalHadShowerModel
GFlashEMShowerModel * theHcalEMShowerModel
GFlashEMShowerModel * theEcalEMShowerModel

Member Data Documentation

GFlashEMShowerModel* ParametrisedEMPhysics::theEcalEMShowerModel
private
GFlashHadronShowerModel* ParametrisedEMPhysics::theEcalHadShowerModel
private
ElectronLimiter* ParametrisedEMPhysics::theElectronLimiter
private

Definition at line 37 of file ParametrisedEMPhysics.h.

Referenced by ConstructProcess().

GFlashEMShowerModel* ParametrisedEMPhysics::theHcalEMShowerModel
private
GFlashHadronShowerModel* ParametrisedEMPhysics::theHcalHadShowerModel
private
edm::ParameterSet ParametrisedEMPhysics::theParSet
private

Definition at line 30 of file ParametrisedEMPhysics.h.

Referenced by ConstructProcess().

ElectronLimiter* ParametrisedEMPhysics::thePositronLimiter
private

Definition at line 38 of file ParametrisedEMPhysics.h.

Referenced by ConstructProcess().