CMS 3D CMS Logo

CMSEmStandardPhysicsEMH.cc
Go to the documentation of this file.
3 
4 #include "G4SystemOfUnits.hh"
5 #include "G4ParticleDefinition.hh"
6 #include "G4EmParameters.hh"
7 #include "G4EmBuilder.hh"
8 
9 #include "G4MscStepLimitType.hh"
10 
11 #include "G4hIonisation.hh"
12 #include "G4hMultipleScattering.hh"
13 #include "G4ionIonisation.hh"
14 
15 #include "G4ParticleTable.hh"
16 #include "G4Gamma.hh"
17 #include "G4Electron.hh"
18 #include "G4Positron.hh"
19 #include "G4GenericIon.hh"
20 
21 #include "G4PhysicsListHelper.hh"
22 #include "G4BuilderType.hh"
23 #include "G4ProcessManager.hh"
24 
25 #include "G4HepEmProcess.hh"
26 
27 #include <string>
28 
30  : G4VPhysicsConstructor("CMSEmStandard_emh") {
31  SetVerboseLevel(ver);
32  G4EmParameters* param = G4EmParameters::Instance();
33  param->SetDefaults();
34  param->SetVerbose(ver);
35  param->SetApplyCuts(true);
36  param->SetStepFunction(0.8, 1 * CLHEP::mm);
37  param->SetMscRangeFactor(0.2);
38  param->SetMscStepLimitType(fUseSafety);
39  param->SetFluo(false);
40  SetPhysicsType(bElectromagnetic);
41  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
42  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
43  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
44  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
45  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
46  fStepLimitType = fUseSafety;
47  if (msc == "UseSafetyPlus") {
48  fStepLimitType = fUseSafetyPlus;
49  }
50  if (msc == "Minimal") {
51  fStepLimitType = fMinimal;
52  }
53  double tcut = p.getParameter<double>("G4TrackingCut") * CLHEP::MeV;
54  param->SetLowestElectronEnergy(tcut);
55  param->SetLowestMuHadEnergy(tcut);
56 }
57 
59 
61  // minimal set of particles for EM physics
62  G4EmBuilder::ConstructMinimalEmSet();
63 }
64 
66  if (verboseLevel > 0) {
67  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct EM Processes";
68  }
69 
70  // This EM builder takes default models of Geant4 10 EMV.
71  // Multiple scattering by WentzelVI for all particles except:
72  // a) e+e- below 100 MeV for which the Urban model is used
73  // b) ions for which Urban model is used
74  G4EmBuilder::PrepareEMPhysics();
75 
76  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
77  // processes used by several particles
78  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
79  G4NuclearStopping* pnuc(nullptr);
80 
81  G4HepEmProcess* hepEmProcess = new G4HepEmProcess();
82  G4Electron::Electron()->GetProcessManager()->AddProcess(hepEmProcess, -1, -1, 1);
83  G4Positron::Positron()->GetProcessManager()->AddProcess(hepEmProcess, -1, -1, 1);
84  G4Gamma::Gamma()->GetProcessManager()->AddProcess(hepEmProcess, -1, -1, 1);
85 
86  // generic ion
87  G4ParticleDefinition* particle = G4GenericIon::GenericIon();
88  G4ionIonisation* ionIoni = new G4ionIonisation();
89  ph->RegisterProcess(hmsc, particle);
90  ph->RegisterProcess(ionIoni, particle);
91 
92  // muons, hadrons ions
93  G4EmBuilder::ConstructCharged(hmsc, pnuc);
94 }
G4MscStepLimitType fStepLimitType
Log< level::Info, true > LogVerbatim
CMSEmStandardPhysicsEMH(G4int ver, const edm::ParameterSet &p)