CMS 3D CMS Logo

CMSEmStandardPhysicsTrackingManager.cc
Go to the documentation of this file.
1 #include "G4Version.hh"
2 #if G4VERSION_NUMBER >= 1100
3 
6 
7 #include "G4CoulombScattering.hh"
8 #include "G4UrbanMscModel.hh"
9 #include "G4WentzelVIModel.hh"
10 #include "G4eBremsstrahlung.hh"
11 #include "G4eCoulombScatteringModel.hh"
12 #include "G4eIonisation.hh"
13 #include "G4eMultipleScattering.hh"
14 #include "G4eplusAnnihilation.hh"
15 
16 #include "G4ElectroVDNuclearModel.hh"
17 #include "G4ElectronNuclearProcess.hh"
18 #include "G4PositronNuclearProcess.hh"
19 
20 #include "G4ComptonScattering.hh"
21 #include "G4GammaConversion.hh"
22 #include "G4PhotoElectricEffect.hh"
23 
24 #include "G4CascadeInterface.hh"
25 #include "G4CrossSectionDataSetRegistry.hh"
26 #include "G4ExcitedStringDecay.hh"
27 #include "G4GammaNuclearXS.hh"
28 #include "G4GammaParticipants.hh"
29 #include "G4GeneratorPrecompoundInterface.hh"
30 #include "G4HadronInelasticProcess.hh"
31 #include "G4HadronicParameters.hh"
32 #include "G4LowEGammaNuclearModel.hh"
33 #include "G4PhotoNuclearCrossSection.hh"
34 #include "G4QGSMFragmentation.hh"
35 #include "G4QGSModel.hh"
36 #include "G4TheoFSGenerator.hh"
37 
38 #include "G4GammaGeneralProcess.hh"
39 #include "G4LossTableManager.hh"
40 
41 #include "G4EmParameters.hh"
42 #include "G4SystemOfUnits.hh"
43 
44 #include "G4Electron.hh"
45 #include "G4Gamma.hh"
46 #include "G4Positron.hh"
47 
48 #include "G4RegionStore.hh"
49 #include "G4Region.hh"
50 #include <string>
51 
52 CMSEmStandardPhysicsTrackingManager *CMSEmStandardPhysicsTrackingManager::masterTrackingManager = nullptr;
53 
54 CMSEmStandardPhysicsTrackingManager::CMSEmStandardPhysicsTrackingManager(const edm::ParameterSet &p) {
55  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
56  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
57  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
58  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
59  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
60  fStepLimitType = fUseSafety;
61  if (msc == "UseSafetyPlus") {
62  fStepLimitType = fUseSafetyPlus;
63  } else if (msc == "Minimal") {
64  fStepLimitType = fMinimal;
65  }
66 
67  G4EmParameters *param = G4EmParameters::Instance();
68  G4double highEnergyLimit = param->MscEnergyLimit();
69 
70  const G4Region *aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
71  const G4Region *bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
72 
73  // e-
74  {
75  G4eMultipleScattering *msc = new G4eMultipleScattering;
76  G4UrbanMscModel *msc1 = new G4UrbanMscModel;
77  G4WentzelVIModel *msc2 = new G4WentzelVIModel;
78  msc1->SetHighEnergyLimit(highEnergyLimit);
79  msc2->SetLowEnergyLimit(highEnergyLimit);
80  msc->SetEmModel(msc1);
81  msc->SetEmModel(msc2);
82 
83  // e-/e+ msc for HCAL and HGCAL using the Urban model
84  if (nullptr != aRegion || nullptr != bRegion) {
85  G4UrbanMscModel *msc3 = new G4UrbanMscModel();
86  msc3->SetHighEnergyLimit(highEnergyLimit);
87  msc3->SetRangeFactor(fRangeFactor);
88  msc3->SetGeomFactor(fGeomFactor);
89  msc3->SetSafetyFactor(fSafetyFactor);
90  msc3->SetLambdaLimit(fLambdaLimit);
91  msc3->SetStepLimitType(fStepLimitType);
92  msc3->SetLocked(true);
93  if (nullptr != aRegion) {
94  msc->AddEmModel(-1, msc3, aRegion);
95  }
96  if (nullptr != bRegion) {
97  msc->AddEmModel(-1, msc3, bRegion);
98  }
99  }
100 
101  electron.msc = msc;
102 
103  electron.ioni = new G4eIonisation;
104  electron.brems = new G4eBremsstrahlung;
105 
106  G4CoulombScattering *ss = new G4CoulombScattering;
107  G4eCoulombScatteringModel *ssm = new G4eCoulombScatteringModel;
108  ssm->SetLowEnergyLimit(highEnergyLimit);
109  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
110  ss->SetEmModel(ssm);
111  ss->SetMinKinEnergy(highEnergyLimit);
112  electron.ss = ss;
113  }
114 
115  // e+
116  {
117  G4eMultipleScattering *msc = new G4eMultipleScattering;
118  G4UrbanMscModel *msc1 = new G4UrbanMscModel;
119  G4WentzelVIModel *msc2 = new G4WentzelVIModel;
120  msc1->SetHighEnergyLimit(highEnergyLimit);
121  msc2->SetLowEnergyLimit(highEnergyLimit);
122  msc->SetEmModel(msc1);
123  msc->SetEmModel(msc2);
124 
125  // e-/e+ msc for HCAL and HGCAL using the Urban model
126  if (nullptr != aRegion || nullptr != bRegion) {
127  G4UrbanMscModel *msc3 = new G4UrbanMscModel();
128  msc3->SetHighEnergyLimit(highEnergyLimit);
129  msc3->SetRangeFactor(fRangeFactor);
130  msc3->SetGeomFactor(fGeomFactor);
131  msc3->SetSafetyFactor(fSafetyFactor);
132  msc3->SetLambdaLimit(fLambdaLimit);
133  msc3->SetStepLimitType(fStepLimitType);
134  msc3->SetLocked(true);
135  if (nullptr != aRegion) {
136  msc->AddEmModel(-1, msc3, aRegion);
137  }
138  if (nullptr != bRegion) {
139  msc->AddEmModel(-1, msc3, bRegion);
140  }
141  }
142 
143  positron.msc = msc;
144 
145  positron.ioni = new G4eIonisation;
146  positron.brems = new G4eBremsstrahlung;
147  positron.annihilation = new G4eplusAnnihilation;
148 
149  G4CoulombScattering *ss = new G4CoulombScattering;
150  G4eCoulombScatteringModel *ssm = new G4eCoulombScatteringModel;
151  ssm->SetLowEnergyLimit(highEnergyLimit);
152  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
153  ss->SetEmModel(ssm);
154  ss->SetMinKinEnergy(highEnergyLimit);
155  positron.ss = ss;
156  }
157 
158  // gamma
159  {
160  gammaProc = new G4GammaGeneralProcess;
161 
162  gammaProc->AddEmProcess(new G4PhotoElectricEffect);
163  gammaProc->AddEmProcess(new G4ComptonScattering);
164  gammaProc->AddEmProcess(new G4GammaConversion);
165 
166  G4HadronInelasticProcess *nuc = new G4HadronInelasticProcess("photonNuclear", G4Gamma::Definition());
167  auto xsreg = G4CrossSectionDataSetRegistry::Instance();
168  G4VCrossSectionDataSet *xs = nullptr;
169  bool useGammaNuclearXS = true;
170  if (useGammaNuclearXS) {
171  xs = xsreg->GetCrossSectionDataSet("GammaNuclearXS");
172  if (nullptr == xs)
173  xs = new G4GammaNuclearXS;
174  } else {
175  xs = xsreg->GetCrossSectionDataSet("PhotoNuclearXS");
176  if (nullptr == xs)
177  xs = new G4PhotoNuclearCrossSection;
178  }
179  nuc->AddDataSet(xs);
180 
181  G4QGSModel<G4GammaParticipants> *theStringModel = new G4QGSModel<G4GammaParticipants>;
182  G4QGSMFragmentation *theFrag = new G4QGSMFragmentation;
183  G4ExcitedStringDecay *theStringDecay = new G4ExcitedStringDecay(theFrag);
184  theStringModel->SetFragmentationModel(theStringDecay);
185 
186  G4GeneratorPrecompoundInterface *theCascade = new G4GeneratorPrecompoundInterface;
187 
188  G4TheoFSGenerator *theModel = new G4TheoFSGenerator;
189  theModel->SetTransport(theCascade);
190  theModel->SetHighEnergyGenerator(theStringModel);
191 
192  G4HadronicParameters *param = G4HadronicParameters::Instance();
193 
194  G4CascadeInterface *cascade = new G4CascadeInterface;
195 
196  // added low-energy model LEND disabled
197  G4double gnLowEnergyLimit = 200 * CLHEP::MeV;
198  if (gnLowEnergyLimit > 0.0) {
199  G4LowEGammaNuclearModel *lemod = new G4LowEGammaNuclearModel;
200  lemod->SetMaxEnergy(gnLowEnergyLimit);
201  nuc->RegisterMe(lemod);
202  cascade->SetMinEnergy(gnLowEnergyLimit - CLHEP::MeV);
203  }
204  cascade->SetMaxEnergy(param->GetMaxEnergyTransitionFTF_Cascade());
205  nuc->RegisterMe(cascade);
206  theModel->SetMinEnergy(param->GetMinEnergyTransitionFTF_Cascade());
207  theModel->SetMaxEnergy(param->GetMaxEnergy());
208  nuc->RegisterMe(theModel);
209 
210  gammaProc->AddHadProcess(nuc);
211  G4LossTableManager::Instance()->SetGammaGeneralProcess(gammaProc);
212  }
213 
214  // Create lepto-nuclear processes last, they access cross section data from
215  // the gamma nuclear process!
216  G4ElectroVDNuclearModel *eModel = new G4ElectroVDNuclearModel;
217 
218  {
219  G4ElectronNuclearProcess *nuc = new G4ElectronNuclearProcess;
220  nuc->RegisterMe(eModel);
221  electron.nuc = nuc;
222  }
223  {
224  G4PositronNuclearProcess *nuc = new G4PositronNuclearProcess;
225  nuc->RegisterMe(eModel);
226  positron.nuc = nuc;
227  }
228 
229  if (masterTrackingManager == nullptr) {
230  masterTrackingManager = this;
231  } else {
232  electron.msc->SetMasterProcess(masterTrackingManager->electron.msc);
233  electron.ss->SetMasterProcess(masterTrackingManager->electron.ss);
234  electron.ioni->SetMasterProcess(masterTrackingManager->electron.ioni);
235  electron.brems->SetMasterProcess(masterTrackingManager->electron.brems);
236  electron.nuc->SetMasterProcess(masterTrackingManager->electron.nuc);
237 
238  positron.msc->SetMasterProcess(masterTrackingManager->positron.msc);
239  positron.ss->SetMasterProcess(masterTrackingManager->positron.ss);
240  positron.ioni->SetMasterProcess(masterTrackingManager->positron.ioni);
241  positron.brems->SetMasterProcess(masterTrackingManager->positron.brems);
242  positron.annihilation->SetMasterProcess(masterTrackingManager->positron.annihilation);
243  positron.nuc->SetMasterProcess(masterTrackingManager->positron.nuc);
244 
245  gammaProc->SetMasterProcess(masterTrackingManager->gammaProc);
246  }
247 }
248 
249 CMSEmStandardPhysicsTrackingManager::~CMSEmStandardPhysicsTrackingManager() {
250  if (masterTrackingManager == this) {
251  masterTrackingManager = nullptr;
252  }
253 }
254 
255 void CMSEmStandardPhysicsTrackingManager::BuildPhysicsTable(const G4ParticleDefinition &part) {
256  if (&part == G4Electron::Definition()) {
257  electron.msc->BuildPhysicsTable(part);
258  electron.ioni->BuildPhysicsTable(part);
259  electron.brems->BuildPhysicsTable(part);
260  electron.ss->BuildPhysicsTable(part);
261  electron.nuc->BuildPhysicsTable(part);
262  } else if (&part == G4Positron::Definition()) {
263  positron.msc->BuildPhysicsTable(part);
264  positron.ioni->BuildPhysicsTable(part);
265  positron.brems->BuildPhysicsTable(part);
266  positron.annihilation->BuildPhysicsTable(part);
267  positron.ss->BuildPhysicsTable(part);
268  positron.nuc->BuildPhysicsTable(part);
269  } else if (&part == G4Gamma::Definition()) {
270  gammaProc->BuildPhysicsTable(part);
271  }
272 }
273 
274 void CMSEmStandardPhysicsTrackingManager::PreparePhysicsTable(const G4ParticleDefinition &part) {
275  if (&part == G4Electron::Definition()) {
276  electron.msc->PreparePhysicsTable(part);
277  electron.ioni->PreparePhysicsTable(part);
278  electron.brems->PreparePhysicsTable(part);
279  electron.ss->PreparePhysicsTable(part);
280  electron.nuc->PreparePhysicsTable(part);
281  } else if (&part == G4Positron::Definition()) {
282  positron.msc->PreparePhysicsTable(part);
283  positron.ioni->PreparePhysicsTable(part);
284  positron.brems->PreparePhysicsTable(part);
285  positron.annihilation->PreparePhysicsTable(part);
286  positron.ss->PreparePhysicsTable(part);
287  positron.nuc->PreparePhysicsTable(part);
288  } else if (&part == G4Gamma::Definition()) {
289  gammaProc->PreparePhysicsTable(part);
290  }
291 }
292 
293 void CMSEmStandardPhysicsTrackingManager::TrackElectron(G4Track *aTrack) {
294  class ElectronPhysics final : public TrackingManagerHelper::Physics {
295  public:
296  ElectronPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
297 
298  void StartTracking(G4Track *aTrack) override {
299  auto &electron = fMgr.electron;
300 
301  electron.msc->StartTracking(aTrack);
302  electron.ioni->StartTracking(aTrack);
303  electron.brems->StartTracking(aTrack);
304  electron.ss->StartTracking(aTrack);
305  electron.nuc->StartTracking(aTrack);
306 
307  fPreviousStepLength = 0;
308  }
309  void EndTracking() override {
310  auto &electron = fMgr.electron;
311 
312  electron.msc->EndTracking();
313  electron.ioni->EndTracking();
314  electron.brems->EndTracking();
315  electron.ss->EndTracking();
316  electron.nuc->EndTracking();
317  }
318 
319  G4double GetPhysicalInteractionLength(const G4Track &track) override {
320  auto &electron = fMgr.electron;
321  G4double physIntLength, proposedSafety = DBL_MAX;
322  G4ForceCondition condition;
323  G4GPILSelection selection;
324 
325  fProposedStep = DBL_MAX;
326  fSelected = -1;
327 
328  physIntLength = electron.nuc->PostStepGPIL(track, fPreviousStepLength, &condition);
329  if (physIntLength < fProposedStep) {
330  fProposedStep = physIntLength;
331  fSelected = 3;
332  }
333 
334  physIntLength = electron.ss->PostStepGPIL(track, fPreviousStepLength, &condition);
335  if (physIntLength < fProposedStep) {
336  fProposedStep = physIntLength;
337  fSelected = 0;
338  }
339 
340  physIntLength = electron.brems->PostStepGPIL(track, fPreviousStepLength, &condition);
341  if (physIntLength < fProposedStep) {
342  fProposedStep = physIntLength;
343  fSelected = 1;
344  }
345 
346  physIntLength = electron.ioni->PostStepGPIL(track, fPreviousStepLength, &condition);
347  if (physIntLength < fProposedStep) {
348  fProposedStep = physIntLength;
349  fSelected = 2;
350  }
351 
352  physIntLength =
353  electron.ioni->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
354  if (physIntLength < fProposedStep) {
355  fProposedStep = physIntLength;
356  fSelected = -1;
357  }
358 
359  physIntLength =
360  electron.msc->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
361  if (physIntLength < fProposedStep) {
362  fProposedStep = physIntLength;
363  // Check if MSC actually wants to win, in most cases it only limits the
364  // step size.
365  if (selection == CandidateForSelection) {
366  fSelected = -1;
367  }
368  }
369 
370  return fProposedStep;
371  }
372 
373  void AlongStepDoIt(G4Track &track, G4Step &step, G4TrackVector &) override {
374  if (step.GetStepLength() == fProposedStep) {
375  step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
376  } else {
377  // Remember that the step was limited by geometry.
378  fSelected = -1;
379  }
380  auto &electron = fMgr.electron;
381  G4VParticleChange *particleChange;
382 
383  particleChange = electron.msc->AlongStepDoIt(track, step);
384  particleChange->UpdateStepForAlongStep(&step);
385  track.SetTrackStatus(particleChange->GetTrackStatus());
386  particleChange->Clear();
387 
388  particleChange = electron.ioni->AlongStepDoIt(track, step);
389  particleChange->UpdateStepForAlongStep(&step);
390  track.SetTrackStatus(particleChange->GetTrackStatus());
391  particleChange->Clear();
392 
393  fPreviousStepLength = step.GetStepLength();
394  }
395 
396  void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
397  if (fSelected < 0) {
398  return;
399  }
400  step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
401 
402  auto &electron = fMgr.electron;
403  G4VProcess *process = nullptr;
404  G4VParticleChange *particleChange = nullptr;
405 
406  switch (fSelected) {
407  case 3:
408  process = electron.nuc;
409  particleChange = electron.nuc->PostStepDoIt(track, step);
410  break;
411  case 0:
412  process = electron.ss;
413  particleChange = electron.ss->PostStepDoIt(track, step);
414  break;
415  case 1:
416  process = electron.brems;
417  particleChange = electron.brems->PostStepDoIt(track, step);
418  break;
419  case 2:
420  process = electron.ioni;
421  particleChange = electron.ioni->PostStepDoIt(track, step);
422  break;
423  }
424 
425  particleChange->UpdateStepForPostStep(&step);
426  step.UpdateTrack();
427 
428  int numSecondaries = particleChange->GetNumberOfSecondaries();
429  for (int i = 0; i < numSecondaries; i++) {
430  G4Track *secondary = particleChange->GetSecondary(i);
431  secondary->SetParentID(track.GetTrackID());
432  secondary->SetCreatorProcess(process);
433  secondaries.push_back(secondary);
434  }
435 
436  track.SetTrackStatus(particleChange->GetTrackStatus());
437  particleChange->Clear();
438  }
439 
440  private:
441  CMSEmStandardPhysicsTrackingManager &fMgr;
442  G4double fPreviousStepLength;
443  G4double fProposedStep;
444  G4int fSelected;
445  };
446 
447  ElectronPhysics physics(*this);
449 }
450 
451 void CMSEmStandardPhysicsTrackingManager::TrackPositron(G4Track *aTrack) {
452  class PositronPhysics final : public TrackingManagerHelper::Physics {
453  public:
454  PositronPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
455 
456  void StartTracking(G4Track *aTrack) override {
457  auto &positron = fMgr.positron;
458 
459  positron.msc->StartTracking(aTrack);
460  positron.ioni->StartTracking(aTrack);
461  positron.brems->StartTracking(aTrack);
462  positron.annihilation->StartTracking(aTrack);
463  positron.ss->StartTracking(aTrack);
464  positron.nuc->StartTracking(aTrack);
465 
466  fPreviousStepLength = 0;
467  }
468  void EndTracking() override {
469  auto &positron = fMgr.positron;
470 
471  positron.msc->EndTracking();
472  positron.ioni->EndTracking();
473  positron.brems->EndTracking();
474  positron.annihilation->EndTracking();
475  positron.ss->EndTracking();
476  positron.nuc->EndTracking();
477  }
478 
479  G4double GetPhysicalInteractionLength(const G4Track &track) override {
480  auto &positron = fMgr.positron;
481  G4double physIntLength, proposedSafety = DBL_MAX;
482  G4ForceCondition condition;
483  G4GPILSelection selection;
484 
485  fProposedStep = DBL_MAX;
486  fSelected = -1;
487 
488  physIntLength = positron.nuc->PostStepGPIL(track, fPreviousStepLength, &condition);
489  if (physIntLength < fProposedStep) {
490  fProposedStep = physIntLength;
491  fSelected = 4;
492  }
493 
494  physIntLength = positron.ss->PostStepGPIL(track, fPreviousStepLength, &condition);
495  if (physIntLength < fProposedStep) {
496  fProposedStep = physIntLength;
497  fSelected = 0;
498  }
499 
500  physIntLength = positron.annihilation->PostStepGPIL(track, fPreviousStepLength, &condition);
501  if (physIntLength < fProposedStep) {
502  fProposedStep = physIntLength;
503  fSelected = 1;
504  }
505 
506  physIntLength = positron.brems->PostStepGPIL(track, fPreviousStepLength, &condition);
507  if (physIntLength < fProposedStep) {
508  fProposedStep = physIntLength;
509  fSelected = 2;
510  }
511 
512  physIntLength = positron.ioni->PostStepGPIL(track, fPreviousStepLength, &condition);
513  if (physIntLength < fProposedStep) {
514  fProposedStep = physIntLength;
515  fSelected = 3;
516  }
517 
518  physIntLength =
519  positron.ioni->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
520  if (physIntLength < fProposedStep) {
521  fProposedStep = physIntLength;
522  fSelected = -1;
523  }
524 
525  physIntLength =
526  positron.msc->AlongStepGPIL(track, fPreviousStepLength, fProposedStep, proposedSafety, &selection);
527  if (physIntLength < fProposedStep) {
528  fProposedStep = physIntLength;
529  // Check if MSC actually wants to win, in most cases it only limits the
530  // step size.
531  if (selection == CandidateForSelection) {
532  fSelected = -1;
533  }
534  }
535 
536  return fProposedStep;
537  }
538 
539  void AlongStepDoIt(G4Track &track, G4Step &step, G4TrackVector &) override {
540  if (step.GetStepLength() == fProposedStep) {
541  step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
542  } else {
543  // Remember that the step was limited by geometry.
544  fSelected = -1;
545  }
546  auto &positron = fMgr.positron;
547  G4VParticleChange *particleChange;
548 
549  particleChange = positron.msc->AlongStepDoIt(track, step);
550  particleChange->UpdateStepForAlongStep(&step);
551  track.SetTrackStatus(particleChange->GetTrackStatus());
552  particleChange->Clear();
553 
554  particleChange = positron.ioni->AlongStepDoIt(track, step);
555  particleChange->UpdateStepForAlongStep(&step);
556  track.SetTrackStatus(particleChange->GetTrackStatus());
557  particleChange->Clear();
558 
559  fPreviousStepLength = step.GetStepLength();
560  }
561 
562  void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
563  if (fSelected < 0) {
564  return;
565  }
566  step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
567 
568  auto &positron = fMgr.positron;
569  G4VProcess *process;
570  G4VParticleChange *particleChange = nullptr;
571 
572  switch (fSelected) {
573  case 4:
574  process = positron.nuc;
575  particleChange = positron.nuc->PostStepDoIt(track, step);
576  break;
577  case 0:
578  process = positron.ss;
579  particleChange = positron.ss->PostStepDoIt(track, step);
580  break;
581  case 1:
582  process = positron.annihilation;
583  particleChange = positron.annihilation->PostStepDoIt(track, step);
584  break;
585  case 2:
586  process = positron.brems;
587  particleChange = positron.brems->PostStepDoIt(track, step);
588  break;
589  case 3:
590  process = positron.ioni;
591  particleChange = positron.ioni->PostStepDoIt(track, step);
592  break;
593  }
594 
595  particleChange->UpdateStepForPostStep(&step);
596  step.UpdateTrack();
597 
598  int numSecondaries = particleChange->GetNumberOfSecondaries();
599  for (int i = 0; i < numSecondaries; i++) {
600  G4Track *secondary = particleChange->GetSecondary(i);
601  secondary->SetParentID(track.GetTrackID());
602  secondary->SetCreatorProcess(process);
603  secondaries.push_back(secondary);
604  }
605 
606  track.SetTrackStatus(particleChange->GetTrackStatus());
607  particleChange->Clear();
608  }
609 
610  G4bool HasAtRestProcesses() override { return true; }
611 
612  void AtRestDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
613  auto &positron = fMgr.positron;
614  // Annihilate the positron at rest.
615  G4VParticleChange *particleChange = positron.annihilation->AtRestDoIt(track, step);
616  particleChange->UpdateStepForAtRest(&step);
617  step.UpdateTrack();
618 
619  int numSecondaries = particleChange->GetNumberOfSecondaries();
620  for (int i = 0; i < numSecondaries; i++) {
621  G4Track *secondary = particleChange->GetSecondary(i);
622  secondary->SetParentID(track.GetTrackID());
623  secondary->SetCreatorProcess(positron.annihilation);
624  secondaries.push_back(secondary);
625  }
626 
627  track.SetTrackStatus(particleChange->GetTrackStatus());
628  particleChange->Clear();
629  }
630 
631  private:
632  CMSEmStandardPhysicsTrackingManager &fMgr;
633  G4double fPreviousStepLength;
634  G4double fProposedStep;
635  G4int fSelected;
636  };
637 
638  PositronPhysics physics(*this);
640 }
641 
642 void CMSEmStandardPhysicsTrackingManager::TrackGamma(G4Track *aTrack) {
643  class GammaPhysics final : public TrackingManagerHelper::Physics {
644  public:
645  GammaPhysics(CMSEmStandardPhysicsTrackingManager &mgr) : fMgr(mgr) {}
646 
647  void StartTracking(G4Track *aTrack) override {
648  fMgr.gammaProc->StartTracking(aTrack);
649 
650  fPreviousStepLength = 0;
651  }
652  void EndTracking() override { fMgr.gammaProc->EndTracking(); }
653 
654  G4double GetPhysicalInteractionLength(const G4Track &track) override {
655  G4double physIntLength;
656  G4ForceCondition condition;
657 
658  fProposedStep = DBL_MAX;
659  fSelected = -1;
660 
661  physIntLength = fMgr.gammaProc->PostStepGPIL(track, fPreviousStepLength, &condition);
662  if (physIntLength < fProposedStep) {
663  fProposedStep = physIntLength;
664  fSelected = 0;
665  }
666 
667  return fProposedStep;
668  }
669 
670  void AlongStepDoIt(G4Track &, G4Step &step, G4TrackVector &) override {
671  if (step.GetStepLength() == fProposedStep) {
672  step.GetPostStepPoint()->SetStepStatus(fAlongStepDoItProc);
673  } else {
674  // Remember that the step was limited by geometry.
675  fSelected = -1;
676  }
677  fPreviousStepLength = step.GetStepLength();
678  }
679 
680  void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries) override {
681  if (fSelected < 0) {
682  return;
683  }
684  step.GetPostStepPoint()->SetStepStatus(fPostStepDoItProc);
685 
686  G4VProcess *process = fMgr.gammaProc;
687  G4VParticleChange *particleChange = fMgr.gammaProc->PostStepDoIt(track, step);
688 
689  particleChange->UpdateStepForPostStep(&step);
690  step.UpdateTrack();
691 
692  int numSecondaries = particleChange->GetNumberOfSecondaries();
693  for (int i = 0; i < numSecondaries; i++) {
694  G4Track *secondary = particleChange->GetSecondary(i);
695  secondary->SetParentID(track.GetTrackID());
696  secondary->SetCreatorProcess(process);
697  secondaries.push_back(secondary);
698  }
699 
700  track.SetTrackStatus(particleChange->GetTrackStatus());
701  particleChange->Clear();
702  }
703 
704  private:
705  CMSEmStandardPhysicsTrackingManager &fMgr;
706  G4double fPreviousStepLength;
707  G4double fProposedStep;
708  G4int fSelected;
709  };
710 
711  GammaPhysics physics(*this);
713 }
714 
715 void CMSEmStandardPhysicsTrackingManager::HandOverOneTrack(G4Track *aTrack) {
716  const G4ParticleDefinition *part = aTrack->GetParticleDefinition();
717 
718  if (part == G4Electron::Definition()) {
719  TrackElectron(aTrack);
720  } else if (part == G4Positron::Definition()) {
721  TrackPositron(aTrack);
722  } else if (part == G4Gamma::Definition()) {
723  TrackGamma(aTrack);
724  }
725 
726  aTrack->SetTrackStatus(fStopAndKill);
727  delete aTrack;
728 }
729 
730 #endif
static void TrackNeutralParticle(G4Track *aTrack, PhysicsImpl &physics)
virtual void StartTracking(G4Track *)
selection
main part
Definition: corrVsCorr.py:100
static void TrackChargedParticle(G4Track *aTrack, PhysicsImpl &physics)
virtual void PostStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries)=0
virtual void AtRestDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries)
part
Definition: HCALResponse.h:20
virtual void AlongStepDoIt(G4Track &track, G4Step &step, G4TrackVector &secondaries)=0
step
Definition: StallMonitor.cc:98
virtual G4double GetPhysicalInteractionLength(const G4Track &track)=0