CMS 3D CMS Logo

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

#include <CMSEmStandardPhysicsTrackingManager.h>

Inheritance diagram for CMSEmStandardPhysicsTrackingManager:

Public Member Functions

void BuildPhysicsTable (const G4ParticleDefinition &) override
 
 CMSEmStandardPhysicsTrackingManager (const edm::ParameterSet &p)
 
void HandOverOneTrack (G4Track *aTrack) override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
 ~CMSEmStandardPhysicsTrackingManager () override
 

Private Member Functions

void TrackElectron (G4Track *aTrack)
 
void TrackGamma (G4Track *aTrack)
 
void TrackPositron (G4Track *aTrack)
 

Private Attributes

struct {
   G4eBremsstrahlung *   brems
 
   G4eIonisation *   ioni
 
   G4eMultipleScattering *   msc
 
   G4ElectronNuclearProcess *   nuc
 
   G4CoulombScattering *   ss
 
electron
 
G4double fGeomFactor
 
G4double fLambdaLimit
 
G4double fRangeFactor
 
G4double fSafetyFactor
 
G4MscStepLimitType fStepLimitType
 
G4GammaGeneralProcess * gammaProc
 
struct {
   G4eplusAnnihilation *   annihilation
 
   G4eBremsstrahlung *   brems
 
   G4eIonisation *   ioni
 
   G4eMultipleScattering *   msc
 
   G4PositronNuclearProcess *   nuc
 
