CMS 3D CMS Logo

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

#include <CMSSQNeutronAnnih.h>

Inheritance diagram for CMSSQNeutronAnnih:

Public Member Functions

G4HadFinalState * ApplyYourself (const G4HadProjectile &aTrack, G4Nucleus &targetNucleus) override
 
 CMSSQNeutronAnnih (double mass)
 
G4double momDistr (G4double x_in)
 
 ~CMSSQNeutronAnnih () override
 

Private Attributes

G4ParticleDefinition * theAntiL
 
G4ParticleDefinition * theK0S
 
G4ParticleDefinition * theProton
 
G4ParticleDefinition * theSQ
 

Detailed Description

Definition at line 13 of file CMSSQNeutronAnnih.h.

Constructor & Destructor Documentation

◆ CMSSQNeutronAnnih()

CMSSQNeutronAnnih::CMSSQNeutronAnnih ( double  mass)

Definition at line 16 of file CMSSQNeutronAnnih.cc.

References EgHLTOffHistBins_cfi::mass, CMSSQ::SQ(), theAntiL, theK0S, theProton, and theSQ.

16  : G4HadronicInteraction("SexaQuark-neutron annihilation") {
17  SetMinEnergy(0.0 * GeV);
18  SetMaxEnergy(100. * TeV);
19 
20  theSQ = CMSSQ::SQ(mass);
21  theK0S = G4KaonZeroShort::KaonZeroShort();
22  theAntiL = G4AntiLambda::AntiLambda();
23  theProton = G4Proton::
24  Proton(); //proton only used when the particle which the sexaquark hits is a deutereon and the neutron dissapears, so what stays behind is a proton
25 }
G4ParticleDefinition * theProton
G4ParticleDefinition * theAntiL
G4ParticleDefinition * theK0S
static CMSSQ * SQ(double mass)
Definition: CMSSQ.cc:41
G4ParticleDefinition * theSQ

◆ ~CMSSQNeutronAnnih()

CMSSQNeutronAnnih::~CMSSQNeutronAnnih ( )
override

Definition at line 27 of file CMSSQNeutronAnnih.cc.

27 {}

Member Function Documentation

◆ ApplyYourself()

G4HadFinalState * CMSSQNeutronAnnih::ApplyYourself ( const G4HadProjectile &  aTrack,
G4Nucleus &  targetNucleus 
)
override

Definition at line 101 of file CMSSQNeutronAnnih.cc.

References A, momDistr(), AlCaHLTBitMon_ParallelJobs::p, SiStripOfflineCRack_cfg::p2, Pi, funct::pow(), mathSSE::sqrt(), theAntiL, theK0S, theProton, and beamSpotPI::Z.

