Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4FTFParticipants.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 //
27 // $Id$
28 // GEANT4 tag $Name: $
29 //
30 // ------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // ---------------- G4FTFParticipants----------------
34 // by Gunter Folger, June 1998.
35 // class finding colliding particles in FTFPartonStringModel
36 // Changed in a part by V. Uzhinsky in oder to put in correcpondence
37 // with original FRITIOF mode. November - December 2006.
38 // Ajusted for (anti) nucleus - nucleus interactions by V. Uzhinsky.
39 // (February 2011)
40 // ------------------------------------------------------------
41 
42 #include <utility>
43 #include <vector>
44 #include <algorithm>
45 
46 #include "G4FTFParticipants.hh"
47 #include "G4ios.hh"
48 #include "Randomize.hh"
49 #include "G4SystemOfUnits.hh"
50 #include "G4FTFParameters.hh" // Uzhi 29.03.08
52 #include "G4VSplitableHadron.hh"
53 
55  theProjectileNucleus(0),
56  currentInteraction(-1)
57 {
58 }
59 
61  , theProjectileNucleus(0), currentInteraction(-1) //A.R. 14-Aug-2012 Coverity fix.
62 {
63  G4Exception("G4FTFParticipants::G4FTFParticipants()","HAD_FTF_001",
64  FatalException," Must not use copy ctor()");
65 }
66 
67 
69 {
70  if ( theProjectileNucleus != NULL ) delete theProjectileNucleus;
71 }
72 
73 //-------------------------------------------------------------------------
74 
76 {
78 
79  theProjectileNucleus = aNucleus;
80 }
81 
83 {
84  return theProjectileNucleus;
85 }
86 
88 {
90  theProjectileNucleus->Init(theA, theZ);
92 }
93 //-------------------------------------------------------------------------
95  G4FTFParameters *theParameters)
96 {
97 //G4cout<<"Participants::GetList"<<G4endl;
98 //G4cout<<"thePrimary "<<thePrimary.GetMomentum()<<G4endl;
99  StartLoop(); // reset Loop over Interactions
100 
101  for(unsigned int i=0; i<theInteractions.size(); i++) delete theInteractions[i];
102  theInteractions.clear();
103 
104  G4double deltaxy=2 * fermi; // Extra nuclear radius
105 //G4cout<<"theProjectileNucleus "<<theProjectileNucleus<<G4endl;
106  if(theProjectileNucleus == 0)
107  { // Hadron-nucleus or anti-baryon-nucleus interactions
108 //G4cout<<"Hadron-nucleus or anti-baryon-nucleus interactions"<<G4endl;
109 
110  G4double impactX(0.), impactY(0.);
111 
112  G4VSplitableHadron * primarySplitable=new G4DiffractiveSplitableHadron(thePrimary);
113 //G4cout<<"Prim in Part "<<primarySplitable->Get4Momentum()<<G4endl;
114  G4double xyradius;
115  xyradius =theNucleus->GetOuterRadius() + deltaxy; // Impact parameter sampling
116 
117 // G4bool nucleusNeedsShift = true; // Uzhi 20 July 2009
118 
119  do
120  {
121  std::pair<G4double, G4double> theImpactParameter;
122  theImpactParameter = theNucleus->ChooseImpactXandY(xyradius);
123  impactX = theImpactParameter.first;
124  impactY = theImpactParameter.second;
125 
126  G4ThreeVector thePosition(impactX, impactY, -DBL_MAX);
127  primarySplitable->SetPosition(thePosition);
128 
130  G4Nucleon * nucleon;
131 
132 G4int TrN(0);
133  while ( (nucleon=theNucleus->GetNextNucleon()) )
134  {
135  G4double impact2= sqr(impactX - nucleon->GetPosition().x()) +
136  sqr(impactY - nucleon->GetPosition().y());
137 
138  if ( theParameters->GetProbabilityOfInteraction(impact2/fermi/fermi)
139  > G4UniformRand() )
140  {
141  primarySplitable->SetStatus(1); // It takes part in the interaction
142 
143  G4VSplitableHadron * targetSplitable=0;
144  if ( ! nucleon->AreYouHit() )
145  {
146  targetSplitable= new G4DiffractiveSplitableHadron(*nucleon);
147  nucleon->Hit(targetSplitable);
148  nucleon->SetBindingEnergy(3.*nucleon->GetBindingEnergy());
149 //G4cout<<" Part nucl "<<TrN<<" "<<nucleon->Get4Momentum()<<G4endl;
150 //G4cout<<" Part nucl "<<G4endl;
151  targetSplitable->SetStatus(1); // It takes part in the interaction
152  }
153  G4InteractionContent * aInteraction =
154  new G4InteractionContent(primarySplitable);
155  aInteraction->SetTarget(targetSplitable);
156  aInteraction->SetTargetNucleon(nucleon); // Uzhi 16.07.09
157  aInteraction->SetStatus(1); // Uzhi Feb26
158  theInteractions.push_back(aInteraction);
159  }
160 TrN++;
161  }
162  } while ( theInteractions.size() == 0 );
163 
164  //if ( theInteractions.size() == 0 ) delete primarySplitable; //A.R. 14-Aug-2012 Coverity fix
165 
166 // G4cout << "Number of Hit nucleons " << theInteractions.size()
167 // << "\t" << impactX/fermi << "\t"<<impactY/fermi
168 // << "\t" << std::sqrt(sqr(impactX)+sqr(impactY))/fermi <<G4endl;
169 
170  return;
171  } // end of if(theProjectileNucleus == 0)
172 
173 //-------------------------------------------------------------------
174 // Projectile and target are nuclei
175 //-------------------------------------------------------------------
176 //VU G4VSplitableHadron * primarySplitable=new G4DiffractiveSplitableHadron(thePrimary);
177 //G4cout<<"Prim in Part "<<primarySplitable->Get4Momentum()<<G4endl;
178 //G4cout<<"Projectile and target are nuclei"<<G4endl;
179 //G4cout<<thePrimary.GetMomentum()<<G4endl;
180 //G4cout<<"Part Pr Tr "<<theProjectileNucleus<<" "<<theNucleus<<G4endl;
181 
182 
183  G4double xyradius;
184  xyradius =theProjectileNucleus->GetOuterRadius() + // Impact parameter sampling
185  theNucleus->GetOuterRadius() + deltaxy;
186 
187  G4double impactX(0.), impactY(0.);
188 
189  do
190  {
191 //G4cout<<"New interaction list"<<G4endl;
192  std::pair<G4double, G4double> theImpactParameter;
193  theImpactParameter = theNucleus->ChooseImpactXandY(xyradius);
194  impactX = theImpactParameter.first;
195  impactY = theImpactParameter.second;
196 //G4cout<<"B "<<std::sqrt(sqr(impactX)+sqr(impactY))/fermi<<G4endl;
197 
198  G4ThreeVector thePosition(impactX, impactY, -DBL_MAX);
199 //VU primarySplitable->SetPosition(thePosition);
200 
202  G4Nucleon * ProjectileNucleon;
203 G4int PrNuclN(0);
204 
205  while ( (ProjectileNucleon=theProjectileNucleus->GetNextNucleon()) )
206  {
207  G4VSplitableHadron * ProjectileSplitable=0;
208 //G4cout<<G4endl<<"Prj N mom "<<ProjectileNucleon->Get4Momentum()<<"-------------"<<G4endl;
210  G4Nucleon * TargetNucleon;
211 
212 G4int TrNuclN(0);
213  while ( (TargetNucleon=theNucleus->GetNextNucleon()) )
214  {
215 //G4cout<<"Trg N mom "<<TargetNucleon->Get4Momentum()<<G4endl;
216  G4double impact2=
217  sqr(impactX+ProjectileNucleon->GetPosition().x()-TargetNucleon->GetPosition().x())+
218  sqr(impactY+ProjectileNucleon->GetPosition().y()-TargetNucleon->GetPosition().y());
219 
220  G4VSplitableHadron * TargetSplitable=0;
221 
222  if ( theParameters->GetProbabilityOfInteraction(impact2/fermi/fermi)
223  > G4UniformRand() )
224  { // An Interaction has happend!
225 //G4cout<<"An Interaction has happend"<<G4endl;
226 //G4cout<<"PrN TrN "<<PrNuclN<<" "<<TrNuclN<<" "<<ProjectileNucleon->GetPosition().z()/fermi<<" "<<TargetNucleon->GetPosition().z()/fermi<<" "<<ProjectileNucleon->GetPosition().z()/fermi + TargetNucleon->GetPosition().z()/fermi <<G4endl;
227 
228  if ( ! ProjectileNucleon->AreYouHit() )
229  { // Projectile nucleon was not involved until now.
230  ProjectileSplitable= new G4DiffractiveSplitableHadron(*ProjectileNucleon);
231  ProjectileNucleon->Hit(ProjectileSplitable);
232  ProjectileNucleon->SetBindingEnergy(3.*ProjectileNucleon->GetBindingEnergy());
233  ProjectileSplitable->SetStatus(1); // It takes part in the interaction
234  }
235  else
236  { // Projectile nucleon was involved before.
237  ProjectileSplitable=ProjectileNucleon->GetSplitableHadron();
238  } // End of if ( ! Projectileucleon->AreYouHit() )
239 
240  if ( ! TargetNucleon->AreYouHit() )
241  { // Target nucleon was not involved until now
242  TargetSplitable= new G4DiffractiveSplitableHadron(*TargetNucleon);
243  TargetNucleon->Hit(TargetSplitable);
244  TargetNucleon->SetBindingEnergy(3.*ProjectileNucleon->GetBindingEnergy());
245  TargetSplitable->SetStatus(1); // It takes part in the interaction
246  }
247  else
248  { // Target nucleon was involved before.
249  TargetSplitable=TargetNucleon->GetSplitableHadron();
250  } // End of if ( ! TargetNeucleon->AreYouHit() )
251 
252  G4InteractionContent * anInteraction =
253  new G4InteractionContent(ProjectileSplitable);
254  anInteraction->SetTarget(TargetSplitable);
255  anInteraction->SetTargetNucleon(TargetNucleon);
256  anInteraction->SetStatus(1); // Uzhi Feb26
257 // anInteraction->SetInteractionTime(ProjectileNucleon->GetPosition().z()+
258 // TargetNucleon->GetPosition().z());
259 //G4cout<<"Z's pr tr "<<ProjectileNucleon->GetPosition().z()/fermi<<" "<<TargetNucleon->GetPosition().z()/fermi<<" "<<ProjectileNucleon->GetPosition().z()/fermi + TargetNucleon->GetPosition().z()/fermi <<G4endl;
260  theInteractions.push_back(anInteraction);
261 //G4cout<<"Ppr tr "<<ProjectileSplitable<<" "<<TargetSplitable<<G4endl;
262  } // End of An Interaction has happend!
263 TrNuclN++;
264  } // End of while ( (TargetNucleon=theNucleus->GetNextNucleon()) )
265 PrNuclN++;
266  } // End of while ( (ProjectileNucleon=theProjectileNucleus->GetNextNucleon()) )
267  } while ( theInteractions.size() == 0 ); // end of while ( theInteractions.size() == 0 )
268 
269 //std::sort(theInteractions.begin(),theInteractions.end()); // ????
270 
271 // G4cout << "Number of primary collisions " << theInteractions.size()
272 // << "\t" << impactX/fermi << "\t"<<impactY/fermi
273 // << "\t" << std::sqrt(sqr(impactX)+sqr(impactY))/fermi <<G4endl;
274 //G4int Uzhi; G4cin >> Uzhi;
275  return;
276 }
277 //--------------------------------------------------------------
278 
279 // Implementation (private) methods