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