Geant4  10.03.p03
 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: G4ElasticHNScattering.cc 100828 2016-11-02 15:25:59Z 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 
52 #include "G4SampleResonance.hh"
53 
54 #include "G4Exp.hh"
55 #include "G4Log.hh"
56 
57 //============================================================================
58 
60 
61 
62 //============================================================================
63 
66  G4FTFParameters* theParameters ) const {
67  projectile->IncrementCollisionCount( 1 );
68  target->IncrementCollisionCount( 1 );
69 
71 
72  // Projectile parameters
73  G4LorentzVector Pprojectile = projectile->Get4Momentum();
74  if ( Pprojectile.z() < 0.0 ) return false;
75  G4bool PutOnMassShell( false );
76  G4double M0projectile = Pprojectile.mag();
77  //if ( M0projectile < projectile->GetDefinition()->GetPDGMass() ) {
78 
79  G4double MminProjectile=BrW.GetMinimumMass(projectile->GetDefinition());
80 
81  if ( M0projectile < MminProjectile ) {
82  PutOnMassShell = true;
83  M0projectile = projectile->GetDefinition()->GetPDGMass();
84  }
85  G4double M0projectile2 = M0projectile * M0projectile;
86  G4double AveragePt2 = theParameters->GetAvaragePt2ofElasticScattering();
87 
88  // Target parameters
89  G4LorentzVector Ptarget = target->Get4Momentum();
90  G4double M0target = Ptarget.mag();
91  //if ( M0target < target->GetDefinition()->GetPDGMass() ) {
92 
93  G4double MminTarget=BrW.GetMinimumMass(target->GetDefinition());
94 
95  if ( M0target < MminTarget ) {
96  PutOnMassShell = true;
97  M0target = target->GetDefinition()->GetPDGMass();
98  }
99  G4double M0target2 = M0target * M0target;
100 
101  // Transform momenta to cms and then rotate parallel to z axis;
102  G4LorentzVector Psum;
103  Psum = Pprojectile + Ptarget;
104  G4LorentzRotation toCms( -1*Psum.boostVector() );
105  G4LorentzVector Ptmp = toCms*Pprojectile;
106  if ( Ptmp.pz() <= 0.0 ) return false;
107  //"String" moving backwards in CMS, abort collision !
108  //G4cout << " abort Collision! " << G4endl;
109  toCms.rotateZ( -1*Ptmp.phi() );
110  toCms.rotateY( -1*Ptmp.theta() );
111  G4LorentzRotation toLab( toCms.inverse() );
112  Pprojectile.transform( toCms );
113  Ptarget.transform( toCms );
114 
115  // Putting on mass-on-shell, if needed
116  G4double PZcms2, PZcms;
117  G4double S = Psum.mag2();
118  G4double SqrtS = std::sqrt( S );
119  if ( SqrtS < M0projectile + M0target ) return false;
120 
121  PZcms2 = ( S*S + sqr( M0projectile2 ) + sqr( M0target2 )
122  - 2*S*M0projectile2 - 2*S*M0target2 - 2*M0projectile2*M0target2 ) / 4.0 / S;
123 
124  if ( PZcms2 < 0.0 ) { // It can be in an interaction with off-shell nuclear nucleon
125  if ( M0projectile > projectile->GetDefinition()->GetPDGMass() ) {
126  // An attempt to de-excite the projectile
127  // It is assumed that the target is in the ground state
128  M0projectile = projectile->GetDefinition()->GetPDGMass();
129  M0projectile2 = M0projectile * M0projectile;
130  PZcms2= ( S*S + sqr( M0projectile2 ) + sqr( M0target2 )
131  - 2*S*M0projectile2 - 2*S*M0target2 - 2*M0projectile2*M0target2 ) / 4.0 / S;
132  if ( PZcms2 < 0.0 ) { return false; } // Nonsuccesful attempt to de-excitate the projectile
133  } else {
134  return false; // The projectile was not excited, but the energy was too low to put
135  // the target nucleon on mass-shell
136  }
137  }
138 
139  PZcms = std::sqrt( PZcms2 );
140 
141  if ( PutOnMassShell ) {
142  if ( Pprojectile.z() > 0.0 ) {
143  Pprojectile.setPz( PZcms );
144  Ptarget.setPz( -PZcms );
145  } else {
146  Pprojectile.setPz( -PZcms );
147  Ptarget.setPz( PZcms );
148  };
149  Pprojectile.setE( std::sqrt( M0projectile2 + Pprojectile.x() * Pprojectile.x() +
150  Pprojectile.y() * Pprojectile.y() + PZcms2 ) );
151  Ptarget.setE( std::sqrt( M0target2 + Ptarget.x() * Ptarget.x() + Ptarget.y() * Ptarget.y() +
152  PZcms2 ) );
153  }
154 
155  G4double maxPtSquare = PZcms2;
156 
157  // Now we can calculate the transferred Pt
158  G4double Pt2;
159  G4double ProjMassT2, ProjMassT;
160  G4double TargMassT2, TargMassT;
161  G4LorentzVector Qmomentum;
162 
163  const G4int maxNumberOfLoops = 1000;
164  G4int loopCounter = 0;
165  do {
166  Qmomentum = G4LorentzVector( GaussianPt( AveragePt2, maxPtSquare ), 0.0 );
167  Pt2 = G4ThreeVector( Qmomentum.vect() ).mag2();
168  ProjMassT2 = M0projectile2 + Pt2;
169  ProjMassT = std::sqrt( ProjMassT2 );
170  TargMassT2 = M0target2 + Pt2;
171  TargMassT = std::sqrt( TargMassT2 );
172  } while ( ( SqrtS < ProjMassT + TargMassT ) &&
173  ++loopCounter < maxNumberOfLoops ); /* Loop checking, 10.08.2015, A.Ribon */
174  if ( loopCounter >= maxNumberOfLoops ) {
175  return false;
176  }
177 
178  PZcms2 = ( S*S + sqr( ProjMassT2 ) + sqr( TargMassT2 )
179  - 2.0*S*ProjMassT2 - 2.0*S*TargMassT2 - 2.0*ProjMassT2*TargMassT2 ) / 4.0 / S;
180 
181  if ( PZcms2 < 0.0 ) { PZcms2 = 0.0; }; // to avoid the exactness problem
182  PZcms = std::sqrt( PZcms2 );
183  Pprojectile.setPz( PZcms );
184  Ptarget.setPz( -PZcms );
185  Pprojectile += Qmomentum;
186  Ptarget -= Qmomentum;
187 
188  // Transform back and update SplitableHadron Participant.
189  Pprojectile.transform( toLab );
190  Ptarget.transform( toLab );
191 
192  // Calculation of the creation time
193  projectile->SetTimeOfCreation( target->GetTimeOfCreation() );
194  projectile->SetPosition( target->GetPosition() );
195 
196  // Creation time and position of target nucleon were determined at
197  // ReggeonCascade() of G4FTFModel
198 
199  projectile->Set4Momentum( Pprojectile );
200  target->Set4Momentum( Ptarget );
201 
202  //projectile->IncrementCollisionCount( 1 );
203  //target->IncrementCollisionCount( 1 );
204 
205  return true;
206 }
207 
208 
209 //============================================================================
210 
211 G4ThreeVector G4ElasticHNScattering::GaussianPt( G4double AveragePt2,
212  G4double maxPtSquare ) const {
213  // @@ this method is used in FTFModel as well. Should go somewhere common!
214  G4double Pt2( 0.0 );
215  if ( AveragePt2 <= 0.0 ) {
216  Pt2 = 0.0;
217  } else {
218  Pt2 = -AveragePt2 * G4Log( 1.0 + G4UniformRand() *
219  ( G4Exp( -maxPtSquare/AveragePt2 ) -1.0 ) );
220  }
221  G4double Pt = std::sqrt( Pt2 );
222  G4double phi = G4UniformRand() * twopi;
223  return G4ThreeVector( Pt * std::cos( phi ), Pt * std::sin( phi ), 0.0 );
224 }
225 
226 
227 //============================================================================
228 
230  throw G4HadronicException( __FILE__, __LINE__,
231  "G4ElasticHNScattering copy contructor not meant to be called" );
232 }
233 
234 
235 //============================================================================
236 
238 
239 
240 //============================================================================
241 
242 const G4ElasticHNScattering & G4ElasticHNScattering::operator=( const G4ElasticHNScattering& ) {
243  throw G4HadronicException( __FILE__, __LINE__,
244  "G4ElasticHNScattering = operator not meant to be called" );
245 }
246 
247 
248 //============================================================================
249 
250 int G4ElasticHNScattering::operator==( const G4ElasticHNScattering& ) const {
251  throw G4HadronicException( __FILE__, __LINE__,
252  "G4ElasticHNScattering == operator not meant to be called" );
253 }
254 
255 
256 //============================================================================
257 
258 int G4ElasticHNScattering::operator!=( const G4ElasticHNScattering& ) const {
259  throw G4HadronicException( __FILE__, __LINE__,
260  "G4ElasticHNScattering != operator not meant to be called" );
261 }
262 
virtual G4bool ElasticScattering(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner, G4FTFParameters *theParameters) const
Hep3Vector boostVector() const
const XML_Char * target
Definition: expat.h:268
double S(double temp)
CLHEP::Hep3Vector G4ThreeVector
void SetTimeOfCreation(G4double aTime)
G4double GetAvaragePt2ofElasticScattering()
HepLorentzVector & rotateZ(double)
int G4int
Definition: G4Types.hh:78
static constexpr double twopi
Definition: G4SIunits.hh:76
const G4ParticleDefinition * GetDefinition() const
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
double mag() const
bool G4bool
Definition: G4Types.hh:79
HepLorentzRotation & transform(const HepBoost &b)
void IncrementCollisionCount(G4int aCount)
const G4LorentzVector & Get4Momentum() const
G4double GetMinimumMass(const G4ParticleDefinition *p) const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
G4double GetPDGMass() const
void SetPosition(const G4ThreeVector &aPosition)
double mag2() const
HepLorentzVector & rotateY(double)
const G4ThreeVector & GetPosition() const
HepLorentzVector & transform(const HepRotation &)
void Set4Momentum(const G4LorentzVector &a4Momentum)
T sqr(const T &x)
Definition: templates.hh:145
double G4double
Definition: G4Types.hh:76
CLHEP::HepLorentzVector G4LorentzVector