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  // EM parameters common for any EM physics configuration
68  G4EmParameters* param = G4EmParameters::Instance();
69  G4int verb = theParSet.getUntrackedParameter<int>("Verbosity", 0);
70  param->SetVerbose(verb);
71 
72  G4double bremth = theParSet.getParameter<double>("G4BremsstrahlungThreshold") * CLHEP::GeV;
73  param->SetBremsstrahlungTh(bremth);
74  G4double mubrth = theParSet.getParameter<double>("G4MuonBremsstrahlungThreshold") * CLHEP::GeV;
75  param->SetMuHadBremsstrahlungTh(mubrth);
76 
77  bool genp = theParSet.getParameter<bool>("G4GammaGeneralProcess");
78  param->SetGeneralProcessActive(genp);
79 
80  bool pe = p.getParameter<bool>("PhotoeffectBelowKShell");
81  int type = p.getParameter<int>("G4TransportWithMSC");
82 
83 #if G4VERSION_NUMBER >= 1110
84  param->SetPhotoeffectBelowKShell(pe);
85  G4TransportationWithMscType trtype = G4TransportationWithMscType::fDisabled;
86  if (type == 1) {
87  trtype = G4TransportationWithMscType::fEnabled;
88  } else if (type == 2) {
89  trtype = G4TransportationWithMscType::fMultipleSteps;
90  }
91  param->SetTransportationWithMsc(trtype);
92 #endif
93 
94  bool mudat = theParSet.getParameter<bool>("ReadMuonData");
95  param->SetRetrieveMuDataFromFile(mudat);
96 
97  bool fluo = theParSet.getParameter<bool>("FlagFluo");
98  param->SetFluo(fluo);
99 
100  bool modifyT = theParSet.getParameter<bool>("ModifyTransportation");
101  double th1 = theParSet.getUntrackedParameter<double>("ThresholdWarningEnergy") * MeV;
102  double th2 = theParSet.getUntrackedParameter<double>("ThresholdImportantEnergy") * MeV;
103  int nt = theParSet.getUntrackedParameter<int>("ThresholdTrials");
104 
105  edm::LogVerbatim("SimG4CoreApplication")
106  << "### ParametrisedEMPhysics parameters:"
107  << "\n verbosity= " << verb << "\n fluoFlag: " << fluo << "\n fluo below K-shell: " << pe
108  << "\n transportation with msc: " << type << "\n modifyTransport: " << modifyT << " Ntrials= " << nt
109  << "\n ThWarning(MeV)= " << th1 / CLHEP::MeV << "\n ThException(MeV)= " << th2 / CLHEP::MeV
110  << "\n read muon data: " << mudat << "\n bremsstrahlung threshold Eth(GeV)= " << bremth / CLHEP::GeV;
111 
112  // Russian roulette and tracking cut for e+-
113  double energyLim = theParSet.getParameter<double>("RusRoElectronEnergyLimit") * MeV;
114  if (energyLim > 0.0) {
115  G4double rrfact[NREG] = {1.0};
116 
117  rrfact[0] = theParSet.getParameter<double>("RusRoEcalElectron");
118  rrfact[1] = theParSet.getParameter<double>("RusRoHcalElectron");
119  rrfact[2] = theParSet.getParameter<double>("RusRoMuonIronElectron");
120  rrfact[3] = theParSet.getParameter<double>("RusRoPreShowerElectron");
121  rrfact[4] = theParSet.getParameter<double>("RusRoCastorElectron");
122  rrfact[5] = theParSet.getParameter<double>("RusRoWorldElectron");
123  for (int i = 0; i < NREG; ++i) {
124  if (rrfact[i] < 1.0) {
125  param->ActivateSecondaryBiasing("eIoni", rname[i], rrfact[i], energyLim);
126  param->ActivateSecondaryBiasing("hIoni", rname[i], rrfact[i], energyLim);
127  edm::LogVerbatim("SimG4CoreApplication")
128  << "ParametrisedEMPhysics: Russian Roulette"
129  << " for e- Prob= " << rrfact[i] << " Elimit(MeV)= " << energyLim / CLHEP::MeV << " inside " << rname[i];
130  }
131  }
132  }
133 }
134 
136  delete m_tpmod;
137  m_tpmod = nullptr;
138 }
139 
141  // minimal set of particles for EM physics
142  G4EmBuilder::ConstructMinimalEmSet();
143 }
144 
146  edm::LogVerbatim("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess() started";
147 
148  // GFlash part
149  bool gem = theParSet.getParameter<bool>("GflashEcal");
150  bool lowEnergyGem = theParSet.getParameter<bool>("LowEnergyGflashEcal");
151  bool ghad = theParSet.getParameter<bool>("GflashHcal");
152  bool gemHad = theParSet.getParameter<bool>("GflashEcalHad");
153  bool ghadHad = theParSet.getParameter<bool>("GflashHcalHad");
154 
155  if (gem || ghad || lowEnergyGem || gemHad || ghadHad) {
156  if (!m_tpmod) {
157  m_tpmod = new TLSmod;
158  }
159  edm::LogVerbatim("SimG4CoreApplication")
160  << "ParametrisedEMPhysics: GFlash Construct for e+-: " << gem << " " << ghad << " " << lowEnergyGem
161  << " for hadrons: " << gemHad << " " << ghadHad;
162 
163  m_tpmod->theFastSimulationManagerProcess = std::make_unique<G4FastSimulationManagerProcess>();
164 
165  if (gem || ghad) {
166  G4Electron::Electron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
167  G4Positron::Positron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
168  } else if (lowEnergyGem) {
169  G4Electron::Electron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
170  G4Positron::Positron()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
171  }
172 
173  if (gemHad || ghadHad) {
174  G4Proton::Proton()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
175  G4AntiProton::AntiProton()->GetProcessManager()->AddDiscreteProcess(
177  G4PionPlus::PionPlus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
178  G4PionMinus::PionMinus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
179  G4KaonPlus::KaonPlus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
180  G4KaonMinus::KaonMinus()->GetProcessManager()->AddDiscreteProcess(m_tpmod->theFastSimulationManagerProcess.get());
181  }
182 
183  if (gem || gemHad || lowEnergyGem) {
184  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("EcalRegion", false);
185 
186  if (!aRegion) {
187  edm::LogWarning("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess: "
188  << "EcalRegion is not defined, GFlash will not be enabled for ECAL!";
189 
190  } else {
191  if (gem) {
192  //Electromagnetic Shower Model for ECAL
194  std::make_unique<GFlashEMShowerModel>("GflashEcalEMShowerModel", aRegion, theParSet);
195  } else if (lowEnergyGem) {
196  //Low energy electromagnetic Shower Model for ECAL
198  std::make_unique<LowEnergyFastSimModel>("LowEnergyFastSimModel", aRegion, theParSet);
199  }
200 
201  if (gemHad) {
202  //Electromagnetic Shower Model for ECAL
204  std::make_unique<GFlashHadronShowerModel>("GflashEcalHadShowerModel", aRegion, theParSet);
205  }
206  }
207  }
208  if (ghad || ghadHad) {
209  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
210  if (!aRegion) {
211  edm::LogWarning("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess: "
212  << "HcalRegion is not defined, GFlash will not be enabled for HCAL!";
213 
214  } else {
215  if (ghad) {
216  //Electromagnetic Shower Model for HCAL
218  std::make_unique<GFlashEMShowerModel>("GflashHcalEMShowerModel", aRegion, theParSet);
219  }
220  if (ghadHad) {
221  //Electromagnetic Shower Model for ECAL
223  std::make_unique<GFlashHadronShowerModel>("GflashHcalHadShowerModel", aRegion, theParSet);
224  }
225  }
226  }
227  }
228 
229  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
230 
231  // Step limiters for e+-
232  bool eLimiter = theParSet.getParameter<bool>("ElectronStepLimit");
233  bool rLimiter = theParSet.getParameter<bool>("ElectronRangeTest");
234  bool pLimiter = theParSet.getParameter<bool>("PositronStepLimit");
235  // Step limiters for hadrons
236  bool pTCut = theParSet.getParameter<bool>("ProtonRegionLimit");
237  bool piTCut = theParSet.getParameter<bool>("PionRegionLimit");
238 
239  std::vector<std::string> regnames = theParSet.getParameter<std::vector<std::string> >("LimitsPerRegion");
240  std::vector<double> limitsE = theParSet.getParameter<std::vector<double> >("EnergyLimitsE");
241  std::vector<double> limitsH = theParSet.getParameter<std::vector<double> >("EnergyLimitsH");
242  std::vector<double> facE = theParSet.getParameter<std::vector<double> >("EnergyFactorsE");
243  std::vector<double> rmsE = theParSet.getParameter<std::vector<double> >("EnergyRMSE");
244  int nlimits = regnames.size();
245  int nlimitsH = 0;
246  std::vector<const G4Region*> reg;
247  std::vector<G4double> rlimE;
248  std::vector<G4double> rlimH;
249  std::vector<G4double> factE;
250  std::vector<G4double> rmsvE;
251  if (0 < nlimits) {
252  G4RegionStore* store = G4RegionStore::GetInstance();
253  for (int i = 0; i < nlimits; ++i) {
254  // apply limiter for whole CMS
255  if (regnames[i] == "all") {
256  reg.clear();
257  rlimE.clear();
258  rlimH.clear();
259  factE.clear();
260  rmsvE.clear();
261  reg.emplace_back(nullptr);
262  rlimE.emplace_back(limitsE[i] * CLHEP::MeV);
263  rlimH.emplace_back(limitsH[i] * CLHEP::MeV);
264  factE.emplace_back(facE[i]);
265  rmsvE.emplace_back(rmsE[i]);
266  nlimitsH = (limitsH[i] > 0) ? 1 : 0;
267  break;
268  }
269  const G4Region* r = store->GetRegion(regnames[i], false);
270  // apply for concrete G4Region
271  if (r && (limitsE[i] > 0.0 || limitsH[i] > 0.0)) {
272  reg.emplace_back(r);
273  rlimE.emplace_back(limitsE[i] * CLHEP::MeV);
274  rlimH.emplace_back(limitsH[i] * CLHEP::MeV);
275  factE.emplace_back(facE[i]);
276  rmsvE.emplace_back(rmsE[i]);
277  if (limitsH[i] > 0) {
278  ++nlimitsH;
279  }
280  }
281  }
282  nlimits = reg.size();
283  }
284 
285  if (eLimiter || rLimiter || 0 < nlimits) {
287  elim->SetRangeCheckFlag(rLimiter);
288  elim->SetFieldCheckFlag(eLimiter);
289  elim->SetTrackingCutPerRegion(reg, rlimE, factE, rmsvE);
290  ph->RegisterProcess(elim, G4Electron::Electron());
291  }
292 
293  if (pLimiter || 0 < nlimits) {
294  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4Positron::Positron());
295  plim->SetFieldCheckFlag(pLimiter);
296  plim->SetTrackingCutPerRegion(reg, rlimE, factE, rmsvE);
297  ph->RegisterProcess(plim, G4Positron::Positron());
298  }
299  if (0 < nlimits && 0 < nlimitsH) {
300  if (pTCut) {
301  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4Proton::Proton());
302  plim->SetFieldCheckFlag(pLimiter);
303  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
304  ph->RegisterProcess(plim, G4Proton::Proton());
305  }
306  if (piTCut) {
307  ElectronLimiter* plim = new ElectronLimiter(theParSet, G4PionPlus::PionPlus());
308  plim->SetFieldCheckFlag(pLimiter);
309  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
310  ph->RegisterProcess(plim, G4PionPlus::PionPlus());
311  plim = new ElectronLimiter(theParSet, G4PionMinus::PionMinus());
312  plim->SetFieldCheckFlag(pLimiter);
313  plim->SetTrackingCutPerRegion(reg, rlimH, factE, rmsvE);
314  ph->RegisterProcess(plim, G4PionMinus::PionMinus());
315  }
316  }
317  // enable fluorescence
318  bool fluo = theParSet.getParameter<bool>("FlagFluo");
319  if (fluo && !G4LossTableManager::Instance()->AtomDeexcitation()) {
320  G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
321  G4LossTableManager::Instance()->SetAtomDeexcitation(de);
322  }
323  // change parameters of transportation
324  bool modifyT = theParSet.getParameter<bool>("ModifyTransportation");
325  if (modifyT) {
326  double th1 = theParSet.getUntrackedParameter<double>("ThresholdWarningEnergy") * MeV;
327  double th2 = theParSet.getUntrackedParameter<double>("ThresholdImportantEnergy") * MeV;
328  int nt = theParSet.getUntrackedParameter<int>("ThresholdTrials");
330  }
331  edm::LogVerbatim("SimG4CoreApplication") << "ParametrisedEMPhysics::ConstructProcess() is done";
332 }
333 
334 void ParametrisedEMPhysics::ModifyTransportation(const G4ParticleDefinition* part, int ntry, double th1, double th2) {
335  G4ProcessManager* man = part->GetProcessManager();
336  G4Transportation* trans = (G4Transportation*)((*(man->GetProcessList()))[0]);
337  if (trans) {
338  trans->SetThresholdWarningEnergy(th1);
339  trans->SetThresholdImportantEnergy(th2);
340  trans->SetThresholdTrials(ntry);
341  edm::LogVerbatim("SimG4CoreApplication")
342  << "ParametrisedEMPhysics: printout level changed for " << part->GetParticleName();
343  }
344 }
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