101  {
102  theParticleChange.Clear();
103  const G4HadProjectile* aParticle = &aTrack;
104  G4double ekin = aParticle->GetKineticEnergy();
105 
106  G4int A = targetNucleus.GetA_asInt();
107  G4int Z = targetNucleus.GetZ_asInt();
108 
109  G4double m_K0S = G4KaonZeroShort::KaonZeroShort()->GetPDGMass();
110  G4double m_L = G4AntiLambda::AntiLambda()->GetPDGMass();
111 
112  //G4double plab = aParticle->GetTotalMomentum();
113 
114  // edm::LogVerbatim("CMSSWNeutronAnnih") << "CMSSQNeutronAnnih: Incident particle p (GeV), total Energy (GeV), particle name, eta ="
115  // << plab/GeV << " "
116  // << aParticle->GetTotalEnergy()/GeV << " "
117  // << aParticle->GetDefinition()->GetParticleName() << " "
118  // << aParticle->Get4Momentum();
119 
120  // Scattered particle referred to axis of incident particle
121  //const G4ParticleDefinition* theParticle = aParticle->GetDefinition();
122 
123  //G4int projPDG = theParticle->GetPDGEncoding();
124  // edm::LogVerbatim("CMSSWNeutronAnnih") << "CMSSQNeutronAnnih: for " << theParticle->GetParticleName()
125  // << " PDGcode= " << projPDG << " on nucleus Z= " << Z
126  // << " A= " << A << " N= " << N;
127 
128  G4LorentzVector lv1 = aParticle->Get4Momentum();
129  edm::LogVerbatim("CMSSWNeutronAnnih") << "The neutron Fermi momentum (mag, x, y, z) "
130  << targetNucleus.GetFermiMomentum().mag() / MeV << " "
131  << targetNucleus.GetFermiMomentum().x() / MeV << " "
132  << targetNucleus.GetFermiMomentum().y() / MeV << " "
133  << targetNucleus.GetFermiMomentum().z() / MeV;
134 
135  //calculate fermi momentum
136 
137  G4double k_neutron = momDistr(G4UniformRand());
138  G4double momentum_neutron = 0.1973 * GeV * k_neutron;
139 
140  G4double theta_neutron = TMath::ACos(2 * G4UniformRand() - 1);
141  G4double phi_neutron = 2. * TMath::Pi() * G4UniformRand();
142 
143  G4double p_neutron_x = momentum_neutron * TMath::Sin(theta_neutron) * TMath::Cos(phi_neutron);
144  G4double p_neutron_y = momentum_neutron * TMath::Sin(theta_neutron) * TMath::Sin(phi_neutron);
145  G4double p_neutron_z = momentum_neutron * TMath::Cos(theta_neutron);
146 
147  //G4LorentzVector lv0(targetNucleus.GetFermiMomentum(), sqrt( pow(G4Neutron::Neutron()->GetPDGMass(),2) + targetNucleus.GetFermiMomentum().mag2() ) );
148  G4LorentzVector lv0(p_neutron_x,
149  p_neutron_y,
150  p_neutron_z,
151  sqrt(pow(G4Neutron::Neutron()->GetPDGMass(), 2) + momentum_neutron * momentum_neutron));
152 
153  //const G4Nucleus* aNucleus = &targetNucleus;
154  G4double BENeutronInNucleus = 0;
155  if (A != 0)
156  BENeutronInNucleus = G4NucleiProperties::GetBindingEnergy(A, Z) / (A);
157 
158  edm::LogVerbatim("CMSSWNeutronAnnih") << "BE of nucleon in the nucleus (GeV): " << BENeutronInNucleus / GeV;
159 
160  G4LorentzVector lvBE(0, 0, 0, BENeutronInNucleus / GeV);
161  G4LorentzVector lv = lv0 + lv1 - lvBE;
162 
163  // kinematiacally impossible ?
164  G4double etot = lv0.e() + lv1.e() - lvBE.e();
165  if (etot < theK0S->GetPDGMass() + theAntiL->GetPDGMass()) {
166  theParticleChange.SetEnergyChange(ekin);
167  theParticleChange.SetMomentumChange(aTrack.Get4Momentum().vect().unit());
168  return &theParticleChange;
169  }
170 
171  float newIonMass = targetNucleus.AtomicMass(A - 1, Z) * 931.5 * MeV;
172  ;
173  G4LorentzVector nlvIon(0, 0, 0, newIonMass);
174 
175  G4double theta_KS0_star = TMath::ACos(2 * G4UniformRand() - 1);
176  G4double phi_KS0_star = 2. * TMath::Pi() * G4UniformRand();
177 
178  G4double p_K0S_star_x = TMath::Sin(theta_KS0_star) * TMath::Cos(phi_KS0_star);
179  G4double p_K0S_star_y = TMath::Sin(theta_KS0_star) * TMath::Sin(phi_KS0_star);
180  G4double p_K0S_star_z = TMath::Cos(theta_KS0_star);
181 
182  G4ThreeVector p(p_K0S_star_x, p_K0S_star_y, p_K0S_star_z);
183  double m0 = lv.m();
184  double m0_2 = m0 * m0;
185  double m1_2 = m_K0S * m_K0S;
186  double m2_2 = m_L * m_L;
187 
188  p *= 0.5 / m0 * sqrt(m0_2 * m0_2 + m1_2 * m1_2 + m2_2 * m2_2 - 2 * m0_2 * m1_2 - 2 * m0_2 * m2_2 - 2 * m1_2 * m2_2);
189  double p2 = p.mag2();
190 
191  G4LorentzVector nlvK0S(p, sqrt(p2 + m1_2));
192  G4LorentzVector nlvAntiL(-p, sqrt(p2 + m2_2));
193 
194  // Boost out of the rest frame.
195  nlvK0S.boost(lv.boostVector());
196  nlvAntiL.boost(lv.boostVector());
197 
198  // now move to implement the interaction
199  theParticleChange.SetStatusChange(stopAndKill);
200  //theParticleChange.SetEnergyChange(ekin); // was 0.0
201 
202  G4DynamicParticle* aSec1 = new G4DynamicParticle(theK0S, nlvK0S);
203  theParticleChange.AddSecondary(aSec1);
204  G4DynamicParticle* aSec2 = new G4DynamicParticle(theAntiL, nlvAntiL);
205  theParticleChange.AddSecondary(aSec2);
206 
207  const G4ParticleDefinition* theRemainingNucleusDef = theProton;
208  if (A != 1)
209  theRemainingNucleusDef = G4IonTable::GetIonTable()->GetIon(Z, A - 1);
210  G4DynamicParticle* aSec3 = new G4DynamicParticle(theRemainingNucleusDef, nlvIon);
211  theParticleChange.AddSecondary(aSec3);
212 
213  // return as is; we don't care about what happens to the nucleus
214  return &theParticleChange;
215 }
const double Pi
Log< level::Info, true > LogVerbatim
G4ParticleDefinition * theProton
G4double momDistr(G4double x_in)
G4ParticleDefinition * theAntiL
G4ParticleDefinition * theK0S
T sqrt(T t)
Definition: SSEVec.h:23
Definition: APVGainStruct.h:7
Power< A, B >::type pow(const A &a, const B &b)
Definition: Power.h:29

