CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
ParametrisedEMPhysics Class Reference

#include <ParametrisedEMPhysics.h>

Inheritance diagram for ParametrisedEMPhysics:

Public Member Functions

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

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 41 of file ParametrisedEMPhysics.cc.

References edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), GeV, theEcalEMShowerModel, theEcalHadShowerModel, theHcalEMShowerModel, theHcalHadShowerModel, and theParSet.

43  : G4VPhysicsConstructor(name), theParSet(p)
44 {
45  theEcalEMShowerModel = nullptr;
46  theEcalHadShowerModel = nullptr;
47  theHcalEMShowerModel = nullptr;
48  theHcalHadShowerModel = nullptr;
49 
50  // bremsstrahlung threshold and EM verbosity
51  G4EmParameters* param = G4EmParameters::Instance();
52  G4int verb = theParSet.getUntrackedParameter<int>("Verbosity",0);
53  param->SetVerbose(verb);
54 
55  G4double bremth = theParSet.getParameter<double>("G4BremsstrahlungThreshold")*GeV;
56  param->SetBremsstrahlungTh(bremth);
57 
58  bool fluo = theParSet.getParameter<bool>("FlagFluo");
59  param->SetFluo(fluo);
60 
61  edm::LogInfo("SimG4CoreApplication")
62  << "ParametrisedEMPhysics::ConstructProcess: bremsstrahlung threshold Eth= "
63  << bremth/GeV << " GeV"
64  << "\n verbosity= " << verb
65  << " fluoFlag: " << fluo;
66 }
T getParameter(std::string const &) const
T getUntrackedParameter(std::string const &, T const &) const
const double GeV
Definition: MathUtil.h:16
edm::ParameterSet theParSet
GFlashHadronShowerModel * theEcalHadShowerModel
GFlashHadronShowerModel * theHcalHadShowerModel
GFlashEMShowerModel * theHcalEMShowerModel
GFlashEMShowerModel * theEcalEMShowerModel
ParametrisedEMPhysics::~ParametrisedEMPhysics ( )
override

Definition at line 68 of file ParametrisedEMPhysics.cc.

References theEcalEMShowerModel, theEcalHadShowerModel, theHcalEMShowerModel, and theHcalHadShowerModel.

68  {
69  delete theEcalEMShowerModel;
70  delete theEcalHadShowerModel;
71  delete theHcalEMShowerModel;
72  delete theHcalHadShowerModel;
73 }
GFlashHadronShowerModel * theEcalHadShowerModel
GFlashHadronShowerModel * theHcalHadShowerModel
GFlashEMShowerModel * theHcalEMShowerModel
GFlashEMShowerModel * theEcalEMShowerModel

Member Function Documentation

void ParametrisedEMPhysics::ConstructParticle ( )
override

Definition at line 75 of file ParametrisedEMPhysics.cc.

76 {
77  G4LeptonConstructor pLeptonConstructor;
78  pLeptonConstructor.ConstructParticle();
79 
80  G4MesonConstructor pMesonConstructor;
81  pMesonConstructor.ConstructParticle();
82 
83  G4BaryonConstructor pBaryonConstructor;
84  pBaryonConstructor.ConstructParticle();
85 
86  G4ShortLivedConstructor pShortLivedConstructor;
87  pShortLivedConstructor.ConstructParticle();
88 
89  G4IonConstructor pConstructor;
90  pConstructor.ConstructParticle();
91 }
void ParametrisedEMPhysics::ConstructProcess ( )
override

Definition at line 93 of file ParametrisedEMPhysics.cc.

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

