source: trunk/source/processes/hadronic/models/parton_string/qgsm/src/G4SingleDiffractiveExcitation.cc@ 1199

Last change on this file since 1199 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 9.4 KB
Line 
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,v 1.1 2007/05/25 07:30:47 gunter Exp $
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
38
39#include "globals.hh"
40#include "Randomize.hh"
41
42#include "G4SingleDiffractiveExcitation.hh"
43#include "G4LorentzRotation.hh"
44#include "G4ThreeVector.hh"
45#include "G4ParticleDefinition.hh"
46#include "G4VSplitableHadron.hh"
47#include "G4ExcitedString.hh"
48//#include "G4ios.hh"
49
50G4SingleDiffractiveExcitation::G4SingleDiffractiveExcitation(G4double sigmaPt, G4double minextraMass,G4double x0mass)
51:
52widthOfPtSquare(-2*sqr(sigmaPt)) , minExtraMass(minextraMass),
53minmass(x0mass)
54{
55}
56
57G4bool G4SingleDiffractiveExcitation::
58 ExciteParticipants(G4VSplitableHadron *projectile, G4VSplitableHadron *target) const
59{
60
61 G4LorentzVector Pprojectile=projectile->Get4Momentum();
62 G4double Mprojectile2=sqr(projectile->GetDefinition()->GetPDGMass() + minExtraMass);
63
64 G4LorentzVector Ptarget=target->Get4Momentum();
65 G4double Mtarget2=sqr(target->GetDefinition()->GetPDGMass() + minExtraMass);
66// G4cout << "E proj, target :" << Pprojectile.e() << ", " <<
67// 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
104
105 G4LorentzRotation toLab(toCms.inverse());
106
107 Pprojectile.transform(toCms);
108 Ptarget.transform(toCms);
109
110 G4LorentzVector Qmomentum;
111 G4int whilecount=0;
112 do {
113// Generate pt
114
115 G4double maxPtSquare=sqr(Ptarget.pz());
116 if (whilecount++ >= 500 && (whilecount%100)==0)
117// G4cout << "G4SingleDiffractiveExcitation::ExciteParticipants possibly looping"
118// << ", loop count/ maxPtSquare : "
119// << whilecount << " / " << maxPtSquare << G4endl;
120 if (whilecount > 1000 )
121 {
122 Qmomentum=G4LorentzVector(0.,0.,0.,0.);
123// G4cout << "G4SingleDiffractiveExcitation::ExciteParticipants: Aborting loop!" << G4endl;
124 return false; // Ignore this interaction
125 }
126 Qmomentum=G4LorentzVector(GaussianPt(widthOfPtSquare,maxPtSquare),0);
127
128
129// Momentum transfer
130 G4double Xmin = minmass / ( Pprojectile.e() + Ptarget.e() );
131 G4double Xmax=1.;
132 G4double Xplus =ChooseX(Xmin,Xmax);
133 G4double Xminus=ChooseX(Xmin,Xmax);
134
135 G4double pt2=G4ThreeVector(Qmomentum.vect()).mag2();
136 G4double Qplus =-1 * pt2 / Xminus/Ptarget.minus();
137 G4double Qminus= pt2 / Xplus /Pprojectile.plus();
138
139 if ( KeepProjectile )
140 {
141 Qminus = (sqr(projectile->GetDefinition()->GetPDGMass()) + pt2 )
142 / (Pprojectile.plus() + Qplus )
143 - Pprojectile.minus();
144 } else
145 {
146 Qplus = Ptarget.plus()
147 - (sqr(target->GetDefinition()->GetPDGMass()) + pt2 )
148 / (Ptarget.minus() - Qminus );
149 }
150
151 Qmomentum.setPz( (Qplus-Qminus)/2 );
152 Qmomentum.setE( (Qplus+Qminus)/2 );
153
154// G4cout << "Qplus / Qminus " << Qplus << " / " << Qminus<<G4endl;
155// G4cout << "pt2 " << pt2 << G4endl;
156// G4cout << "Qmomentum " << Qmomentum << G4endl;
157// G4cout << " Masses (P/T) : " << (Pprojectile+Qmomentum).mag() <<
158// " / " << (Ptarget-Qmomentum).mag() << G4endl;
159
160 } while ( (Ptarget-Qmomentum).mag2() <= Mtarget2
161 || (Pprojectile+Qmomentum).mag2() <= Mprojectile2
162 || (Ptarget-Qmomentum).e() < 0.
163 || (Pprojectile+Qmomentum).e() < 0. );
164
165
166// G4double Ecms=Pprojectile.e() + Ptarget.e();
167
168 Pprojectile += Qmomentum;
169
170 Ptarget -= Qmomentum;
171
172// G4cout << "Pprojectile.e() : " << Pprojectile.e() << G4endl;
173// G4cout << "Ptarget.e() : " << Ptarget.e() << G4endl;
174
175// G4cout << "end event_______________________________________________"<<G4endl;
176//
177
178
179// G4cout << "Pprojectile with Q : " << Pprojectile << G4endl;
180// G4cout << "Ptarget with Q : " << Ptarget << G4endl;
181// G4cout << "Projectile back: " << toLab * Pprojectile << G4endl;
182// G4cout << "Target back: " << toLab * Ptarget << G4endl;
183
184// Transform back and update SplitableHadron Participant.
185 Pprojectile.transform(toLab);
186 Ptarget.transform(toLab);
187
188// G4cout << "G4SingleDiffractiveExcitation- Target mass " << Ptarget.mag() << G4endl;
189// G4cout << "G4SingleDiffractiveExcitation- Projectile mass " << Pprojectile.mag() << G4endl;
190
191 target->Set4Momentum(Ptarget);
192 projectile->Set4Momentum(Pprojectile);
193
194
195 return true;
196}
197
198
199
200
201// --------- private methods ----------------------
202
203G4double G4SingleDiffractiveExcitation::ChooseX(G4double Xmin, G4double Xmax) const
204{
205// choose an x between Xmin and Xmax with P(x) ~ 1/x
206
207// to be improved...
208
209 G4double range=Xmax-Xmin;
210
211 if ( Xmin <= 0. || range <=0. )
212 {
213 G4cout << " Xmin, range : " << Xmin << " , " << range << G4endl;
214 throw G4HadronicException(__FILE__, __LINE__, "G4SingleDiffractiveExcitation::ChooseX : Invalid arguments ");
215 }
216
217 G4double x;
218 do {
219 x=Xmin + G4UniformRand() * range;
220 } while ( Xmin/x < G4UniformRand() );
221
222// cout << "DiffractiveX "<<x<<G4endl;
223 return x;
224}
225
226G4ThreeVector G4SingleDiffractiveExcitation::GaussianPt(G4double widthSquare, G4double maxPtSquare) const
227{ // @@ this method is used in FTFModel as well. Should go somewhere common!
228
229 G4double pt2;
230
231 do {
232 pt2=widthSquare * std::log( G4UniformRand() );
233 } while ( pt2 > maxPtSquare);
234
235 pt2=std::sqrt(pt2);
236
237 G4double phi=G4UniformRand() * twopi;
238
239 return G4ThreeVector (pt2*std::cos(phi), pt2*std::sin(phi), 0.);
240}
241
242G4SingleDiffractiveExcitation::G4SingleDiffractiveExcitation(const G4SingleDiffractiveExcitation &)
243: G4QGSDiffractiveExcitation(),
244widthOfPtSquare(0) , minExtraMass(0),
245minmass(0)
246{
247 throw G4HadronicException(__FILE__, __LINE__, "G4SingleDiffractiveExcitation copy contructor not meant to be called");
248}
249
250
251G4SingleDiffractiveExcitation::~G4SingleDiffractiveExcitation()
252{
253}
254
255
256const G4SingleDiffractiveExcitation & G4SingleDiffractiveExcitation::operator=(const G4SingleDiffractiveExcitation &)
257{
258 throw G4HadronicException(__FILE__, __LINE__, "G4SingleDiffractiveExcitation = operator meant to be called");
259 return *this;
260}
261
262
263int G4SingleDiffractiveExcitation::operator==(const G4SingleDiffractiveExcitation &) const
264{
265 throw G4HadronicException(__FILE__, __LINE__, "G4SingleDiffractiveExcitation == operator meant to be called");
266 return false;
267}
268
269int G4SingleDiffractiveExcitation::operator!=(const G4SingleDiffractiveExcitation &) const
270{
271 throw G4HadronicException(__FILE__, __LINE__, "G4SingleDiffractiveExcitation != operator meant to be called");
272 return true;
273}
274
275
276
277
278
Note: See TracBrowser for help on using the repository browser.