Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4SingleDiffractiveExcitation.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: G4SingleDiffractiveExcitation.cc 100828 2016-11-02 15:25:59Z gcosmo $
28 // ------------------------------------------------------------
29 // GEANT 4 class implemetation file
30 //
31 // ---------------- G4SingleDiffractiveExcitation --------------
32 // by Gunter Folger, October 1998.
33 // diffractive Excitation used by strings models
34 // Take a projectile and a target
35 // excite the projectile and target
36 // ------------------------------------------------------------
37 
39 #include "globals.hh"
40 #include "G4PhysicalConstants.hh"
41 #include "G4SystemOfUnits.hh"
42 #include "Randomize.hh"
43 #include "G4LorentzRotation.hh"
44 #include "G4ThreeVector.hh"
45 #include "G4ParticleDefinition.hh"
46 #include "G4VSplitableHadron.hh"
47 #include "G4ExcitedString.hh"
48 
49 #include "G4Log.hh"
50 
51 
53 : widthOfPtSquare(-2*sqr(sigmaPt)) , minExtraMass(minextraMass), minmass(x0mass)
54 {}
55 
57 {}
58 
61 {
62  G4LorentzVector Pprojectile=projectile->Get4Momentum();
63  G4double Mprojectile2=sqr(projectile->GetDefinition()->GetPDGMass() + minExtraMass);
64 
65  G4LorentzVector Ptarget=target->Get4Momentum();
66  G4double Mtarget2=sqr(target->GetDefinition()->GetPDGMass() + minExtraMass);
67  //G4cout << "E proj, target :" << Pprojectile.e() << ", " << Ptarget.e() << G4endl;
68 
69  G4bool KeepProjectile= G4UniformRand() > 0.5;
70 
71  // reset the min.mass of the non diffractive particle to its value, ( minus a bit for rounding...)
72  if ( KeepProjectile )
73  {
74  //cout << " Projectile fix" << G4endl;
75  Mprojectile2 = sqr(projectile->GetDefinition()->GetPDGMass() * (1-perCent) );
76  } else {
77  //cout << " Target fix" << G4endl;
78  Mtarget2=sqr(target->GetDefinition()->GetPDGMass() * (1-perCent) );
79  }
80 
81  // Transform momenta to cms and then rotate parallel to z axis;
82 
83  G4LorentzVector Psum;
84  Psum=Pprojectile+Ptarget;
85 
86  G4LorentzRotation toCms(-1*Psum.boostVector());
87 
88  G4LorentzVector Ptmp=toCms*Pprojectile;
89 
90  if ( Ptmp.pz() <= 0. )
91  {
92  // "String" moving backwards in CMS, abort collision !!
93  //G4cout << " abort Collision!! " << G4endl;
94  return false;
95  }
96 
97  toCms.rotateZ(-1*Ptmp.phi());
98  toCms.rotateY(-1*Ptmp.theta());
99 
100  //G4cout << "Pprojectile be4 boost " << Pprojectile << G4endl;
101  //G4cout << "Ptarget be4 boost : " << Ptarget << G4endl;
102 
103  G4LorentzRotation toLab(toCms.inverse());
104 
105  Pprojectile.transform(toCms);
106  Ptarget.transform(toCms);
107 
108  G4LorentzVector Qmomentum;
109  G4int whilecount=0;
110  do {
111  // Generate pt
112 
113  G4double maxPtSquare=sqr(Ptarget.pz());
114  if (whilecount++ >= 500 && (whilecount%100)==0)
115  //G4cout << "G4SingleDiffractiveExcitation::ExciteParticipants possibly looping"
116  // << ", loop count/ maxPtSquare : "
117  // << whilecount << " / " << maxPtSquare << G4endl;
118 
119  if (whilecount > 1000 )
120  {
121  Qmomentum=G4LorentzVector(0.,0.,0.,0.);
122  //G4cout << "G4SingleDiffractiveExcitation::ExciteParticipants: Aborting loop!" << G4endl;
123  return false; // Ignore this interaction
124  }
125  Qmomentum=G4LorentzVector(GaussianPt(widthOfPtSquare,maxPtSquare),0);
126 
127  // Momentum transfer
128  G4double Xmin = minmass / ( Pprojectile.e() + Ptarget.e() );
129  G4double Xmax=1.;
130  G4double Xplus =ChooseX(Xmin,Xmax);
131  G4double Xminus=ChooseX(Xmin,Xmax);
132 
133  G4double pt2=G4ThreeVector(Qmomentum.vect()).mag2();
134  G4double Qplus =-1 * pt2 / Xminus/Ptarget.minus();
135  G4double Qminus= pt2 / Xplus /Pprojectile.plus();
136 
137  if ( KeepProjectile )
138  {
139  Qminus = (sqr(projectile->GetDefinition()->GetPDGMass()) + pt2 )
140  / (Pprojectile.plus() + Qplus ) - Pprojectile.minus();
141  } else {
142  Qplus = Ptarget.plus() - (sqr(target->GetDefinition()->GetPDGMass()) + pt2 )
143  / (Ptarget.minus() - Qminus );
144  }
145 
146  Qmomentum.setPz( (Qplus-Qminus)/2 );
147  Qmomentum.setE( (Qplus+Qminus)/2 );
148 
149  //G4cout << "Qplus / Qminus " << Qplus << " / " << Qminus<<G4endl;
150  //G4cout << "pt2 " << pt2 << G4endl;
151  //G4cout << "Qmomentum " << Qmomentum << G4endl;
152  //G4cout << " Masses (P/T) : " << (Pprojectile+Qmomentum).mag() <<
153  // " / " << (Ptarget-Qmomentum).mag() << G4endl;
154 
155  } while ( (Ptarget-Qmomentum).mag2() <= Mtarget2 /* Loop checking, 26.10.2015, A.Ribon */
156  || (Pprojectile+Qmomentum).mag2() <= Mprojectile2
157  || (Ptarget-Qmomentum).e() < 0.
158  || (Pprojectile+Qmomentum).e() < 0. );
159 
160  //G4double Ecms=Pprojectile.e() + Ptarget.e();
161 
162  Pprojectile += Qmomentum;
163  Ptarget -= Qmomentum;
164 
165  //G4cout << "Pprojectile.e() : " << Pprojectile.e() << G4endl;
166  //G4cout << "Ptarget.e() : " << Ptarget.e() << G4endl;
167  //G4cout << "end event_______________________________________________"<<G4endl;
168  //G4cout << "Pprojectile with Q : " << Pprojectile << G4endl;
169  //G4cout << "Ptarget with Q : " << Ptarget << G4endl;
170  //G4cout << "Projectile back: " << toLab * Pprojectile << G4endl;
171  //G4cout << "Target back: " << toLab * Ptarget << G4endl;
172 
173  // Transform back and update SplitableHadron Participant.
174  Pprojectile.transform(toLab);
175  Ptarget.transform(toLab);
176 
177  //G4cout << "G4SingleDiffractiveExcitation- Target mass " << Ptarget.mag() << G4endl;
178  //G4cout << "G4SingleDiffractiveExcitation- Projectile mass " << Pprojectile.mag() << G4endl;
179 
180  target->Set4Momentum(Ptarget);
181  projectile->Set4Momentum(Pprojectile);
182 
183  return true;
184 }
185 
186 
187 // --------- private methods ----------------------
188 
189 G4double G4SingleDiffractiveExcitation::ChooseX(G4double Xmin, G4double Xmax) const
190 {
191  // choose an x between Xmin and Xmax with P(x) ~ 1/x
192  // to be improved...
193 
194  G4double range=Xmax-Xmin;
195 
196  if ( Xmin <= 0. || range <=0. )
197  {
198  G4cout << " Xmin, range : " << Xmin << " , " << range << G4endl;
199  throw G4HadronicException(__FILE__, __LINE__, "G4SingleDiffractiveExcitation::ChooseX : Invalid arguments ");
200  }
201 
202  G4double x;
203  do {
204  x=Xmin + G4UniformRand() * range;
205  } while ( Xmin/x < G4UniformRand() ); /* Loop checking, 26.10.2015, A.Ribon */
206 
207  //cout << "DiffractiveX "<<x<<G4endl;
208  return x;
209 }
210 
211 G4ThreeVector G4SingleDiffractiveExcitation::GaussianPt(G4double widthSquare, G4double maxPtSquare) const
212 {
213  // @@ this method is used in FTFModel as well. Should go somewhere common!
214 
215  G4double pt2;
216 
217  const G4int maxNumberOfLoops = 1000;
218  G4int loopCounter = -1;
219  do {
220  pt2=widthSquare * G4Log( G4UniformRand() );
221  } while ( ( pt2 > maxPtSquare) && ++loopCounter < maxNumberOfLoops ); /* Loop checking, 26.10.2015, A.Ribon */
222  if ( loopCounter >= maxNumberOfLoops ) pt2 = 0.0;
223 
224  pt2=std::sqrt(pt2);
225 
226  G4double phi=G4UniformRand() * twopi;
227 
228  return G4ThreeVector (pt2*std::cos(phi), pt2*std::sin(phi), 0.);
229 }
230 
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
Hep3Vector boostVector() const
const XML_Char * target
Definition: expat.h:268
CLHEP::Hep3Vector G4ThreeVector
static constexpr double perCent
Definition: G4SIunits.hh:332
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
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
HepLorentzRotation & transform(const HepBoost &b)
const G4LorentzVector & Get4Momentum() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double GetPDGMass() const
G4SingleDiffractiveExcitation(G4double sigmaPt=0.6 *CLHEP::GeV, G4double minExtraMass=250 *CLHEP::MeV, G4double x0mass=250 *CLHEP::MeV)
HepLorentzVector & rotateY(double)
#define G4endl
Definition: G4ios.hh:61
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