CMS 3D CMS Logo

CMSSIMPInelasticProcess.cc
Go to the documentation of this file.
1 
5 
6 #include "G4Types.hh"
7 #include "G4SystemOfUnits.hh"
8 #include "G4HadProjectile.hh"
9 #include "G4ElementVector.hh"
10 #include "G4Track.hh"
11 #include "G4Step.hh"
12 #include "G4Element.hh"
13 #include "G4ParticleChange.hh"
14 #include "G4NucleiProperties.hh"
15 #include "G4Nucleus.hh"
16 
17 #include "G4HadronicException.hh"
18 #include "G4HadronicProcessStore.hh"
19 #include "G4HadronicInteraction.hh"
20 
21 #include "G4HadronInelasticDataSet.hh"
22 #include "G4ParticleDefinition.hh"
23 
26  : G4HadronicProcess(processName, fHadronic) {
27  AddDataSet(new CMSSIMPInelasticXS());
29 }
30 
32 
33 G4bool CMSSIMPInelasticProcess::IsApplicable(const G4ParticleDefinition& aP) {
34  return theParticle->GetParticleType() == aP.GetParticleType();
35 }
36 
37 G4VParticleChange* CMSSIMPInelasticProcess::PostStepDoIt(const G4Track& aTrack, const G4Step&) {
38  // if primary is not Alive then do nothing
39  theTotalResult->Clear();
40  theTotalResult->Initialize(aTrack);
41  theTotalResult->ProposeWeight(aTrack.GetWeight());
42  if (aTrack.GetTrackStatus() != fAlive) {
43  return theTotalResult;
44  }
45 
46  // Find cross section at end of step and check if <= 0
47  //
48  G4DynamicParticle* aParticle = const_cast<G4DynamicParticle*>(aTrack.GetDynamicParticle());
49 
50  // change this SIMP particle in a neutron
51  aParticle->SetPDGcode(2112);
52  aParticle->SetDefinition(G4Neutron::Neutron());
53 
54  G4Nucleus* target = GetTargetNucleusPointer();
55  const G4Material* aMaterial = aTrack.GetMaterial();
56  const G4Element* anElement = GetCrossSectionDataStore()->SampleZandA(aParticle, aMaterial, *target);
57 
58  // Next check for illegal track status
59  //
60  if (aTrack.GetTrackStatus() != fAlive && aTrack.GetTrackStatus() != fSuspend) {
61  if (aTrack.GetTrackStatus() == fStopAndKill || aTrack.GetTrackStatus() == fKillTrackAndSecondaries ||
62  aTrack.GetTrackStatus() == fPostponeToNextEvent) {
63  G4ExceptionDescription ed;
64  ed << "CMSSIMPInelasticProcess: track in unusable state - " << aTrack.GetTrackStatus() << G4endl;
65  ed << "CMSSIMPInelasticProcess: returning unchanged track " << G4endl;
66  DumpState(aTrack, "PostStepDoIt", ed);
67  G4Exception("CMSSIMPInelasticProcess::PostStepDoIt", "had004", JustWarning, ed);
68  }
69  // No warning for fStopButAlive which is a legal status here
70  return theTotalResult;
71  }
72 
73  // Initialize the hadronic projectile from the track
74  thePro.Initialise(aTrack);
75  G4HadronicInteraction* anInteraction = GetHadronicInteractionList()[0];
76 
77  G4HadFinalState* result = nullptr;
78  G4int reentryCount = 0;
79 
80  do {
81  try {
82  // Call the interaction
83  result = anInteraction->ApplyYourself(thePro, *target);
84  ++reentryCount;
85  } catch (G4HadronicException& aR) {
86  G4ExceptionDescription ed;
87  aR.Report(ed);
88  ed << "Call for " << anInteraction->GetModelName() << G4endl;
89  ed << "Target element " << anElement->GetName() << " Z= " << target->GetZ_asInt()
90  << " A= " << target->GetA_asInt() << G4endl;
91  DumpState(aTrack, "ApplyYourself", ed);
92  ed << " ApplyYourself failed" << G4endl;
93  G4Exception("CMSSIMPInelasticProcess::PostStepDoIt", "had006", FatalException, ed);
94  }
95 
96  // Check the result for catastrophic energy non-conservation
97  result = CheckResult(thePro, *target, result);
98  if (reentryCount > 100) {
99  G4ExceptionDescription ed;
100  ed << "Call for " << anInteraction->GetModelName() << G4endl;
101  ed << "Target element " << anElement->GetName() << " Z= " << target->GetZ_asInt()
102  << " A= " << target->GetA_asInt() << G4endl;
103  DumpState(aTrack, "ApplyYourself", ed);
104  ed << " ApplyYourself does not completed after 100 attempts" << G4endl;
105  G4Exception("CMSSIMPInelasticProcess::PostStepDoIt", "had006", FatalException, ed);
106  }
107  } while (!result);
108  // Check whether kaon0 or anti_kaon0 are present between the secondaries:
109  // if this is the case, transform them into either kaon0S or kaon0L,
110  // with equal, 50% probability, keeping their dynamical masses (and
111  // the other kinematical properties).
112  // When this happens - very rarely - a "JustWarning" exception is thrown.
113  G4int nSec = result->GetNumberOfSecondaries();
114  if (nSec > 0) {
115  for (G4int i = 0; i < nSec; ++i) {
116  G4DynamicParticle* dynamicParticle = result->GetSecondary(i)->GetParticle();
117  const G4ParticleDefinition* particleDefinition = dynamicParticle->GetParticleDefinition();
118  if (particleDefinition == G4KaonZero::Definition() || particleDefinition == G4AntiKaonZero::Definition()) {
119  G4ParticleDefinition* newPart;
120  if (G4UniformRand() > 0.5) {
121  newPart = G4KaonZeroShort::Definition();
122  } else {
123  newPart = G4KaonZeroLong::Definition();
124  }
125  dynamicParticle->SetDefinition(newPart);
126  G4ExceptionDescription ed;
127  ed << " Hadronic model " << anInteraction->GetModelName() << G4endl;
128  ed << " created " << particleDefinition->GetParticleName() << G4endl;
129  ed << " -> forced to be " << newPart->GetParticleName() << G4endl;
130  G4Exception("G4HadronicProcess::PostStepDoIt", "had007", JustWarning, ed);
131  }
132  }
133  }
134 
135  result->SetTrafoToLab(thePro.GetTrafoToLab());
136 
137  ClearNumberOfInteractionLengthLeft();
138 
139  FillResult(result, aTrack);
140 
141  return theTotalResult;
142 }
mps_fire.i
i
Definition: mps_fire.py:428
CMSSIMPInelasticXS
Definition: CMSSIMPInelasticXS.h:9
CMSSIMPInelasticProcess::PostStepDoIt
G4VParticleChange * PostStepDoIt(const G4Track &aTrack, const G4Step &aStep) override
Definition: CMSSIMPInelasticProcess.cc:37
CMSSIMPInelasticXS.h
CMSSIMPInelasticProcess::theParticle
G4ParticleDefinition * theParticle
Definition: CMSSIMPInelasticProcess.h:24
CMSSIMP.h
CMSSIMP::SIMP
static CMSSIMP * SIMP()
Definition: CMSSIMP.cc:39
CMSSIMPInelasticProcess::CMSSIMPInelasticProcess
CMSSIMPInelasticProcess(const G4String &processName="SIMPInelastic")
Definition: CMSSIMPInelasticProcess.cc:25
CMSSIMPInelasticProcess::~CMSSIMPInelasticProcess
~CMSSIMPInelasticProcess() override
Definition: CMSSIMPInelasticProcess.cc:31
CMSSIMPInelasticProcess::IsApplicable
G4bool IsApplicable(const G4ParticleDefinition &aParticleType) override
Definition: CMSSIMPInelasticProcess.cc:33
CMSSIMPInelasticProcess.h
SimL1EmulatorRepack_CalouGT_cff.processName
processName
Definition: SimL1EmulatorRepack_CalouGT_cff.py:17
filterCSVwithJSON.target
target
Definition: filterCSVwithJSON.py:32
mps_fire.result
result
Definition: mps_fire.py:311