   G4CoulombScattering *   ss
 
positron
 

Static Private Attributes

static CMSEmStandardPhysicsTrackingManagermasterTrackingManager = nullptr
 

Detailed Description

Definition at line 20 of file CMSEmStandardPhysicsTrackingManager.h.

Constructor & Destructor Documentation

◆ CMSEmStandardPhysicsTrackingManager()

CMSEmStandardPhysicsTrackingManager::CMSEmStandardPhysicsTrackingManager ( const edm::ParameterSet p)

Definition at line 51 of file CMSEmStandardPhysicsTrackingManager.cc.

References annihilation, brems, electron, fGeomFactor, fLambdaLimit, fRangeFactor, fSafetyFactor, fStepLimitType, g4SimHits_cfi::G4GammaGeneralProcess, gammaProc, ioni, masterTrackingManager, msc, nuc, AlCaHLTBitMon_ParallelJobs::p, positron, ss, and AlCaHLTBitMon_QueryRunRegistry::string.

51  {
52  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
53  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
54  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
55  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
56  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
57  fStepLimitType = fUseSafety;
58  if (msc == "UseSafetyPlus") {
59  fStepLimitType = fUseSafetyPlus;
60  } else if (msc == "Minimal") {
61  fStepLimitType = fMinimal;
62  }
63 
64  G4EmParameters *param = G4EmParameters::Instance();
65  G4double highEnergyLimit = param->MscEnergyLimit();
66 
67  const G4Region *aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
68  const G4Region *bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
69 
70  // e-
71  {
72  G4eMultipleScattering *msc = new G4eMultipleScattering;
73  G4UrbanMscModel *msc1 = new G4UrbanMscModel;
74  G4WentzelVIModel *msc2 = new G4WentzelVIModel;
75  msc1->SetHighEnergyLimit(highEnergyLimit);
76  msc2->SetLowEnergyLimit(highEnergyLimit);
77  msc->SetEmModel(msc1);
78  msc->SetEmModel(msc2);
79 
80  // e-/e+ msc for HCAL and HGCAL using the Urban model
81  if (nullptr != aRegion || nullptr != bRegion) {
82  G4UrbanMscModel *msc3 = new G4UrbanMscModel();
83  msc3->SetHighEnergyLimit(highEnergyLimit);
84  msc3->SetRangeFactor(fRangeFactor);
85  msc3->SetGeomFactor(fGeomFactor);
86  msc3->SetSafetyFactor(fSafetyFactor);
87  msc3->SetLambdaLimit(fLambdaLimit);
88  msc3->SetStepLimitType(fStepLimitType);
89  msc3->SetLocked(true);
90  if (nullptr != aRegion) {
91  msc->AddEmModel(-1, msc3, aRegion);
92  }
93  if (nullptr != bRegion) {
94  msc->AddEmModel(-1, msc3, bRegion);
95  }
96  }
97 
98  electron.msc = msc;
99 
100  electron.ioni = new G4eIonisation;
101  electron.brems = new G4eBremsstrahlung;
102 
103  G4CoulombScattering *ss = new G4CoulombScattering;
104  G4eCoulombScatteringModel *ssm = new G4eCoulombScatteringModel;
105  ssm->SetLowEnergyLimit(highEnergyLimit);
106  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
107  ss->SetEmModel(ssm);
108  ss->SetMinKinEnergy(highEnergyLimit);
109  electron.ss = ss;
110  }
111 
112  // e+
113  {
114  G4eMultipleScattering *msc = new G4eMultipleScattering;
115  G4UrbanMscModel *msc1 = new G4UrbanMscModel;
116  G4WentzelVIModel *msc2 = new G4WentzelVIModel;
117  msc1->SetHighEnergyLimit(highEnergyLimit);
118  msc2->SetLowEnergyLimit(highEnergyLimit);
119  msc->SetEmModel(msc1);
120  msc->SetEmModel(msc2);
121 
122  // e-/e+ msc for HCAL and HGCAL using the Urban model
123  if (nullptr != aRegion || nullptr != bRegion) {
124  G4UrbanMscModel *msc3 = new G4UrbanMscModel();
125  msc3->SetHighEnergyLimit(highEnergyLimit);
126  msc3->SetRangeFactor(fRangeFactor);
127  msc3->SetGeomFactor(fGeomFactor);
128  msc3->SetSafetyFactor(fSafetyFactor);
129  msc3->SetLambdaLimit(fLambdaLimit);
130  msc3->SetStepLimitType(fStepLimitType);
131  msc3->SetLocked(true);
132  if (nullptr != aRegion) {
133  msc->AddEmModel(-1, msc3, aRegion);
134  }
135  if (nullptr != bRegion) {
136  msc->AddEmModel(-1, msc3, bRegion);
137  }
138  }
139 
140  positron.msc = msc;
141 
142  positron.ioni = new G4eIonisation;
143  positron.brems = new G4eBremsstrahlung;
144  positron.annihilation = new G4eplusAnnihilation;
145 
146  G4CoulombScattering *ss = new G4CoulombScattering;
147  G4eCoulombScatteringModel *ssm = new G4eCoulombScatteringModel;
148  ssm->SetLowEnergyLimit(highEnergyLimit);
149  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
150  ss->SetEmModel(ssm);
151  ss->SetMinKinEnergy(highEnergyLimit);
152  positron.ss = ss;
153  }
154 
155  // gamma
156  {
158 
159  gammaProc->AddEmProcess(new G4PhotoElectricEffect);
160  gammaProc->AddEmProcess(new G4ComptonScattering);
161  gammaProc->AddEmProcess(new G4GammaConversion);
162 
163  G4HadronInelasticProcess *nuc = new G4HadronInelasticProcess("photonNuclear", G4Gamma::Definition());
164  auto xsreg = G4CrossSectionDataSetRegistry::Instance();
165  G4VCrossSectionDataSet *xs = nullptr;
166  bool useGammaNuclearXS = true;
167  if (useGammaNuclearXS) {
168  xs = xsreg->GetCrossSectionDataSet("GammaNuclearXS");
169  if (nullptr == xs)
170  xs = new G4GammaNuclearXS;
171  } else {
172  xs = xsreg->GetCrossSectionDataSet("PhotoNuclearXS");
173  if (nullptr == xs)
174  xs = new G4PhotoNuclearCrossSection;
175  }
176  nuc->AddDataSet(xs);
177 
178  G4QGSModel<G4GammaParticipants> *theStringModel = new G4QGSModel<G4GammaParticipants>;
179  G4QGSMFragmentation *theFrag = new G4QGSMFragmentation;
180  G4ExcitedStringDecay *theStringDecay = new G4ExcitedStringDecay(theFrag);
181  theStringModel->SetFragmentationModel(theStringDecay);
182 
183  G4GeneratorPrecompoundInterface *theCascade = new G4GeneratorPrecompoundInterface;
184 
185  G4TheoFSGenerator *theModel = new G4TheoFSGenerator;
186  theModel->SetTransport(theCascade);
187  theModel->SetHighEnergyGenerator(theStringModel);
188 
189  G4HadronicParameters *param = G4HadronicParameters::Instance();
190 
191  G4CascadeInterface *cascade = new G4CascadeInterface;
192 
193  // added low-energy model LEND disabled
194  G4double gnLowEnergyLimit = 200 * CLHEP::MeV;
195  if (gnLowEnergyLimit > 0.0) {
196  G4LowEGammaNuclearModel *lemod = new G4LowEGammaNuclearModel;
197  lemod->SetMaxEnergy(gnLowEnergyLimit);
198  nuc->RegisterMe(lemod);
199  cascade->SetMinEnergy(gnLowEnergyLimit - CLHEP::MeV);
200  }
201  cascade->SetMaxEnergy(param->GetMaxEnergyTransitionFTF_Cascade());
202  nuc->RegisterMe(cascade);
203  theModel->SetMinEnergy(param->GetMinEnergyTransitionFTF_Cascade());
204  theModel->SetMaxEnergy(param->GetMaxEnergy());
205  nuc->RegisterMe(theModel);
206 
207  gammaProc->AddHadProcess(nuc);
208  G4LossTableManager::Instance()->SetGammaGeneralProcess(gammaProc);
209  }
210 
211  // Create lepto-nuclear processes last, they access cross section data from
212  // the gamma nuclear process!
213  G4ElectroVDNuclearModel *eModel = new G4ElectroVDNuclearModel;
214 
215  {
216  G4ElectronNuclearProcess *nuc = new G4ElectronNuclearProcess;
217  nuc->RegisterMe(eModel);
218  electron.nuc = nuc;
219  }
220  {
221  G4PositronNuclearProcess *nuc = new G4PositronNuclearProcess;
222  nuc->RegisterMe(eModel);
223  positron.nuc = nuc;
224  }
225 
226  if (masterTrackingManager == nullptr) {
227  masterTrackingManager = this;
228  } else {
229  electron.msc->SetMasterProcess(masterTrackingManager->electron.msc);
230  electron.ss->SetMasterProcess(masterTrackingManager->electron.ss);
231  electron.ioni->SetMasterProcess(masterTrackingManager->electron.ioni);
232  electron.brems->SetMasterProcess(masterTrackingManager->electron.brems);
233  electron.nuc->SetMasterProcess(masterTrackingManager->electron.nuc);
234 
235  positron.msc->SetMasterProcess(masterTrackingManager->positron.msc);
236  positron.ss->SetMasterProcess(masterTrackingManager->positron.ss);
237  positron.ioni->SetMasterProcess(masterTrackingManager->positron.ioni);
238  positron.brems->SetMasterProcess(masterTrackingManager->positron.brems);
239  positron.annihilation->SetMasterProcess(masterTrackingManager->positron.annihilation);
240  positron.nuc->SetMasterProcess(masterTrackingManager->positron.nuc);
241 
242  gammaProc->SetMasterProcess(masterTrackingManager->gammaProc);
243  }
244 }
struct CMSEmStandardPhysicsTrackingManager::@1036 positron
struct CMSEmStandardPhysicsTrackingManager::@1035 electron
static CMSEmStandardPhysicsTrackingManager * masterTrackingManager

◆ ~CMSEmStandardPhysicsTrackingManager()

CMSEmStandardPhysicsTrackingManager::~CMSEmStandardPhysicsTrackingManager ( )
override

Definition at line 246 of file CMSEmStandardPhysicsTrackingManager.cc.

References masterTrackingManager.

246  {
247  if (masterTrackingManager == this) {
248  masterTrackingManager = nullptr;
249  }
250 }
static CMSEmStandardPhysicsTrackingManager * masterTrackingManager

Member Function Documentation

◆ BuildPhysicsTable()

void CMSEmStandardPhysicsTrackingManager::BuildPhysicsTable ( const G4ParticleDefinition &  part)
override

Definition at line 252 of file CMSEmStandardPhysicsTrackingManager.cc.

References electron, gammaProc, and positron.

252  {
253  if (&part == G4Electron::Definition()) {
254  electron.msc->BuildPhysicsTable(part);
255  electron.ioni->BuildPhysicsTable(part);
256  electron.brems->BuildPhysicsTable(part);
257  electron.ss->BuildPhysicsTable(part);
258  electron.nuc->BuildPhysicsTable(part);
259  } else if (&part == G4Positron::Definition()) {
260  positron.msc->BuildPhysicsTable(part);
261  positron.ioni->BuildPhysicsTable(part);
262  positron.brems->BuildPhysicsTable(part);
263  positron.annihilation->BuildPhysicsTable(part);
264  positron.ss->BuildPhysicsTable(part);
265  positron.nuc->BuildPhysicsTable(part);
266  } else if (&part == G4Gamma::Definition()) {
267  gammaProc->BuildPhysicsTable(part);
268  }
269 }
struct CMSEmStandardPhysicsTrackingManager::@1036 positron
struct CMSEmStandardPhysicsTrackingManager::@1035 electron
part
Definition: HCALResponse.h:20

◆ HandOverOneTrack()

void CMSEmStandardPhysicsTrackingManager::HandOverOneTrack ( G4Track *  aTrack)
override

Definition at line 712 of file CMSEmStandardPhysicsTrackingManager.cc.

References TrackElectron(), TrackGamma(), and TrackPositron().

712  {
713  const G4ParticleDefinition *part = aTrack->GetParticleDefinition();
714 
715  if (part == G4Electron::Definition()) {
716  TrackElectron(aTrack);
717  } else if (part == G4Positron::Definition()) {
718  TrackPositron(aTrack);
719  } else if (part == G4Gamma::Definition()) {
720  TrackGamma(aTrack);
721  }
722 
723  aTrack->SetTrackStatus(fStopAndKill);
724  delete aTrack;
725 }
part
Definition: HCALResponse.h:20

◆ PreparePhysicsTable()

void CMSEmStandardPhysicsTrackingManager::PreparePhysicsTable ( const G4ParticleDefinition &  part)
override

Definition at line 271 of file CMSEmStandardPhysicsTrackingManager.cc.

References electron, gammaProc, and positron.

271  {
272  if (&part == G4Electron::Definition()) {
273  electron.msc->PreparePhysicsTable(part);
274  electron.ioni->PreparePhysicsTable(part);
275  electron.brems->PreparePhysicsTable(part);
276  electron.ss->PreparePhysicsTable(part);
277  electron.nuc->PreparePhysicsTable(part);
278  } else if (&part == G4Positron::Definition()) {
279  positron.msc->PreparePhysicsTable(part);
280  positron.ioni->PreparePhysicsTable(part);
281  positron.brems->PreparePhysicsTable(part);
282  positron.annihilation->PreparePhysicsTable(part);
283  positron.ss->PreparePhysicsTable(part);
284  positron.nuc->PreparePhysicsTable(part);
285  } else if (&part == G4Gamma::Definition()) {
286  gammaProc->PreparePhysicsTable(part);
287  }
288 }
struct CMSEmStandardPhysicsTrackingManager::@1036 positron
struct CMSEmStandardPhysicsTrackingManager::@1035 electron
part
Definition: HCALResponse.h:20

◆ TrackElectron()

void CMSEmStandardPhysicsTrackingManager::TrackElectron ( G4Track *  aTrack)
private

Definition at line 290 of file CMSEmStandardPhysicsTrackingManager.cc.

References electron, mps_fire::i, EcalCondDBWriter_cfi::physics, LaserDQM_cfg::process, corrVsCorr::selection, HLT_2023v12_cff::track, and TrackingManagerHelper::TrackChargedParticle().

Referenced by HandOverOneTrack().

290  {
291  class ElectronPhysics final : public TrackingManagerHelper::Physics {
292  public:
293  ElectronPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
294 
295  void StartTracking(G4Track *aTrack) override {
296  auto &electron = fMgr.electron;
297 
298  electron.msc->StartTracking(aTrack);
299  electron.ioni->StartTracking(aTrack);
300  electron.brems->StartTracking(aTrack);
301  electron.ss->StartTracking(aTrack);
302  electron.nuc->StartTracking(aTrack);
303 
304  fPreviousStepLength = 0;
305  }
306  void EndTracking() override {
307  auto &electron = fMgr.electron;
308 
309  electron.msc->EndTracking();
310  electron.ioni->EndTracking();
311  electron.brems->EndTracking();
312  electron.ss->EndTracking();
313  electron.nuc->EndTracking();
314  }
315 
316  G4double GetPhysicalInteractionLength(const G4Track &track) override {
317  auto &electron = fMgr.electron;
318  G4double physIntLength, proposedSafety = DBL_MAX;
319  G4ForceCondition condition;
320  G4GPILSelection selection;
321 
322  fProposedStep = DBL_MAX;
323  fSelected = -1;
324 
325  physIntLength = electron.nuc->PostStepGPIL(track, fPreviousStepLength, &condition);
326  if (physIntLength < fProposedStep) {
327  fProposedStep = physIntLength;
328  fSelected = 3;
329  }
330 
331  physIntLength = electron.ss->PostStepGPIL(track, fPreviousStepLength, &condition);
332  if (physIntLength < fProposedStep) {
333  fProposedStep = physIntLength;
334  fSelected = 0;
335  }
336 
337  physIntLength = electron.brems->PostStepGPIL(track, fPreviousStepLength, &condition);
338  if (physIntLength < fProposedStep) {
339  fProposedStep = physIntLength;
340  fSelected = 1;
341  }
342 
343  physIntLength = electron.ioni->PostStepGPIL(track, fPreviousStepLength, &condition);
344  if (physIntLength < fProposedStep) {
345  fProposedStep = physIntLength;
346  fSelected = 2;
347  }
348 
349  physIntLength =
350  electron.ioni->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
351  if (physIntLength < fProposedStep) {
352  fProposedStep = physIntLength;
353  fSelected = -1;
354  }
355 
356  physIntLength =
357  electron.msc->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
358  if (physIntLength < fProposedStep) {
359  fProposedStep = physIntLength;
360  // Check if MSC actually wants to win, in most cases it only limits the
361  // step size.
362  if (selection == CandidateForSelection) {
363  fSelected = -1;
364  }
365  }
366 
367  return fProposedStep;
368  }
369 
370  void AlongStepDoIt(G4Track &track, G4Step &step, G4TrackVector &) override {
371  if (step.GetStepLength() == fProposedStep) {
372  step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
373  } else {
374  // Remember that the step was limited by geometry.
375  fSelected = -1;
376  }
377  auto &electron = fMgr.electron;
378  G4VParticleChange *particleChange;
379 
380  particleChange = electron.msc->AlongStepDoIt(track, step);
381  particleChange->UpdateStepForAlongStep(&step);
382  track.SetTrackStatus(particleChange->GetTrackStatus());
383  particleChange->Clear();
384 
385  particleChange = electron.ioni->AlongStepDoIt(track, step);
386  particleChange->UpdateStepForAlongStep(&step);
387  track.SetTrackStatus(particleChange->GetTrackStatus());
388  particleChange->Clear();
389 
390  fPreviousStepLength = step.GetStepLength();
391  }
392 
393  void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
394  if (fSelected < 0) {
395  return;
396  }
397  step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
398 
399  auto &electron = fMgr.electron;
400  G4VProcess *process = nullptr;
401  G4VParticleChange *particleChange = nullptr;
402 
403  switch (fSelected) {
404  case 3:
405  process = electron.nuc;
406  particleChange = electron.nuc->PostStepDoIt(track, step);
407  break;
408  case 0:
409  process = electron.ss;
410  particleChange = electron.ss->PostStepDoIt(track, step);
411  break;
412  case 1:
413  process = electron.brems;
414  particleChange = electron.brems->PostStepDoIt(track, step);
415  break;
416  case 2:
417  process = electron.ioni;
418  particleChange = electron.ioni->PostStepDoIt(track, step);
419  break;
420  }
421 
422  particleChange->UpdateStepForPostStep(&step);
423  step.UpdateTrack();
424 
425  int numSecondaries = particleChange->GetNumberOfSecondaries();
426  for (int i = 0; i < numSecondaries; i++) {
427  G4Track *secondary = particleChange->GetSecondary(i);
428  secondary->SetParentID(track.GetTrackID());
429  secondary->SetCreatorProcess(process);
430  secondaries.push_back(secondary);
431  }
432 
433  track.SetTrackStatus(particleChange->GetTrackStatus());
434  particleChange->Clear();
435  }
436 
437  private:
439  G4double fPreviousStepLength;
440  G4double fProposedStep;
441  G4int fSelected;
442  };
443 
444  ElectronPhysics physics(*this);
446 }
selection
main part
Definition: corrVsCorr.py:100
struct CMSEmStandardPhysicsTrackingManager::@1035 electron
static void TrackChargedParticle(G4Track *aTrack, PhysicsImpl &physics)
step
Definition: StallMonitor.cc:98

◆ TrackGamma()

void CMSEmStandardPhysicsTrackingManager::TrackGamma ( G4Track *  aTrack)
private

Definition at line 639 of file CMSEmStandardPhysicsTrackingManager.cc.

References mps_fire::i, EcalCondDBWriter_cfi::physics, LaserDQM_cfg::process, HLT_2023v12_cff::track, and TrackingManagerHelper::TrackNeutralParticle().

Referenced by HandOverOneTrack().

639  {
640  class GammaPhysics final : public TrackingManagerHelper::Physics {
641  public:
642  GammaPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
643 
644  void StartTracking(G4Track *aTrack) override {
645  fMgr.gammaProc->StartTracking(aTrack);
646 
647  fPreviousStepLength = 0;
648  }
649  void EndTracking() override { fMgr.gammaProc->EndTracking(); }
650 
651  G4double GetPhysicalInteractionLength(const G4Track &track) override {
652  G4double physIntLength;
653  G4ForceCondition condition;
654 
655  fProposedStep = DBL_MAX;
656  fSelected = -1;
657 
658  physIntLength = fMgr.gammaProc->PostStepGPIL(track, fPreviousStepLength, &condition);
659  if (physIntLength < fProposedStep) {
660  fProposedStep = physIntLength;
661  fSelected = 0;
662  }
663 
664  return fProposedStep;
665  }
666 
667  void AlongStepDoIt(G4Track &, G4Step &step, G4TrackVector &) override {
668  if (step.GetStepLength() == fProposedStep) {
669  step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
670  } else {
671  // Remember that the step was limited by geometry.
672  fSelected = -1;
673  }
674  fPreviousStepLength = step.GetStepLength();
675  }
676 
677  void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
678  if (fSelected < 0) {
679  return;
680  }
681  step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
682 
683  G4VProcess *process = fMgr.gammaProc;
684  G4VParticleChange *particleChange = fMgr.gammaProc->PostStepDoIt(track, step);
685 
686  particleChange->UpdateStepForPostStep(&step);
687  step.UpdateTrack();
688 
689  int numSecondaries = particleChange->GetNumberOfSecondaries();
690  for (int i = 0; i < numSecondaries; i++) {
691  G4Track *secondary = particleChange->GetSecondary(i);
692  secondary->SetParentID(track.GetTrackID());
693  secondary->SetCreatorProcess(process);
694  secondaries.push_back(secondary);
695  }
696 
697  track.SetTrackStatus(particleChange->GetTrackStatus());
698  particleChange->Clear();
699  }
700 
701  private:
703  G4double fPreviousStepLength;
704  G4double fProposedStep;
705  G4int fSelected;
706  };
707 
708  GammaPhysics physics(*this);
710 }
static void TrackNeutralParticle(G4Track *aTrack, PhysicsImpl &physics)
step
Definition: StallMonitor.cc:98

