CMS 3D CMS Logo

ParametrisedEMPhysics.cc
Go to the documentation of this file.
1 //
2 // Joanna Weng 08.2005
3 // Physics process for Gflash parameterisation
4 // modified by Soon Yung Jun, Dongwook Jang
5 // V.Ivanchenko rename the class, cleanup, and move
6 // to SimG4Core/Application - 2012/08/14
7 
14 
15 #include "G4FastSimulationManagerProcess.hh"
16 #include "G4ProcessManager.hh"
17 #include "G4HadronicParameters.hh"
18 #include "G4EmBuilder.hh"
19 
20 #include "G4RegionStore.hh"
21 #include "G4Electron.hh"
22 #include "G4Positron.hh"
23 #include "G4MuonMinus.hh"
24 #include "G4MuonPlus.hh"
25 #include "G4PionMinus.hh"
26 #include "G4PionPlus.hh"
27 #include "G4KaonMinus.hh"
28 #include "G4KaonPlus.hh"
29 #include "G4Proton.hh"
30 #include "G4AntiProton.hh"
31 
32 #include "G4EmParameters.hh"
33 #include "G4LossTableManager.hh"
34 #include "G4PhysicsListHelper.hh"
35 #include "G4ProcessManager.hh"
36 #include "G4SystemOfUnits.hh"
37 #include "G4Transportation.hh"
38 #include "G4UAtomicDeexcitation.hh"
39 #include "G4Version.hh"
40 #include <memory>
41 
42 #include <string>
43 #include <vector>
44 
45 const G4int NREG = 7;
46 const G4String rname[NREG] = {"EcalRegion",
47  "HcalRegion",
48  "HGcalRegion",
49  "MuonIron",
50  "PreshowerRegion",
51  "CastorRegion",
52  "DefaultRegionForTheWorld"};
53 
55  std::unique_ptr<GFlashEMShowerModel> theEcalEMShowerModel;
56  std::unique_ptr<LowEnergyFastSimModel> theLowEnergyFastSimModel;
57  std::unique_ptr<GFlashEMShowerModel> theHcalEMShowerModel;
58  std::unique_ptr<GFlashHadronShowerModel> theEcalHadShowerModel;
59  std::unique_ptr<GFlashHadronShowerModel> theHcalHadShowerModel;
60  std::unique_ptr<G4FastSimulationManagerProcess> theFastSimulationManagerProcess;
61 };
62 
64 
66  : G4VPhysicsConstructor(name), theParSet(p) {
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 #if G4VERSION_NUMBER >= 1110
80  bool genn = theParSet.getParameter<bool>("G4NeutronGeneralProcess");
81  G4HadronicParameters::Instance()->SetEnableNeutronGeneralProcess(genn);
82 
83  if (theParSet.getParameter<bool>("G4TransportWithMSC"))
84  param->SetTransportationWithMsc(G4TransportationWithMscType::fEnabled);
85 #endif
86 
87  bool mudat = theParSet.getParameter<bool>("ReadMuonData");
88  param->SetRetrieveMuDataFromFile(mudat);
89 
90  bool fluo = theParSet.getParameter<bool>("FlagFluo");
91  param->SetFluo(fluo);
92 
93  bool modifyT = theParSet.getParameter<bool>("ModifyTransportation");
94  double th1 = theParSet.getUntrackedParameter<double>("ThresholdWarningEnergy") * MeV;
95  double th2 = theParSet.getUntrackedParameter<double>("ThresholdImportantEnergy") * MeV;
96  int nt = theParSet.getUntrackedParameter<int>("ThresholdTrials");
97 
98  edm::LogVerbatim("SimG4CoreApplication")
99  << "ParametrisedEMPhysics::ConstructProcess: bremsstrahlung threshold Eth= " << bremth / GeV << " GeV"
100  << "\n verbosity= " << verb << " fluoFlag: " << fluo
101  << " modifyTransport: " << modifyT << " Ntrials= " << nt
102  << "\n ThWarning(MeV)= " << th1 << " ThException(MeV)= " << th2;
103 
104  // Russian roulette and tracking cut for e+-
105  double energyLim = theParSet.getParameter<double>("RusRoElectronEnergyLimit") * MeV;
106  if (energyLim > 0.0) {
107  G4double rrfact[NREG] = {1.0};
108 
109  rrfact[0] = theParSet.getParameter<double>("RusRoEcalElectron");
110  rrfact[1] = theParSet.getParameter<double>("RusRoHcalElectron");
111  rrfact[2] = theParSet.getParameter<double>("RusRoMuonIronElectron");
112  rrfact[3] = theParSet.getParameter<double>("RusRoPreShowerElectron");
113  rrfact[4] = theParSet.getParameter<double>("RusRoCastorElectron");
114  rrfact[5] = theParSet.getParameter<double>("RusRoWorldElectron");
115  for (int i = 0; i < NREG; ++i) {
116  if (rrfact[i] < 1.0) {
117  param->ActivateSecondaryBiasing("eIoni", rname[i], rrfact[i], energyLim);
118  param->ActivateSecondaryBiasing("hIoni", rname[i], rrfact[i], energyLim);
119  edm::LogVerbatim("SimG4CoreApplication")
120  << "ParametrisedEMPhysics: Russian Roulette"
121  << " for e- Prob= " << rrfact[i] << " Elimit(MeV)= " << energyLim / CLHEP::MeV << " inside " << rname[i];
122  }
123  }
124  }
125 }
126 
128  delete m_tpmod;
129  m_tpmod = nullptr;
130 }
131 
133  // minimal set of particles for EM physics
134  G4EmBuilder::ConstructMinimalEmSet();
135 }
136 
138  edm::LogVerbatim("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess() started";
139 
140  // GFlash part
141  bool gem = theParSet.getParameter<bool>("GflashEcal");
142  bool lowEnergyGem = theParSet.getParameter<bool>("LowEnergyGflashEcal");
143  bool ghad = theParSet.getParameter<bool>("GflashHcal");
144  bool gemHad = theParSet.getParameter<bool>("GflashEcalHad");
145  bool ghadHad = theParSet.getParameter<bool>("GflashHcalHad");
146 
147  if (gem || ghad || lowEnergyGem || gemHad || ghadHad) {
148  if (!m_tpmod) {
149  m_tpmod = new TLSmod;
150  }
151  edm::LogVerbatim("SimG4CoreApplication")
152  << "ParametrisedEMPhysics: GFlash Construct for e+-: " << gem << " " << ghad << " " << lowEnergyGem
153  << " for hadrons: " << gemHad << " " << ghadHad;
154 
155  m_tpmod->theFastSimulationManagerProcess = std::make_unique<G4FastSimulationManagerProcess>();
156 
157  if (gem || ghad) {
158  G4Electron::Electron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
159  G4Positron::Positron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
160  } else if (lowEnergyGem) {
161  G4Electron::Electron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
162  G4Positron::Positron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
163  }
164 
165  if (gemHad || ghadHad) {
166  G4Proton::Proton()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
167  G4AntiProton::AntiProton()->GetProcessManager()->AddDiscreteProcess(
169  G4PionPlus::PionPlus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
170  G4PionMinus::PionMinus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
171  G4KaonPlus::KaonPlus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
172  G4KaonMinus::KaonMinus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
173  }
174 
175  if (gem || gemHad || lowEnergyGem) {
176  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("EcalRegion", false);
177 
178  if (!aRegion) {
179  edm::LogWarning("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess: "
180  << "EcalRegion is not defined, GFlash will not be enabled for ECAL!";
181 
182  } else {
183  if (gem) {
184  //Electromagnetic Shower Model for ECAL
186  std::make_unique<GFlashEMShowerModel>("GflashEcalEMShowerModel", aRegion, theParSet);
187  } else if (lowEnergyGem) {
188  //Low energy electromagnetic Shower Model for ECAL
190  std::make_unique<LowEnergyFastSimModel>("LowEnergyFastSimModel", aRegion, theParSet);
191  }
192 
193  if (gemHad) {
194  //Electromagnetic Shower Model for ECAL
196  std::make_unique<GFlashHadronShowerModel>("GflashEcalHadShowerModel", aRegion, theParSet);
197  }
198  }
199  }
200  if (ghad || ghadHad) {
201  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
202  if (!aRegion) {
203  edm::LogWarning("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess: "
204  << "HcalRegion is not defined, GFlash will not be enabled for HCAL!";
205 
206  } else {
207  if (ghad) {
208  //Electromagnetic Shower Model for HCAL
210  std::make_unique<GFlashEMShowerModel>("GflashHcalEMShowerModel", aRegion, theParSet);
211  }
212  if (ghadHad) {
213  //Electromagnetic Shower Model for ECAL
215  std::make_unique<GFlashHadronShowerModel>("GflashHcalHadShowerModel", aRegion, theParSet);
216  }
217  }
218  }
219  }
220 
221  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
222 
223  // Step limiters for e+-
224  bool eLimiter = theParSet.getParameter<bool>("ElectronStepLimit");
225  bool rLimiter = theParSet.getParameter<bool>("ElectronRangeTest");
226  bool pLimiter = theParSet.getParameter<bool>("PositronStepLimit");
227  // Step limiters for hadrons
228  bool pTCut = theParSet.getParameter<bool>("ProtonRegionLimit");
229  bool piTCut = theParSet.getParameter<bool>("PionRegionLimit");
230 
231  std::vector<std::string> regnames = theParSet.getParameter<std::vector<std::string> >("LimitsPerRegion");
232  std::vector<double> limitsE = theParSet.getParameter<std::vector<double> >("EnergyLimitsE");
233  std::vector<double> limitsH = theParSet.getParameter<std::vector<double> >("EnergyLimitsH");
234  std::vector<double> facE = theParSet.getParameter<std::vector<double> >("EnergyFactorsE");
235  std::vector<double> rmsE = theParSet.getParameter<std::vector<double> >("EnergyRMSE");
236  int nlimits = regnames.size();
237  int nlimitsH = 0;
238  std::vector<const G4Region*> reg;
239  std::vector<G4double> rlimE;
240  std::vector<G4double> rlimH;
241  std::vector<G4double> factE;
242  std::vector<G4double> rmsvE;
243  if (0 < nlimits) {
244  G4RegionStore* store = G4RegionStore::GetInstance();
245  for (int i = 0; i < nlimits; ++i) {
246  // apply limiter for whole CMS
247  if (regnames[i] == "all") {
248  reg.clear();
249  rlimE.clear();
250  rlimH.clear();
251  factE.clear();
252  rmsvE.clear();
253  reg.emplace_back(nullptr);
254  rlimE.emplace_back(limitsE[i] * CLHEP::MeV);
255  rlimH.emplace_back(limitsH[i] * CLHEP::MeV);
256  factE.emplace_back(facE[i]);
257  rmsvE.emplace_back(rmsE[i]);
258  nlimitsH = (limitsH[i] > 0) ? 1 : 0;
259  break;
260  }
261  const G4Region* r = store->GetRegion(regnames[i], false);
262  // apply for concrete G4Region
263  if (r && (limitsE[i] > 0.0 || limitsH[i] > 0.0)) {
264  reg.emplace_back(r);
265  rlimE.emplace_back(limitsE[i] * CLHEP::MeV);
266  rlimH.emplace_back(limitsH[i] * CLHEP::MeV);
267  factE.emplace_back(facE[i]);
268  rmsvE.emplace_back(rmsE[i]);
269  if (limitsH[i] > 0) {
270  ++nlimitsH;
271  }
272  }
273  }
274  nlimits = reg.size();
275  }
276 
277  if (eLimiter || rLimiter || 0 < nlimits) {
279  elim->SetRangeCheckFlag(rLimiter);
280  elim->SetFieldCheckFlag(eLimiter);
281  elim->SetTrackingCutPerRegion(reg, rlimE, factE, rmsvE);
282  ph->RegisterProcess(elim, G4Electron::Electron());
283  }
284 
285  if (pLimiter || 0 < nlimits) {
286  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4Positron::Positron());
287  plim->SetFieldCheckFlag(pLimiter);
288  plim->SetTrackingCutPerRegion(reg, rlimE, factE, rmsvE);
289  ph->RegisterProcess(plim, G4Positron::Positron());
290  }
291  if (0 < nlimits && 0 < nlimitsH) {
292  if (pTCut) {
293  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4Proton::Proton());
294  plim->SetFieldCheckFlag(pLimiter);
295  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
296  ph->RegisterProcess(plim, G4Proton::Proton());
297  }
298  if (piTCut) {
299  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4PionPlus::PionPlus());
300  plim->SetFieldCheckFlag(pLimiter);
301  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
302  ph->RegisterProcess(plim, G4PionPlus::PionPlus());
303  plim = new ElectronLimiter(theParSet, G4PionMinus::PionMinus());
304  plim->SetFieldCheckFlag(pLimiter);
305  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
306  ph->RegisterProcess(plim, G4PionMinus::PionMinus());
307  }
308  }
309  // enable fluorescence
310  bool fluo = theParSet.getParameter<bool>("FlagFluo");
311  if (fluo && !G4LossTableManager::Instance()->AtomDeexcitation()) {
312  G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
313  G4LossTableManager::Instance()->SetAtomDeexcitation(de);
314  }
315  // change parameters of transportation
316  bool modifyT = theParSet.getParameter<bool>("ModifyTransportation");
317  if (modifyT) {
318  double th1 = theParSet.getUntrackedParameter<double>("ThresholdWarningEnergy") * MeV;
319  double th2 = theParSet.getUntrackedParameter<double>("ThresholdImportantEnergy") * MeV;
320  int nt = theParSet.getUntrackedParameter<int>("ThresholdTrials");
322  }
323  edm::LogVerbatim("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess() is done";
324 }
325 
326 void ParametrisedEMPhysics::ModifyTransportation(const G4ParticleDefinition* part, int ntry, double th1, double th2) {
327  G4ProcessManager* man = part->GetProcessManager();
328  G4Transportation* trans = (G4Transportation*)((*(man->GetProcessList()))[0]);
329  if (trans) {
330  trans->SetThresholdWarningEnergy(th1);
331  trans->SetThresholdImportantEnergy(th2);
332  trans->SetThresholdTrials(ntry);
333  edm::LogVerbatim("SimG4CoreApplication")
334  << "ParametrisedEMPhysics: printout level changed for " << part->GetParticleName();
335  }
336 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
genp
produce generated paricles in acceptance #
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
part
Definition: HCALResponse.h:20
ParametrisedEMPhysics(const std::string &name, const edm::ParameterSet &p)
const G4String rname[NREG]
std::unique_ptr< GFlashHadronShowerModel > theEcalHadShowerModel
const G4int NREG
Log< level::Warning, false > LogWarning