Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ElasticHNScattering.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 // ------------------------------------------------------------
29 // GEANT 4 class implemetation file
30 //
31 // ---------------- G4ElasticHNScattering --------------
32 // by V. Uzhinsky, March 2008.
33 // elastic scattering used by Fritiof model
34 // Take a projectile and a target
35 // scatter the projectile and target
36 // ---------------------------------------------------------------------
37 
38 
39 #include "globals.hh"
40 #include "Randomize.hh"
41 #include "G4PhysicalConstants.hh"
42 
43 #include "G4ElasticHNScattering.hh"
44 #include "G4LorentzRotation.hh"
45 #include "G4ThreeVector.hh"
46 #include "G4ParticleDefinition.hh"
47 #include "G4VSplitableHadron.hh"
48 #include "G4ExcitedString.hh"
49 #include "G4FTFParameters.hh"
50 //#include "G4ios.hh"
51 
53 {
54 }
55 
59  G4FTFParameters *theParameters) const
60 {
61 // -------------------- Projectile parameters -----------------------------------
62  G4LorentzVector Pprojectile=projectile->Get4Momentum();
63 
64  if(Pprojectile.z() < 0.)
65  {
66  target->SetStatus(2);
67  return false;
68  }
69 
70  G4bool PutOnMassShell(false);
71 
72  G4double M0projectile = Pprojectile.mag();
73  if(M0projectile < projectile->GetDefinition()->GetPDGMass())
74  {
75  PutOnMassShell=true;
76  M0projectile=projectile->GetDefinition()->GetPDGMass();
77  }
78 
79  G4double Mprojectile2 = M0projectile * M0projectile;
80 
81  G4double AveragePt2=theParameters->GetAvaragePt2ofElasticScattering();
82 
83 // -------------------- Target parameters ----------------------------------------------
84 
85  G4LorentzVector Ptarget=target->Get4Momentum();
86 
87  G4double M0target = Ptarget.mag();
88 
89  if(M0target < target->GetDefinition()->GetPDGMass())
90  {
91  PutOnMassShell=true;
92  M0target=target->GetDefinition()->GetPDGMass();
93  }
94 
95  G4double Mtarget2 = M0target * M0target;
96 
97 // Transform momenta to cms and then rotate parallel to z axis;
98 
99  G4LorentzVector Psum;
100  Psum=Pprojectile+Ptarget;
101 
102  G4LorentzRotation toCms(-1*Psum.boostVector());
103 
104  G4LorentzVector Ptmp=toCms*Pprojectile;
105 
106  if ( Ptmp.pz() <= 0. )
107  {
108  // "String" moving backwards in CMS, abort collision !!
109  //G4cout << " abort Collision!! " << G4endl;
110  target->SetStatus(2);
111  return false;
112  }
113 
114  toCms.rotateZ(-1*Ptmp.phi());
115  toCms.rotateY(-1*Ptmp.theta());
116 
117  G4LorentzRotation toLab(toCms.inverse());
118 
119  Pprojectile.transform(toCms);
120  Ptarget.transform(toCms);
121 
122 // ---------------------- Putting on mass-on-shell, if needed ------------------------
123  G4double PZcms2, PZcms;
124 
125  G4double S=Psum.mag2();
126 // G4double SqrtS=std::sqrt(S);
127 
128  PZcms2=(S*S+Mprojectile2*Mprojectile2+Mtarget2*Mtarget2-
129  2*S*Mprojectile2-2*S*Mtarget2-2*Mprojectile2*Mtarget2)/4./S;
130 
131  if(PZcms2 < 0.)
132  { // It can be in an interaction with off-shell nuclear nucleon
133  if(M0projectile > projectile->GetDefinition()->GetPDGMass())
134  { // An attempt to de-excite the projectile
135  // It is assumed that the target is in the ground state
136  M0projectile = projectile->GetDefinition()->GetPDGMass();
137  Mprojectile2=M0projectile*M0projectile;
138  PZcms2=(S*S+Mprojectile2*Mprojectile2+Mtarget2*Mtarget2-
139  2*S*Mprojectile2 - 2*S*Mtarget2 - 2*Mprojectile2*Mtarget2)
140  /4./S;
141 
142  if(PZcms2 < 0.){ return false;} // Non succesful attempt after the de-excitation
143  }
144  else // if(M0projectile > projectile->GetDefinition()->GetPDGMass())
145  {
146  target->SetStatus(2);
147  return false; // The projectile was not excited,
148  // but the energy was too low to put
149  // the target nucleon on mass-shell
150  } // end of if(M0projectile > projectile->GetDefinition()->GetPDGMass())
151  } // end of if(PZcms2 < 0.)
152 
153  PZcms = std::sqrt(PZcms2);
154 
155  if(PutOnMassShell)
156  {
157  if(Pprojectile.z() > 0.)
158  {
159  Pprojectile.setPz( PZcms);
160  Ptarget.setPz( -PZcms);
161  }
162  else // if(Pprojectile.z() > 0.)
163  {
164  Pprojectile.setPz(-PZcms);
165  Ptarget.setPz( PZcms);
166  };
167 
168  Pprojectile.setE(std::sqrt(Mprojectile2+
169  Pprojectile.x()*Pprojectile.x()+
170  Pprojectile.y()*Pprojectile.y()+
171  PZcms2));
172  Ptarget.setE(std::sqrt( Mtarget2 +
173  Ptarget.x()*Ptarget.x()+
174  Ptarget.y()*Ptarget.y()+
175  PZcms2));
176  } // end of if(PutOnMassShell)
177 
178  G4double maxPtSquare = PZcms2;
179 
180 // ------ Now we can calculate the transfered Pt --------------------------
181  G4double Pt2;
182  G4double ProjMassT2; //, ProjMassT;
183  G4double TargMassT2; //, TargMassT;
184 
185  G4LorentzVector Qmomentum;
186  Qmomentum=G4LorentzVector(GaussianPt(AveragePt2,maxPtSquare),0);
187 
188  Pt2=G4ThreeVector(Qmomentum.vect()).mag2();
189 
190  ProjMassT2=Mprojectile2+Pt2;
191 // ProjMassT =std::sqrt(ProjMassT2);
192 
193  TargMassT2=Mtarget2+Pt2;
194 // TargMassT =std::sqrt(TargMassT2);
195 
196  PZcms2=(S*S+ProjMassT2*ProjMassT2+
197  TargMassT2*TargMassT2-
198  2.*S*ProjMassT2-2.*S*TargMassT2-
199  2.*ProjMassT2*TargMassT2)/4./S;
200  if(PZcms2 < 0 ) {PZcms2=0;};// to avoid the exactness problem
201  PZcms =std::sqrt(PZcms2);
202 
203  Pprojectile.setPz( PZcms);
204  Ptarget.setPz( -PZcms);
205 
206  Pprojectile += Qmomentum;
207  Ptarget -= Qmomentum;
208 
209 // Transform back and update SplitableHadron Participant.
210  Pprojectile.transform(toLab);
211  Ptarget.transform(toLab);
212 /* // Maybe it will be needed for an exact calculations--------------------
213  G4double TargetMomentum=std::sqrt(Ptarget.x()*Ptarget.x()+
214  Ptarget.y()*Ptarget.y()+
215  Ptarget.z()*Ptarget.z());
216 */
217 
218 // Calculation of the creation time ---------------------
219  projectile->SetTimeOfCreation(target->GetTimeOfCreation());
220  projectile->SetPosition(target->GetPosition());
221 // Creation time and position of target nucleon were determined at
222 // ReggeonCascade() of G4FTFModel
223 // ------------------------------------------------------
224 
225  projectile->Set4Momentum(Pprojectile);
226  target->Set4Momentum(Ptarget);
227 
228  projectile->IncrementCollisionCount(1);
229  target->IncrementCollisionCount(1);
230 
231  return true;
232 }
233 
234 
235 // --------- private methods ----------------------
236 
237 G4ThreeVector G4ElasticHNScattering::GaussianPt(G4double AveragePt2, G4double maxPtSquare) const
238 { // @@ this method is used in FTFModel as well. Should go somewhere common!
239 
240  G4double Pt2(0.);
241  if(AveragePt2 <= 0.) {Pt2=0.;}
242  else
243  {
244  Pt2 = -AveragePt2 * std::log(1. + G4UniformRand() *
245  (std::exp(-maxPtSquare/AveragePt2)-1.));
246  }
247  G4double Pt=std::sqrt(Pt2);
248  G4double phi=G4UniformRand() * twopi;
249 
250  return G4ThreeVector (Pt*std::cos(phi), Pt*std::sin(phi), 0.);
251 }
252 
254 {
255  throw G4HadronicException(__FILE__, __LINE__, "G4ElasticHNScattering copy contructor not meant to be called");
256 }
257 
258 
260 {
261 }
262 
263 
264 const G4ElasticHNScattering & G4ElasticHNScattering::operator=(const G4ElasticHNScattering &)
265 {
266  throw G4HadronicException(__FILE__, __LINE__, "G4ElasticHNScattering = operator not meant to be called");
267  //return *this; //A.R. 25-Jul-2012 : fix Coverity
268 }
269 
270 
271 int G4ElasticHNScattering::operator==(const G4ElasticHNScattering &) const
272 {
273  throw G4HadronicException(__FILE__, __LINE__, "G4ElasticHNScattering == operator not meant to be called");
274  //return false; //A.R. 25-Jul-2012 : fix Coverity
275 }
276 
277 int G4ElasticHNScattering::operator!=(const G4ElasticHNScattering &) const
278 {
279  throw G4HadronicException(__FILE__, __LINE__, "G4ElasticHNScattering != operator not meant to be called");
280  //return true; //A.R. 25-Jul-2012 : fix Coverity
281 }