◆ TrackPositron()

void CMSEmStandardPhysicsTrackingManager::TrackPositron ( G4Track *  aTrack)
private

Definition at line 448 of file CMSEmStandardPhysicsTrackingManager.cc.

References mps_fire::i, EcalCondDBWriter_cfi::physics, positron, LaserDQM_cfg::process, corrVsCorr::selection, HLT_2023v12_cff::track, and TrackingManagerHelper::TrackChargedParticle().

Referenced by HandOverOneTrack().

448  {
449  class PositronPhysics final : public TrackingManagerHelper::Physics {
450  public:
451  PositronPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
452 
453  void StartTracking(G4Track *aTrack) override {
454  auto &positron = fMgr.positron;
455 
456  positron.msc->StartTracking(aTrack);
457  positron.ioni->StartTracking(aTrack);
458  positron.brems->StartTracking(aTrack);
459  positron.annihilation->StartTracking(aTrack);
460  positron.ss->StartTracking(aTrack);
461  positron.nuc->StartTracking(aTrack);
462 
463  fPreviousStepLength = 0;
464  }
465  void EndTracking() override {
466  auto &positron = fMgr.positron;
467 
468  positron.msc->EndTracking();
469  positron.ioni->EndTracking();
470  positron.brems->EndTracking();
471  positron.annihilation->EndTracking();
472  positron.ss->EndTracking();
473  positron.nuc->EndTracking();
474  }
475 
476  G4double GetPhysicalInteractionLength(const G4Track &track) override {
477  auto &positron = fMgr.positron;
478  G4double physIntLength, proposedSafety = DBL_MAX;
479  G4ForceCondition condition;
480  G4GPILSelection selection;
481 
482  fProposedStep = DBL_MAX;
483  fSelected = -1;
484 
485  physIntLength = positron.nuc->PostStepGPIL(track, fPreviousStepLength, &condition);
486  if (physIntLength < fProposedStep) {
487  fProposedStep = physIntLength;
488  fSelected = 4;
489  }
490 
491  physIntLength = positron.ss->PostStepGPIL(track, fPreviousStepLength, &condition);
492  if (physIntLength < fProposedStep) {
493  fProposedStep = physIntLength;
494  fSelected = 0;
495  }
496 
497  physIntLength = positron.annihilation->PostStepGPIL(track, fPreviousStepLength, &condition);
498  if (physIntLength < fProposedStep) {
499  fProposedStep = physIntLength;
500  fSelected = 1;
501  }
502 
503  physIntLength = positron.brems->PostStepGPIL(track, fPreviousStepLength, &condition);
504  if (physIntLength < fProposedStep) {
505  fProposedStep = physIntLength;
506  fSelected = 2;
507  }
508 
509  physIntLength = positron.ioni->PostStepGPIL(track, fPreviousStepLength, &condition);
510  if (physIntLength < fProposedStep) {
511  fProposedStep = physIntLength;
512  fSelected = 3;
513  }
514 
515  physIntLength =
516  positron.ioni->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
517  if (physIntLength < fProposedStep) {
518  fProposedStep = physIntLength;
519  fSelected = -1;
520  }
521 
522  physIntLength =
523  positron.msc->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
524  if (physIntLength < fProposedStep) {
525  fProposedStep = physIntLength;
526  // Check if MSC actually wants to win, in most cases it only limits the
527  // step size.
528  if (selection == CandidateForSelection) {
529  fSelected = -1;
530  }
531  }
532 
533  return fProposedStep;
534  }
535 
536  void AlongStepDoIt(G4Track &track, G4Step &step, G4TrackVector &) override {
537  if (step.GetStepLength() == fProposedStep) {
538  step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
539  } else {
540  // Remember that the step was limited by geometry.
541  fSelected = -1;
542  }
543  auto &positron = fMgr.positron;
544  G4VParticleChange *particleChange;
545 
546  particleChange = positron.msc->AlongStepDoIt(track, step);
547  particleChange->UpdateStepForAlongStep(&step);
548  track.SetTrackStatus(particleChange->GetTrackStatus());
549  particleChange->Clear();
550 
551  particleChange = positron.ioni->AlongStepDoIt(track, step);
552  particleChange->UpdateStepForAlongStep(&step);
553  track.SetTrackStatus(particleChange->GetTrackStatus());
554  particleChange->Clear();
555 
556  fPreviousStepLength = step.GetStepLength();
557  }
558 
559  void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
560  if (fSelected < 0) {
561  return;
562  }
563  step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
564 
565  auto &positron = fMgr.positron;
566  G4VProcess *process;
567  G4VParticleChange *particleChange = nullptr;
568 
569  switch (fSelected) {
570  case 4:
571  process = positron.nuc;
572  particleChange = positron.nuc->PostStepDoIt(track, step);
573  break;
574  case 0:
575  process = positron.ss;
576  particleChange = positron.ss->PostStepDoIt(track, step);
577  break;
578  case 1:
579  process = positron.annihilation;
580  particleChange = positron.annihilation->PostStepDoIt(track, step);
581  break;
582  case 2:
583  process = positron.brems;
584  particleChange = positron.brems->PostStepDoIt(track, step);
585  break;
586  case 3:
587  process = positron.ioni;
588  particleChange = positron.ioni->PostStepDoIt(track, step);
589  break;
590  }
591 
592  particleChange->UpdateStepForPostStep(&step);
593  step.UpdateTrack();
594 
595  int numSecondaries = particleChange->GetNumberOfSecondaries();
596  for (int i = 0; i < numSecondaries; i++) {
597  G4Track *secondary = particleChange->GetSecondary(i);
598  secondary->SetParentID(track.GetTrackID());
599  secondary->SetCreatorProcess(process);
600  secondaries.push_back(secondary);
601  }
602 
603  track.SetTrackStatus(particleChange->GetTrackStatus());
604  particleChange->Clear();
605  }
606 
607  G4bool HasAtRestProcesses() override { return true; }
608 
609  void AtRestDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
610  auto &positron = fMgr.positron;
611  // Annihilate the positron at rest.
612  G4VParticleChange *particleChange = positron.annihilation->AtRestDoIt(track, step);
613  particleChange->UpdateStepForAtRest(&step);
614  step.UpdateTrack();
615 
616  int numSecondaries = particleChange->GetNumberOfSecondaries();
617  for (int i = 0; i < numSecondaries; i++) {
618  G4Track *secondary = particleChange->GetSecondary(i);
619  secondary->SetParentID(track.GetTrackID());
620  secondary->SetCreatorProcess(positron.annihilation);
621  secondaries.push_back(secondary);
622  }
623 
624  track.SetTrackStatus(particleChange->GetTrackStatus());
625  particleChange->Clear();
626  }
627 
628  private:
630  G4double fPreviousStepLength;
631  G4double fProposedStep;
632  G4int fSelected;
633  };
634 
635  PositronPhysics physics(*this);
637 }
struct CMSEmStandardPhysicsTrackingManager::@1036 positron
selection
main part
Definition: corrVsCorr.py:100
static void TrackChargedParticle(G4Track *aTrack, PhysicsImpl &physics)
step
Definition: StallMonitor.cc:98

