[819] | 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 | // |
---|
[1055] | 27 | // $Id: G4QString.cc,v 1.5 2009/02/23 09:49:24 mkossov Exp $ |
---|
| 28 | // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ |
---|
[819] | 29 | // |
---|
| 30 | // ------------------------------------------------------------ |
---|
| 31 | // GEANT 4 class implementation file |
---|
| 32 | // |
---|
| 33 | // ---------------- G4QString ---------------- |
---|
| 34 | // by Mikhail Kossov Oct, 2006 |
---|
| 35 | // class for excited string used by Parton String Models |
---|
| 36 | // For comparison mirror member functions are taken from G4 classes: |
---|
| 37 | // G4FragmentingString |
---|
| 38 | // G4ExcitedStringDecay |
---|
[1055] | 39 | // --------------------------------------------------------------------- |
---|
| 40 | // Short description: If partons from the G4QPartonPair are close in |
---|
| 41 | // rapidity, they create Quasmons, but if they are far in the rapidity |
---|
| 42 | // space, they can not interact directly. Say the bottom parton (quark) |
---|
| 43 | // has rapidity 0, and the top parton (antiquark) has rapidity 8, then |
---|
| 44 | // the top quark splits in two by radiating gluon, and each part has |
---|
| 45 | // rapidity 4, then the gluon splits in quark-antiquark pair (rapidity |
---|
| 46 | // 2 each), and then the quark gadiates anothe gluon and reachs rapidity |
---|
| 47 | // 1. Now it can interact with the bottom antiquark, creating a Quasmon |
---|
| 48 | // or a hadron. The intermediate partons is the string ladder. |
---|
| 49 | // --------------------------------------------------------------------- |
---|
[819] | 50 | |
---|
| 51 | //#define debug |
---|
| 52 | |
---|
| 53 | #include "G4QString.hh" |
---|
| 54 | #include <algorithm> |
---|
| 55 | // Static parameters definition |
---|
| 56 | G4double G4QString::MassCut=350.*MeV; // minimum mass cut for the string |
---|
| 57 | G4double G4QString::ClusterMass=150.*MeV; // minimum cluster mass for fragmentation |
---|
| 58 | G4double G4QString::SigmaQT=0.5*GeV; // quarkTransverseMomentum distribution parameter |
---|
| 59 | G4double G4QString::DiquarkSuppress=0.1; // is Diquark suppression parameter |
---|
| 60 | G4double G4QString::DiquarkBreakProb=0.1; // is Diquark breaking probability |
---|
| 61 | G4double G4QString::SmoothParam=0.9; // QGS model parameter |
---|
| 62 | G4double G4QString::StrangeSuppress=0.435;// Strangeness suppression (u:d:s=1:1:0.3 ?M.K.) |
---|
| 63 | G4double G4QString::widthOfPtSquare=-0.72*GeV*GeV; // pt -width2 forStringExcitation |
---|
| 64 | G4int G4QString::StringLoopInterrupt=999; // String fragmentation LOOP limit |
---|
| 65 | G4int G4QString::ClusterLoopInterrupt=500;// Cluster fragmentation LOOP limit |
---|
| 66 | |
---|
| 67 | G4QString::G4QString() : theDirection(0),thePosition(G4ThreeVector(0.,0.,0.)), |
---|
| 68 | hadronizer(new G4QHadronBuilder){} |
---|
| 69 | |
---|
| 70 | G4QString::G4QString(G4QParton* Color, G4QParton* AntiColor, G4int Direction) : |
---|
| 71 | hadronizer(new G4QHadronBuilder) |
---|
| 72 | { |
---|
| 73 | ExciteString(Color, AntiColor, Direction); |
---|
| 74 | } |
---|
| 75 | |
---|
| 76 | G4QString::G4QString(G4QParton* QCol,G4QParton* Gluon,G4QParton* QAntiCol,G4int Direction): |
---|
[1055] | 77 | theDirection(Direction),thePosition(QCol->GetPosition()),hadronizer(new G4QHadronBuilder) |
---|
[819] | 78 | { |
---|
| 79 | thePartons.push_back(QCol); |
---|
| 80 | thePartons.push_back(Gluon); |
---|
| 81 | thePartons.push_back(QAntiCol); |
---|
| 82 | } |
---|
| 83 | |
---|
| 84 | G4QString::G4QString(const G4QString &right) : theDirection(right.GetDirection()), |
---|
| 85 | thePosition(right.GetPosition()), hadronizer(new G4QHadronBuilder){} |
---|
| 86 | |
---|
| 87 | G4QString::~G4QString() |
---|
| 88 | {if(thePartons.size())std::for_each(thePartons.begin(),thePartons.end(),DeleteQParton());} |
---|
| 89 | |
---|
| 90 | G4LorentzVector G4QString::Get4Momentum() const |
---|
| 91 | { |
---|
[1055] | 92 | G4LorentzVector momentum(0.,0.,0.,0.); |
---|
| 93 | for(unsigned i=0; i<thePartons.size(); i++) momentum += thePartons[i]->Get4Momentum(); |
---|
| 94 | return momentum; |
---|
[819] | 95 | } |
---|
| 96 | |
---|
| 97 | void G4QString::LorentzRotate(const G4LorentzRotation & rotation) |
---|
| 98 | { |
---|
[1055] | 99 | for(unsigned i=0; i<thePartons.size(); i++) |
---|
| 100 | thePartons[i]->Set4Momentum(rotation*thePartons[i]->Get4Momentum()); |
---|
[819] | 101 | } |
---|
| 102 | |
---|
| 103 | void G4QString::InsertParton(G4QParton* aParton, const G4QParton* addafter) |
---|
| 104 | { |
---|
| 105 | G4QPartonVector::iterator insert_index; // Begin by default (? M.K.) |
---|
[1055] | 106 | |
---|
| 107 | if(addafter != NULL) |
---|
| 108 | { |
---|
| 109 | insert_index=std::find(thePartons.begin(), thePartons.end(), addafter); |
---|
| 110 | if (insert_index == thePartons.end()) // No object addafter in thePartons |
---|
| 111 | { |
---|
| 112 | G4cerr<<"***G4QString::InsertParton: Address Parton is not found"<<G4endl; |
---|
[819] | 113 | G4Exception("G4QString::InsertParton:","72",FatalException,"NoAddressParton"); |
---|
[1055] | 114 | } |
---|
| 115 | } |
---|
| 116 | thePartons.insert(insert_index+1, aParton); |
---|
[819] | 117 | } |
---|
| 118 | |
---|
| 119 | G4LorentzRotation G4QString::TransformToCenterOfMass() |
---|
| 120 | { |
---|
[1055] | 121 | G4LorentzRotation toCms(-1*Get4Momentum().boostVector()); |
---|
| 122 | for(unsigned i=0; i<thePartons.size(); i++) |
---|
| 123 | thePartons[i]->Set4Momentum(toCms*thePartons[i]->Get4Momentum()); |
---|
| 124 | return toCms; |
---|
[819] | 125 | } |
---|
| 126 | |
---|
| 127 | G4LorentzRotation G4QString::TransformToAlignedCms() |
---|
| 128 | { |
---|
[1055] | 129 | G4LorentzVector momentum=Get4Momentum(); |
---|
| 130 | G4LorentzRotation toAlignedCms(-1*momentum.boostVector()); |
---|
| 131 | momentum= toAlignedCms* thePartons[0]->Get4Momentum(); |
---|
| 132 | toAlignedCms.rotateZ(-1*momentum.phi()); |
---|
| 133 | toAlignedCms.rotateY(-1*momentum.theta()); |
---|
| 134 | for(unsigned index=0; index<thePartons.size(); index++) |
---|
| 135 | { |
---|
| 136 | momentum=toAlignedCms * thePartons[index]->Get4Momentum(); |
---|
| 137 | thePartons[index]->Set4Momentum(momentum); |
---|
| 138 | } |
---|
| 139 | return toAlignedCms; |
---|
[819] | 140 | } |
---|
| 141 | |
---|
| 142 | void G4QString::Boost(G4ThreeVector& Velocity) |
---|
| 143 | { |
---|
| 144 | for(unsigned cParton=0; cParton<thePartons.size(); cParton++) |
---|
| 145 | { |
---|
| 146 | G4LorentzVector Mom = thePartons[cParton]->Get4Momentum(); |
---|
| 147 | Mom.boost(Velocity); |
---|
| 148 | thePartons[cParton]->Set4Momentum(Mom); |
---|
| 149 | } |
---|
| 150 | } |
---|
| 151 | |
---|
| 152 | G4QParton* G4QString::GetColorParton(void) const |
---|
| 153 | { |
---|
| 154 | G4QParton* start = *(thePartons.begin()); |
---|
| 155 | G4QParton* end = *(thePartons.end()-1); |
---|
| 156 | G4int Encoding = start->GetPDGCode(); |
---|
| 157 | if (Encoding<-1000 || (Encoding < 1000 && Encoding > 0)) return start; |
---|
| 158 | return end; |
---|
| 159 | } |
---|
| 160 | |
---|
| 161 | G4QParton* G4QString::GetAntiColorParton(void) const |
---|
| 162 | { |
---|
| 163 | G4QParton* start = *(thePartons.begin()); |
---|
| 164 | G4QParton* end = *(thePartons.end()-1); |
---|
| 165 | G4int Encoding = start->GetPDGCode(); |
---|
| 166 | if (Encoding < -1000 || (Encoding < 1000 && Encoding > 0)) return end; |
---|
| 167 | return start; |
---|
| 168 | } |
---|
| 169 | |
---|
| 170 | // Fill parameters |
---|
| 171 | void G4QString::SetParameters(G4double mCut,G4double clustM, G4double sigQT,G4double DQSup, |
---|
| 172 | G4double DQBU, G4double smPar, G4double SSup, G4double SigPt, G4int SLmax, G4int CLmax) |
---|
| 173 | {// ============================================================================= |
---|
| 174 | MassCut = mCut; // minimum mass cut for the string |
---|
| 175 | ClusterMass = clustM; // minimum cluster mass for the fragmentation |
---|
| 176 | SigmaQT = sigQT; // quark transverse momentum distribution parameter |
---|
| 177 | DiquarkSuppress = DQSup; // is Diquark suppression parameter |
---|
| 178 | DiquarkBreakProb= DQBU; // is Diquark breaking probability |
---|
| 179 | SmoothParam = smPar; // QGS model parameter |
---|
| 180 | StrangeSuppress = SSup; // Strangeness suppression parameter |
---|
[1055] | 181 | widthOfPtSquare = -2*SigPt*SigPt; // width^2 of pt for string excitation |
---|
[819] | 182 | StringLoopInterrupt = SLmax; // String fragmentation LOOP limit |
---|
| 183 | ClusterLoopInterrupt= CLmax; // Cluster fragmentation LOOP limit |
---|
| 184 | } |
---|
| 185 | |
---|
| 186 | // Pt distribution @@ one can use 1/(1+A*Pt^2)^B |
---|
| 187 | G4ThreeVector G4QString::GaussianPt(G4double widthSquare, G4double maxPtSquare) const |
---|
| 188 | { |
---|
[1055] | 189 | G4double pt2; do{pt2=widthSquare*std::log(G4UniformRand());} while (pt2>maxPtSquare); |
---|
| 190 | pt2=std::sqrt(pt2); |
---|
| 191 | G4double phi=G4UniformRand()*twopi; |
---|
| 192 | return G4ThreeVector(pt2*std::cos(phi),pt2*std::sin(phi),0.); |
---|
[819] | 193 | } // End of GaussianPt |
---|
| 194 | |
---|
| 195 | // Diffractively excite the string |
---|
| 196 | void G4QString::DiffString(G4QHadron* hadron, G4bool isProjectile) |
---|
| 197 | { |
---|
[1055] | 198 | hadron->SplitUp(); |
---|
| 199 | G4QParton* start = hadron->GetNextParton(); |
---|
| 200 | if( start==NULL) |
---|
| 201 | { |
---|
[819] | 202 | G4cerr<<"***G4QString::DiffString: No start parton found, nothing is done"<<G4endl; |
---|
[1055] | 203 | return; |
---|
| 204 | } |
---|
| 205 | G4QParton* end = hadron->GetNextParton(); |
---|
| 206 | if( end==NULL) |
---|
| 207 | { |
---|
[819] | 208 | G4cerr<<"***G4QString::DiffString: No end parton found, nothing is done"<<G4endl; |
---|
[1055] | 209 | return; |
---|
| 210 | } |
---|
| 211 | if(isProjectile) ExciteString(end, start, 1); // 1 = Projectile |
---|
| 212 | else ExciteString(start, end,-1); // -1 = Target |
---|
| 213 | SetPosition(hadron->GetPosition()); |
---|
[819] | 214 | // momenta of string ends |
---|
[1055] | 215 | G4double ptSquared= hadron->Get4Momentum().perp2(); |
---|
[819] | 216 | G4double hmins=hadron->Get4Momentum().minus(); |
---|
| 217 | G4double hplus=hadron->Get4Momentum().plus(); |
---|
[1055] | 218 | G4double transMassSquared=hplus*hmins; |
---|
| 219 | G4double maxMomentum = std::sqrt(transMassSquared) - std::sqrt(ptSquared); |
---|
| 220 | G4double maxAvailMomentumSquared = maxMomentum * maxMomentum; |
---|
| 221 | G4ThreeVector pt=GaussianPt(widthOfPtSquare,maxAvailMomentumSquared); |
---|
| 222 | G4LorentzVector Pstart(G4LorentzVector(pt,0.)); |
---|
| 223 | G4LorentzVector Pend(hadron->Get4Momentum().px(), hadron->Get4Momentum().py(), 0.); |
---|
| 224 | Pend-=Pstart; |
---|
| 225 | G4double tm1=hmins+(Pend.perp2()-Pstart.perp2())/hplus; |
---|
| 226 | G4double tm2=std::sqrt( std::max(0., tm1*tm1-4*Pend.perp2()*hmins/hplus ) ); |
---|
| 227 | G4int Sign= isProjectile ? TARGET : PROJECTILE; |
---|
| 228 | G4double endMinus = 0.5 * (tm1 + Sign*tm2); |
---|
| 229 | G4double startMinus= hmins - endMinus; |
---|
| 230 | G4double startPlus = Pstart.perp2() / startMinus; |
---|
| 231 | G4double endPlus = hplus - startPlus; |
---|
| 232 | Pstart.setPz(0.5*(startPlus - startMinus)); |
---|
| 233 | Pstart.setE (0.5*(startPlus + startMinus)); |
---|
| 234 | Pend.setPz (0.5*(endPlus - endMinus)); |
---|
| 235 | Pend.setE (0.5*(endPlus + endMinus)); |
---|
| 236 | start->Set4Momentum(Pstart); |
---|
| 237 | end->Set4Momentum(Pend); |
---|
[819] | 238 | #ifdef debug |
---|
[1055] | 239 | G4cout<<"G4QString::DiffString: StartQ="<<start->GetPDGCode()<<start->Get4Momentum()<<"(" |
---|
[819] | 240 | <<start->Get4Momentum().mag()<<"), EndQ="<<end->GetPDGCode()<<end ->Get4Momentum() |
---|
| 241 | <<"("<<end->Get4Momentum().mag()<<"), sumOfEnds="<<Pstart+Pend<<", H4M(original)=" |
---|
| 242 | <<hadron->Get4Momentum()<<G4endl; |
---|
| 243 | #endif |
---|
| 244 | } // End of DiffString (The string is excited diffractively) |
---|
| 245 | |
---|
| 246 | // Excite the string by two partons |
---|
| 247 | void G4QString::ExciteString(G4QParton* Color, G4QParton* AntiColor, G4int Direction) |
---|
| 248 | { |
---|
[1055] | 249 | theDirection = Direction; |
---|
[819] | 250 | thePosition = Color->GetPosition(); |
---|
| 251 | thePartons.push_back(Color); |
---|
| 252 | thePartons.push_back(AntiColor); |
---|
| 253 | } // End of ExciteString |
---|
| 254 | |
---|
| 255 | // LUND Longitudinal fragmentation |
---|
| 256 | G4double G4QString::GetLundLightConeZ(G4double zmin, G4double zmax, G4int , // @@ ? M.K. |
---|
| 257 | G4QHadron* pHadron, G4double Px, G4double Py) |
---|
| 258 | { |
---|
| 259 | static const G4double alund = 0.7/GeV/GeV; |
---|
| 260 | // If blund get restored, you MUST adapt the calculation of zOfMaxyf. |
---|
| 261 | //static const G4double blund = 1; |
---|
| 262 | G4double z, yf; |
---|
| 263 | G4double Mt2 = Px*Px + Py*Py + pHadron->GetMass2(); |
---|
| 264 | G4double zOfMaxyf=alund*Mt2/(alund*Mt2+1.); |
---|
| 265 | G4double maxYf=(1.-zOfMaxyf)/zOfMaxyf * std::exp(-alund*Mt2/zOfMaxyf); |
---|
| 266 | do |
---|
| 267 | { |
---|
| 268 | z = zmin + G4UniformRand()*(zmax-zmin); |
---|
| 269 | // yf = std::pow(1. - z, blund)/z*std::exp(-alund*Mt2/z); |
---|
[1055] | 270 | yf = (1-z)/z * std::exp(-alund*Mt2/z); |
---|
[819] | 271 | } while (G4UniformRand()*maxYf>yf); |
---|
| 272 | return z; |
---|
| 273 | } // End of GetLundLightConeZ |
---|
| 274 | |
---|
| 275 | |
---|
| 276 | // QGSM Longitudinal fragmentation |
---|
| 277 | G4double G4QString::GetQGSMLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding, |
---|
[1055] | 278 | G4QHadron* , G4double, G4double) // @@ M.K. |
---|
[819] | 279 | { |
---|
| 280 | static const G4double arho=0.5; |
---|
| 281 | static const G4double aphi=0.; |
---|
| 282 | static const G4double an=-0.5; |
---|
| 283 | static const G4double ala=-0.75; |
---|
| 284 | static const G4double aksi=-1.; |
---|
| 285 | static const G4double alft=0.5; |
---|
| 286 | G4double z; |
---|
| 287 | G4double theA(0), d2, yf; |
---|
| 288 | G4int absCode = std::abs(PartonEncoding); |
---|
| 289 | // @@ Crazy algorithm is used for simple power low... |
---|
| 290 | if (absCode < 10) // Quarks (@@ 9 can be a gluon) |
---|
| 291 | { |
---|
| 292 | if(absCode == 1 || absCode == 2) theA = arho; |
---|
| 293 | else if(absCode == 3) theA = aphi; |
---|
| 294 | else |
---|
| 295 | { |
---|
| 296 | G4cerr<<"***G4QString::GetQGSMLightConeZ: CHIPS is SU(3), quakCode="<<absCode<<G4endl; |
---|
| 297 | G4Exception("G4QString::GetQGSMLightConeZ:","72",FatalException,"WrongQuark"); |
---|
| 298 | } |
---|
[1055] | 299 | do |
---|
[819] | 300 | { |
---|
| 301 | z = zmin + G4UniformRand()*(zmax - zmin); |
---|
| 302 | d2 = alft - theA; |
---|
| 303 | yf = std::pow( 1.-z, d2); |
---|
| 304 | } |
---|
| 305 | while (G4UniformRand()>yf); |
---|
| 306 | } |
---|
| 307 | else // Di-quarks (@@ Crazy codes are not checked) |
---|
| 308 | { |
---|
| 309 | if (absCode==3101||absCode==3103||absCode==3201||absCode==3203) d2=alft-ala-ala+arho; |
---|
| 310 | else if(absCode==1103||absCode==2101||absCode==2203||absCode==2103) d2=alft-an-an+arho; |
---|
| 311 | else d2=alft-aksi-aksi+arho; |
---|
| 312 | do |
---|
| 313 | { |
---|
| 314 | z = zmin + G4UniformRand()*(zmax - zmin); |
---|
| 315 | yf = std::pow(1.-z, d2); |
---|
| 316 | } |
---|
| 317 | while (G4UniformRand()>yf); |
---|
| 318 | } |
---|
| 319 | return z; |
---|
| 320 | } // End of GetQGSMLightConeZ |
---|
| 321 | |
---|
| 322 | // Diffractively excite the string |
---|
| 323 | G4QHadronVector* G4QString::FragmentString(G4bool QL) |
---|
| 324 | { |
---|
| 325 | // Can no longer modify Parameters for Fragmentation. |
---|
[1055] | 326 | // PastInitPhase=true; // Now static |
---|
| 327 | |
---|
| 328 | // check if string has enough mass to fragment... |
---|
| 329 | G4QHadronVector* LeftVector=LightFragmentationTest(); |
---|
| 330 | if(LeftVector) return LeftVector; |
---|
| 331 | |
---|
| 332 | LeftVector = new G4QHadronVector; |
---|
| 333 | G4QHadronVector* RightVector = new G4QHadronVector; |
---|
[819] | 334 | |
---|
| 335 | // this should work but it's only a semi deep copy. |
---|
[1055] | 336 | // %GF G4ExcitedString theStringInCMS(theString); |
---|
[819] | 337 | G4QString* theStringInCMS=CPExcited(); // must be deleted |
---|
[1055] | 338 | G4LorentzRotation toCms=theStringInCMS->TransformToAlignedCms(); |
---|
| 339 | G4bool success=false; |
---|
| 340 | G4bool inner_sucess=true; |
---|
| 341 | G4int attempt=0; |
---|
| 342 | while (!success && attempt++<StringLoopInterrupt) //@@ It's Loop with break |
---|
| 343 | { |
---|
| 344 | G4QString* currentString = new G4QString(*theStringInCMS); |
---|
| 345 | if(LeftVector->size()) |
---|
[819] | 346 | { |
---|
| 347 | std::for_each(LeftVector->begin(), LeftVector->end(), DeleteQHadron()); |
---|
[1055] | 348 | LeftVector->clear(); |
---|
[819] | 349 | } |
---|
[1055] | 350 | if(RightVector->size()) |
---|
| 351 | { |
---|
[819] | 352 | std::for_each(RightVector->begin(), RightVector->end(), DeleteQHadron()); |
---|
[1055] | 353 | RightVector->clear(); |
---|
| 354 | } |
---|
| 355 | inner_sucess=true; // set false on failure.. |
---|
| 356 | while (!StopFragmentation()) |
---|
| 357 | { // Split current string into hadron + new string |
---|
| 358 | G4QString* newString=0; // used as output from SplitUp... |
---|
| 359 | G4QHadron* Hadron=Splitup(QL); |
---|
| 360 | if(Hadron && IsFragmentable()) |
---|
| 361 | { |
---|
| 362 | if(currentString->GetDecayDirection()>0) LeftVector->push_back(Hadron); |
---|
| 363 | else RightVector->push_back(Hadron); |
---|
| 364 | delete currentString; |
---|
| 365 | currentString=newString; |
---|
| 366 | } |
---|
[819] | 367 | else |
---|
| 368 | { |
---|
[1055] | 369 | // abandon ... start from the beginning |
---|
| 370 | if (newString) delete newString; |
---|
| 371 | if (Hadron) delete Hadron; |
---|
| 372 | inner_sucess=false; |
---|
| 373 | break; |
---|
| 374 | } |
---|
| 375 | } |
---|
| 376 | // Split current string into 2 final Hadrons |
---|
| 377 | if( inner_sucess && SplitLast(currentString, LeftVector, RightVector) ) success=true; |
---|
| 378 | delete currentString; |
---|
| 379 | } |
---|
| 380 | delete theStringInCMS; |
---|
| 381 | if (!success) |
---|
| 382 | { |
---|
| 383 | if(RightVector->size()) |
---|
| 384 | { |
---|
[819] | 385 | std::for_each(RightVector->begin(), RightVector->end(), DeleteQHadron()); |
---|
[1055] | 386 | RightVector->clear(); |
---|
| 387 | } |
---|
| 388 | delete RightVector; |
---|
| 389 | if(LeftVector->size()) |
---|
[819] | 390 | { |
---|
| 391 | std::for_each(LeftVector->begin(), LeftVector->end(), DeleteQHadron()); |
---|
[1055] | 392 | LeftVector->clear(); |
---|
[819] | 393 | } |
---|
[1055] | 394 | return LeftVector; |
---|
| 395 | } |
---|
| 396 | // Join Left and Right Vectors into LeftVector in correct order. |
---|
| 397 | while(!RightVector->empty()) |
---|
| 398 | { |
---|
| 399 | LeftVector->push_back(RightVector->back()); |
---|
| 400 | RightVector->erase(RightVector->end()-1); |
---|
| 401 | } |
---|
| 402 | delete RightVector; |
---|
| 403 | CalculateHadronTimePosition(Get4Momentum().mag(), LeftVector); |
---|
| 404 | G4LorentzRotation toObserverFrame(toCms.inverse()); |
---|
| 405 | for(unsigned C1 = 0; C1 < LeftVector->size(); C1++) |
---|
| 406 | { |
---|
| 407 | G4QHadron* Hadron = (*LeftVector)[C1]; |
---|
| 408 | G4LorentzVector Momentum = Hadron->Get4Momentum(); |
---|
| 409 | Momentum = toObserverFrame*Momentum; |
---|
| 410 | Hadron->Set4Momentum(Momentum); |
---|
| 411 | G4LorentzVector Coordinate(Hadron->GetPosition(), Hadron->GetFormationTime()); |
---|
| 412 | Momentum = toObserverFrame*Coordinate; |
---|
| 413 | Hadron->SetFormationTime(Momentum.e()); |
---|
| 414 | Hadron->SetPosition(GetPosition()+Momentum.vect()); |
---|
| 415 | } |
---|
| 416 | return LeftVector; |
---|
[819] | 417 | } // End of FragmentLundString |
---|
| 418 | |
---|
| 419 | // Creates a string, using only the end-partons of the string |
---|
| 420 | G4QString* G4QString::CPExcited() |
---|
| 421 | { |
---|
[1055] | 422 | G4QParton* LeftParton = new G4QParton(GetLeftParton()); |
---|
| 423 | G4QParton* RightParton= new G4QParton(GetRightParton()); |
---|
| 424 | return new G4QString(LeftParton,RightParton,GetDirection()); |
---|
[819] | 425 | } // End of CPExcited |
---|
| 426 | |
---|
| 427 | // Simple decay of the string |
---|
| 428 | G4QHadronVector* G4QString::LightFragmentationTest() |
---|
| 429 | { |
---|
| 430 | // Check string decay threshold |
---|
[1055] | 431 | |
---|
| 432 | G4QHadronVector* result=0; // return 0 when string exceeds the mass cut |
---|
| 433 | |
---|
| 434 | G4QHadronPair hadrons((G4QHadron*)0, (G4QHadron*)0); // pair of hadrons |
---|
| 435 | if(sqr(FragmentationMass(0,&hadrons)+MassCut)<Mass2()) return 0; //Par MassCut |
---|
| 436 | |
---|
| 437 | result = new G4QHadronVector; |
---|
[819] | 438 | |
---|
[1055] | 439 | if(hadrons.second==0) // Second hadron exists |
---|
| 440 | { |
---|
| 441 | // Substitute string by a light hadron, Note that Energy is not conserved here! @@ |
---|
| 442 | G4ThreeVector Mom3 = Get4Momentum().vect(); |
---|
| 443 | G4LorentzVector Mom(Mom3, std::sqrt(Mom3.mag2() + hadrons.first->GetMass2()) ); |
---|
[819] | 444 | result->push_back(new G4QHadron(hadrons.first, 0, GetPosition(), Mom)); |
---|
[1055] | 445 | } |
---|
[819] | 446 | else |
---|
[1055] | 447 | { |
---|
| 448 | //... string was qq--qqbar type: Build two stable hadrons, |
---|
| 449 | G4LorentzVector Mom1, Mom2; |
---|
| 450 | Sample4Momentum(&Mom1, hadrons.first->GetMass(), &Mom2, hadrons.second->GetMass(), |
---|
| 451 | Get4Momentum().mag()); |
---|
| 452 | result->push_back(new G4QHadron(hadrons.first, 0, GetPosition(), Mom1)); |
---|
| 453 | result->push_back(new G4QHadron(hadrons.second,0, GetPosition(), Mom2)); |
---|
[819] | 454 | G4ThreeVector Velocity = Get4Momentum().boostVector(); |
---|
| 455 | G4int L=result->size(); if(L) for(G4int i=0; i<L; i++) (*result)[i]->Boost(Velocity); |
---|
[1055] | 456 | } |
---|
| 457 | return result; |
---|
[819] | 458 | } // End of LightFragmentationTest |
---|
| 459 | |
---|
| 460 | // Calculate Fragmentation Mass (if pdefs # 0, returns two hadrons) |
---|
[1055] | 461 | G4double G4QString::FragmentationMass(G4QHcreate build, G4QHadronPair* pdefs) |
---|
[819] | 462 | { |
---|
| 463 | G4double mass; |
---|
| 464 | // Example how to use an interface to different member functions |
---|
[1055] | 465 | if(build==0) build=&G4QHadronBuilder::BuildLowSpin; // @@ Build S Hadrons? |
---|
[819] | 466 | G4QHadron* Hadron1 = 0; // @@ Not initialized |
---|
| 467 | G4QHadron* Hadron2 = 0; |
---|
| 468 | if(!FourQuarkString() ) |
---|
| 469 | { |
---|
| 470 | // spin 0 meson or spin 1/2 barion will be built |
---|
| 471 | Hadron1 = (hadronizer->*build)(GetLeftParton(), GetRightParton()); |
---|
| 472 | mass = Hadron1->GetMass(); |
---|
| 473 | } |
---|
| 474 | else |
---|
| 475 | { |
---|
| 476 | // string is qq--qqbar: Build two stable hadrons, with extra uubar or ddbar quark pair |
---|
[1055] | 477 | G4int iflc = (G4UniformRand() < 0.5)? 1 : 2; |
---|
| 478 | if (GetLeftParton()->GetPDGCode() < 0) iflc = -iflc; |
---|
| 479 | //... theSpin = 4; spin 3/2 baryons will be built |
---|
| 480 | Hadron1 = (hadronizer->*build)(GetLeftParton(),CreateParton(iflc)); |
---|
| 481 | Hadron2 = (hadronizer->*build)(GetRightParton(),CreateParton(-iflc)); |
---|
[819] | 482 | mass = Hadron1->GetMass() + Hadron2->GetMass(); |
---|
| 483 | } |
---|
[1055] | 484 | if(pdefs) // need to return hadrons as well.... |
---|
| 485 | { |
---|
| 486 | pdefs->first = Hadron1; |
---|
| 487 | pdefs->second = Hadron2; |
---|
| 488 | } |
---|
[819] | 489 | return mass; |
---|
| 490 | } // End of FragmentationMass |
---|
| 491 | |
---|
| 492 | // Checks that the string is qq-(qq)bar string |
---|
| 493 | G4bool G4QString::FourQuarkString() const |
---|
| 494 | { |
---|
| 495 | return GetLeftParton()->GetParticleSubType() == "di_quark" |
---|
| 496 | && GetRightParton()->GetParticleSubType()== "di_quark"; |
---|
| 497 | } // End of FourQuarkString |
---|
| 498 | |
---|
| 499 | void G4QString::CalculateHadronTimePosition(G4double StringMass, G4QHadronVector* Hadrons) |
---|
| 500 | { |
---|
| 501 | // `yo-yo` formation time |
---|
| 502 | static const G4double kappa = 1.0 * GeV/fermi; |
---|
| 503 | static const G4double dkappa = kappa+kappa; |
---|
| 504 | for(unsigned c1 = 0; c1 < Hadrons->size(); c1++) |
---|
| 505 | { |
---|
| 506 | G4double SumPz = 0; |
---|
| 507 | G4double SumE = 0; |
---|
| 508 | for(unsigned c2 = 0; c2 < c1; c2++) |
---|
| 509 | { |
---|
| 510 | G4LorentzVector hc2M=(*Hadrons)[c2]->Get4Momentum(); |
---|
| 511 | SumPz += hc2M.pz(); |
---|
| 512 | SumE += hc2M.e(); |
---|
| 513 | } |
---|
| 514 | G4QHadron* hc1=(*Hadrons)[c1]; |
---|
| 515 | G4LorentzVector hc1M=hc1->Get4Momentum(); |
---|
| 516 | G4double HadronE = hc1M.e(); |
---|
| 517 | G4double HadronPz= hc1M.pz(); |
---|
| 518 | hc1->SetFormationTime((StringMass-SumPz-SumPz+HadronE-HadronPz)/dkappa); |
---|
| 519 | hc1->SetPosition(G4ThreeVector(0,0,(StringMass-SumE-SumE-HadronE+HadronPz)/dkappa)); |
---|
| 520 | } |
---|
| 521 | } // End of CalculateHadronTimePosition |
---|
| 522 | |
---|
| 523 | void G4QString::SetLeftPartonStable() |
---|
| 524 | { |
---|
| 525 | theStableParton=GetLeftParton(); |
---|
| 526 | theDecayParton=GetRightParton(); |
---|
| 527 | decaying=Right; |
---|
| 528 | } |
---|
| 529 | |
---|
| 530 | void G4QString::SetRightPartonStable() |
---|
| 531 | { |
---|
| 532 | theStableParton=GetRightParton(); |
---|
| 533 | theDecayParton=GetLeftParton(); |
---|
| 534 | decaying=Left; |
---|
| 535 | } |
---|
| 536 | |
---|
| 537 | G4int G4QString::GetDecayDirection() const |
---|
| 538 | { |
---|
[1055] | 539 | if (decaying == Left ) return +1; |
---|
| 540 | else if (decaying == Right) return -1; |
---|
| 541 | else |
---|
| 542 | { |
---|
| 543 | G4cerr<<"***G4QString::GetDecayDirection: wrong DecayDirection="<<decaying<<G4endl; |
---|
[819] | 544 | G4Exception("G4QString::GetDecayDirection:","72",FatalException,"WrongDecayDirection"); |
---|
[1055] | 545 | } |
---|
| 546 | return 0; |
---|
[819] | 547 | } |
---|
| 548 | |
---|
| 549 | G4ThreeVector G4QString::StablePt() |
---|
| 550 | { |
---|
[1055] | 551 | if (decaying == Left ) return Ptright; |
---|
| 552 | else if (decaying == Right ) return Ptleft; |
---|
| 553 | else |
---|
| 554 | { |
---|
| 555 | G4cerr<<"***G4QString::StablePt: wrong DecayDirection="<<decaying<<G4endl; |
---|
[819] | 556 | G4Exception("G4QString::StablePt:","72",FatalException,"WrongDecayDirection"); |
---|
[1055] | 557 | } |
---|
| 558 | return G4ThreeVector(); |
---|
[819] | 559 | } |
---|
| 560 | |
---|
| 561 | G4ThreeVector G4QString::DecayPt() |
---|
| 562 | { |
---|
[1055] | 563 | if (decaying == Left ) return Ptleft; |
---|
| 564 | else if (decaying == Right ) return Ptright; |
---|
| 565 | else |
---|
| 566 | { |
---|
| 567 | G4cerr<<"***G4QString::DecayPt: wrong DecayDirection="<<decaying<<G4endl; |
---|
[819] | 568 | G4Exception("G4QString::DecayPt:","72",FatalException,"WrongDecayDirection"); |
---|
[1055] | 569 | } |
---|
| 570 | return G4ThreeVector(); |
---|
[819] | 571 | } |
---|
| 572 | |
---|
| 573 | G4double G4QString::LightConeDecay() |
---|
| 574 | { |
---|
[1055] | 575 | if (decaying == Left ) return Pplus; |
---|
| 576 | else if (decaying == Right ) return Pminus; |
---|
| 577 | else |
---|
| 578 | { |
---|
| 579 | G4cerr<<"***G4QString::LightConeDecay: wrong DecayDirection="<<decaying<<G4endl; |
---|
[819] | 580 | G4Exception("G4QString::LightConeDecay:","72",FatalException,"WrongDecayDirection"); |
---|
[1055] | 581 | } |
---|
| 582 | return 0; |
---|
[819] | 583 | } |
---|
| 584 | |
---|
| 585 | G4LorentzVector G4QString::GetFragmentation4Mom() const |
---|
| 586 | { |
---|
| 587 | G4LorentzVector momentum(Ptleft+Ptright,0.5*(Pplus+Pminus)); |
---|
[1055] | 588 | momentum.setPz(0.5*(Pplus-Pminus)); |
---|
| 589 | return momentum; |
---|
[819] | 590 | } |
---|
| 591 | |
---|
| 592 | // Random choice of string end to use it for creating the hadron (decay) |
---|
| 593 | G4QHadron* G4QString::Splitup(G4bool QL) |
---|
| 594 | { |
---|
| 595 | SideOfDecay = (G4UniformRand() < 0.5)? 1: -1; |
---|
| 596 | if(SideOfDecay<0) SetLeftPartonStable(); // Decay Right parton |
---|
| 597 | else SetRightPartonStable(); // Decay Left parton |
---|
| 598 | G4QParton* newStringEnd; |
---|
| 599 | G4QHadron* Hadron; |
---|
| 600 | if(DecayIsQuark()) Hadron=QuarkSplitup(GetDecayParton(), newStringEnd); // MF1 |
---|
| 601 | else Hadron= DiQuarkSplitup(GetDecayParton(), newStringEnd); // MF2 |
---|
| 602 | // create new String from old, ie. keep Left and Right order, but replace decay |
---|
| 603 | G4LorentzVector* HadronMomentum=SplitEandP(Hadron, QL); // MF3 |
---|
| 604 | if(HadronMomentum) |
---|
| 605 | { |
---|
[1055] | 606 | G4ThreeVector Pos(0.,0.,0.); |
---|
| 607 | Hadron->Set4Momentum(*HadronMomentum); |
---|
| 608 | UpdateString(newStringEnd, HadronMomentum); |
---|
| 609 | delete HadronMomentum; |
---|
[819] | 610 | } |
---|
| 611 | return Hadron; |
---|
| 612 | } // End of Splitup |
---|
| 613 | |
---|
| 614 | void G4QString::UpdateString(G4QParton* decay, const G4LorentzVector* mom) |
---|
| 615 | { |
---|
[1055] | 616 | decaying=None; |
---|
| 617 | if(decaying == Left) |
---|
| 618 | { |
---|
[819] | 619 | G4QParton* theFirst = thePartons[0]; |
---|
[1055] | 620 | delete theFirst; |
---|
| 621 | theFirst = decay; |
---|
| 622 | Ptleft -= mom->vect(); |
---|
| 623 | Ptleft.setZ(0.); |
---|
| 624 | } |
---|
[819] | 625 | else if (decaying == Right) |
---|
[1055] | 626 | { |
---|
[819] | 627 | G4QParton* theLast = thePartons[thePartons.size()-1]; |
---|
| 628 | delete theLast; |
---|
[1055] | 629 | theLast = decay; |
---|
| 630 | Ptright -= mom->vect(); |
---|
| 631 | Ptright.setZ(0.); |
---|
| 632 | } |
---|
[819] | 633 | else |
---|
[1055] | 634 | { |
---|
| 635 | G4cerr<<"***G4QString::UpdateString: wrong oldDecay="<<decaying<<G4endl; |
---|
[819] | 636 | G4Exception("G4QString::UpdateString","72",FatalException,"WrongDecayDirection"); |
---|
[1055] | 637 | } |
---|
| 638 | Pplus -= mom->e() + mom->pz(); |
---|
| 639 | Pminus -= mom->e() - mom->pz(); |
---|
[819] | 640 | } // End of UpdateString |
---|
| 641 | |
---|
| 642 | // QL=true for QGSM and QL=false for Lund fragmentation |
---|
| 643 | G4LorentzVector* G4QString::SplitEandP(G4QHadron* pHadron, G4bool QL) |
---|
| 644 | { |
---|
| 645 | G4double HadronMass = pHadron->GetMass(); |
---|
| 646 | // calculate and assign hadron transverse momentum component HadronPx andHadronPy |
---|
| 647 | G4ThreeVector HadronPt = SampleQuarkPt() + DecayPt(); |
---|
| 648 | HadronPt.setZ(0.); |
---|
| 649 | //... sample z to define hadron longitudinal momentum and energy |
---|
| 650 | //... but first check the available phase space |
---|
| 651 | G4double DecayQuarkMass2 = sqr(GetDecayParton()->GetMass()); // Mass of quark? M.K. |
---|
| 652 | G4double HadronMass2T = HadronMass*HadronMass + HadronPt.mag2(); |
---|
[1055] | 653 | if (DecayQuarkMass2 + HadronMass2T >= SmoothParam*Mass2() ) return 0; // restart! |
---|
[819] | 654 | //... then compute allowed z region z_min <= z <= z_max |
---|
| 655 | G4double zMin = HadronMass2T/Mass2(); |
---|
| 656 | G4double zMax = 1. - DecayQuarkMass2/Mass2(); |
---|
[1055] | 657 | if (zMin >= zMax) return 0; // have to start all over! |
---|
[819] | 658 | G4double z=0; |
---|
| 659 | if(QL) z = GetQGSMLightConeZ(zMin, zMax, GetDecayParton()->GetPDGCode(), pHadron, |
---|
| 660 | HadronPt.x(), HadronPt.y()); |
---|
| 661 | else z = GetLundLightConeZ(zMin, zMax, GetDecayParton()->GetPDGCode(), pHadron, |
---|
| 662 | HadronPt.x(), HadronPt.y()); |
---|
| 663 | //... now compute hadron longitudinal momentum and energy |
---|
| 664 | // longitudinal hadron momentum component HadronPz |
---|
| 665 | G4double zl= z*LightConeDecay(); |
---|
| 666 | G4double HadronE = (zl + HadronMass2T/zl)/2; |
---|
| 667 | HadronPt.setZ( GetDecayDirection() * (zl - HadronE) ); |
---|
| 668 | G4LorentzVector* a4Momentum= new G4LorentzVector(HadronPt,HadronE); |
---|
| 669 | |
---|
| 670 | return a4Momentum; |
---|
| 671 | } |
---|
| 672 | |
---|
| 673 | G4ThreeVector G4QString::SampleQuarkPt() |
---|
| 674 | { |
---|
| 675 | G4double Pt = SigmaQT * std::sqrt( -std::log(G4UniformRand())); |
---|
| 676 | G4double phi = twopi*G4UniformRand(); |
---|
| 677 | return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0); |
---|
| 678 | } |
---|
| 679 | |
---|
| 680 | void G4QString::Sample4Momentum(G4LorentzVector* Mom, G4double Mass, |
---|
| 681 | G4LorentzVector* AntiMom, G4double AntiMass, G4double InitialMass) |
---|
| 682 | { |
---|
| 683 | G4double m2 = Mass*Mass; |
---|
| 684 | G4double am2= AntiMass*AntiMass; |
---|
| 685 | G4double dub=InitialMass*InitialMass - m2 - am2; |
---|
| 686 | G4double r_val = dub - 4*m2*am2; |
---|
| 687 | G4double Pabs = (r_val > 0.) ? std::sqrt(r_val)/(InitialMass*InitialMass) : 0; |
---|
| 688 | //... sample unit vector |
---|
| 689 | G4double r = G4UniformRand(); // @@ G4RandomDirection() |
---|
| 690 | G4double pz = 1. - r - r; // cos(theta) |
---|
| 691 | G4double st = std::sqrt(1. - pz * pz) * Pabs; |
---|
| 692 | G4double phi= twopi*G4UniformRand(); |
---|
| 693 | G4double px = st*std::cos(phi); |
---|
| 694 | G4double py = st*std::sin(phi); |
---|
| 695 | pz *= Pabs; |
---|
| 696 | G4double p2=Pabs*Pabs; |
---|
| 697 | Mom->setPx(px); Mom->setPy(py); Mom->setPz(pz); |
---|
| 698 | Mom->setE(std::sqrt(p2 + Mass*Mass)); |
---|
| 699 | AntiMom->setPx(-px); AntiMom->setPy(-py); AntiMom->setPz(-pz); |
---|
| 700 | AntiMom->setE (std::sqrt(Pabs*Pabs + AntiMass*AntiMass)); |
---|
| 701 | } // End of Sample4Momentum |
---|
| 702 | |
---|
| 703 | G4bool G4QString::SplitLast(G4QString* string, G4QHadronVector* LeftVector, |
---|
| 704 | G4QHadronVector* RightVector) |
---|
| 705 | { |
---|
| 706 | //... perform last cluster decay |
---|
| 707 | G4ThreeVector ClusterVel =string->Get4Momentum().boostVector(); |
---|
| 708 | G4double ResidualMass =string->Mass(); |
---|
| 709 | G4double ClusterMassCut = ClusterMass; |
---|
| 710 | G4int cClusterInterrupt = 0; |
---|
| 711 | G4QHadron* LeftHadron; |
---|
| 712 | G4QHadron* RightHadron; |
---|
| 713 | do |
---|
| 714 | { |
---|
| 715 | if(cClusterInterrupt++ >= ClusterLoopInterrupt) return false; |
---|
[1055] | 716 | G4QParton* quark = 0; |
---|
| 717 | string->SetLeftPartonStable(); // to query quark contents.. |
---|
| 718 | if(string->DecayIsQuark() && string->StableIsQuark()) // There're quarks on clusterEnds |
---|
| 719 | LeftHadron= QuarkSplitup(string->GetLeftParton(), quark); |
---|
| 720 | else |
---|
[819] | 721 | { |
---|
[1055] | 722 | //... there is a Diquark on cluster ends |
---|
| 723 | G4int IsParticle; |
---|
| 724 | if(string->StableIsQuark())IsParticle=(string->GetLeftParton()->GetPDGCode()>0)?-1:1; |
---|
| 725 | else IsParticle=(string->GetLeftParton()->GetPDGCode()>0)?1:-1; |
---|
[819] | 726 | G4QPartonPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted |
---|
| 727 | quark = QuarkPair.GetParton2(); |
---|
| 728 | LeftHadron=hadronizer->Build(QuarkPair.GetParton1(), string->GetLeftParton()); |
---|
[1055] | 729 | } |
---|
[819] | 730 | RightHadron = hadronizer->Build(string->GetRightParton(), quark); |
---|
| 731 | //... repeat procedure, if mass of cluster is too low to produce hadrons |
---|
| 732 | //... ClusterMassCut = 0.15*GeV model parameter |
---|
[1055] | 733 | if ( quark->GetParticleSubType()== "quark" ) ClusterMassCut = 0.; |
---|
| 734 | else ClusterMassCut = ClusterMass; |
---|
[819] | 735 | } while(ResidualMass <= LeftHadron->GetMass() + RightHadron->GetMass() + ClusterMassCut); |
---|
| 736 | //... compute hadron momenta and energies |
---|
| 737 | G4LorentzVector LeftMom, RightMom; |
---|
| 738 | G4ThreeVector Pos; |
---|
| 739 | Sample4Momentum(&LeftMom,LeftHadron->GetMass(),&RightMom,RightHadron->GetMass(), |
---|
| 740 | ResidualMass); |
---|
| 741 | LeftMom.boost(ClusterVel); |
---|
| 742 | RightMom.boost(ClusterVel); |
---|
| 743 | LeftVector->push_back(new G4QHadron(LeftHadron, 0, Pos, LeftMom)); |
---|
| 744 | RightVector->push_back(new G4QHadron(RightHadron, 0, Pos, RightMom)); |
---|
| 745 | |
---|
| 746 | return true; |
---|
| 747 | } // End of SplitLast |
---|
| 748 | |
---|
[1055] | 749 | G4QHadron* G4QString::QuarkSplitup(G4QParton* decay, G4QParton *&created) |
---|
[819] | 750 | { |
---|
[1055] | 751 | G4int IsParticle=(decay->GetPDGCode()>0) ? -1 : +1; // a quark needs antiquark or diquark |
---|
[819] | 752 | G4QPartonPair QuarkPair = CreatePartonPair(IsParticle); |
---|
| 753 | created = QuarkPair.GetParton2(); |
---|
| 754 | return hadronizer->Build(QuarkPair.GetParton1(), decay); |
---|
| 755 | } // End of QuarkSplitup |
---|
| 756 | |
---|
| 757 | // |
---|
| 758 | G4QHadron* G4QString::DiQuarkSplitup(G4QParton* decay, G4QParton* &created) |
---|
| 759 | { |
---|
| 760 | //... can Diquark break or not? |
---|
| 761 | if (G4UniformRand() < DiquarkBreakProb ) |
---|
| 762 | { |
---|
| 763 | //... Diquark break |
---|
| 764 | G4int stableQuarkEncoding = decay->GetPDGCode()/1000; |
---|
| 765 | G4int decayQuarkEncoding = (decay->GetPDGCode()/100)%10; |
---|
| 766 | if (G4UniformRand() < 0.5) |
---|
| 767 | { |
---|
| 768 | G4int Swap = stableQuarkEncoding; |
---|
| 769 | stableQuarkEncoding = decayQuarkEncoding; |
---|
| 770 | decayQuarkEncoding = Swap; |
---|
| 771 | } |
---|
| 772 | G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1; |
---|
[1055] | 773 | // if we have a quark, we need antiquark) |
---|
[819] | 774 | G4QPartonPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted |
---|
| 775 | //... Build new Diquark |
---|
| 776 | G4int QuarkEncoding=QuarkPair.GetParton2()->GetPDGCode(); |
---|
| 777 | G4int i10 = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding)); |
---|
| 778 | G4int i20 = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding)); |
---|
| 779 | G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3; |
---|
| 780 | G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin); |
---|
| 781 | created = CreateParton(NewDecayEncoding); |
---|
| 782 | G4QParton* decayQuark=CreateParton(decayQuarkEncoding); |
---|
| 783 | return hadronizer->Build(QuarkPair.GetParton1(), decayQuark); |
---|
| 784 | } |
---|
| 785 | else |
---|
| 786 | { |
---|
| 787 | //... Diquark does not break |
---|
| 788 | G4int IsParticle=(decay->GetPDGCode()>0) ? +1 : -1; |
---|
[1055] | 789 | // if we have a diquark, we need quark) |
---|
[819] | 790 | G4QPartonPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted |
---|
| 791 | created = QuarkPair.GetParton2(); |
---|
| 792 | return hadronizer->Build(QuarkPair.GetParton1(), decay); |
---|
| 793 | } |
---|
| 794 | } // End of DiQuarkSplitup |
---|
| 795 | |
---|
| 796 | G4QPartonPair G4QString::CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks) |
---|
| 797 | { |
---|
| 798 | // NeedParticle = {+1 for Particle, -1 for AntiParticle} |
---|
| 799 | if(AllowDiquarks && G4UniformRand() < DiquarkSuppress) |
---|
| 800 | { |
---|
| 801 | // Create a Diquark - AntiDiquark pair , first in pair is anti to IsParticle |
---|
| 802 | G4int q1 = SampleQuarkFlavor(); |
---|
| 803 | G4int q2 = SampleQuarkFlavor(); |
---|
| 804 | G4int spin = (q1 != q2 && G4UniformRand() <= 0.5)? 1 : 3; |
---|
| 805 | // Convention: quark with higher PDG number is first |
---|
| 806 | G4int PDGcode = (std::max(q1,q2) * 1000 + std::min(q1,q2) * 100 + spin) * NeedParticle; |
---|
| 807 | return G4QPartonPair(CreateParton(-PDGcode),CreateParton(PDGcode)); |
---|
| 808 | } |
---|
| 809 | else |
---|
| 810 | { |
---|
| 811 | // Create a Quark - AntiQuark pair, first in pair IsParticle |
---|
| 812 | G4int PDGcode=SampleQuarkFlavor()*NeedParticle; |
---|
| 813 | return G4QPartonPair(CreateParton(PDGcode),CreateParton(-PDGcode)); |
---|
| 814 | } |
---|
| 815 | } // End of CreatePartonPair |
---|