◆ momDistr()

G4double CMSSQNeutronAnnih::momDistr ( G4double  x_in)

Definition at line 30 of file CMSSQNeutronAnnih.cc.

References mps_fire::i, and x.

Referenced by ApplyYourself().

30  {
31  const int n_entries = 50;
32 
33  G4double CDF_k[n_entries] = {0, 0.1, 0.2, 0.3, 0.4, 0.5, 0.6, 0.7, 0.8, 0.9, 1, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6,
34  1.7, 1.8, 1.9, 2, 2.1, 2.2, 2.3, 2.4, 2.5, 2.6, 2.7, 2.8, 2.9, 3, 3.1, 3.2, 3.3,
35  3.4, 3.5, 3.6, 3.7, 3.8, 3.9, 4., 4.1, 4.2, 4.3, 4.4, 4.5, 4.6, 4.7, 4.8, 4.9};
36 
37  G4double x[n_entries] = {0,
38  0.0038033182,
39  0.0187291764,
40  0.0510409777,
41  0.1048223609,
42  0.1807862863,
43  0.2756514534,
44  0.3825832103,
45  0.4926859745,
46  0.5970673837,
47  0.6887542272,
48  0.7637748784,
49  0.8212490273,
50  0.8627259608,
51  0.8911605331,
52  0.9099115186,
53  0.9220525854,
54  0.9300190818,
55  0.9355376091,
56  0.9397242185,
57  0.9432387722,
58  0.946438928,
59  0.9495023924,
60  0.9525032995,
61  0.9554669848,
62  0.9583936672,
63  0.9612770117,
64  0.9641067202,
65  0.9668727859,
66  0.9695676121,
67  0.9721815799,
68  0.9747092981,
69  0.9771426396,
70  0.9794740235,
71  0.9816956807,
72  0.9838003583,
73  0.9857816165,
74  0.9876331761,
75  0.9893513365,
76  0.9909333198,
77  0.992378513,
78  0.9936885054,
79  0.9948665964,
80  0.9959179448,
81  0.9968491104,
82  0.9976680755,
83  0.9983832508,
84  0.9990041784,
85  0.9995400073,
86  1};
87 
88  //now interpolate the above points for x_in
89  if (x_in <= 0.0)
90  return 0.;
91  if (x_in >= 1.0)
92  return CDF_k[n_entries - 1];
93  for (int i = 1; i < n_entries; i++) {
94  if (x[i] >= x_in) {
95  return (CDF_k[i] - CDF_k[i - 1]) / (x[i] - x[i - 1]) * (x_in - x[i - 1]) + CDF_k[i - 1];
96  }
97  }
98  return 0.0;
99 }

Member Data Documentation

◆ theAntiL

G4ParticleDefinition* CMSSQNeutronAnnih::theAntiL
private

Definition at line 26 of file CMSSQNeutronAnnih.h.

Referenced by ApplyYourself(), and CMSSQNeutronAnnih().

◆ theK0S

G4ParticleDefinition* CMSSQNeutronAnnih::theK0S
private

Definition at line 25 of file CMSSQNeutronAnnih.h.

Referenced by ApplyYourself(), and CMSSQNeutronAnnih().

◆ theProton

G4ParticleDefinition* CMSSQNeutronAnnih::theProton
private

Definition at line 27 of file CMSSQNeutronAnnih.h.

Referenced by ApplyYourself(), and CMSSQNeutronAnnih().

◆ theSQ

G4ParticleDefinition* CMSSQNeutronAnnih::theSQ
private

Definition at line 24 of file CMSSQNeutronAnnih.h.

Referenced by CMSSQNeutronAnnih().