CMS 3D CMS Logo

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

#include <ParametrisedEMPhysics.h>

Inheritance diagram for ParametrisedEMPhysics:

Classes

struct  TLSmod
 

Public Member Functions

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

Private Member Functions

void ModifyTransportation (const G4ParticleDefinition *, int ntry, double th1, double th2)
 

Private Attributes

edm::ParameterSet theParSet
 

Static Private Attributes

static G4ThreadLocal TLSmodm_tpmod = nullptr
 

Detailed Description

Definition at line 16 of file ParametrisedEMPhysics.h.

Constructor & Destructor Documentation

◆ ParametrisedEMPhysics()

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

Definition at line 64 of file ParametrisedEMPhysics.cc.

References EgammaValidation_cff::genp, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), mps_fire::i, NREG, nt, AlCaHLTBitMon_ParallelJobs::p, rname, and theParSet.

65  : G4VPhysicsConstructor(name), theParSet(p) {
66  // EM parameters common for any EM physics configuration
67  G4EmParameters* param = G4EmParameters::Instance();
68  G4int verb = theParSet.getUntrackedParameter<int>("Verbosity", 0);
69  param->SetVerbose(verb);
70 
71  G4double bremth = theParSet.getParameter<double>("G4BremsstrahlungThreshold") * CLHEP::GeV;
72  param->SetBremsstrahlungTh(bremth);
73  G4double mubrth = theParSet.getParameter<double>("G4MuonBremsstrahlungThreshold") * CLHEP::GeV;
74  param->SetMuHadBremsstrahlungTh(mubrth);
75 
76  bool genp = theParSet.getParameter<bool>("G4GammaGeneralProcess");
77  param->SetGeneralProcessActive(genp);
78 
79  bool pe = p.getParameter<bool>("PhotoeffectBelowKShell");
80  param->SetPhotoeffectBelowKShell(pe);
81 
82  int type = p.getParameter<int>("G4TransportWithMSC");
83  G4TransportationWithMscType trtype = G4TransportationWithMscType::fDisabled;
84  if (type == 1) {
85  trtype = G4TransportationWithMscType::fEnabled;
86  } else if (type == 2) {
87  trtype = G4TransportationWithMscType::fMultipleSteps;
88  }
89  param->SetTransportationWithMsc(trtype);
90 
91  bool mudat = theParSet.getParameter<bool>("ReadMuonData");
92  param->SetRetrieveMuDataFromFile(mudat);
93 
94  bool fluo = theParSet.getParameter<bool>("FlagFluo");
95  param->SetFluo(fluo);
96 
97  bool modifyT = theParSet.getParameter<bool>("ModifyTransportation");
98  double th1 = theParSet.getUntrackedParameter<double>("ThresholdWarningEnergy") * CLHEP::MeV;
99  double th2 = theParSet.getUntrackedParameter<double>("ThresholdImportantEnergy") * CLHEP::MeV;
100  int nt = theParSet.getUntrackedParameter<int>("ThresholdTrials");
101 
102  edm::LogVerbatim("SimG4CoreApplication")
103  << "### ParametrisedEMPhysics parameters:"
104  << "\n verbosity= " << verb << "\n fluoFlag: " << fluo << "\n fluo below K-shell: " << pe
105  << "\n transportation with msc: " << type << "\n modifyTransport: " << modifyT << " Ntrials= " << nt
106  << "\n ThWarning(MeV)= " << th1 / CLHEP::MeV << "\n ThException(MeV)= " << th2 / CLHEP::MeV
107  << "\n read muon data: " << mudat << "\n bremsstrahlung threshold Eth(GeV)= " << bremth / CLHEP::GeV;
108 
109  // Russian roulette and tracking cut for e+-
110  double energyLim = theParSet.getParameter<double>("RusRoElectronEnergyLimit") * CLHEP::MeV;
111  if (energyLim > 0.0) {
112  G4double rrfact[NREG] = {1.0};
113 
114  rrfact[0] = theParSet.getParameter<double>("RusRoEcalElectron");
115  rrfact[1] = theParSet.getParameter<double>("RusRoHcalElectron");
116  rrfact[2] = theParSet.getParameter<double>("RusRoMuonIronElectron");
117  rrfact[3] = theParSet.getParameter<double>("RusRoPreShowerElectron");
118  rrfact[4] = theParSet.getParameter<double>("RusRoCastorElectron");
119  rrfact[5] = theParSet.getParameter<double>("RusRoWorldElectron");
120  for (int i = 0; i < NREG; ++i) {
121  if (rrfact[i] < 1.0) {
122  param->ActivateSecondaryBiasing("eIoni", rname[i], rrfact[i], energyLim);
123  param->ActivateSecondaryBiasing("hIoni", rname[i], rrfact[i], energyLim);
124  edm::LogVerbatim("SimG4CoreApplication")
125  << "ParametrisedEMPhysics: Russian Roulette"
126  << " for e- Prob= " << rrfact[i] << " Elimit(MeV)= " << energyLim / CLHEP::MeV << " inside " << rname[i];
127  }
128  }
129  }
130 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
genp
produce generated paricles in acceptance #
edm::ParameterSet theParSet
T getUntrackedParameter(std::string const &, T const &) const
int nt
Definition: AMPTWrapper.h:42
const G4String rname[NREG]
const G4int NREG

◆ ~ParametrisedEMPhysics()

ParametrisedEMPhysics::~ParametrisedEMPhysics ( )
override

Definition at line 132 of file ParametrisedEMPhysics.cc.

References m_tpmod.

132  {
133  delete m_tpmod;
134  m_tpmod = nullptr;
135 }
static G4ThreadLocal TLSmod * m_tpmod

Member Function Documentation

◆ ConstructParticle()

void ParametrisedEMPhysics::ConstructParticle ( )
override

Definition at line 137 of file ParametrisedEMPhysics.cc.

137  {
138  // minimal set of particles for EM physics
139  G4EmBuilder::ConstructMinimalEmSet();
140 }

◆ ConstructProcess()

void ParametrisedEMPhysics::ConstructProcess ( )
override

Definition at line 142 of file ParametrisedEMPhysics.cc.

References nanoDQM_cfi::Electron, mixOne_premix_on_sim_cfi::gem, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), mps_fire::i, m_tpmod, ModifyTransportation(), nt, alignCSCRings::r, ElectronLimiter::SetFieldCheckFlag(), ElectronLimiter::SetRangeCheckFlag(), ElectronLimiter::SetTrackingCutPerRegion(), ParametrisedEMPhysics::TLSmod::theEcalEMShowerModel, ParametrisedEMPhysics::TLSmod::theEcalHadShowerModel, ParametrisedEMPhysics::TLSmod::theFastSimulationManagerProcess, ParametrisedEMPhysics::TLSmod::theHcalEMShowerModel, ParametrisedEMPhysics::TLSmod::theHcalHadShowerModel, ParametrisedEMPhysics::TLSmod::theLowEnergyFastSimModel, and theParSet.

