Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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$
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"
41#include "G4SystemOfUnits.hh"
42#include "Randomize.hh"
43#include "G4LorentzRotation.hh"
44#include "G4ThreeVector.hh"
46#include "G4VSplitableHadron.hh"
47#include "G4ExcitedString.hh"
48
50:
51widthOfPtSquare(-2*sqr(sigmaPt)) , minExtraMass(minextraMass),
52minmass(x0mass)
53{}
54
56{}
57
60{
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() << ", " <<
68 // Ptarget.e() << G4endl;
69
70 G4bool KeepProjectile= G4UniformRand() > 0.5;
71
72 // reset the min.mass of the non diffractive particle to its value, ( minus a bit for rounding...)
73 if ( KeepProjectile )
74 {
75 // cout << " Projectile fix" << G4endl;
76 Mprojectile2 = sqr(projectile->GetDefinition()->GetPDGMass() * (1-perCent) );
77 } else {
78 // cout << " Target fix" << G4endl;
79 Mtarget2=sqr(target->GetDefinition()->GetPDGMass() * (1-perCent) );
80 }
81
82 // Transform momenta to cms and then rotate parallel to z axis;
83
84 G4LorentzVector Psum;
85 Psum=Pprojectile+Ptarget;
86
87 G4LorentzRotation toCms(-1*Psum.boostVector());
88
89 G4LorentzVector Ptmp=toCms*Pprojectile;
90
91 if ( Ptmp.pz() <= 0. )
92 {
93 // "String" moving backwards in CMS, abort collision !!
94 // G4cout << " abort Collision!! " << G4endl;
95 return false;
96 }
97
98 toCms.rotateZ(-1*Ptmp.phi());
99 toCms.rotateY(-1*Ptmp.theta());
100
101 // G4cout << "Pprojectile be4 boost " << Pprojectile << G4endl;
102 // G4cout << "Ptarget be4 boost : " << Ptarget << G4endl;
103
104
105
106 G4LorentzRotation toLab(toCms.inverse());
107
108 Pprojectile.transform(toCms);
109 Ptarget.transform(toCms);
110
111 G4LorentzVector Qmomentum;
112 G4int whilecount=0;
113 do {
114 // Generate pt
115
116 G4double maxPtSquare=sqr(Ptarget.pz());
117 if (whilecount++ >= 500 && (whilecount%100)==0)
118 // G4cout << "G4SingleDiffractiveExcitation::ExciteParticipants possibly looping"
119 // << ", loop count/ maxPtSquare : "
120 // << whilecount << " / " << maxPtSquare << G4endl;
121 if (whilecount > 1000 )
122 {
123 Qmomentum=G4LorentzVector(0.,0.,0.,0.);
124 // G4cout << "G4SingleDiffractiveExcitation::ExciteParticipants: Aborting loop!" << G4endl;
125 return false; // Ignore this interaction
126 }
127 Qmomentum=G4LorentzVector(GaussianPt(widthOfPtSquare,maxPtSquare),0);
128
129
130 // Momentum transfer
131 G4double Xmin = minmass / ( Pprojectile.e() + Ptarget.e() );
132 G4double Xmax=1.;
133 G4double Xplus =ChooseX(Xmin,Xmax);
134 G4double Xminus=ChooseX(Xmin,Xmax);
135
136 G4double pt2=G4ThreeVector(Qmomentum.vect()).mag2();
137 G4double Qplus =-1 * pt2 / Xminus/Ptarget.minus();
138 G4double Qminus= pt2 / Xplus /Pprojectile.plus();
139
140 if ( KeepProjectile )
141 {
142 Qminus = (sqr(projectile->GetDefinition()->GetPDGMass()) + pt2 )
143 / (Pprojectile.plus() + Qplus )
144 - Pprojectile.minus();
145 } else
146 {
147 Qplus = Ptarget.plus()
148 - (sqr(target->GetDefinition()->GetPDGMass()) + pt2 )
149 / (Ptarget.minus() - Qminus );
150 }
151
152 Qmomentum.setPz( (Qplus-Qminus)/2 );
153 Qmomentum.setE( (Qplus+Qminus)/2 );
154
155 // G4cout << "Qplus / Qminus " << Qplus << " / " << Qminus<<G4endl;
156 // G4cout << "pt2 " << pt2 << G4endl;
157 // G4cout << "Qmomentum " << Qmomentum << G4endl;
158 // G4cout << " Masses (P/T) : " << (Pprojectile+Qmomentum).mag() <<
159 // " / " << (Ptarget-Qmomentum).mag() << G4endl;
160
161 } while ( (Ptarget-Qmomentum).mag2() <= Mtarget2
162 || (Pprojectile+Qmomentum).mag2() <= Mprojectile2
163 || (Ptarget-Qmomentum).e() < 0.
164 || (Pprojectile+Qmomentum).e() < 0. );
165
166
167 // G4double Ecms=Pprojectile.e() + Ptarget.e();
168
169 Pprojectile += Qmomentum;
170
171 Ptarget -= Qmomentum;
172
173 // G4cout << "Pprojectile.e() : " << Pprojectile.e() << G4endl;
174 // G4cout << "Ptarget.e() : " << Ptarget.e() << G4endl;
175
176 // G4cout << "end event_______________________________________________"<<G4endl;
177 //
178
179
180 // G4cout << "Pprojectile with Q : " << Pprojectile << G4endl;
181 // G4cout << "Ptarget with Q : " << Ptarget << G4endl;
182 // G4cout << "Projectile back: " << toLab * Pprojectile << G4endl;
183 // G4cout << "Target back: " << toLab * Ptarget << G4endl;
184
185 // Transform back and update SplitableHadron Participant.
186 Pprojectile.transform(toLab);
187 Ptarget.transform(toLab);
188
189 // G4cout << "G4SingleDiffractiveExcitation- Target mass " << Ptarget.mag() << G4endl;
190 // G4cout << "G4SingleDiffractiveExcitation- Projectile mass " << Pprojectile.mag() << G4endl;
191
192 target->Set4Momentum(Ptarget);
193 projectile->Set4Momentum(Pprojectile);
194
195
196 return true;
197}
198
199
200
201
202// --------- private methods ----------------------
203
204G4double G4SingleDiffractiveExcitation::ChooseX(G4double Xmin, G4double Xmax) const
205{
206 // choose an x between Xmin and Xmax with P(x) ~ 1/x
207
208 // to be improved...
209
210 G4double range=Xmax-Xmin;
211
212 if ( Xmin <= 0. || range <=0. )
213 {
214 G4cout << " Xmin, range : " << Xmin << " , " << range << G4endl;
215 throw G4HadronicException(__FILE__, __LINE__, "G4SingleDiffractiveExcitation::ChooseX : Invalid arguments ");
216 }
217
218 G4double x;
219 do {
220 x=Xmin + G4UniformRand() * range;
221 } while ( Xmin/x < G4UniformRand() );
222
223 // cout << "DiffractiveX "<<x<<G4endl;
224 return x;
225}
226
227G4ThreeVector G4SingleDiffractiveExcitation::GaussianPt(G4double widthSquare, G4double maxPtSquare) const
228{ // @@ this method is used in FTFModel as well. Should go somewhere common!
229
230 G4double pt2;
231
232 do {
233 pt2=widthSquare * std::log( G4UniformRand() );
234 } while ( pt2 > maxPtSquare);
235
236 pt2=std::sqrt(pt2);
237
238 G4double phi=G4UniformRand() * twopi;
239
240 return G4ThreeVector (pt2*std::cos(phi), pt2*std::sin(phi), 0.);
241}
242
243
244
245
246
CLHEP::HepLorentzVector G4LorentzVector
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
#define G4UniformRand()
Definition: Randomize.hh:53
double mag2() const
HepLorentzRotation & rotateY(double delta)
HepLorentzRotation & rotateZ(double delta)
HepLorentzRotation inverse() const
double theta() const
Hep3Vector boostVector() const
Hep3Vector vect() const
double minus() const
HepLorentzVector & transform(const HepRotation &)
G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner) const
G4SingleDiffractiveExcitation(G4double sigmaPt=0.6 *CLHEP::GeV, G4double minExtraMass=250 *CLHEP::MeV, G4double x0mass=250 *CLHEP::MeV)
G4ParticleDefinition * GetDefinition() const
void Set4Momentum(const G4LorentzVector &a4Momentum)
const G4LorentzVector & Get4Momentum() const
T sqr(const T &x)
Definition: templates.hh:145