93  {
94 
95  // GFlash part
96  bool gem = theParSet.getParameter<bool>("GflashEcal");
97  bool ghad = theParSet.getParameter<bool>("GflashHcal");
98  bool gemHad = theParSet.getParameter<bool>("GflashEcalHad");
99  bool ghadHad = theParSet.getParameter<bool>("GflashHcalHad");
100 
101  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
102  if(gem || ghad || gemHad || ghadHad) {
103  edm::LogInfo("SimG4CoreApplication")
104  << "ParametrisedEMPhysics: GFlash Construct for e+-: "
105  << gem << " " << ghad << " for hadrons: " << gemHad << " " << ghadHad;
106 
107  G4FastSimulationManagerProcess * theFastSimulationManagerProcess =
108  new G4FastSimulationManagerProcess();
109 
110  if(gem || ghad) {
111  ph->RegisterProcess(theFastSimulationManagerProcess, G4Electron::Electron());
112  ph->RegisterProcess(theFastSimulationManagerProcess, G4Positron::Positron());
113  }
114  if(gemHad || ghadHad) {
115  ph->RegisterProcess(theFastSimulationManagerProcess, G4Proton::Proton());
116  ph->RegisterProcess(theFastSimulationManagerProcess, G4AntiProton::AntiProton());
117  ph->RegisterProcess(theFastSimulationManagerProcess, G4PionPlus::PionPlus());
118  ph->RegisterProcess(theFastSimulationManagerProcess, G4PionMinus::PionMinus());
119  ph->RegisterProcess(theFastSimulationManagerProcess, G4KaonPlus::KaonPlus());
120  ph->RegisterProcess(theFastSimulationManagerProcess, G4KaonMinus::KaonMinus());
121  }
122 
123  if(gem || gemHad) {
124  G4Region* aRegion =
125  G4RegionStore::GetInstance()->GetRegion("EcalRegion");
126 
127  if(!aRegion){
128  edm::LogInfo("SimG4CoreApplication")
129  << "ParametrisedEMPhysics::ConstructProcess: "
130  << "EcalRegion is not defined, GFlash will not be enabled for ECAL!";
131 
132  } else {
133  if(gem) {
134 
135  //Electromagnetic Shower Model for ECAL
137  new GFlashEMShowerModel("GflashEcalEMShowerModel",aRegion,theParSet);
138  //std::cout << "GFlash is defined for EcalRegion" << std::endl;
139  }
140  if(gemHad) {
141 
142  //Electromagnetic Shower Model for ECAL
144  new GFlashHadronShowerModel("GflashEcalHadShowerModel",aRegion,theParSet);
145  //std::cout << "GFlash is defined for EcalRegion" << std::endl;
146  }
147  }
148  }
149  if(ghad || ghadHad) {
150  G4Region* aRegion =
151  G4RegionStore::GetInstance()->GetRegion("HcalRegion");
152  if(!aRegion) {
153  edm::LogInfo("SimG4CoreApplication")
154  << "ParametrisedEMPhysics::ConstructProcess: "
155  << "HcalRegion is not defined, GFlash will not be enabled for HCAL!";
156 
157  } else {
158  if(ghad) {
159 
160  //Electromagnetic Shower Model for HCAL
162  new GFlashEMShowerModel("GflashHcalEMShowerModel",aRegion,theParSet);
163  //std::cout << "GFlash is defined for HcalRegion" << std::endl;
164  }
165  if(ghadHad) {
166 
167  //Electromagnetic Shower Model for ECAL
169  new GFlashHadronShowerModel("GflashHcalHadShowerModel",aRegion,theParSet);
170  //std::cout << "GFlash is defined for EcalRegion" << std::endl;
171  }
172  }
173  }
174  }
175  // Russian roulette and tracking cut for e+-
176  const G4int NREG = 6;
177  const G4String rname[NREG] = {"EcalRegion", "HcalRegion", "MuonIron",
178  "PreshowerRegion","CastorRegion",
179  "DefaultRegionForTheWorld"};
180  G4double rrfact[NREG] = { 1.0 };
181 
182  double energyLim =
183  theParSet.getParameter<double>("RusRoElectronEnergyLimit")*MeV;
184  if(energyLim > 0.0) {
185  G4EmProcessOptions opt;
186  rrfact[0] = theParSet.getParameter<double>("RusRoEcalElectron");
187  rrfact[1] = theParSet.getParameter<double>("RusRoHcalElectron");
188  rrfact[2] = theParSet.getParameter<double>("RusRoMuonIronElectron");
189  rrfact[3] = theParSet.getParameter<double>("RusRoPreShowerElectron");
190  rrfact[4] = theParSet.getParameter<double>("RusRoCastorElectron");
191  rrfact[5] = theParSet.getParameter<double>("RusRoWorldElectron");
192  for(int i=0; i<NREG; ++i) {
193  if(rrfact[i] < 1.0) {
194  opt.ActivateSecondaryBiasing("eIoni",rname[i],rrfact[i],energyLim);
195  opt.ActivateSecondaryBiasing("hIoni",rname[i],rrfact[i],energyLim);
196  edm::LogInfo("SimG4CoreApplication")
197  << "ParametrisedEMPhysics: Russian Roulette"
198  << " for e- Prob= " << rrfact[i]
199  << " Elimit(MeV)= " << energyLim/CLHEP::MeV
200  << " inside " << rname[i];
201  }
202  }
203  }
204 
205  // Step limiters for e+-
206  bool eLimiter = theParSet.getParameter<bool>("ElectronStepLimit");
207  bool rLimiter = theParSet.getParameter<bool>("ElectronRangeTest");
208  bool pLimiter = theParSet.getParameter<bool>("PositronStepLimit");
209 
210  if(eLimiter || rLimiter) {
214  ph->RegisterProcess(theElectronLimiter, G4Electron::Electron());
215  }
216 
217  if(pLimiter){
220  ph->RegisterProcess(theElectronLimiter, G4Positron::Positron());
221  }
222  // enable fluorescence
223  bool fluo = theParSet.getParameter<bool>("FlagFluo");
224  if(fluo && !G4LossTableManager::Instance()->AtomDeexcitation()) {
225  G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
226  G4LossTableManager::Instance()->SetAtomDeexcitation(de);
227  }
228 }
ElectronLimiter * theElectronLimiter
T getParameter(std::string const &) const
void SetRangeCheckFlag(G4bool)
edm::ParameterSet theParSet
ElectronLimiter * thePositronLimiter
A class for AMC data.
Definition: AMC13Event.h:6
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(), and ParametrisedEMPhysics().

ElectronLimiter* ParametrisedEMPhysics::thePositronLimiter
private

Definition at line 38 of file ParametrisedEMPhysics.h.

Referenced by ConstructProcess().