Member Data Documentation

◆ annihilation

G4eplusAnnihilation* CMSEmStandardPhysicsTrackingManager::annihilation

◆ brems

G4eBremsstrahlung* CMSEmStandardPhysicsTrackingManager::brems

◆ electron

struct { ... } CMSEmStandardPhysicsTrackingManager::electron

◆ fGeomFactor

G4double CMSEmStandardPhysicsTrackingManager::fGeomFactor
private

◆ fLambdaLimit

G4double CMSEmStandardPhysicsTrackingManager::fLambdaLimit
private

◆ fRangeFactor

G4double CMSEmStandardPhysicsTrackingManager::fRangeFactor
private

◆ fSafetyFactor

G4double CMSEmStandardPhysicsTrackingManager::fSafetyFactor
private

◆ fStepLimitType

G4MscStepLimitType CMSEmStandardPhysicsTrackingManager::fStepLimitType
private

◆ gammaProc

G4GammaGeneralProcess* CMSEmStandardPhysicsTrackingManager::gammaProc
private

◆ ioni

G4eIonisation* CMSEmStandardPhysicsTrackingManager::ioni

◆ masterTrackingManager

CMSEmStandardPhysicsTrackingManager * CMSEmStandardPhysicsTrackingManager::masterTrackingManager = nullptr
staticprivate

◆ msc

G4eMultipleScattering* CMSEmStandardPhysicsTrackingManager::msc

◆ nuc [1/2]

G4ElectronNuclearProcess* CMSEmStandardPhysicsTrackingManager::nuc

◆ nuc [2/2]

G4PositronNuclearProcess* CMSEmStandardPhysicsTrackingManager::nuc

Definition at line 56 of file CMSEmStandardPhysicsTrackingManager.h.

◆ positron

struct { ... } CMSEmStandardPhysicsTrackingManager::positron

◆ ss

G4CoulombScattering* CMSEmStandardPhysicsTrackingManager::ss