source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/body/src/G4QParton.cc @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 11.1 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: G4QParton.cc,v 1.9 2009/07/17 16:54:57 mkossov Exp $
28// GEANT4 tag $Name: hadr-chips-V09-03-08 $
29//
30// ------------------------------------------------------------
31//      GEANT 4 class implementation file
32//
33//      ---------------- G4QParton ----------------
34//             by Mikhail Kossov, Oct 2006.
35//   class for Parton (inside a string) used by Parton String Models
36//   For comparison mirror member functions are taken from G4 class:
37//   G4Parton
38// ------------------------------------------------------------------------
39// Short description: The Quark-Gluon String consists of the partons, which
40// are quarks and some times gluons.
41// ------------------------------------------------------------------------
42
43//#define debug
44
45#include "G4QParton.hh"
46
47G4QParton::G4QParton() // By default creates only quarks (not di-quarks)
48{
49  // CHIPS is working only with u, d, and s quarks (SU(3)xSU(3)) (no gluons! M.K.)
50  // Random Flavor/Colour/Spin definition for default constructor (with .3 s-suppresion)
51  PGGCode=(G4int)(2.3*G4UniformRand())+1; //@@ Additional parameter of s/u (M.K.)
52  theType=1;
53#ifdef debug
54  G4cout<<"....G4QParton::DefConstructer: PDG = "<<PGGCode<<", Type="<<theType<<G4endl;
55#endif
56  // random colour (1,2,3)=(R,G,B) for quarks and (-1,-2,-3)=(aR,aG,aB) for anti-quarks
57  theColour = (G4int)(3*G4UniformRand())+1;
58  if(theColour>3) theColour = 3;                   // Should never happend
59  theSpinZ = (G4int)(2*G4UniformRand()) - 0.5;
60  QCont = G4QContent(0,0,0,0,0,0); 
61  // Default definition (initialization)
62  theX = 0.;
63  thePosition=G4ThreeVector(0.,0.,0.);
64  theMomentum=G4LorentzVector(0.,0.,0.,0.);
65}
66
67G4QParton::G4QParton(G4int PDG)
68{
69  SetPDGCode(PDG);
70  // Default definition (initialization)
71  theX = 0.;
72  thePosition=G4ThreeVector(0.,0.,0.);
73  theMomentum=G4LorentzVector(0.,0.,0.,0.);
74}
75
76G4QParton::G4QParton(const G4QParton &right)
77{
78  PGGCode     = right.PGGCode;
79  QCont       = right.QCont;
80  theType     = right.theType;
81  theMomentum = right.theMomentum;
82  thePosition = right.thePosition;
83  theX        = right.theX;
84  theColour   = right.theColour;
85  theSpinZ    = right.theSpinZ;
86#ifdef debug
87  G4cout<<"G4QParton::RCopyConstructer: PDG="<<PGGCode<<", Col="<<theColour<<", Sz="
88        <<theSpinZ<<G4endl;
89#endif
90}
91
92G4QParton::G4QParton(const G4QParton* right)
93{
94  PGGCode       = right->PGGCode;
95  QCont         = right->QCont;
96  theType       = right->theType;
97  theMomentum   = right->theMomentum;
98  thePosition   = right->thePosition;
99  theX          = right->theX;
100  theColour     = right->theColour;
101  theSpinZ      = right->theSpinZ;
102#ifdef debug
103  G4cout<<"G4QParton::PCopyConstructer: PDG="<<PGGCode<<", Col="<<theColour<<", Sz="
104        <<theSpinZ<<G4endl;
105#endif
106}
107
108const G4QParton& G4QParton::operator=(const G4QParton &right)
109{
110  if(this != &right)                          // Beware of self assignment
111  {
112    PGGCode     = right.GetPDGCode();
113    QCont       = right.QCont;
114    theType     = right.GetType();
115    theMomentum = right.Get4Momentum();
116    thePosition = right.GetPosition();
117    theX        = right.theX;
118    theColour   = right.theColour;
119    theSpinZ    = right.theSpinZ; 
120#ifdef debug
121    G4cout<<"G4QParton::=Constructer: PDG="<<PGGCode<<", Col="<<theColour<<", Sz="
122          <<theSpinZ<<G4endl;
123#endif
124  }
125  return *this;
126}
127
128G4QParton::~G4QParton() {}
129
130// Redefine the parton nature without changing x, 4Mom, Pos etc.
131void G4QParton::SetPDGCode(G4int PDG)
132{
133  PGGCode=PDG;
134  G4int aPDG=std::abs(PDG);
135  if(aPDG < 3304 && aPDG > 1100 && aPDG%100 < 4) // di-quark
136  {
137    theType=2;
138    G4int cPDG=aPDG/100;
139    if(PDG>0)
140    {
141      if     (cPDG==11) QCont=G4QContent(2,0,0,0,0,0);   // dd
142      else if(cPDG==21) QCont=G4QContent(1,1,0,0,0,0);   // ud
143      else if(cPDG==22) QCont=G4QContent(0,2,0,0,0,0);   // uu
144      else if(cPDG==31) QCont=G4QContent(1,0,1,0,0,0);   // sd
145      else if(cPDG==32) QCont=G4QContent(0,1,1,0,0,0);   // su
146      else if(cPDG==33) QCont=G4QContent(0,0,2,0,0,0);   // ss
147      else
148      {
149        G4cerr<<"***G4QParton::SetPDGCode: bad di-quark PDG="<<PDG<<G4endl;
150        G4Exception("G4QParton::SetPDGCode:","72",FatalException,"Not SU(3) DiQuark");
151      }
152    }
153    else
154    {
155      if     (cPDG==11) QCont=G4QContent(0,0,0,2,0,0);   // anti-dd
156      else if(cPDG==21) QCont=G4QContent(0,0,0,1,1,0);   // anti-ud
157      else if(cPDG==22) QCont=G4QContent(0,0,0,0,2,0);   // anti-uu
158      else if(cPDG==31) QCont=G4QContent(0,0,0,1,0,1);   // anti-sd
159      else if(cPDG==32) QCont=G4QContent(0,0,0,0,1,1);   // anti-su
160      else if(cPDG==33) QCont=G4QContent(0,0,0,0,0,2);   // anti-ss
161      else
162      {
163        G4cerr<<"***G4QParton::SetPDGCode: bad anti-di-quark PDG="<<PDG<<G4endl;
164        G4Exception("G4QParton::SetPDGCode:","72",FatalException,"Not SU(3) AntiDiQuark");
165      }
166    }
167  }
168  else if(aPDG && aPDG<4)                        // quark
169  {
170    theType=1;
171    if(PDG>0)
172    {
173      if     (PDG==1) QCont=G4QContent(1,0,0,0,0,0);   // d
174      else if(PDG==2) QCont=G4QContent(0,1,0,0,0,0);   // u
175      else if(PDG==3) QCont=G4QContent(0,0,1,0,0,0);   // s
176      else
177      {
178        G4cerr<<"***G4QParton::SetPDGCode: bad quark PDG="<<PDG<<G4endl;
179        G4Exception("G4QParton::SetPDGCode:","72",FatalException,"Not SU(3) Quark");
180      }
181    }
182    else
183    {
184      if     (PDG==-1) QCont=G4QContent(0,0,0,1,0,0);  // anti-d
185      else if(PDG==-2) QCont=G4QContent(0,0,0,0,1,0);  // anti-u
186      else if(PDG==-3) QCont=G4QContent(0,0,0,0,0,1);  // anti-s
187      else
188      {
189        G4cerr<<"***G4QParton::SetPDGCode: bad anti-quark PDG="<<PDG<<G4endl;
190        G4Exception("G4QParton::SetPDGCode:","72",FatalException,"Not SU(3) Anti-Quark");
191      }
192    }
193  }
194  else if(aPDG==9 || aPDG==21)                   // gluon
195  {
196    theType=0;
197    QCont=G4QContent(0,0,0,0,0,0);
198  }
199  else
200  {
201    G4cerr<<"***G4QParton::SetPDGCode: wrong gluon/quark/diquark PDG="<<PDG<<G4endl;
202    G4Exception("G4QParton::SetPDGCode:","72",FatalException,"WrongPartonPDG");
203  }
204#ifdef debug
205  G4cout<<"....G4QParton::SetPDGCode: PDG = "<<PDG<<", Type="<<theType<<G4endl;
206#endif
207  //
208  // colour by random in (1,2,3)=(R,G,B) for quarks and
209  //                  in (-1,-2,-3)=(Rbar,Gbar,Bbar) for anti-quarks:
210  G4int RGB=(G4int)(3*G4UniformRand())+1;
211  if(theType==1)
212  {
213    if(PDG>0) theColour = RGB;
214    else      theColour =-RGB;
215  }
216  // colour by random in (-1,-2,-3)=(Rbar,Gbar,Bbar)=(GB,RB,RG) for di-quarks and
217  //                  in (1,2,3)=(R,G,B)=(GB,RB,RG) for anti-di-quarks:
218  else if(theType==2)
219  {
220    if(PDG>0) theColour =-RGB;
221    else      theColour = RGB;
222  }
223  // ColourByRandom (-11,-12,-13,-21,...,-33)=(RRbar,RGbar,RBbar,...,BBbar) for gluons
224  else theColour = -(RGB*10 + (G4int)(3*G4UniformRand())+1);
225  //
226  // spin-z choosen at random from PDG-encoded spin:
227  //
228  G4double            dPDGSpin=1.;        // Quark 2S
229  if     (theType==0) dPDGSpin=2.;        // Gluon 2S
230  else if(theType==2) dPDGSpin=aPDG%10-1; // Di-quark 2S
231  theSpinZ = (G4int)((dPDGSpin+1)*G4UniformRand())-dPDGSpin/2;
232}
233
234// QGS x+/x- logic of the Energy and Pz calculation
235void G4QParton::DefineMomentumInZ(G4double aLightConeMomentum, G4bool aDirection)
236{
237  G4LorentzVector a4Momentum = Get4Momentum();
238  aLightConeMomentum*=theX;
239  G4double TransverseMass2 = sqr(a4Momentum.px()) + sqr(a4Momentum.py());
240  a4Momentum.setPz(0.5*(aLightConeMomentum - TransverseMass2/aLightConeMomentum) *
241                                                                      (aDirection? 1: -1));
242  a4Momentum.setE( 0.5*(aLightConeMomentum + TransverseMass2/aLightConeMomentum));
243  Set4Momentum(a4Momentum);
244} 
245
246// Reduce DiQ-aDiQ to Q-aQ (true if succeeded). General function of the QPartons operations
247G4bool G4QParton::ReduceDiQADiQ(G4QParton* d1, G4QParton* d2)
248{
249  G4bool result=false;
250  G4int sPDG=d1->GetPDGCode();
251  G4int nPDG=d2->GetPDGCode();
252#ifdef debug
253  G4cout<<"G4QParton::ReduceDiQADiQ: **Called** LPDG="<<sPDG<<", RPDG="<<nPDG<<G4endl;
254#endif
255  G4int        qPDG=sPDG;
256  if(qPDG<-99) qPDG=(-qPDG)/100;
257  else         qPDG/=100;
258  G4int        dPDG=nPDG;
259  if(dPDG<-99) dPDG=(-dPDG)/100;
260  else         dPDG/=100;
261  G4int L1=qPDG/10;
262  G4int L2=qPDG%10;
263  G4int R1=dPDG/10;
264  G4int R2=dPDG%10;
265  if(L1==R1 || L1==R2 || L2==R1 || L2==R2) // Annihilation condition
266  {
267    if     (L1==R1)
268    {
269      if(sPDG>0) sPDG=L2;
270      else       sPDG=-L2;
271      if(nPDG>0) nPDG=R2;
272      else       nPDG=-R2;
273#ifdef debug
274      G4cout<<"G4QParton::ReDiQADiQ:L2="<<L2<<",R2="<<R2<<",L="<<sPDG<<",R="<<nPDG<<G4endl;
275#endif
276    }
277    else if(L1==R2)
278    {
279      if(sPDG>0) sPDG=L2;
280      else       sPDG=-L2;
281      if(nPDG>0) nPDG=R1;
282      else       nPDG=-R1;
283#ifdef debug
284      G4cout<<"G4QParton::ReDiQADiQ:L2="<<L2<<",R1="<<R1<<",L="<<sPDG<<",R="<<nPDG<<G4endl;
285#endif
286    }
287    else if(L2==R1)
288    {
289      if(sPDG>0) sPDG=L1;
290      else       sPDG=-L1;
291      if(nPDG>0) nPDG=R2;
292      else       nPDG=-R2;
293#ifdef debug
294      G4cout<<"G4QParton::ReDiQADiQ:L1="<<L1<<",R2="<<R2<<",L="<<sPDG<<",R="<<nPDG<<G4endl;
295#endif
296    }
297    else //(L2==R2)
298    {
299      if(sPDG>0) sPDG=L1;
300      else       sPDG=-L1;
301      if(nPDG>0) nPDG=R1;
302      else       nPDG=-R1;
303#ifdef debug
304      G4cout<<"G4QParton::ReDiQADiQ:L1="<<L1<<",R1="<<R1<<",L="<<sPDG<<",R="<<nPDG<<G4endl;
305#endif
306    }
307    d1->SetPDGCode(sPDG);             // Reset the left quark
308    d2->SetPDGCode(nPDG);            // Reset the right quark
309    result=true;
310#ifdef debug
311    G4cout<<"G4QParton::ReduceDiQADiQ:AfterReduction,L="<<sPDG<<",R="<<nPDG<<G4endl;
312#endif
313  }
314#ifdef debug
315  else G4cout<<"-Warning-G4QParton::ReduceDiQADiQ:DQ-aDQ reduction to Q-aQ Failed"<<G4endl;
316#endif
317  return result;
318}
Note: See TracBrowser for help on using the repository browser.