CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
CMSEmStandardPhysics Class Reference

#include <CMSEmStandardPhysics.h>

Inheritance diagram for CMSEmStandardPhysics:

Public Member Functions

 CMSEmStandardPhysics (G4int ver, const edm::ParameterSet &p)
 
void ConstructParticle () override
 
void ConstructProcess () override
 
 ~CMSEmStandardPhysics () override
 

Private Attributes

G4double fGeomFactor
 
G4double fLambdaLimit
 
G4double fRangeFactor
 
G4double fSafetyFactor
 
G4MscStepLimitType fStepLimitType
 

Detailed Description

Definition at line 17 of file CMSEmStandardPhysics.h.

Constructor & Destructor Documentation

◆ CMSEmStandardPhysics()

CMSEmStandardPhysics::CMSEmStandardPhysics ( G4int  ver,
const edm::ParameterSet p 
)

Definition at line 45 of file CMSEmStandardPhysics.cc.

46  : G4VPhysicsConstructor("CMSEmStandard_emm") {
47  SetVerboseLevel(ver);
48  G4EmParameters* param = G4EmParameters::Instance();
49  param->SetDefaults();
50  param->SetVerbose(ver);
51  param->SetApplyCuts(true);
52  param->SetStepFunction(0.8, 1 * CLHEP::mm);
53  param->SetMscRangeFactor(0.2);
54  param->SetMscStepLimitType(fMinimal);
55  param->SetFluo(false);
56  SetPhysicsType(bElectromagnetic);
57  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
58  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
59  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
60  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
61  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
62  fStepLimitType = fUseSafety;
63  if (msc == "UseSafetyPlus") {
64  fStepLimitType = fUseSafetyPlus;
65  }
66  if (msc == "Minimal") {
67  fStepLimitType = fMinimal;
68  }
69 }

References fGeomFactor, fLambdaLimit, fRangeFactor, fSafetyFactor, fStepLimitType, AlCaHLTBitMon_ParallelJobs::p, and AlCaHLTBitMon_QueryRunRegistry::string.

◆ ~CMSEmStandardPhysics()

CMSEmStandardPhysics::~CMSEmStandardPhysics ( )
override

Definition at line 71 of file CMSEmStandardPhysics.cc.

71 {}

Member Function Documentation

◆ ConstructParticle()

void CMSEmStandardPhysics::ConstructParticle ( )
override

Definition at line 73 of file CMSEmStandardPhysics.cc.

73  {
74  // minimal set of particles for EM physics
75  G4EmBuilder::ConstructMinimalEmSet();
76 }

◆ ConstructProcess()

void CMSEmStandardPhysics::ConstructProcess ( )
override

Definition at line 78 of file CMSEmStandardPhysics.cc.

