test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ToyModelHadronicProcess.cc
Go to the documentation of this file.
1 //Geant4 include
2 #include "G4ProcessManager.hh"
3 #include "G4ParticleTable.hh"
4 #include "G4Track.hh"
5 #include "G4FermiPhaseSpaceDecay.hh"
6 //Our includes
7 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
10 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh"
12 
13 #include <vector>
14 
15 using namespace CLHEP;
16 
17 ToyModelHadronicProcess::ToyModelHadronicProcess(HadronicProcessHelper * aHelper,
18  const G4String& processName) :
19  G4VDiscreteProcess(processName), m_verboseLevel(0), m_helper(aHelper), m_detachCloud(true)
20 {
21  m_decay = new G4FermiPhaseSpaceDecay();
22 }
23 
24 G4bool ToyModelHadronicProcess::IsApplicable(const G4ParticleDefinition& aP)
25 {
26  return m_helper->applicabilityTester(aP);
27 }
28 
29 G4double ToyModelHadronicProcess::GetMicroscopicCrossSection(const G4DynamicParticle *particle,
30  const G4Element *element,
31  G4double /*temperature*/)
32 {
33  //Get the cross section for this particle/element combination from the ProcessHelper
34  G4double inclusiveCrossSection = m_helper->inclusiveCrossSection(particle,element);
35 
36  if(m_verboseLevel >= 3) {
37  G4cout << "ToyModelHadronicProcess: Return cross section " << inclusiveCrossSection
38  << G4endl;
39  }
40  return inclusiveCrossSection;
41 }
42 
43 G4double ToyModelHadronicProcess::GetMeanFreePath(const G4Track& aTrack, G4double,
44  G4ForceCondition*)
45 {
46  G4Material *aMaterial = aTrack.GetMaterial();
47  const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
48  G4double sigma = 0.0;
49 
50  G4int nElements = aMaterial->GetNumberOfElements();
51 
52  const G4double *theAtomicNumDensityVector =
53  aMaterial->GetAtomicNumDensityVector();
54  G4double aTemp = aMaterial->GetTemperature();
55 
56  for( G4int i=0; i<nElements; ++i )
57  {
58  G4double xSection =
59  GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[i], aTemp);
60  sigma += theAtomicNumDensityVector[i] * xSection;
61  }
62  G4double x = DBL_MAX;
63  if(sigma > 0.0) { x = 1./sigma; }
64  return x;
65 }
66 
67 G4VParticleChange* ToyModelHadronicProcess::PostStepDoIt(const G4Track& track,
68  const G4Step& /* step*/)
69 {
70  const G4TouchableHandle thisTouchable(track.GetTouchableHandle());
71 
72  // A little setting up
73  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
74  m_particleChange.Initialize(track);
75 
76  //This will contain RHadron Def + RHad momentum
77  const G4DynamicParticle* incidentRHadron = track.GetDynamicParticle();
78  const G4ThreeVector aPosition = track.GetPosition();
79 
80  CustomParticle* CustomIncident = dynamic_cast<CustomParticle*>(incidentRHadron->GetDefinition());
81 
82  //This will contain Cloud Def + scaled momentum
83  G4DynamicParticle* cloudParticle = new G4DynamicParticle();
84 
85  //set cloud definition
86  if(CustomIncident==0) {
87  G4cout << "ToyModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!"
88  << G4endl;
89  } else {
90  cloudParticle->SetDefinition(CustomIncident->GetCloud());
91  }
92 
93  //compute scaled momentum
94  double scale = cloudParticle->GetDefinition()->GetPDGMass()
95  /incidentRHadron->GetDefinition()->GetPDGMass();
96  G4LorentzVector cloudMomentum(incidentRHadron->GetMomentum()*scale,
97  cloudParticle->GetDefinition()->GetPDGMass());
98  G4LorentzVector gluinoMomentum(incidentRHadron->GetMomentum()*(1.-scale),
99  CustomIncident->GetSpectator()->GetPDGMass());
100 
101  cloudParticle->Set4Momentum(cloudMomentum);
102 
103  const G4DynamicParticle *incidentParticle;
104  if(! m_detachCloud) { incidentParticle = incidentRHadron; }
105  else { incidentParticle = cloudParticle; }
106 
107  if(m_verboseLevel >= 3)
108  {
109  G4cout << "ToyModelHadronicProcess::PostStepDoIt After scaling " << G4endl;
110  G4cout << " RHadron "<<incidentRHadron->GetDefinition()->GetParticleName()
111  <<" pdgmass= " << incidentRHadron->GetDefinition()->GetPDGMass() / GeV
112  << " P= " << incidentRHadron->Get4Momentum() / GeV<<" mass= "
113  << incidentRHadron->Get4Momentum().m() / GeV <<G4endl;
114  G4cout << " Cloud "<<cloudParticle->GetDefinition()->GetParticleName()
115  <<" pdgmass= " << cloudParticle->GetDefinition()->GetPDGMass() / GeV
116  << " P= " << cloudParticle->Get4Momentum() / GeV<<" mass= "
117  << cloudParticle->Get4Momentum().m() / GeV <<G4endl;
118  G4cout << " Incident pdgmass= "
119  << incidentParticle->GetDefinition()->GetPDGMass() / GeV
120  << " P= " << incidentParticle->Get4Momentum() / GeV<<" mass= "
121  << incidentParticle->Get4Momentum().m() / GeV <<G4endl;
122  }
123 
124  const G4ThreeVector position = track.GetPosition();
125  std::vector<G4ParticleDefinition*> newParticles; // this will contain clouds
126  std::vector<G4ParticleDefinition*> newParticlesRHadron; // this will contain r-hadron
127 
128  //Get the final state particles and target
129  G4ParticleDefinition* targetParticle;
130  HadronicProcessHelper::ReactionProduct reactionProduct =
131  m_helper->finalState(incidentRHadron,track.GetMaterial(),targetParticle);
132 
133  // Fill a list of the new particles to create including the incident if it survives
134  for(HadronicProcessHelper::ReactionProduct::iterator it = reactionProduct.begin();
135  it != reactionProduct.end() ;
136  ++it)
137  {
138  G4ParticleDefinition * productDefinition =theParticleTable->FindParticle(*it);
139  newParticlesRHadron.push_back(productDefinition);
140  CustomParticle* cProd = dynamic_cast<CustomParticle*>(productDefinition);
141 
142  if( cProd!=0 && m_detachCloud)
143  productDefinition=cProd->GetCloud();
144 
145  newParticles.push_back(productDefinition);
146  }
147 
148  int numberOfSecondaries = reactionProduct.size();
149 
150  if(m_verboseLevel >= 2) {
151  G4cout << "ToyModelHadronicProcess::PostStepDoIt N secondaries: "
152  << numberOfSecondaries << G4endl;
153  }
154 
155  // Getting 4-momenta
156  G4LorentzVector sum4Momentum = incidentParticle->Get4Momentum();
157  G4LorentzVector target4Momentum(0,0,0,targetParticle->GetPDGMass());
158  sum4Momentum += target4Momentum;
159  G4ThreeVector cmBoost = sum4Momentum.boostVector();
160  G4double M = sum4Momentum.m();
161 
162  if(m_verboseLevel >= 2) {
163  G4cout << "ToyModelHadronicProcess::PostStepDoIt Kinematics in GeV: " << G4endl
164  << " Boost = " << cmBoost / GeV<< G4endl
165  << " 4P Lab = " << sum4Momentum / GeV << G4endl
166  << " Ecm = " << M / GeV << G4endl;
167  }
168  std::vector<G4double> m;
169 
170  //Fill the masses
171  for(int ip=0; ip<numberOfSecondaries; ++ip)
172  {
173  m.push_back(newParticles[ip]->GetPDGMass());
174  }
175 
176  std::vector<G4LorentzVector*>* fourMomenta = m_decay->Decay(M, m);
177 
178  if(fourMomenta) {
179 
180  //incident particle has to be killed
181  m_particleChange.SetNumberOfSecondaries(numberOfSecondaries);
182  m_particleChange.ProposeTrackStatus(fStopAndKill);
183  if(m_verboseLevel >= 3) {
184  G4cout << "Incident does not survive: stopAndKill + set num secondaries to "
185  << numberOfSecondaries << G4endl;
186  }
187 
188  for (int ip=0; ip <numberOfSecondaries; ++ip) {
189  // transform to Lab
190  (*fourMomenta)[ip]->boost(cmBoost);
191 
192  // added cloud to R-hadron
193  if (newParticlesRHadron[ip]==incidentRHadron->GetDefinition()) {
194  if(m_detachCloud) {
195  if(m_verboseLevel >= 3) {
196  G4cout << "ToyModelHadronicProcess::PostStepDoIt Add gluino momentum "
197  << *((*fourMomenta)[ip])/GeV <<"(m="<< (*fourMomenta)[ip]->m()/ GeV
198  <<") + " <<gluinoMomentum/GeV<<G4endl;
199  }
200  G4LorentzVector p4_new = gluinoMomentum + *((*fourMomenta)[ip]);
201  G4ThreeVector momentum = p4_new.vect();
202  double rhMass = newParticlesRHadron[ip]->GetPDGMass();
203  (*fourMomenta)[ip]->setVectM(momentum,rhMass);
204  if(m_verboseLevel >= 3) {
205  G4cout << " = " << *((*fourMomenta)[ip])/GeV <<"(m="<< (*fourMomenta)[ip]->m() / GeV
206  <<") vs. "<<rhMass/GeV
207  << G4endl;
208  }
209  }
210  }
211  //Create new dynamic particle
212  G4DynamicParticle* productDynParticle = new G4DynamicParticle();
213  //Set the pDef to the dynParte
214  productDynParticle->SetDefinition(newParticlesRHadron[ip]);
215  //Set the 4-vector to dynPart
216  productDynParticle->Set4Momentum(*((*fourMomenta)[ip]));
217 
218  if(m_verboseLevel >= 3) {
219  G4cout << "ToyModelHadronicProcess::PostStepDoIt: " << *((*fourMomenta)[ip])/GeV
220  << "(m="<< (*fourMomenta)[ip]->m()/ GeV<<") "
221  << newParticlesRHadron[ip]->GetParticleName()
222  << G4endl;
223  }
224  //Create a G4Track
225  G4Track* productTrack = new G4Track(productDynParticle,
226  track.GetGlobalTime(),
227  position);
228  productTrack->SetTouchableHandle(thisTouchable);
229  m_particleChange.AddSecondary(productTrack);
230  delete (*fourMomenta)[ip];
231  }
232  delete fourMomenta;
233  }
234  //clear interaction length
235  ClearNumberOfInteractionLengthLeft();
236  return &m_particleChange;
237 }
238 
int i
Definition: DBlmapReader.cc:9
const double GeV
Definition: MathUtil.h:16
G4ParticleDefinition * GetCloud()
T x() const
Cartesian x coordinate.
G4ParticleDefinition * GetSpectator()
static int position[264][3]
Definition: ReadPGInfo.cc:509