| 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: G4VLongitudinalStringDecay.cc,v 1.8 2007/04/24 14:55:23 gunter Exp $
|
|---|
| 28 | // GEANT4 tag $Name: geant4-09-01-patch-02 $
|
|---|
| 29 | //
|
|---|
| 30 | // -----------------------------------------------------------------------------
|
|---|
| 31 | // GEANT 4 class implementation file
|
|---|
| 32 | //
|
|---|
| 33 | // History: first implementation, Maxim Komogorov, 1-Jul-1998
|
|---|
| 34 | // redesign Gunter Folger, August/September 2001
|
|---|
| 35 | // -----------------------------------------------------------------------------
|
|---|
| 36 | #include "G4ios.hh"
|
|---|
| 37 | #include "Randomize.hh"
|
|---|
| 38 | #include "G4VLongitudinalStringDecay.hh"
|
|---|
| 39 | #include "G4FragmentingString.hh"
|
|---|
| 40 |
|
|---|
| 41 | #include "G4ParticleDefinition.hh"
|
|---|
| 42 | #include "G4ParticleTypes.hh"
|
|---|
| 43 | #include "G4ParticleChange.hh"
|
|---|
| 44 | #include "G4VShortLivedParticle.hh"
|
|---|
| 45 | #include "G4ShortLivedConstructor.hh"
|
|---|
| 46 | #include "G4ParticleTable.hh"
|
|---|
| 47 | #include "G4ShortLivedTable.hh"
|
|---|
| 48 | #include "G4PhaseSpaceDecayChannel.hh"
|
|---|
| 49 | #include "G4VDecayChannel.hh"
|
|---|
| 50 | #include "G4DecayTable.hh"
|
|---|
| 51 |
|
|---|
| 52 | #include "G4DiQuarks.hh"
|
|---|
| 53 | #include "G4Quarks.hh"
|
|---|
| 54 | #include "G4Gluons.hh"
|
|---|
| 55 |
|
|---|
| 56 | //------------------------debug switches
|
|---|
| 57 | //#define DEBUG_LightFragmentationTest 1
|
|---|
| 58 |
|
|---|
| 59 |
|
|---|
| 60 | //********************************************************************************
|
|---|
| 61 | // Constructors
|
|---|
| 62 |
|
|---|
| 63 | G4VLongitudinalStringDecay::G4VLongitudinalStringDecay()
|
|---|
| 64 | {
|
|---|
| 65 | MassCut = 0.35*GeV;
|
|---|
| 66 | ClusterMass = 0.15*GeV;
|
|---|
| 67 |
|
|---|
| 68 | SmoothParam = 0.9;
|
|---|
| 69 | StringLoopInterrupt = 1000;
|
|---|
| 70 | ClusterLoopInterrupt = 500;
|
|---|
| 71 |
|
|---|
| 72 | // Changable Parameters below.
|
|---|
| 73 |
|
|---|
| 74 | SigmaQT = 0.5 * GeV;
|
|---|
| 75 |
|
|---|
| 76 | StrangeSuppress = 0.44; // 27 % strange quarks produced, ie. u:d:s=1:1:0.27
|
|---|
| 77 | DiquarkSuppress = 0.1;
|
|---|
| 78 | DiquarkBreakProb = 0.1;
|
|---|
| 79 |
|
|---|
| 80 | //... pspin_meson is probability to create vector meson
|
|---|
| 81 | pspin_meson = 0.5;
|
|---|
| 82 |
|
|---|
| 83 | //... pspin_barion is probability to create 3/2 barion
|
|---|
| 84 | pspin_barion = 0.5;
|
|---|
| 85 |
|
|---|
| 86 | //... vectorMesonMix[] is quark mixing parameters for vector mesons (Variable spin = 3)
|
|---|
| 87 | vectorMesonMix.resize(6);
|
|---|
| 88 | vectorMesonMix[0] = 0.5;
|
|---|
| 89 | vectorMesonMix[1] = 0.0;
|
|---|
| 90 | vectorMesonMix[2] = 0.5;
|
|---|
| 91 | vectorMesonMix[3] = 0.0;
|
|---|
| 92 | vectorMesonMix[4] = 1.0;
|
|---|
| 93 | vectorMesonMix[5] = 1.0;
|
|---|
| 94 |
|
|---|
| 95 | //... scalarMesonMix[] is quark mixing parameters for scalar mesons (Variable spin=1)
|
|---|
| 96 | scalarMesonMix.resize(6);
|
|---|
| 97 | scalarMesonMix[0] = 0.5;
|
|---|
| 98 | scalarMesonMix[1] = 0.25;
|
|---|
| 99 | scalarMesonMix[2] = 0.5;
|
|---|
| 100 | scalarMesonMix[3] = 0.25;
|
|---|
| 101 | scalarMesonMix[4] = 1.0;
|
|---|
| 102 | scalarMesonMix[5] = 0.5;
|
|---|
| 103 |
|
|---|
| 104 | // Parameters may be changed until the first fragmentation starts
|
|---|
| 105 | PastInitPhase=false;
|
|---|
| 106 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 107 | scalarMesonMix,vectorMesonMix);
|
|---|
| 108 | }
|
|---|
| 109 |
|
|---|
| 110 |
|
|---|
| 111 | G4VLongitudinalStringDecay::~G4VLongitudinalStringDecay()
|
|---|
| 112 | {
|
|---|
| 113 | delete hadronizer;
|
|---|
| 114 | }
|
|---|
| 115 |
|
|---|
| 116 | //=============================================================================================-------------
|
|---|
| 117 |
|
|---|
| 118 | // Operators
|
|---|
| 119 |
|
|---|
| 120 | //const & G4VLongitudinalStringDecay::operator=(const G4VLongitudinalStringDecay &)
|
|---|
| 121 | // {
|
|---|
| 122 | // }
|
|---|
| 123 |
|
|---|
| 124 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 125 |
|
|---|
| 126 | int G4VLongitudinalStringDecay::operator==(const G4VLongitudinalStringDecay &) const
|
|---|
| 127 | {
|
|---|
| 128 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator== forbidden");
|
|---|
| 129 | return false;
|
|---|
| 130 | }
|
|---|
| 131 |
|
|---|
| 132 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 133 |
|
|---|
| 134 | int G4VLongitudinalStringDecay::operator!=(const G4VLongitudinalStringDecay &) const
|
|---|
| 135 | {
|
|---|
| 136 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator!= forbidden");
|
|---|
| 137 | return true;
|
|---|
| 138 | }
|
|---|
| 139 |
|
|---|
| 140 | //==========================================================================================================
|
|---|
| 141 |
|
|---|
| 142 | G4int G4VLongitudinalStringDecay::SampleQuarkFlavor(void)
|
|---|
| 143 | {
|
|---|
| 144 | return (1 + (int)(G4UniformRand()/StrangeSuppress));
|
|---|
| 145 | }
|
|---|
| 146 |
|
|---|
| 147 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 148 |
|
|---|
| 149 | G4VLongitudinalStringDecay::pDefPair G4VLongitudinalStringDecay::CreatePartonPair(G4int NeedParticle,G4bool AllowDiquarks)
|
|---|
| 150 | {
|
|---|
| 151 | // NeedParticle = +1 for Particle, -1 for Antiparticle
|
|---|
| 152 |
|
|---|
| 153 | if ( AllowDiquarks && G4UniformRand() < DiquarkSuppress )
|
|---|
| 154 | {
|
|---|
| 155 | // Create a Diquark - AntiDiquark pair , first in pair is anti to IsParticle
|
|---|
| 156 | G4int q1 = SampleQuarkFlavor();
|
|---|
| 157 | G4int q2 = SampleQuarkFlavor();
|
|---|
| 158 | G4int spin = (q1 != q2 && G4UniformRand() <= 0.5)? 1 : 3;
|
|---|
| 159 | // convention: quark with higher PDG number is first
|
|---|
| 160 | G4int PDGcode = (std::max(q1,q2) * 1000 + std::min(q1,q2) * 100 + spin) * NeedParticle;
|
|---|
| 161 | return pDefPair (FindParticle(-PDGcode),FindParticle(PDGcode));
|
|---|
| 162 |
|
|---|
| 163 |
|
|---|
| 164 | } else {
|
|---|
| 165 | // Create a Quark - AntiQuark pair, first in pair IsParticle
|
|---|
| 166 | G4int PDGcode=SampleQuarkFlavor()*NeedParticle;
|
|---|
| 167 | return pDefPair (FindParticle(PDGcode),FindParticle(-PDGcode));
|
|---|
| 168 | }
|
|---|
| 169 |
|
|---|
| 170 | }
|
|---|
| 171 |
|
|---|
| 172 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 173 |
|
|---|
| 174 | // G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt()
|
|---|
| 175 | // {
|
|---|
| 176 | // G4double width_param= 2.0 * GeV*GeV;
|
|---|
| 177 | // G4double R = G4UniformRand();
|
|---|
| 178 | // G4double Pt = std::sqrt(width_param*R/(1-R));
|
|---|
| 179 | // G4double phi = 2.*pi*G4UniformRand();
|
|---|
| 180 | // return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0);
|
|---|
| 181 | // }
|
|---|
| 182 | G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt()
|
|---|
| 183 | {
|
|---|
| 184 | G4double Pt = -std::log(G4UniformRand());
|
|---|
| 185 | Pt = SigmaQT * std::sqrt(Pt);
|
|---|
| 186 | G4double phi = 2.*pi*G4UniformRand();
|
|---|
| 187 | return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0);
|
|---|
| 188 | }
|
|---|
| 189 |
|
|---|
| 190 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 191 |
|
|---|
| 192 | void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons)
|
|---|
| 193 | {
|
|---|
| 194 | // `yo-yo` formation time
|
|---|
| 195 | const G4double kappa = 1.0 * GeV/fermi;
|
|---|
| 196 | for(size_t c1 = 0; c1 < Hadrons->size(); c1++)
|
|---|
| 197 | {
|
|---|
| 198 | G4double SumPz = 0;
|
|---|
| 199 | G4double SumE = 0;
|
|---|
| 200 | for(size_t c2 = 0; c2 < c1; c2++)
|
|---|
| 201 | {
|
|---|
| 202 | SumPz += Hadrons->operator[](c2)->Get4Momentum().pz();
|
|---|
| 203 | SumE += Hadrons->operator[](c2)->Get4Momentum().e();
|
|---|
| 204 | }
|
|---|
| 205 | G4double HadronE = Hadrons->operator[](c1)->Get4Momentum().e();
|
|---|
| 206 | G4double HadronPz = Hadrons->operator[](c1)->Get4Momentum().pz();
|
|---|
| 207 | Hadrons->operator[](c1)->SetFormationTime((theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa));
|
|---|
| 208 | G4ThreeVector aPosition(0, 0, (theInitialStringMass - 2.*SumE - HadronE + HadronPz)/(2.*kappa));
|
|---|
| 209 | Hadrons->operator[](c1)->SetPosition(aPosition);
|
|---|
| 210 | }
|
|---|
| 211 | }
|
|---|
| 212 |
|
|---|
| 213 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 214 |
|
|---|
| 215 | /*
|
|---|
| 216 | void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons)
|
|---|
| 217 | {
|
|---|
| 218 | // 'constituent' formation time
|
|---|
| 219 | const G4double kappa = 1.0 * GeV/fermi;
|
|---|
| 220 | for(G4int c1 = 0; c1 < Hadrons->length(); c1++)
|
|---|
| 221 | {
|
|---|
| 222 | G4double SumPz = 0;
|
|---|
| 223 | G4double SumE = 0;
|
|---|
| 224 | for(G4int c2 = 0; c2 <= c1; c2++)
|
|---|
| 225 | {
|
|---|
| 226 | SumPz += Hadrons->at(c2)->Get4Momentum().pz();
|
|---|
| 227 | SumE += Hadrons->at(c2)->Get4Momentum().e();
|
|---|
| 228 | }
|
|---|
| 229 | Hadrons->at(c1)->SetFormationTime((theInitialStringMass - 2.*SumPz)/(2.*kappa));
|
|---|
| 230 | G4ThreeVector aPosition(0, 0, (theInitialStringMass - 2.*SumE)/(2.*kappa));
|
|---|
| 231 | Hadrons->at(c1)->SetPosition(aPosition);
|
|---|
| 232 | }
|
|---|
| 233 | c1 = Hadrons->length()-1;
|
|---|
| 234 | Hadrons->at(c1)->SetFormationTime(Hadrons->at(c1-1)->GetFormationTime());
|
|---|
| 235 | Hadrons->at(c1)->SetPosition(Hadrons->at(c1-1)->GetPosition());
|
|---|
| 236 | }
|
|---|
| 237 | */
|
|---|
| 238 |
|
|---|
| 239 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 240 |
|
|---|
| 241 | G4ParticleDefinition *
|
|---|
| 242 | G4VLongitudinalStringDecay::QuarkSplitup(G4ParticleDefinition*
|
|---|
| 243 | decay, G4ParticleDefinition *&created)
|
|---|
| 244 | {
|
|---|
| 245 | G4int IsParticle=(decay->GetPDGEncoding()>0) ? -1 : +1; // if we have a quark, we need antiquark (or diquark)
|
|---|
| 246 | pDefPair QuarkPair = CreatePartonPair(IsParticle);
|
|---|
| 247 | created = QuarkPair.second;
|
|---|
| 248 | return hadronizer->Build(QuarkPair.first, decay);
|
|---|
| 249 |
|
|---|
| 250 | }
|
|---|
| 251 |
|
|---|
| 252 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 253 |
|
|---|
| 254 | G4ParticleDefinition *G4VLongitudinalStringDecay::DiQuarkSplitup(
|
|---|
| 255 | G4ParticleDefinition* decay,
|
|---|
| 256 | G4ParticleDefinition *&created)
|
|---|
| 257 | {
|
|---|
| 258 | //... can Diquark break or not?
|
|---|
| 259 | if (G4UniformRand() < DiquarkBreakProb ){
|
|---|
| 260 | //... Diquark break
|
|---|
| 261 |
|
|---|
| 262 | G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000;
|
|---|
| 263 | G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10;
|
|---|
| 264 | if (G4UniformRand() < 0.5)
|
|---|
| 265 | {
|
|---|
| 266 | G4int Swap = stableQuarkEncoding;
|
|---|
| 267 | stableQuarkEncoding = decayQuarkEncoding;
|
|---|
| 268 | decayQuarkEncoding = Swap;
|
|---|
| 269 | }
|
|---|
| 270 |
|
|---|
| 271 | G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1;
|
|---|
| 272 | // if we have a quark, we need antiquark)
|
|---|
| 273 | pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
|
|---|
| 274 | //... Build new Diquark
|
|---|
| 275 | G4int QuarkEncoding=QuarkPair.second->GetPDGEncoding();
|
|---|
| 276 | G4int i10 = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
|
|---|
| 277 | G4int i20 = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
|
|---|
| 278 | G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3;
|
|---|
| 279 | G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin);
|
|---|
| 280 | created = FindParticle(NewDecayEncoding);
|
|---|
| 281 | G4ParticleDefinition * decayQuark=FindParticle(decayQuarkEncoding);
|
|---|
| 282 | G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decayQuark);
|
|---|
| 283 | return had;
|
|---|
| 284 | // return hadronizer->Build(QuarkPair.first, decayQuark);
|
|---|
| 285 |
|
|---|
| 286 | } else {
|
|---|
| 287 | //... Diquark does not break
|
|---|
| 288 |
|
|---|
| 289 | G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1;
|
|---|
| 290 | // if we have a diquark, we need quark)
|
|---|
| 291 | pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
|
|---|
| 292 | created = QuarkPair.second;
|
|---|
| 293 |
|
|---|
| 294 | G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
|
|---|
| 295 | return had;
|
|---|
| 296 | // return G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
|
|---|
| 297 | }
|
|---|
| 298 | }
|
|---|
| 299 |
|
|---|
| 300 | //-----------------------------------------------------------------------------------------
|
|---|
| 301 |
|
|---|
| 302 | G4KineticTrack * G4VLongitudinalStringDecay::Splitup(
|
|---|
| 303 | G4FragmentingString *string,
|
|---|
| 304 | G4FragmentingString *&newString)
|
|---|
| 305 | {
|
|---|
| 306 |
|
|---|
| 307 | //... random choice of string end to use for creating the hadron (decay)
|
|---|
| 308 | SideOfDecay = (G4UniformRand() < 0.5)? 1: -1;
|
|---|
| 309 | if (SideOfDecay < 0)
|
|---|
| 310 | {
|
|---|
| 311 | string->SetLeftPartonStable();
|
|---|
| 312 | } else
|
|---|
| 313 | {
|
|---|
| 314 | string->SetRightPartonStable();
|
|---|
| 315 | }
|
|---|
| 316 |
|
|---|
| 317 | G4ParticleDefinition *newStringEnd;
|
|---|
| 318 | G4ParticleDefinition * HadronDefinition;
|
|---|
| 319 | if (string->DecayIsQuark())
|
|---|
| 320 | {
|
|---|
| 321 | HadronDefinition= QuarkSplitup(string->GetDecayParton(), newStringEnd);
|
|---|
| 322 | } else {
|
|---|
| 323 | HadronDefinition= DiQuarkSplitup(string->GetDecayParton(), newStringEnd);
|
|---|
| 324 | }
|
|---|
| 325 | // create new String from old, ie. keep Left and Right order, but replace decay
|
|---|
| 326 | G4LorentzVector* HadronMomentum=SplitEandP(HadronDefinition, string);
|
|---|
| 327 |
|
|---|
| 328 | G4KineticTrack * Hadron =0;
|
|---|
| 329 | if ( HadronMomentum != 0 ) {
|
|---|
| 330 |
|
|---|
| 331 | G4ThreeVector Pos;
|
|---|
| 332 | Hadron = new G4KineticTrack(HadronDefinition, 0,Pos, *HadronMomentum);
|
|---|
| 333 |
|
|---|
| 334 | newString=new G4FragmentingString(*string,newStringEnd,
|
|---|
| 335 | HadronMomentum);
|
|---|
| 336 |
|
|---|
| 337 | delete HadronMomentum;
|
|---|
| 338 | }
|
|---|
| 339 | return Hadron;
|
|---|
| 340 | }
|
|---|
| 341 |
|
|---|
| 342 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 343 |
|
|---|
| 344 | G4ExcitedString *G4VLongitudinalStringDecay::CPExcited(const G4ExcitedString & in)
|
|---|
| 345 | {
|
|---|
| 346 | G4Parton *Left=new G4Parton(*in.GetLeftParton());
|
|---|
| 347 | G4Parton *Right=new G4Parton(*in.GetRightParton());
|
|---|
| 348 | return new G4ExcitedString(Left,Right,in.GetDirection());
|
|---|
| 349 | }
|
|---|
| 350 |
|
|---|
| 351 | G4double G4VLongitudinalStringDecay::FragmentationMass(
|
|---|
| 352 | const G4FragmentingString *
|
|---|
| 353 | const string,
|
|---|
| 354 | Pcreate build,
|
|---|
| 355 | pDefPair * pdefs)
|
|---|
| 356 | {
|
|---|
| 357 |
|
|---|
| 358 | G4double mass;
|
|---|
| 359 | static G4bool NeedInit(true);
|
|---|
| 360 | static std::vector<double> nomix;
|
|---|
| 361 | static G4HadronBuilder * minMassHadronizer;
|
|---|
| 362 | if ( NeedInit )
|
|---|
| 363 | {
|
|---|
| 364 | NeedInit = false;
|
|---|
| 365 | nomix.resize(6);
|
|---|
| 366 | for ( G4int i=0; i<6 ; i++ ) nomix[i]=0;
|
|---|
| 367 | // minMassHadronizer=new G4HadronBuilder(pspin_meson,pspin_barion,nomix,nomix);
|
|---|
| 368 | minMassHadronizer=hadronizer;
|
|---|
| 369 | }
|
|---|
| 370 |
|
|---|
| 371 | if ( build==0 ) build=&G4HadronBuilder::BuildLowSpin;
|
|---|
| 372 |
|
|---|
| 373 | G4ParticleDefinition *Hadron1, *Hadron2=0;
|
|---|
| 374 |
|
|---|
| 375 | if (!string->FourQuarkString() )
|
|---|
| 376 | {
|
|---|
| 377 | // spin 0 meson or spin 1/2 barion will be built
|
|---|
| 378 |
|
|---|
| 379 | Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
|
|---|
| 380 | string->GetRightParton());
|
|---|
| 381 | mass= (Hadron1)->GetPDGMass();
|
|---|
| 382 | } else
|
|---|
| 383 | {
|
|---|
| 384 | //... string is qq--qqbar: Build two stable hadrons,
|
|---|
| 385 | //... with extra uubar or ddbar quark pair
|
|---|
| 386 | G4int iflc = (G4UniformRand() < 0.5)? 1 : 2;
|
|---|
| 387 | if (string->GetLeftParton()->GetPDGEncoding() < 0) iflc = -iflc;
|
|---|
| 388 |
|
|---|
| 389 | //... theSpin = 4; spin 3/2 baryons will be built
|
|---|
| 390 | Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),FindParticle(iflc));
|
|---|
| 391 | Hadron2 =(minMassHadronizer->*build)(string->GetRightParton(),FindParticle(-iflc));
|
|---|
| 392 | mass = (Hadron1)->GetPDGMass() + (Hadron2)->GetPDGMass();
|
|---|
| 393 | }
|
|---|
| 394 |
|
|---|
| 395 | if ( pdefs != 0 )
|
|---|
| 396 | { // need to return hadrons as well....
|
|---|
| 397 | pdefs->first = Hadron1;
|
|---|
| 398 | pdefs->second = Hadron2;
|
|---|
| 399 | }
|
|---|
| 400 |
|
|---|
| 401 | return mass;
|
|---|
| 402 | }
|
|---|
| 403 |
|
|---|
| 404 | G4KineticTrackVector* G4VLongitudinalStringDecay::LightFragmentationTest(const
|
|---|
| 405 | G4ExcitedString * const string)
|
|---|
| 406 | {
|
|---|
| 407 | // Check string decay threshold
|
|---|
| 408 |
|
|---|
| 409 | G4KineticTrackVector * result=0; // return 0 when string exceeds the mass cut
|
|---|
| 410 |
|
|---|
| 411 | pDefPair hadrons((G4ParticleDefinition *)0,(G4ParticleDefinition *)0);
|
|---|
| 412 | G4FragmentingString aString(*string);
|
|---|
| 413 | if ( sqr(FragmentationMass(&aString,0,&hadrons)+MassCut) < aString.Mass2()) {
|
|---|
| 414 | return 0;
|
|---|
| 415 | }
|
|---|
| 416 |
|
|---|
| 417 | result=new G4KineticTrackVector;
|
|---|
| 418 |
|
|---|
| 419 | if ( hadrons.second ==0 )
|
|---|
| 420 | {
|
|---|
| 421 | // Substitute string by light hadron, Note that Energy is not conserved here!
|
|---|
| 422 |
|
|---|
| 423 | #ifdef DEBUG_LightFragmentationTest
|
|---|
| 424 | G4cout << "VlongSF Warning replacing string by single hadron "
|
|---|
| 425 | << hadrons.first->GetParticleName()
|
|---|
| 426 | << "string .. " << string->Get4Momentum() << " "
|
|---|
| 427 | << string->Get4Momentum().m() << G4endl;
|
|---|
| 428 | #endif
|
|---|
| 429 |
|
|---|
| 430 | G4ThreeVector Mom3 = string->Get4Momentum().vect();
|
|---|
| 431 | G4LorentzVector Mom(Mom3,
|
|---|
| 432 | std::sqrt(Mom3.mag2() + sqr(hadrons.first->GetPDGMass())));
|
|---|
| 433 | result->push_back(new G4KineticTrack(hadrons.first, 0, string->GetPosition(), Mom));
|
|---|
| 434 | } else
|
|---|
| 435 | {
|
|---|
| 436 | //... string was qq--qqbar type: Build two stable hadrons,
|
|---|
| 437 |
|
|---|
| 438 | #ifdef DEBUG_LightFragmentationTest
|
|---|
| 439 | G4cout << "VlongSF Warning replacing qq-qqbar string by TWO hadrons "
|
|---|
| 440 | << hadrons.first->GetParticleName() << " / "
|
|---|
| 441 | << hadrons.second->GetParticleName()
|
|---|
| 442 | << "string .. " << string->Get4Momentum() << " "
|
|---|
| 443 | << string->Get4Momentum().m() << G4endl;
|
|---|
| 444 | #endif
|
|---|
| 445 |
|
|---|
| 446 | G4LorentzVector Mom1, Mom2;
|
|---|
| 447 | Sample4Momentum(&Mom1, hadrons.first->GetPDGMass(),
|
|---|
| 448 | &Mom2,hadrons.second->GetPDGMass(),
|
|---|
| 449 | string->Get4Momentum().mag());
|
|---|
| 450 | result->push_back(new G4KineticTrack(hadrons.first, 0, string->GetPosition(), Mom1));
|
|---|
| 451 | result->push_back(new G4KineticTrack(hadrons.second, 0, string->GetPosition(), Mom2));
|
|---|
| 452 | G4ThreeVector Velocity = string->Get4Momentum().boostVector();
|
|---|
| 453 | result->Boost(Velocity);
|
|---|
| 454 | }
|
|---|
| 455 |
|
|---|
| 456 | return result;
|
|---|
| 457 |
|
|---|
| 458 | }
|
|---|
| 459 |
|
|---|
| 460 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 461 |
|
|---|
| 462 | G4ParticleDefinition* G4VLongitudinalStringDecay::FindParticle(G4int Encoding)
|
|---|
| 463 | {
|
|---|
| 464 | G4ParticleDefinition* ptr = G4ParticleTable::GetParticleTable()->FindParticle(Encoding);
|
|---|
| 465 | if (ptr == NULL)
|
|---|
| 466 | {
|
|---|
| 467 | G4cout << "Particle with encoding "<<Encoding<<" does not exist!!!"<<G4endl;
|
|---|
| 468 | throw G4HadronicException(__FILE__, __LINE__, "Check your particle table");
|
|---|
| 469 | }
|
|---|
| 470 | return ptr;
|
|---|
| 471 | }
|
|---|
| 472 |
|
|---|
| 473 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 474 |
|
|---|
| 475 | void G4VLongitudinalStringDecay::SetSigmaTransverseMomentum(G4double aValue)
|
|---|
| 476 | {
|
|---|
| 477 | if ( PastInitPhase ) {
|
|---|
| 478 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetSigmaTransverseMomentum after FragmentString() not allowed");
|
|---|
| 479 | } else {
|
|---|
| 480 | SigmaQT = aValue;
|
|---|
| 481 | }
|
|---|
| 482 | }
|
|---|
| 483 |
|
|---|
| 484 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 485 |
|
|---|
| 486 | void G4VLongitudinalStringDecay::SetStrangenessSuppression(G4double aValue)
|
|---|
| 487 | {
|
|---|
| 488 | if ( PastInitPhase ) {
|
|---|
| 489 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetStrangenessSuppression after FragmentString() not allowed");
|
|---|
| 490 | } else {
|
|---|
| 491 | StrangeSuppress = aValue;
|
|---|
| 492 | }
|
|---|
| 493 | }
|
|---|
| 494 |
|
|---|
| 495 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 496 |
|
|---|
| 497 | void G4VLongitudinalStringDecay::SetDiquarkSuppression(G4double aValue)
|
|---|
| 498 | {
|
|---|
| 499 | if ( PastInitPhase ) {
|
|---|
| 500 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkSuppression after FragmentString() not allowed");
|
|---|
| 501 | } else {
|
|---|
| 502 | DiquarkSuppress = aValue;
|
|---|
| 503 | }
|
|---|
| 504 | }
|
|---|
| 505 |
|
|---|
| 506 | void G4VLongitudinalStringDecay::SetDiquarkBreakProbability(G4double aValue)
|
|---|
| 507 | {
|
|---|
| 508 | if ( PastInitPhase ) {
|
|---|
| 509 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkBreakProbability after FragmentString() not allowed");
|
|---|
| 510 | } else {
|
|---|
| 511 | DiquarkBreakProb = aValue;
|
|---|
| 512 | }
|
|---|
| 513 | }
|
|---|
| 514 |
|
|---|
| 515 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 516 |
|
|---|
| 517 | void G4VLongitudinalStringDecay::SetVectorMesonProbability(G4double aValue)
|
|---|
| 518 | {
|
|---|
| 519 | if ( PastInitPhase ) {
|
|---|
| 520 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonProbability after FragmentString() not allowed");
|
|---|
| 521 | } else {
|
|---|
| 522 | pspin_meson = aValue;
|
|---|
| 523 | delete hadronizer;
|
|---|
| 524 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 525 | scalarMesonMix,vectorMesonMix);
|
|---|
| 526 | }
|
|---|
| 527 | }
|
|---|
| 528 |
|
|---|
| 529 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 530 |
|
|---|
| 531 | void G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability(G4double aValue)
|
|---|
| 532 | {
|
|---|
| 533 | if ( PastInitPhase ) {
|
|---|
| 534 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability after FragmentString() not allowed");
|
|---|
| 535 | } else {
|
|---|
| 536 | pspin_barion = aValue;
|
|---|
| 537 | delete hadronizer;
|
|---|
| 538 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 539 | scalarMesonMix,vectorMesonMix);
|
|---|
| 540 | }
|
|---|
| 541 | }
|
|---|
| 542 |
|
|---|
| 543 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 544 |
|
|---|
| 545 | void G4VLongitudinalStringDecay::SetScalarMesonMixings(std::vector<G4double> aVector)
|
|---|
| 546 | {
|
|---|
| 547 | if ( PastInitPhase ) {
|
|---|
| 548 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings after FragmentString() not allowed");
|
|---|
| 549 | } else {
|
|---|
| 550 | if ( aVector.size() < 6 )
|
|---|
| 551 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings( argument Vector too small");
|
|---|
| 552 | scalarMesonMix[0] = aVector[0];
|
|---|
| 553 | scalarMesonMix[1] = aVector[1];
|
|---|
| 554 | scalarMesonMix[2] = aVector[2];
|
|---|
| 555 | scalarMesonMix[3] = aVector[3];
|
|---|
| 556 | scalarMesonMix[4] = aVector[4];
|
|---|
| 557 | scalarMesonMix[5] = aVector[5];
|
|---|
| 558 | delete hadronizer;
|
|---|
| 559 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 560 | scalarMesonMix,vectorMesonMix);
|
|---|
| 561 | }
|
|---|
| 562 | }
|
|---|
| 563 |
|
|---|
| 564 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 565 |
|
|---|
| 566 | void G4VLongitudinalStringDecay::SetVectorMesonMixings(std::vector<G4double> aVector)
|
|---|
| 567 | {
|
|---|
| 568 | if ( PastInitPhase ) {
|
|---|
| 569 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings after FragmentString() not allowed");
|
|---|
| 570 | } else {
|
|---|
| 571 | if ( aVector.size() < 6 )
|
|---|
| 572 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings( argument Vector too small");
|
|---|
| 573 | vectorMesonMix[0] = aVector[0];
|
|---|
| 574 | vectorMesonMix[1] = aVector[1];
|
|---|
| 575 | vectorMesonMix[2] = aVector[2];
|
|---|
| 576 | vectorMesonMix[3] = aVector[3];
|
|---|
| 577 | vectorMesonMix[4] = aVector[4];
|
|---|
| 578 | vectorMesonMix[5] = aVector[5];
|
|---|
| 579 | delete hadronizer;
|
|---|
| 580 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 581 | scalarMesonMix,vectorMesonMix);
|
|---|
| 582 |
|
|---|
| 583 | }
|
|---|
| 584 | }
|
|---|
| 585 | //*******************************************************************************************************
|
|---|