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