Geant4  10.00.p01
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: G4ElasticHNScattering.cc 74627 2013-10-17 07:04:38Z gcosmo $
28 //
29 
30 // ------------------------------------------------------------
31 // GEANT 4 class implemetation file
32 //
33 // ---------------- G4ElasticHNScattering --------------
34 // by V. Uzhinsky, March 2008.
35 // elastic scattering used by Fritiof model
36 // Take a projectile and a target
37 // scatter the projectile and target
38 // ---------------------------------------------------------------------
39 
40 #include "globals.hh"
41 #include "Randomize.hh"
42 #include "G4PhysicalConstants.hh"
43 
44 #include "G4ElasticHNScattering.hh"
45 #include "G4LorentzRotation.hh"
46 #include "G4ThreeVector.hh"
47 #include "G4ParticleDefinition.hh"
48 #include "G4VSplitableHadron.hh"
49 #include "G4ExcitedString.hh"
50 #include "G4FTFParameters.hh"
51 //#include "G4ios.hh"
52 
53 
54 //============================================================================
55 
57 
58 
59 //============================================================================
60 
62  G4VSplitableHadron* target,
63  G4FTFParameters* theParameters ) const {
64  projectile->IncrementCollisionCount( 1 );
65  target->IncrementCollisionCount( 1 );
66 
67  // Projectile parameters
68  G4LorentzVector Pprojectile = projectile->Get4Momentum();
69  if ( Pprojectile.z() < 0.0 ) return false;
70  G4bool PutOnMassShell( false );
71  G4double M0projectile = Pprojectile.mag();
72  if ( M0projectile < projectile->GetDefinition()->GetPDGMass() ) {
73  PutOnMassShell = true;
74  M0projectile = projectile->GetDefinition()->GetPDGMass();
75  }
76  G4double M0projectile2 = M0projectile * M0projectile;
77  G4double AveragePt2 = theParameters->GetAvaragePt2ofElasticScattering();
78 
79  // Target parameters
80  G4LorentzVector Ptarget = target->Get4Momentum();
81  G4double M0target = Ptarget.mag();
82  if ( M0target < target->GetDefinition()->GetPDGMass() ) {
83  PutOnMassShell = true;
84  M0target = target->GetDefinition()->GetPDGMass();
85  }
86  G4double M0target2 = M0target * M0target;
87 
88  // Transform momenta to cms and then rotate parallel to z axis;
89  G4LorentzVector Psum;
90  Psum = Pprojectile + Ptarget;
91  G4LorentzRotation toCms( -1*Psum.boostVector() );
92  G4LorentzVector Ptmp = toCms*Pprojectile;
93  if ( Ptmp.pz() <= 0.0 ) return false;
94  // "String" moving backwards in CMS, abort collision !
95  //G4cout << " abort Collision! " << G4endl;
96  toCms.rotateZ( -1*Ptmp.phi() );
97  toCms.rotateY( -1*Ptmp.theta() );
98  G4LorentzRotation toLab( toCms.inverse() );
99  Pprojectile.transform( toCms );
100  Ptarget.transform( toCms );
101 
102  // Putting on mass-on-shell, if needed
103  G4double PZcms2, PZcms;
104  G4double S = Psum.mag2();
105  G4double SqrtS = std::sqrt( S );
106  if ( SqrtS < M0projectile + M0target ) return false;
107 
108  PZcms2 = ( S*S + sqr( M0projectile2 ) + sqr( M0target2 )
109  - 2*S*M0projectile2 - 2*S*M0target2 - 2*M0projectile2*M0target2 ) / 4.0 / S;
110 
111  if ( PZcms2 < 0.0 ) { // It can be in an interaction with off-shell nuclear nucleon
112  if ( M0projectile > projectile->GetDefinition()->GetPDGMass() ) {
113  // An attempt to de-excite the projectile
114  // It is assumed that the target is in the ground state
115  M0projectile = projectile->GetDefinition()->GetPDGMass();
116  M0projectile2 = M0projectile * M0projectile;
117  PZcms2= ( S*S + sqr( M0projectile2 ) + sqr( M0target2 )
118  - 2*S*M0projectile2 - 2*S*M0target2 - 2*M0projectile2*M0target2 ) / 4.0 / S;
119  if ( PZcms2 < 0.0 ) { return false; } // Nonsuccesful attempt to de-excitate the projectile
120  } else {
121  return false; // The projectile was not excited, but the energy was too low to put
122  // the target nucleon on mass-shell
123  }
124  }
125 
126  PZcms = std::sqrt( PZcms2 );
127 
128  if ( PutOnMassShell ) {
129  if ( Pprojectile.z() > 0.0 ) {
130  Pprojectile.setPz( PZcms );
131  Ptarget.setPz( -PZcms );
132  } else {
133  Pprojectile.setPz( -PZcms );
134  Ptarget.setPz( PZcms );
135  };
136  Pprojectile.setE( std::sqrt( M0projectile2 + Pprojectile.x() * Pprojectile.x() +
137  Pprojectile.y() * Pprojectile.y() + PZcms2 ) );
138  Ptarget.setE( std::sqrt( M0target2 + Ptarget.x() * Ptarget.x() + Ptarget.y() * Ptarget.y() +
139  PZcms2 ) );
140  }
141 
142  G4double maxPtSquare = PZcms2;
143 
144  // Now we can calculate the transferred Pt
145  G4double Pt2;
146  G4double ProjMassT2, ProjMassT;
147  G4double TargMassT2, TargMassT;
148  G4LorentzVector Qmomentum;
149 
150  do {
151  Qmomentum = G4LorentzVector( GaussianPt( AveragePt2, maxPtSquare ), 0.0 );
152  Pt2 = G4ThreeVector( Qmomentum.vect() ).mag2();
153  ProjMassT2 = M0projectile2 + Pt2;
154  ProjMassT = std::sqrt( ProjMassT2 );
155  TargMassT2 = M0target2 + Pt2;
156  TargMassT = std::sqrt( TargMassT2 );
157  } while ( SqrtS < ProjMassT + TargMassT );
158 
159  PZcms2 = ( S*S + sqr( ProjMassT2 ) + sqr( TargMassT2 )
160  - 2.0*S*ProjMassT2 - 2.0*S*TargMassT2 - 2.0*ProjMassT2*TargMassT2 ) / 4.0 / S;
161 
162  if ( PZcms2 < 0.0 ) { PZcms2 = 0.0; }; // to avoid the exactness problem
163  PZcms = std::sqrt( PZcms2 );
164  Pprojectile.setPz( PZcms );
165  Ptarget.setPz( -PZcms );
166  Pprojectile += Qmomentum;
167  Ptarget -= Qmomentum;
168 
169  // Transform back and update SplitableHadron Participant.
170  Pprojectile.transform( toLab );
171  Ptarget.transform( toLab );
172 
173  // Calculation of the creation time
174  projectile->SetTimeOfCreation( target->GetTimeOfCreation() );
175  projectile->SetPosition( target->GetPosition() );
176 
177  // Creation time and position of target nucleon were determined at
178  // ReggeonCascade() of G4FTFModel
179 
180  projectile->Set4Momentum( Pprojectile );
181  target->Set4Momentum( Ptarget );
182 
183  //projectile->IncrementCollisionCount( 1 );
184  //target->IncrementCollisionCount( 1 );
185 
186  return true;
187 }
188 
189 
190 //============================================================================
191 
193  G4double maxPtSquare ) const {
194  // @@ this method is used in FTFModel as well. Should go somewhere common!
195  G4double Pt2( 0.0 );
196  if ( AveragePt2 <= 0.0 ) {
197  Pt2 = 0.0;
198  } else {
199  Pt2 = -AveragePt2 * std::log( 1.0 + G4UniformRand() *
200  ( std::exp( -maxPtSquare/AveragePt2 ) -1.0 ) );
201  }
202  G4double Pt = std::sqrt( Pt2 );
203  G4double phi = G4UniformRand() * twopi;
204  return G4ThreeVector( Pt * std::cos( phi ), Pt * std::sin( phi ), 0.0 );
205 }
206 
207 
208 //============================================================================
209 
211  throw G4HadronicException( __FILE__, __LINE__,
212  "G4ElasticHNScattering copy contructor not meant to be called" );
213 }
214 
215 
216 //============================================================================
217 
219 
220 
221 //============================================================================
222 
224  throw G4HadronicException( __FILE__, __LINE__,
225  "G4ElasticHNScattering = operator not meant to be called" );
226 }
227 
228 
229 //============================================================================
230 
232  throw G4HadronicException( __FILE__, __LINE__,
233  "G4ElasticHNScattering == operator not meant to be called" );
234 }
235 
236 
237 //============================================================================
238 
240  throw G4HadronicException( __FILE__, __LINE__,
241  "G4ElasticHNScattering != operator not meant to be called" );
242 }
virtual G4bool ElasticScattering(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner, G4FTFParameters *theParameters) const
CLHEP::Hep3Vector G4ThreeVector
int operator==(const G4ElasticHNScattering &right) const
CLHEP::HepLorentzRotation G4LorentzRotation
void SetTimeOfCreation(G4double aTime)
G4double GetAvaragePt2ofElasticScattering()
G4ParticleDefinition * GetDefinition() const
#define G4UniformRand()
Definition: Randomize.hh:87
const G4ElasticHNScattering & operator=(const G4ElasticHNScattering &right)
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector GaussianPt(G4double AveragePt2, G4double maxPtSquare) const
void IncrementCollisionCount(G4int aCount)
const G4LorentzVector & Get4Momentum() const
G4double GetPDGMass() const
void SetPosition(const G4ThreeVector &aPosition)
const G4ThreeVector & GetPosition() const
int operator!=(const G4ElasticHNScattering &right) const
void Set4Momentum(const G4LorentzVector &a4Momentum)
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
CLHEP::HepLorentzVector G4LorentzVector