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 "G4InelasticInteraction.hh"
6 #include "Randomize.hh"
7 //Our includes
8 #include "SimG4Core/CustomPhysics/interface/ToyModelHadronicProcess.hh"
11 #include "SimG4Core/CustomPhysics/interface/HadronicProcessHelper.hh"
13 
14 ToyModelHadronicProcess::ToyModelHadronicProcess(HadronicProcessHelper * aHelper, const G4String& processName) :
15  G4VDiscreteProcess(processName), m_verboseLevel(0), m_helper(aHelper), m_detachCloud(true)
16 {
17  // Instantiating helper class
18  //m_helper = HadronicProcessHelper::instance();
19  //m_verboseLevel=0;
20  //m_detachCloud=true;
21 
22 }
23 
24 
25 G4bool ToyModelHadronicProcess::IsApplicable(const G4ParticleDefinition& aP)
26 {
27  return m_helper->applicabilityTester(aP);
28 }
29 
30 G4double ToyModelHadronicProcess::GetMicroscopicCrossSection(const G4DynamicParticle *particle,
31  const G4Element *element,
32  G4double /*temperature*/)
33 {
34  //Get the cross section for this particle/element combination from the ProcessHelper
35  G4double inclusiveCrossSection = m_helper->inclusiveCrossSection(particle,element);
36 
37  // Need to provide Set-methods for these in time
38  G4double highestEnergyLimit = 10 * TeV ;
39  G4double lowestEnergyLimit = 1 * eV;
40  G4double particleEnergy = particle->GetKineticEnergy();
41 
42  if (particleEnergy > highestEnergyLimit || particleEnergy < lowestEnergyLimit){
43  if(m_verboseLevel >= 1) std::cout << "ToyModelHadronicProcess: Energy out of bounds [" <<
44  lowestEnergyLimit / MeV << "MeV , " <<
45  highestEnergyLimit / MeV << "MeV ] while it is " << particleEnergy/MeV <<
46  std::endl;
47  return 0;
48  } else {
49  if(m_verboseLevel >= 3) std::cout << "ToyModelHadronicProcess: Return cross section " << inclusiveCrossSection << std::endl;
50  return inclusiveCrossSection;
51  }
52 
53 }
54 
55 
56 G4double ToyModelHadronicProcess::GetMeanFreePath(const G4Track& aTrack, G4double, G4ForceCondition*)
57 {
58  G4Material *aMaterial = aTrack.GetMaterial();
59  const G4DynamicParticle *aParticle = aTrack.GetDynamicParticle();
60  G4double sigma = 0.0;
61 
62  G4int nElements = aMaterial->GetNumberOfElements();
63 
64  const G4double *theAtomicNumDensityVector =
65  aMaterial->GetAtomicNumDensityVector();
66  G4double aTemp = aMaterial->GetTemperature();
67 
68  for( G4int i=0; i<nElements; ++i )
69  {
70  G4double xSection =
71  GetMicroscopicCrossSection( aParticle, (*aMaterial->GetElementVector())[i], aTemp);
72  sigma += theAtomicNumDensityVector[i] * xSection;
73  }
74 
75  return 1.0/sigma;
76 
77 }
78 
79 
80 G4VParticleChange* ToyModelHadronicProcess::PostStepDoIt(const G4Track& track,
81  const G4Step& /* step*/)
82 {
83 
84  const G4TouchableHandle thisTouchable(track.GetTouchableHandle());
85 
86  // A little setting up
87  G4ParticleTable* theParticleTable = G4ParticleTable::GetParticleTable();
88  m_particleChange.Initialize(track);
89  // G4DynamicParticle* incidentRHadron = const_cast<G4DynamicParticle*>(track.GetDynamicParticle()); //This will contain RHadron Def + RHad momentum
90  const G4DynamicParticle* incidentRHadron = track.GetDynamicParticle(); //This will contain RHadron Def + RHad momentum
91 // double E_0 = incidentRHadron->GetKineticEnergy();
92 // const G4int theIncidentPDG = incidentRHadron->GetDefinition()->GetPDGEncoding();
93  const G4ThreeVector aPosition = track.GetPosition();
94 
95 // double gamma = incidentRHadron->GetTotalEnergy()/incidentRHadron->GetDefinition()->GetPDGMass();
96 
97  CustomParticle* CustomIncident = dynamic_cast<CustomParticle*>(incidentRHadron->GetDefinition());
98  G4DynamicParticle* cloudParticle = new G4DynamicParticle(); //This will contain Cloud Def + scaled momentum
99 
100  // G4int rHadronPdg = incidentRHadron->GetDefinition()->GetPDGEncoding();
101 
102  //set cloud definition
103  if(CustomIncident==0)
104  {
105  std::cout << "ToyModelHadronicProcess::PostStepDoIt Definition of particle cloud not available!!" << std::endl;
106  } else {
107  cloudParticle->SetDefinition(CustomIncident->GetCloud());
108  }
109 
110  //compute scaled momentum
111  double scale=cloudParticle->GetDefinition()->GetPDGMass()/incidentRHadron->GetDefinition()->GetPDGMass();
112  G4LorentzVector cloudMomentum;
113  cloudMomentum.setVectM(incidentRHadron->GetMomentum()*scale,cloudParticle->GetDefinition()->GetPDGMass());
114  G4LorentzVector gluinoMomentum;
115  gluinoMomentum.setVectM(incidentRHadron->GetMomentum()*(1.-scale),CustomIncident->GetSpectator()->GetPDGMass());
116 
117  cloudParticle->Set4Momentum(cloudMomentum);
118 
119  const G4DynamicParticle *incidentParticle;
120  if(! m_detachCloud) incidentParticle = incidentRHadron;
121  else incidentParticle= cloudParticle;
122 
123  if(m_verboseLevel >= 3)
124  {
125  std::cout << "ToyModelHadronicProcess::PostStepDoIt After scaling " << std::endl;
126  std::cout << " RHadron "<<incidentRHadron->GetDefinition()->GetParticleName()<<" pdgmass= " << incidentRHadron->GetDefinition()->GetPDGMass() / GeV
127  << " P= " << incidentRHadron->Get4Momentum() / GeV<<" mass= "<< incidentRHadron->Get4Momentum().m() / GeV <<std::endl;
128  std::cout << " Cloud "<<cloudParticle->GetDefinition()->GetParticleName()<<" pdgmass= " << cloudParticle->GetDefinition()->GetPDGMass() / GeV
129  << " P= " << cloudParticle->Get4Momentum() / GeV<<" mass= "<< cloudParticle->Get4Momentum().m() / GeV <<std::endl;
130  std::cout << " Incident pdgmass= " << incidentParticle->GetDefinition()->GetPDGMass() / GeV
131  << " P= " << incidentParticle->Get4Momentum() / GeV<<" mass= "<< incidentParticle->Get4Momentum().m() / GeV <<std::endl;
132  }
133 
134  const G4ThreeVector position = track.GetPosition();
135  const G4int incidentParticlePDG = incidentRHadron->GetDefinition()->GetPDGEncoding();
136  std::vector<G4ParticleDefinition*> newParticles; // this will contain clouds
137  std::vector<G4ParticleDefinition*> newParticlesRHadron; // this will contain r-hadron
138 
139  G4bool incidentSurvives = false;
140 
141  //Get the final state particles and target
142  G4ParticleDefinition* targetParticle;
143  HadronicProcessHelper::ReactionProduct reactionProduct = m_helper->finalState(incidentRHadron,track.GetMaterial(),targetParticle);
144 
145  // Fill a list of the new particles to create
146  //(i.e. reaction products without the incident if it survives)
147 
148  for(HadronicProcessHelper::ReactionProduct::iterator it = reactionProduct.begin();
149  it != reactionProduct.end() ;
150  it++)
151  {
152  G4ParticleDefinition * productDefinition =theParticleTable->FindParticle(*it);
153 
154  if (productDefinition->GetPDGEncoding()==incidentParticlePDG)
155  incidentSurvives = true;
156 
157  newParticlesRHadron.push_back(productDefinition);
158 
159  CustomParticle* cProd = 0;
160  cProd = dynamic_cast<CustomParticle*>(productDefinition);
161 
162  if( cProd!=0 && m_detachCloud)
163  productDefinition=cProd->GetCloud();
164 
165  newParticles.push_back(productDefinition);
166  }
167 
168  int numberOfSecondaries = reactionProduct.size();
169 
170  if(m_verboseLevel >= 2) std::cout << "ToyModelHadronicProcess::PostStepDoIt N secondaries: "
171  << numberOfSecondaries << std::endl;
172 
173 
174  //************ My toy model ********************
175  // 2 -> 2 goes to CM, chooses a random direction and fires the particles off back to back
176  // 2 -> 3 Effectively two two-body decays
177 
178  // Getting fourmomenta
179  const G4LorentzVector incident4Momentum = incidentParticle->Get4Momentum();
180  const G4LorentzVector target4Momentum(0,0,0,targetParticle->GetPDGMass());
181  const G4LorentzVector sum4Momentum = incident4Momentum + target4Momentum;
182  const G4ThreeVector cmBoost = sum4Momentum.boostVector();//The boost from CM to lab
183  const G4LorentzVector cm4Momentum = sum4Momentum.rest4Vector();
184 
185  if(m_verboseLevel >= 2) std::cout << "ToyModelHadronicProcess::PostStepDoIt Kinematics in GeV: " << std::endl
186  << " Boost = " << cmBoost / GeV<< std::endl
187  << " 4P CM = " << cm4Momentum / GeV << std::endl
188  << " 4P Inc = " << incident4Momentum / GeV << std::endl
189  << " 4P Targ = " << target4Momentum / GeV<< std::endl;
190 
191  //Choosing random direction
192  const G4double phi_p = 2*pi*G4UniformRand()-pi ;
193  const G4double theta_p = pi*G4UniformRand() ;
194  const G4ThreeVector randomDirection(sin(theta_p)*cos(phi_p),
195  sin(theta_p)*sin(phi_p),
196  cos(theta_p));
197 
198 
199  std::vector<G4double> m;
200  std::vector<G4LorentzVector> fourMomenta;
201 
202  //Fill the masses
203  for(int ip=0;ip<numberOfSecondaries;ip++)
204  {
205  m.push_back(newParticles[ip]->GetPDGMass());
206  }
207 
208 
209  if (numberOfSecondaries==2){
210  // 2 -> 2
211 
212  //Get the CM energy
213  G4double energy = cm4Momentum.e();
214  G4ThreeVector p[2];
215 
216  //Size of momenta in CM
217 
218 
219  // Energy conservation:
220  G4double cmMomentum = 1/(2*energy)*sqrt(energy*energy*energy*energy + m[0]*m[0]*m[0]*m[0] + m[1]*m[1]*m[1]*m[1]
221  - 2*(m[0]*m[0] + m[1]*m[1])*energy*energy -2*m[0]*m[0]*m[1]*m[1]);
222  p[0] = cmMomentum * randomDirection;
223  p[1] = -p[0];
224 
225  if(m_verboseLevel >= 2) std::cout << "ToyModelHadronicProcess::PostStepDoIt 2->2: " << std::endl
226  << " Pcm(GeV) = " << cmMomentum / GeV << std::endl;
227 
228  for(int ip=0;ip<2;ip++)
229  {
230  //Compute energy
231  G4double e = sqrt(p[ip].mag2() + m[ip]*m[ip]);
232  //Set 4-vectory
233  fourMomenta.push_back(G4LorentzVector(p[ip],e));
234  //Boost back to lab
235  fourMomenta[ip].boost(cmBoost);
236 
237  if(m_verboseLevel >= 2) std::cout << " particle " << ip <<" Plab(GeV) = " << fourMomenta[ip] /GeV << std::endl;
238  }
239 
240  } else if (numberOfSecondaries==3) {
241 
242  // 2 -> 3
243  //Size of momenta in CM
244  for (std::vector<G4double>::iterator it=m.begin();it!=m.end();it++)
245  fourMomenta.push_back(G4LorentzVector(0,0,0,*it));
246 
247  Decay3Body KinCalc;
248  KinCalc.doDecay(cm4Momentum, fourMomenta[0], fourMomenta[1], fourMomenta[2] );
249 
250  //Rotating the plane to a random orientation, and boosting home
251  CLHEP::HepRotation rotation(randomDirection,G4UniformRand()*2*pi);
252  for (std::vector<G4LorentzVector>::iterator it = fourMomenta.begin();
253  it!=fourMomenta.end();
254  it++)
255  {
256  *it *= rotation;
257  (*it).boost(cmBoost);
258  }
259  if(m_verboseLevel >= 3) G4cout<<"Momentum-check: "<<incident4Momentum /GeV<<" GeV vs "
260  << (fourMomenta[0]+fourMomenta[1]+fourMomenta[2])/GeV<<G4endl;
261 
262  }
263 
264 
265  //Now we have the fourMomenta of all the products (coming from 2->2 or 2->3)
266  if(incidentSurvives){
267  //if incident particle survives the number of secondaries is n-1
268  m_particleChange.SetNumberOfSecondaries(numberOfSecondaries-1);
269  if(m_verboseLevel >= 3) std::cout << "Incident survives: set num secondaries to " << numberOfSecondaries-1 << std::endl;
270 
271  } else {
272  //incident particle has to be killed and number of secondaries is n
273  m_particleChange.SetNumberOfSecondaries(numberOfSecondaries);
274  m_particleChange.ProposeTrackStatus(fStopAndKill);
275  if(m_verboseLevel >= 3) std::cout << "Incident does not survive: stopAndKill + set num secondaries to " << numberOfSecondaries << std::endl;
276  }
277  double e_kin;
278 
279  for (int ip=0; ip <numberOfSecondaries;ip++)
280  {
281  if (newParticlesRHadron[ip]==incidentRHadron->GetDefinition()) // does incident paricle survive?
282  {
283  // incidentSurvives = true; //yes! Modify its dynamic properties
284  if(m_detachCloud)
285  {
286  if(m_verboseLevel >= 3)
287  std::cout << "ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
288  fourMomenta[ip]/GeV <<"(m="<< fourMomenta[ip].m()/ GeV<<") + " <<gluinoMomentum/GeV<<std::endl; ;
289  G4LorentzVector p4_new = gluinoMomentum+fourMomenta[ip];
290  G4ThreeVector momentum = p4_new.vect();
291  double rhMass=newParticlesRHadron[ip]->GetPDGMass() ;
292  e_kin = sqrt(momentum.mag()*momentum.mag()+rhMass*rhMass)-rhMass;
293  // fourMomenta[ip]=G4LorentzVector(momentum,sqrt(momentum.mag2()+rhMass*rhMass));
294  fourMomenta[ip].setVectM(momentum,rhMass);
295 
296 // double virt=(p4_new-fourMomenta[ip]).m()/MeV;
297 
298  if(m_verboseLevel >= 3)
299  std::cout << " = " << fourMomenta[ip]/GeV <<"(m="<< fourMomenta[ip].m() / GeV<<") vs. "<<rhMass/GeV
300  << std::endl;
301  }
302  m_particleChange.ProposeMomentumDirection(fourMomenta[ip].vect()/fourMomenta[ip].vect().mag());
303  m_particleChange.ProposeEnergy(fourMomenta[ip].e()-fourMomenta[ip].mag());
304  if(m_verboseLevel >= 3) std::cout << "ToyModelHadronicProcess::PostStepDoIt Propose momentum " << fourMomenta[ip]/GeV << std::endl;
305  } else { //this particle is not the incident one
306  //Create new dynamic particle
307  G4DynamicParticle* productDynParticle = new G4DynamicParticle();
308  //Set the pDef to the dynParte
309  productDynParticle->SetDefinition(newParticlesRHadron[ip]);
310  //Set the 4-vector to dynPart
311  if(newParticlesRHadron[ip]!=newParticles[ip] && m_detachCloud ) // check if it is a cloud, second check is useless
312  {
313  G4LorentzVector p4_new;
314  p4_new.setVectM(gluinoMomentum.vect()+fourMomenta[ip].vect(),productDynParticle->GetDefinition()->GetPDGMass());
315  // productDynParticle->SetMomentum(fourMomenta[ip].vect()+gluinoMomentum.vect());
316  productDynParticle->Set4Momentum(p4_new);
317 
318 // double virt=(gluinoMomentum+fourMomenta[ip]-p4_new).m()/MeV;
319 
320  if(m_verboseLevel >= 3)
321  std::cout << "ToyModelHadronicProcess::PostStepDoIt Add gluino momentum " <<
322  fourMomenta[ip]/GeV
323  <<"(m="<< fourMomenta[ip].m()/ GeV<<") + "
324  <<gluinoMomentum/GeV
325  <<" = " << productDynParticle->Get4Momentum()/GeV
326  <<" => "<<productDynParticle->Get4Momentum().m()/GeV
327  <<" vs. "<<productDynParticle->GetDefinition()->GetPDGMass()/GeV
328  << std::endl;
329  }
330  else
331  {
332  productDynParticle->Set4Momentum(fourMomenta[ip]);
333  }
334  //Create a G4Track
335  G4Track* productTrack = new G4Track(productDynParticle,
336  track.GetGlobalTime(),
337  position);
338  productTrack->SetTouchableHandle(thisTouchable);
339  //Append to the result
340  if(m_verboseLevel >= 3) std::cout << "ToyModelHadronicProcess::PostStepDoIt Add secondary with 4-Momentum " << productDynParticle->Get4Momentum()/GeV << std::endl;
341  m_particleChange.AddSecondary(productTrack);
342  }
343  }
344 
345 
346 
347 
348 
349 
350 
351 
352 
353  //clear interaction length
354  ClearNumberOfInteractionLengthLeft();
355  //return the result
356  return &m_particleChange;
357 
358 }
359 
360 const G4DynamicParticle* ToyModelHadronicProcess::FindRhadron(G4ParticleChange* aParticleChange)
361 {
362  G4int nsec = aParticleChange->GetNumberOfSecondaries();
363  if (nsec==0) return 0;
364  int i = 0;
365  G4bool found = false;
366  while (i!=nsec && !found){
367  if (dynamic_cast<CustomParticle*>(aParticleChange->GetSecondary(i)->GetDynamicParticle()->GetDefinition())!=0) found = true;
368  i++;
369  }
370  i--;
371  if(found) return aParticleChange->GetSecondary(i)->GetDynamicParticle();
372  return 0;
373 }
int i
Definition: DBlmapReader.cc:9
void doDecay(const G4LorentzVector &mother, G4LorentzVector &daughter1, G4LorentzVector &daughter2, G4LorentzVector &daughter3)
Definition: Decay3Body.cc:23
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
G4ParticleDefinition * GetCloud()
static int position[TOTALCHAMBERS][3]
Definition: ReadPGInfo.cc:509
T sqrt(T t)
Definition: SSEVec.h:28
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
G4ParticleDefinition * GetSpectator()
tuple cout
Definition: gather_cfg.py:41
double pi