142  {
143  edm::LogVerbatim("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess() started";
144 
145  // GFlash part
146  bool gem = theParSet.getParameter<bool>("GflashEcal");
147  bool lowEnergyGem = theParSet.getParameter<bool>("LowEnergyGflashEcal");
148  bool ghad = theParSet.getParameter<bool>("GflashHcal");
149  bool gemHad = theParSet.getParameter<bool>("GflashEcalHad");
150  bool ghadHad = theParSet.getParameter<bool>("GflashHcalHad");
151 
152  if (gem || ghad || lowEnergyGem || gemHad || ghadHad) {
153  if (!m_tpmod) {
154  m_tpmod = new TLSmod;
155  }
156  edm::LogVerbatim("SimG4CoreApplication")
157  << "ParametrisedEMPhysics: GFlash Construct for e+-: " << gem << " " << ghad << " " << lowEnergyGem
158  << " for hadrons: " << gemHad << " " << ghadHad;
159 
160  m_tpmod->theFastSimulationManagerProcess = std::make_unique<G4FastSimulationManagerProcess>();
161 
162  if (gem || ghad) {
163  G4Electron::Electron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
164  G4Positron::Positron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
165  } else if (lowEnergyGem) {
166  G4Electron::Electron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
167  G4Positron::Positron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
168  }
169 
170  if (gemHad || ghadHad) {
171  G4Proton::Proton()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
172  G4AntiProton::AntiProton()->GetProcessManager()->AddDiscreteProcess(
174  G4PionPlus::PionPlus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
175  G4PionMinus::PionMinus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
176  G4KaonPlus::KaonPlus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
177  G4KaonMinus::KaonMinus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
178  }
179 
180  if (gem || gemHad || lowEnergyGem) {
181  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("EcalRegion", false);
182 
183  if (!aRegion) {
184  edm::LogWarning("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess: "
185  << "EcalRegion is not defined, GFlash will not be enabled for ECAL!";
186 
187  } else {
188  if (gem) {
189  //Electromagnetic Shower Model for ECAL
191  std::make_unique<GFlashEMShowerModel>("GflashEcalEMShowerModel", aRegion, theParSet);
192  } else if (lowEnergyGem) {
193  //Low energy electromagnetic Shower Model for ECAL
195  std::make_unique<LowEnergyFastSimModel>("LowEnergyFastSimModel", aRegion, theParSet);
196  }
197 
198  if (gemHad) {
199  //Electromagnetic Shower Model for ECAL
201  std::make_unique<GFlashHadronShowerModel>("GflashEcalHadShowerModel", aRegion, theParSet);
202  }
203  }
204  }
205  if (ghad || ghadHad) {
206  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
207  if (!aRegion) {
208  edm::LogWarning("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess: "
209  << "HcalRegion is not defined, GFlash will not be enabled for HCAL!";
210 
211  } else {
212  if (ghad) {
213  //Electromagnetic Shower Model for HCAL
215  std::make_unique<GFlashEMShowerModel>("GflashHcalEMShowerModel", aRegion, theParSet);
216  }
217  if (ghadHad) {
218  //Electromagnetic Shower Model for ECAL
220  std::make_unique<GFlashHadronShowerModel>("GflashHcalHadShowerModel", aRegion, theParSet);
221  }
222  }
223  }
224  }
225 
226  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
227 
228  // Step limiters for e+-
229  bool eLimiter = theParSet.getParameter<bool>("ElectronStepLimit");
230  bool rLimiter = theParSet.getParameter<bool>("ElectronRangeTest");
231  bool pLimiter = theParSet.getParameter<bool>("PositronStepLimit");
232  // Step limiters for hadrons
233  bool pTCut = theParSet.getParameter<bool>("ProtonRegionLimit");
234  bool piTCut = theParSet.getParameter<bool>("PionRegionLimit");
235 
236  std::vector<std::string> regnames = theParSet.getParameter<std::vector<std::string> >("LimitsPerRegion");
237  std::vector<double> limitsE = theParSet.getParameter<std::vector<double> >("EnergyLimitsE");
238  std::vector<double> limitsH = theParSet.getParameter<std::vector<double> >("EnergyLimitsH");
239  std::vector<double> facE = theParSet.getParameter<std::vector<double> >("EnergyFactorsE");
240  std::vector<double> rmsE = theParSet.getParameter<std::vector<double> >("EnergyRMSE");
241  int nlimits = regnames.size();
242  int nlimitsH = 0;
243  std::vector<const G4Region*> reg;
244  std::vector<G4double> rlimE;
245  std::vector<G4double> rlimH;
246  std::vector<G4double> factE;
247  std::vector<G4double> rmsvE;
248  if (0 < nlimits) {
249  G4RegionStore* store = G4RegionStore::GetInstance();
250  for (int i = 0; i < nlimits; ++i) {
251  // apply limiter for whole CMS
252  if (regnames[i] == "all") {
253  reg.clear();
254  rlimE.clear();
255  rlimH.clear();
256  factE.clear();
257  rmsvE.clear();
258  reg.emplace_back(nullptr);
259  rlimE.emplace_back(limitsE[i] * CLHEP::MeV);
260  rlimH.emplace_back(limitsH[i] * CLHEP::MeV);
261  factE.emplace_back(facE[i]);
262  rmsvE.emplace_back(rmsE[i]);
263  nlimitsH = (limitsH[i] > 0) ? 1 : 0;
264  break;
265  }
266  const G4Region* r = store->GetRegion(regnames[i], false);
267  // apply for concrete G4Region
268  if (r && (limitsE[i] > 0.0 || limitsH[i] > 0.0)) {
269  reg.emplace_back(r);
270  rlimE.emplace_back(limitsE[i] * CLHEP::MeV);
271  rlimH.emplace_back(limitsH[i] * CLHEP::MeV);
272  factE.emplace_back(facE[i]);
273  rmsvE.emplace_back(rmsE[i]);
274  if (limitsH[i] > 0) {
275  ++nlimitsH;
276  }
277  }
278  }
279  nlimits = reg.size();
280  }
281 
282  if (eLimiter || rLimiter || 0 < nlimits) {
284  elim->SetRangeCheckFlag(rLimiter);
285  elim->SetFieldCheckFlag(eLimiter);
286  elim->SetTrackingCutPerRegion(reg, rlimE, factE, rmsvE);
287  ph->RegisterProcess(elim, G4Electron::Electron());
288  }
289 
290  if (pLimiter || 0 < nlimits) {
291  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4Positron::Positron());
292  plim->SetFieldCheckFlag(pLimiter);
293  plim->SetTrackingCutPerRegion(reg, rlimE, factE, rmsvE);
294  ph->RegisterProcess(plim, G4Positron::Positron());
295  }
296  if (0 < nlimits && 0 < nlimitsH) {
297  if (pTCut) {
298  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4Proton::Proton());
299  plim->SetFieldCheckFlag(pLimiter);
300  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
301  ph->RegisterProcess(plim, G4Proton::Proton());
302  }
303  if (piTCut) {
304  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4PionPlus::PionPlus());
305  plim->SetFieldCheckFlag(pLimiter);
306  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
307  ph->RegisterProcess(plim, G4PionPlus::PionPlus());
308  plim = new ElectronLimiter(theParSet, G4PionMinus::PionMinus());
309  plim->SetFieldCheckFlag(pLimiter);
310  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
311  ph->RegisterProcess(plim, G4PionMinus::PionMinus());
312  }
313  }
314  // enable fluorescence
315  bool fluo = theParSet.getParameter<bool>("FlagFluo");
316  if (fluo && !G4LossTableManager::Instance()->AtomDeexcitation()) {
317  G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
318  G4LossTableManager::Instance()->SetAtomDeexcitation(de);
319  }
320  // change parameters of transportation
321  bool modifyT = theParSet.getParameter<bool>("ModifyTransportation");
322  if (modifyT) {
323  double th1 = theParSet.getUntrackedParameter<double>("ThresholdWarningEnergy") * CLHEP::MeV;
324  double th2 = theParSet.getUntrackedParameter<double>("ThresholdImportantEnergy") * CLHEP::MeV;
325  int nt = theParSet.getUntrackedParameter<int>("ThresholdTrials");
327  }
328  edm::LogVerbatim("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess() is done";
329 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
void SetRangeCheckFlag(G4bool)
edm::ParameterSet theParSet
void SetTrackingCutPerRegion(std::vector< const G4Region *> &, std::vector< G4double > &, std::vector< G4double > &, std::vector< G4double > &)
static G4ThreadLocal TLSmod * m_tpmod
void ModifyTransportation(const G4ParticleDefinition *, int ntry, double th1, double th2)
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< G4FastSimulationManagerProcess > theFastSimulationManagerProcess
std::unique_ptr< GFlashEMShowerModel > theEcalEMShowerModel
void SetFieldCheckFlag(G4bool)
std::unique_ptr< LowEnergyFastSimModel > theLowEnergyFastSimModel
int nt
Definition: AMPTWrapper.h:42
std::unique_ptr< GFlashEMShowerModel > theHcalEMShowerModel
std::unique_ptr< GFlashHadronShowerModel > theHcalHadShowerModel
std::unique_ptr< GFlashHadronShowerModel > theEcalHadShowerModel
Log< level::Warning, false > LogWarning

◆ ModifyTransportation()

void ParametrisedEMPhysics::ModifyTransportation ( const G4ParticleDefinition *  part,
int  ntry,
double  th1,
double  th2 
)
private

Definition at line 331 of file ParametrisedEMPhysics.cc.

Referenced by ConstructProcess().

331  {
332  G4ProcessManager* man = part->GetProcessManager();
333  G4Transportation* trans = (G4Transportation*)((*(man->GetProcessList()))[0]);
334  if (trans) {
335  trans->SetThresholdWarningEnergy(th1);
336  trans->SetThresholdImportantEnergy(th2);
337  trans->SetThresholdTrials(ntry);
338  edm::LogVerbatim("SimG4CoreApplication")
339  << "ParametrisedEMPhysics: printout level changed for " << part->GetParticleName();
340  }
341 }
Log< level::Info, true > LogVerbatim
part
Definition: HCALResponse.h:20

Member Data Documentation

◆ m_tpmod

G4ThreadLocal ParametrisedEMPhysics::TLSmod * ParametrisedEMPhysics::m_tpmod = nullptr
staticprivate

Definition at line 29 of file ParametrisedEMPhysics.h.

Referenced by ConstructProcess(), and ~ParametrisedEMPhysics().

◆ theParSet

edm::ParameterSet ParametrisedEMPhysics::theParSet
private

Definition at line 27 of file ParametrisedEMPhysics.h.

Referenced by ConstructProcess(), and ParametrisedEMPhysics().