78  {
79  if (verboseLevel > 0) {
80  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct EM Processes";
81  }
82 
83  // This EM builder takes default models of Geant4 10 EMV.
84  // Multiple scattering by WentzelVI for all particles except:
85  // a) e+e- below 100 MeV for which the Urban model is used
86  // b) ions for which Urban model is used
87  G4EmBuilder::PrepareEMPhysics();
88 
89  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
90  // processes used by several particles
91  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
92  G4NuclearStopping* pnuc(nullptr);
93 
94  // high energy limit for e+- scattering models
95  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
96 
97  const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
98  const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
99 
100  // Add gamma EM Processes
101  G4ParticleDefinition* particle = G4Gamma::Gamma();
102 
103  G4PhotoElectricEffect* pee = new G4PhotoElectricEffect();
104 
105  if (G4EmParameters::Instance()->GeneralProcessActive()) {
106  G4GammaGeneralProcess* sp = new G4GammaGeneralProcess();
107  sp->AddEmProcess(pee);
108  sp->AddEmProcess(new G4ComptonScattering());
109  sp->AddEmProcess(new G4GammaConversion());
110  G4LossTableManager::Instance()->SetGammaGeneralProcess(sp);
111  ph->RegisterProcess(sp, particle);
112 
113  } else {
114  ph->RegisterProcess(pee, particle);
115  ph->RegisterProcess(new G4ComptonScattering(), particle);
116  ph->RegisterProcess(new G4GammaConversion(), particle);
117  }
118 
119  // e-
120  particle = G4Electron::Electron();
121 
122  G4eIonisation* eioni = new G4eIonisation();
123 
124  G4eMultipleScattering* msc = new G4eMultipleScattering;
125  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
126  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
127  msc1->SetHighEnergyLimit(highEnergyLimit);
128  msc2->SetLowEnergyLimit(highEnergyLimit);
129  msc->SetEmModel(msc1);
130  msc->SetEmModel(msc2);
131 
132  // e-/e+ msc for HCAL and HGCAL using the Urban model
133  if (nullptr != aRegion || nullptr != bRegion) {
134  G4UrbanMscModel* msc3 = new G4UrbanMscModel();
135  msc3->SetHighEnergyLimit(highEnergyLimit);
136  msc3->SetRangeFactor(fRangeFactor);
137  msc3->SetGeomFactor(fGeomFactor);
138  msc3->SetSafetyFactor(fSafetyFactor);
139  msc3->SetLambdaLimit(fLambdaLimit);
140  msc3->SetStepLimitType(fStepLimitType);
141  msc3->SetLocked(true);
142 
143  if (nullptr != aRegion) {
144  msc->AddEmModel(-1, msc3, aRegion);
145  }
146  if (nullptr != bRegion) {
147  msc->AddEmModel(-1, msc3, bRegion);
148  }
149  }
150 
151  // single scattering
152  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
153  G4CoulombScattering* ss = new G4CoulombScattering();
154  ss->SetEmModel(ssm);
155  ss->SetMinKinEnergy(highEnergyLimit);
156  ssm->SetLowEnergyLimit(highEnergyLimit);
157  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
158 
159  ph->RegisterProcess(msc, particle);
160  ph->RegisterProcess(eioni, particle);
161  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
162  ph->RegisterProcess(ss, particle);
163 
164  // e+
165  particle = G4Positron::Positron();
166  eioni = new G4eIonisation();
167 
168  msc = new G4eMultipleScattering();
169  msc1 = new G4UrbanMscModel();
170  msc2 = new G4WentzelVIModel();
171  msc1->SetHighEnergyLimit(highEnergyLimit);
172  msc2->SetLowEnergyLimit(highEnergyLimit);
173  msc->SetEmModel(msc1);
174  msc->SetEmModel(msc2);
175 
176  // e-/e+ msc for HCAL and HGCAL using the Urban model
177  if (nullptr != aRegion || nullptr != bRegion) {
178  G4UrbanMscModel* msc3 = new G4UrbanMscModel();
179  msc3->SetHighEnergyLimit(highEnergyLimit);
180  msc3->SetRangeFactor(fRangeFactor);
181  msc3->SetGeomFactor(fGeomFactor);
182  msc3->SetSafetyFactor(fSafetyFactor);
183  msc3->SetLambdaLimit(fLambdaLimit);
184  msc3->SetStepLimitType(fStepLimitType);
185  msc3->SetLocked(true);
186 
187  if (nullptr != aRegion) {
188  msc->AddEmModel(-1, msc3, aRegion);
189  }
190  if (nullptr != bRegion) {
191  msc->AddEmModel(-1, msc3, bRegion);
192  }
193  }
194 
195  // single scattering
196  ssm = new G4eCoulombScatteringModel();
197  ss = new G4CoulombScattering();
198  ss->SetEmModel(ssm);
199  ss->SetMinKinEnergy(highEnergyLimit);
200  ssm->SetLowEnergyLimit(highEnergyLimit);
201  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
202 
203  ph->RegisterProcess(msc, particle);
204  ph->RegisterProcess(eioni, particle);
205  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
206  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
207  ph->RegisterProcess(ss, particle);
208 
209  // generic ion
210  particle = G4GenericIon::GenericIon();
211  G4ionIonisation* ionIoni = new G4ionIonisation();
212  ph->RegisterProcess(hmsc, particle);
213  ph->RegisterProcess(ionIoni, particle);
214 
215  // muons, hadrons ions
216  G4EmBuilder::ConstructCharged(hmsc, pnuc);
217 }

References nanoDQM_cff::Electron, fGeomFactor, fLambdaLimit, fRangeFactor, fSafetyFactor, fStepLimitType, contentValuesCheck::ss, and HLT_FULL_cff::verboseLevel.

Member Data Documentation

◆ fGeomFactor

G4double CMSEmStandardPhysics::fGeomFactor
private

Definition at line 27 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fLambdaLimit

G4double CMSEmStandardPhysics::fLambdaLimit
private

Definition at line 29 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fRangeFactor

G4double CMSEmStandardPhysics::fRangeFactor
private

Definition at line 26 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fSafetyFactor

G4double CMSEmStandardPhysics::fSafetyFactor
private

Definition at line 28 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fStepLimitType

G4MscStepLimitType CMSEmStandardPhysics::fStepLimitType
private

Definition at line 30 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

HLT_FULL_cff.verboseLevel
verboseLevel
Definition: HLT_FULL_cff.py:8514
CMSEmStandardPhysics::fLambdaLimit
G4double fLambdaLimit
Definition: CMSEmStandardPhysics.h:29
contentValuesCheck.ss
ss
Definition: contentValuesCheck.py:33
CMSEmStandardPhysics::fSafetyFactor
G4double fSafetyFactor
Definition: CMSEmStandardPhysics.h:28
CMSEmStandardPhysics::fGeomFactor
G4double fGeomFactor
Definition: CMSEmStandardPhysics.h:27
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
CMSEmStandardPhysics::fRangeFactor
G4double fRangeFactor
Definition: CMSEmStandardPhysics.h:26
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
nanoDQM_cff.Electron
Electron
Definition: nanoDQM_cff.py:95
CMSEmStandardPhysics::fStepLimitType
G4MscStepLimitType fStepLimitType
Definition: CMSEmStandardPhysics.h:30
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128