| 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.18 2009/08/03 13:21:25 vuzhinsk Exp $
|
|---|
| 28 | // GEANT4 tag $Name: geant4-09-03 $
|
|---|
| 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 | SigmaQT = 0.5 * GeV;
|
|---|
| 74 |
|
|---|
| 75 | StrangeSuppress = 0.44; // 27 % strange quarks produced, ie. u:d:s=1:1:0.27
|
|---|
| 76 | DiquarkSuppress = 0.07;
|
|---|
| 77 | DiquarkBreakProb = 0.1;
|
|---|
| 78 |
|
|---|
| 79 | //... pspin_meson is probability to create vector meson
|
|---|
| 80 | pspin_meson = 0.5;
|
|---|
| 81 |
|
|---|
| 82 | //... pspin_barion is probability to create 3/2 barion
|
|---|
| 83 | pspin_barion = 0.5;
|
|---|
| 84 |
|
|---|
| 85 | //... vectorMesonMix[] is quark mixing parameters for vector mesons (Variable spin = 3)
|
|---|
| 86 | vectorMesonMix.resize(6);
|
|---|
| 87 | vectorMesonMix[0] = 0.5;
|
|---|
| 88 | vectorMesonMix[1] = 0.0;
|
|---|
| 89 | vectorMesonMix[2] = 0.5;
|
|---|
| 90 | vectorMesonMix[3] = 0.0;
|
|---|
| 91 | vectorMesonMix[4] = 1.0;
|
|---|
| 92 | vectorMesonMix[5] = 1.0;
|
|---|
| 93 |
|
|---|
| 94 | //... scalarMesonMix[] is quark mixing parameters for scalar mesons (Variable spin=1)
|
|---|
| 95 | scalarMesonMix.resize(6);
|
|---|
| 96 | scalarMesonMix[0] = 0.5;
|
|---|
| 97 | scalarMesonMix[1] = 0.25;
|
|---|
| 98 | scalarMesonMix[2] = 0.5;
|
|---|
| 99 | scalarMesonMix[3] = 0.25;
|
|---|
| 100 | scalarMesonMix[4] = 1.0;
|
|---|
| 101 | scalarMesonMix[5] = 0.5;
|
|---|
| 102 |
|
|---|
| 103 | // Parameters may be changed until the first fragmentation starts
|
|---|
| 104 | PastInitPhase=false;
|
|---|
| 105 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 106 | scalarMesonMix,vectorMesonMix);
|
|---|
| 107 | Kappa = 1.0 * GeV/fermi;
|
|---|
| 108 |
|
|---|
| 109 |
|
|---|
| 110 | }
|
|---|
| 111 |
|
|---|
| 112 |
|
|---|
| 113 | G4VLongitudinalStringDecay::~G4VLongitudinalStringDecay()
|
|---|
| 114 | {
|
|---|
| 115 | delete hadronizer;
|
|---|
| 116 | }
|
|---|
| 117 |
|
|---|
| 118 | //=============================================================================
|
|---|
| 119 |
|
|---|
| 120 | // Operators
|
|---|
| 121 |
|
|---|
| 122 | //const & G4VLongitudinalStringDecay::operator=(const G4VLongitudinalStringDecay &)
|
|---|
| 123 | // {
|
|---|
| 124 | // }
|
|---|
| 125 |
|
|---|
| 126 | //-----------------------------------------------------------------------------
|
|---|
| 127 |
|
|---|
| 128 | int G4VLongitudinalStringDecay::operator==(const G4VLongitudinalStringDecay &) const
|
|---|
| 129 | {
|
|---|
| 130 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator== forbidden");
|
|---|
| 131 | return false;
|
|---|
| 132 | }
|
|---|
| 133 |
|
|---|
| 134 | //-------------------------------------------------------------------------------------
|
|---|
| 135 |
|
|---|
| 136 | int G4VLongitudinalStringDecay::operator!=(const G4VLongitudinalStringDecay &) const
|
|---|
| 137 | {
|
|---|
| 138 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator!= forbidden");
|
|---|
| 139 | return true;
|
|---|
| 140 | }
|
|---|
| 141 |
|
|---|
| 142 | //***********************************************************************************
|
|---|
| 143 |
|
|---|
| 144 | // For changing Mass Cut used for selection of very small mass strings
|
|---|
| 145 | void G4VLongitudinalStringDecay::SetMassCut(G4double aValue){MassCut=aValue;}
|
|---|
| 146 |
|
|---|
| 147 | //-----------------------------------------------------------------------------
|
|---|
| 148 |
|
|---|
| 149 | // For handling a string with very low mass
|
|---|
| 150 |
|
|---|
| 151 | G4KineticTrackVector* G4VLongitudinalStringDecay::LightFragmentationTest(const
|
|---|
| 152 | G4ExcitedString * const string)
|
|---|
| 153 | {
|
|---|
| 154 | // Check string decay threshold
|
|---|
| 155 |
|
|---|
| 156 | G4KineticTrackVector * result=0; // return 0 when string exceeds the mass cut
|
|---|
| 157 |
|
|---|
| 158 | pDefPair hadrons((G4ParticleDefinition *)0,(G4ParticleDefinition *)0);
|
|---|
| 159 |
|
|---|
| 160 | G4FragmentingString aString(*string);
|
|---|
| 161 |
|
|---|
| 162 | if ( sqr(FragmentationMass(&aString,0,&hadrons)+MassCut) < aString.Mass2()) {
|
|---|
| 163 | return 0;
|
|---|
| 164 | }
|
|---|
| 165 |
|
|---|
| 166 | // The string mass is very low ---------------------------
|
|---|
| 167 |
|
|---|
| 168 | result=new G4KineticTrackVector;
|
|---|
| 169 |
|
|---|
| 170 | if ( hadrons.second ==0 )
|
|---|
| 171 | {
|
|---|
| 172 | // Substitute string by light hadron, Note that Energy is not conserved here!
|
|---|
| 173 |
|
|---|
| 174 | /*
|
|---|
| 175 | #ifdef DEBUG_LightFragmentationTest
|
|---|
| 176 | G4cout << "VlongSF Warning replacing string by single hadron " <<G4endl;
|
|---|
| 177 | G4cout << hadrons.first->GetParticleName()
|
|---|
| 178 | << "string .. " << string->Get4Momentum() << " "
|
|---|
| 179 | << string->Get4Momentum().m() << G4endl;
|
|---|
| 180 | #endif
|
|---|
| 181 | G4cout << "VlongSF Warning replacing string by single hadron " <<G4endl;
|
|---|
| 182 | G4cout << hadrons.first->GetParticleName() << " string .. " <<G4endl;
|
|---|
| 183 | G4cout << string->Get4Momentum() << " " << string->Get4Momentum().m() << G4endl;
|
|---|
| 184 | */
|
|---|
| 185 | G4ThreeVector Mom3 = string->Get4Momentum().vect();
|
|---|
| 186 | G4LorentzVector Mom(Mom3,
|
|---|
| 187 | std::sqrt(Mom3.mag2() +
|
|---|
| 188 | sqr(hadrons.first->GetPDGMass())));
|
|---|
| 189 | result->push_back(new G4KineticTrack(hadrons.first, 0,
|
|---|
| 190 | string->GetPosition(),
|
|---|
| 191 | Mom));
|
|---|
| 192 | } else
|
|---|
| 193 | {
|
|---|
| 194 | //... string was qq--qqbar type: Build two stable hadrons,
|
|---|
| 195 |
|
|---|
| 196 | #ifdef DEBUG_LightFragmentationTest
|
|---|
| 197 | G4cout << "VlongSF Warning replacing qq-qqbar string by TWO hadrons "
|
|---|
| 198 | << hadrons.first->GetParticleName() << " / "
|
|---|
| 199 | << hadrons.second->GetParticleName()
|
|---|
| 200 | << "string .. " << string->Get4Momentum() << " "
|
|---|
| 201 | << string->Get4Momentum().m() << G4endl;
|
|---|
| 202 | #endif
|
|---|
| 203 |
|
|---|
| 204 | // Uzhi Formation time in the case???
|
|---|
| 205 |
|
|---|
| 206 | G4LorentzVector Mom1, Mom2;
|
|---|
| 207 | Sample4Momentum(&Mom1, hadrons.first->GetPDGMass(),
|
|---|
| 208 | &Mom2,hadrons.second->GetPDGMass(),
|
|---|
| 209 | string->Get4Momentum().mag());
|
|---|
| 210 |
|
|---|
| 211 | result->push_back(new G4KineticTrack(hadrons.first, 0,
|
|---|
| 212 | string->GetPosition(),
|
|---|
| 213 | Mom1));
|
|---|
| 214 | result->push_back(new G4KineticTrack(hadrons.second, 0,
|
|---|
| 215 | string->GetPosition(),
|
|---|
| 216 | Mom2));
|
|---|
| 217 |
|
|---|
| 218 | G4ThreeVector Velocity = string->Get4Momentum().boostVector();
|
|---|
| 219 | result->Boost(Velocity);
|
|---|
| 220 | }
|
|---|
| 221 |
|
|---|
| 222 | return result;
|
|---|
| 223 |
|
|---|
| 224 | }
|
|---|
| 225 |
|
|---|
| 226 | //----------------------------------------------------------------------------------------
|
|---|
| 227 |
|
|---|
| 228 | G4double G4VLongitudinalStringDecay::FragmentationMass(
|
|---|
| 229 | const G4FragmentingString * const string,
|
|---|
| 230 | Pcreate build, pDefPair * pdefs )
|
|---|
| 231 | {
|
|---|
| 232 |
|
|---|
| 233 | G4double mass;
|
|---|
| 234 | static G4bool NeedInit(true);
|
|---|
| 235 | static std::vector<double> nomix;
|
|---|
| 236 | static G4HadronBuilder * minMassHadronizer;
|
|---|
| 237 | if ( NeedInit )
|
|---|
| 238 | {
|
|---|
| 239 | NeedInit = false;
|
|---|
| 240 | nomix.resize(6);
|
|---|
| 241 | for ( G4int i=0; i<6 ; i++ ) nomix[i]=0;
|
|---|
| 242 |
|
|---|
| 243 | // minMassHadronizer=new G4HadronBuilder(pspin_meson,pspin_barion,nomix,nomix);
|
|---|
| 244 | minMassHadronizer=hadronizer;
|
|---|
| 245 | }
|
|---|
| 246 |
|
|---|
| 247 | if ( build==0 ) build=&G4HadronBuilder::BuildLowSpin;
|
|---|
| 248 |
|
|---|
| 249 | G4ParticleDefinition *Hadron1, *Hadron2=0;
|
|---|
| 250 |
|
|---|
| 251 | if (!string->FourQuarkString() )
|
|---|
| 252 | {
|
|---|
| 253 | // spin 0 meson or spin 1/2 barion will be built
|
|---|
| 254 |
|
|---|
| 255 | Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
|
|---|
| 256 | string->GetRightParton());
|
|---|
| 257 | mass= (Hadron1)->GetPDGMass();
|
|---|
| 258 | } else
|
|---|
| 259 | {
|
|---|
| 260 | //... string is qq--qqbar: Build two stable hadrons,
|
|---|
| 261 | //... with extra uubar or ddbar quark pair
|
|---|
| 262 | G4int iflc = (G4UniformRand() < 0.5)? 1 : 2;
|
|---|
| 263 | if (string->GetLeftParton()->GetPDGEncoding() < 0) iflc = -iflc;
|
|---|
| 264 |
|
|---|
| 265 | //... theSpin = 4; spin 3/2 baryons will be built
|
|---|
| 266 | Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
|
|---|
| 267 | FindParticle(iflc) );
|
|---|
| 268 | Hadron2 = (minMassHadronizer->*build)(string->GetRightParton(),
|
|---|
| 269 | FindParticle(-iflc) );
|
|---|
| 270 | mass = (Hadron1)->GetPDGMass() + (Hadron2)->GetPDGMass();
|
|---|
| 271 | }
|
|---|
| 272 |
|
|---|
| 273 | if ( pdefs != 0 )
|
|---|
| 274 | { // need to return hadrons as well....
|
|---|
| 275 | pdefs->first = Hadron1;
|
|---|
| 276 | pdefs->second = Hadron2;
|
|---|
| 277 | }
|
|---|
| 278 |
|
|---|
| 279 | return mass;
|
|---|
| 280 | }
|
|---|
| 281 |
|
|---|
| 282 | //----------------------------------------------------------------------------
|
|---|
| 283 |
|
|---|
| 284 | G4ParticleDefinition* G4VLongitudinalStringDecay::FindParticle(G4int Encoding)
|
|---|
| 285 | {
|
|---|
| 286 | G4ParticleDefinition* ptr = G4ParticleTable::GetParticleTable()->FindParticle(Encoding);
|
|---|
| 287 | if (ptr == NULL)
|
|---|
| 288 | {
|
|---|
| 289 | G4cout << "Particle with encoding "<<Encoding<<" does not exist!!!"<<G4endl;
|
|---|
| 290 | throw G4HadronicException(__FILE__, __LINE__, "Check your particle table");
|
|---|
| 291 | }
|
|---|
| 292 | return ptr;
|
|---|
| 293 | }
|
|---|
| 294 |
|
|---|
| 295 | //-----------------------------------------------------------------------------
|
|---|
| 296 | // virtual void Sample4Momentum(G4LorentzVector* Mom, G4double Mass,
|
|---|
| 297 | // G4LorentzVector* AntiMom, G4double AntiMass,
|
|---|
| 298 | // G4double InitialMass)=0;
|
|---|
| 299 | //-----------------------------------------------------------------------------
|
|---|
| 300 |
|
|---|
| 301 | //*********************************************************************************
|
|---|
| 302 | // For decision on continue or stop string fragmentation
|
|---|
| 303 | // virtual G4bool StopFragmenting(const G4FragmentingString * const string)=0;
|
|---|
| 304 | // virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0;
|
|---|
| 305 |
|
|---|
| 306 | // If a string can not fragment, make last break into 2 hadrons
|
|---|
| 307 | // virtual G4bool SplitLast(G4FragmentingString * string,
|
|---|
| 308 | // G4KineticTrackVector * LeftVector,
|
|---|
| 309 | // G4KineticTrackVector * RightVector)=0;
|
|---|
| 310 | //-----------------------------------------------------------------------------
|
|---|
| 311 | //
|
|---|
| 312 | // If a string fragments, do the following
|
|---|
| 313 | //
|
|---|
| 314 | // For transver of a string to its CMS frame
|
|---|
| 315 | //-----------------------------------------------------------------------------
|
|---|
| 316 |
|
|---|
| 317 | G4ExcitedString *G4VLongitudinalStringDecay::CPExcited(const G4ExcitedString & in)
|
|---|
| 318 | {
|
|---|
| 319 | G4Parton *Left=new G4Parton(*in.GetLeftParton());
|
|---|
| 320 | G4Parton *Right=new G4Parton(*in.GetRightParton());
|
|---|
| 321 | return new G4ExcitedString(Left,Right,in.GetDirection());
|
|---|
| 322 | }
|
|---|
| 323 |
|
|---|
| 324 | //-----------------------------------------------------------------------------
|
|---|
| 325 |
|
|---|
| 326 | G4KineticTrack * G4VLongitudinalStringDecay::Splitup(
|
|---|
| 327 | G4FragmentingString *string,
|
|---|
| 328 | G4FragmentingString *&newString)
|
|---|
| 329 | {
|
|---|
| 330 | //G4cout<<"In G4VLong String Dec######################"<<G4endl;
|
|---|
| 331 |
|
|---|
| 332 | //... random choice of string end to use for creating the hadron (decay)
|
|---|
| 333 | SideOfDecay = (G4UniformRand() < 0.5)? 1: -1;
|
|---|
| 334 | if (SideOfDecay < 0)
|
|---|
| 335 | {
|
|---|
| 336 | string->SetLeftPartonStable();
|
|---|
| 337 | } else
|
|---|
| 338 | {
|
|---|
| 339 | string->SetRightPartonStable();
|
|---|
| 340 | }
|
|---|
| 341 |
|
|---|
| 342 | G4ParticleDefinition *newStringEnd;
|
|---|
| 343 | G4ParticleDefinition * HadronDefinition;
|
|---|
| 344 | if (string->DecayIsQuark())
|
|---|
| 345 | {
|
|---|
| 346 | HadronDefinition= QuarkSplitup(string->GetDecayParton(), newStringEnd);
|
|---|
| 347 | } else {
|
|---|
| 348 | HadronDefinition= DiQuarkSplitup(string->GetDecayParton(), newStringEnd);
|
|---|
| 349 | }
|
|---|
| 350 | // create new String from old, ie. keep Left and Right order, but replace decay
|
|---|
| 351 |
|
|---|
| 352 | newString=new G4FragmentingString(*string,newStringEnd); // To store possible
|
|---|
| 353 | // quark containt of new string
|
|---|
| 354 | G4LorentzVector* HadronMomentum=SplitEandP(HadronDefinition, string, newString);
|
|---|
| 355 |
|
|---|
| 356 | delete newString; newString=0; // Uzhi 20.06.08
|
|---|
| 357 |
|
|---|
| 358 | G4KineticTrack * Hadron =0;
|
|---|
| 359 | if ( HadronMomentum != 0 ) {
|
|---|
| 360 |
|
|---|
| 361 | G4ThreeVector Pos;
|
|---|
| 362 | Hadron = new G4KineticTrack(HadronDefinition, 0,Pos, *HadronMomentum);
|
|---|
| 363 |
|
|---|
| 364 | newString=new G4FragmentingString(*string,newStringEnd,
|
|---|
| 365 | HadronMomentum);
|
|---|
| 366 |
|
|---|
| 367 | //G4cout<<"Out G4VLong String Dec######################"<<G4endl;
|
|---|
| 368 | //G4cout<<"newString 4Mom"<<newString->Get4Momentum()<<G4endl;
|
|---|
| 369 | //G4cout<<"newString Pl "<<newString->LightConePlus()<<G4endl;
|
|---|
| 370 | //G4cout<<"newString Mi "<<newString->LightConeMinus()<<G4endl;
|
|---|
| 371 | //G4cout<<"newString Pts "<<newString->StablePt()<<G4endl;
|
|---|
| 372 | //G4cout<<"newString Ptd "<<newString->DecayPt()<<G4endl;
|
|---|
| 373 | //G4cout<<"newString M2 "<<newString->Mass2()<<G4endl;
|
|---|
| 374 | //G4cout<<"newString M2 "<<newString->Mass()<<G4endl;
|
|---|
| 375 | delete HadronMomentum;
|
|---|
| 376 | }
|
|---|
| 377 | return Hadron;
|
|---|
| 378 | }
|
|---|
| 379 |
|
|---|
| 380 | //--------------------------------------------------------------------------------------
|
|---|
| 381 |
|
|---|
| 382 | G4ParticleDefinition *
|
|---|
| 383 | G4VLongitudinalStringDecay::QuarkSplitup(G4ParticleDefinition*
|
|---|
| 384 | decay, G4ParticleDefinition *&created)
|
|---|
| 385 | {
|
|---|
| 386 | G4int IsParticle=(decay->GetPDGEncoding()>0) ? -1 : +1; // if we have a quark,
|
|---|
| 387 | // we need antiquark
|
|---|
| 388 | // (or diquark)
|
|---|
| 389 | pDefPair QuarkPair = CreatePartonPair(IsParticle);
|
|---|
| 390 | created = QuarkPair.second;
|
|---|
| 391 | return hadronizer->Build(QuarkPair.first, decay);
|
|---|
| 392 |
|
|---|
| 393 | }
|
|---|
| 394 |
|
|---|
| 395 | //-----------------------------------------------------------------------------
|
|---|
| 396 |
|
|---|
| 397 | G4ParticleDefinition *G4VLongitudinalStringDecay::DiQuarkSplitup(
|
|---|
| 398 | G4ParticleDefinition* decay,
|
|---|
| 399 | G4ParticleDefinition *&created)
|
|---|
| 400 | {
|
|---|
| 401 | //... can Diquark break or not?
|
|---|
| 402 | if (G4UniformRand() < DiquarkBreakProb ){
|
|---|
| 403 | //... Diquark break
|
|---|
| 404 |
|
|---|
| 405 | G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000;
|
|---|
| 406 | G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10;
|
|---|
| 407 | if (G4UniformRand() < 0.5)
|
|---|
| 408 | {
|
|---|
| 409 | G4int Swap = stableQuarkEncoding;
|
|---|
| 410 | stableQuarkEncoding = decayQuarkEncoding;
|
|---|
| 411 | decayQuarkEncoding = Swap;
|
|---|
| 412 | }
|
|---|
| 413 |
|
|---|
| 414 | G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1;
|
|---|
| 415 | // if we have a quark, we need antiquark)
|
|---|
| 416 | pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
|
|---|
| 417 | //... Build new Diquark
|
|---|
| 418 | G4int QuarkEncoding=QuarkPair.second->GetPDGEncoding();
|
|---|
| 419 | G4int i10 = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
|
|---|
| 420 | G4int i20 = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
|
|---|
| 421 | G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3;
|
|---|
| 422 | G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin);
|
|---|
| 423 | created = FindParticle(NewDecayEncoding);
|
|---|
| 424 | G4ParticleDefinition * decayQuark=FindParticle(decayQuarkEncoding);
|
|---|
| 425 | G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decayQuark);
|
|---|
| 426 | return had;
|
|---|
| 427 | // return hadronizer->Build(QuarkPair.first, decayQuark);
|
|---|
| 428 |
|
|---|
| 429 | } else {
|
|---|
| 430 | //... Diquark does not break
|
|---|
| 431 |
|
|---|
| 432 | G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1;
|
|---|
| 433 | // if we have a diquark, we need quark)
|
|---|
| 434 | pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
|
|---|
| 435 | created = QuarkPair.second;
|
|---|
| 436 |
|
|---|
| 437 | G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
|
|---|
| 438 | return had;
|
|---|
| 439 | // return G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
|
|---|
| 440 | }
|
|---|
| 441 | }
|
|---|
| 442 |
|
|---|
| 443 | //-----------------------------------------------------------------------------
|
|---|
| 444 |
|
|---|
| 445 | G4int G4VLongitudinalStringDecay::SampleQuarkFlavor(void)
|
|---|
| 446 | {
|
|---|
| 447 | return (1 + (int)(G4UniformRand()/StrangeSuppress));
|
|---|
| 448 | }
|
|---|
| 449 |
|
|---|
| 450 | //-----------------------------------------------------------------------------
|
|---|
| 451 |
|
|---|
| 452 | G4VLongitudinalStringDecay::pDefPair G4VLongitudinalStringDecay::CreatePartonPair(G4int NeedParticle,G4bool AllowDiquarks)
|
|---|
| 453 | {
|
|---|
| 454 | // NeedParticle = +1 for Particle, -1 for Antiparticle
|
|---|
| 455 |
|
|---|
| 456 | if ( AllowDiquarks && G4UniformRand() < DiquarkSuppress )
|
|---|
| 457 | {
|
|---|
| 458 | // Create a Diquark - AntiDiquark pair , first in pair is anti to IsParticle
|
|---|
| 459 | G4int q1 = SampleQuarkFlavor();
|
|---|
| 460 | G4int q2 = SampleQuarkFlavor();
|
|---|
| 461 | G4int spin = (q1 != q2 && G4UniformRand() <= 0.5)? 1 : 3;
|
|---|
| 462 | // convention: quark with higher PDG number is first
|
|---|
| 463 | G4int PDGcode = (std::max(q1,q2) * 1000 + std::min(q1,q2) * 100 + spin) * NeedParticle;
|
|---|
| 464 | return pDefPair (FindParticle(-PDGcode),FindParticle(PDGcode));
|
|---|
| 465 |
|
|---|
| 466 |
|
|---|
| 467 | } else {
|
|---|
| 468 | // Create a Quark - AntiQuark pair, first in pair IsParticle
|
|---|
| 469 | G4int PDGcode=SampleQuarkFlavor()*NeedParticle;
|
|---|
| 470 | return pDefPair (FindParticle(PDGcode),FindParticle(-PDGcode));
|
|---|
| 471 | }
|
|---|
| 472 |
|
|---|
| 473 | }
|
|---|
| 474 |
|
|---|
| 475 | //-----------------------------------------------------------------------------
|
|---|
| 476 |
|
|---|
| 477 | // G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt()
|
|---|
| 478 | // {
|
|---|
| 479 | // G4double width_param= 2.0 * GeV*GeV;
|
|---|
| 480 | // G4double R = G4UniformRand();
|
|---|
| 481 | // G4double Pt = std::sqrt(width_param*R/(1-R));
|
|---|
| 482 | // G4double phi = 2.*pi*G4UniformRand();
|
|---|
| 483 | // return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0);
|
|---|
| 484 | // }
|
|---|
| 485 |
|
|---|
| 486 | G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt(G4double ptMax)
|
|---|
| 487 | {
|
|---|
| 488 | G4double Pt;
|
|---|
| 489 | if ( ptMax < 0 ) {
|
|---|
| 490 | // sample full gaussian
|
|---|
| 491 | Pt = -std::log(G4UniformRand());
|
|---|
| 492 | } else {
|
|---|
| 493 | // sample in limited range
|
|---|
| 494 | Pt = -std::log(CLHEP::RandFlat::shoot(std::exp(-sqr(ptMax)/sqr(SigmaQT)), 1.));
|
|---|
| 495 | }
|
|---|
| 496 | Pt = SigmaQT * std::sqrt(Pt);
|
|---|
| 497 | G4double phi = 2.*pi*G4UniformRand();
|
|---|
| 498 | return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0);
|
|---|
| 499 | }
|
|---|
| 500 |
|
|---|
| 501 | //******************************************************************************
|
|---|
| 502 |
|
|---|
| 503 | void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons)
|
|---|
| 504 | {
|
|---|
| 505 |
|
|---|
| 506 | // `yo-yo` formation time
|
|---|
| 507 | // const G4double kappa = 1.0 * GeV/fermi/4.; // Uzhi String tension 1.06.08
|
|---|
| 508 | G4double kappa = GetStringTensionParameter();
|
|---|
| 509 | //G4cout<<"Kappa "<<kappa<<G4endl; // Uzhi 20.06.08
|
|---|
| 510 | //G4int Uzhi; G4cin>>Uzhi; // Uzhi 20.06.08
|
|---|
| 511 | //G4cout<<"Number of hadrons "<<Hadrons->size()<<G4endl; // Vova
|
|---|
| 512 | for(size_t c1 = 0; c1 < Hadrons->size(); c1++)
|
|---|
| 513 | {
|
|---|
| 514 | G4double SumPz = 0;
|
|---|
| 515 | G4double SumE = 0;
|
|---|
| 516 | for(size_t c2 = 0; c2 < c1; c2++)
|
|---|
| 517 | {
|
|---|
| 518 | SumPz += Hadrons->operator[](c2)->Get4Momentum().pz();
|
|---|
| 519 | SumE += Hadrons->operator[](c2)->Get4Momentum().e();
|
|---|
| 520 | }
|
|---|
| 521 | G4double HadronE = Hadrons->operator[](c1)->Get4Momentum().e();
|
|---|
| 522 | G4double HadronPz = Hadrons->operator[](c1)->Get4Momentum().pz();
|
|---|
| 523 | Hadrons->operator[](c1)->SetFormationTime(
|
|---|
| 524 | (theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa)/c_light); //Uzhi1.07.09c_light
|
|---|
| 525 |
|
|---|
| 526 | G4ThreeVector aPosition(0, 0,
|
|---|
| 527 | (theInitialStringMass - 2.*SumE - HadronE + HadronPz)/(2.*kappa));
|
|---|
| 528 | Hadrons->operator[](c1)->SetPosition(aPosition);
|
|---|
| 529 |
|
|---|
| 530 | /*
|
|---|
| 531 | G4cout<<"kappa "<<kappa<<G4endl;
|
|---|
| 532 | G4cout<<c1<<" Partic time, position Old "<<
|
|---|
| 533 | (theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa)<<' '<<
|
|---|
| 534 | (theInitialStringMass - 2.*SumE - HadronE + HadronPz)/(2.*kappa)<<G4endl;
|
|---|
| 535 |
|
|---|
| 536 | G4cout<<c1<<" Partic time, position New "<<
|
|---|
| 537 | (theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa)/c_light<<' '<<
|
|---|
| 538 | (theInitialStringMass - 2.*SumE - HadronE + HadronPz)/(2.*kappa)<<G4endl;
|
|---|
| 539 |
|
|---|
| 540 | G4cout<<"fermi "<<fermi<<" 1/fermi "<<1./fermi<<' '<<"1*fermi/c "<<1.*fermi/c_light<<G4endl;
|
|---|
| 541 | G4int Uzhi; G4cin>>Uzhi; // Uzhi 20.06.08
|
|---|
| 542 | */
|
|---|
| 543 | }
|
|---|
| 544 | }
|
|---|
| 545 |
|
|---|
| 546 | //-----------------------------------------------------------------------------
|
|---|
| 547 |
|
|---|
| 548 | /*
|
|---|
| 549 | void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons)
|
|---|
| 550 | {
|
|---|
| 551 | // 'constituent' formation time
|
|---|
| 552 | const G4double kappa = 1.0 * GeV/fermi;
|
|---|
| 553 | for(G4int c1 = 0; c1 < Hadrons->length(); c1++)
|
|---|
| 554 | {
|
|---|
| 555 | G4double SumPz = 0;
|
|---|
| 556 | G4double SumE = 0;
|
|---|
| 557 | for(G4int c2 = 0; c2 <= c1; c2++)
|
|---|
| 558 | {
|
|---|
| 559 | SumPz += Hadrons->at(c2)->Get4Momentum().pz();
|
|---|
| 560 | SumE += Hadrons->at(c2)->Get4Momentum().e();
|
|---|
| 561 | }
|
|---|
| 562 | Hadrons->at(c1)->SetFormationTime((theInitialStringMass - 2.*SumPz)/(2.*kappa));
|
|---|
| 563 | G4ThreeVector aPosition(0, 0, (theInitialStringMass - 2.*SumE)/(2.*kappa));
|
|---|
| 564 | Hadrons->at(c1)->SetPosition(aPosition);
|
|---|
| 565 | }
|
|---|
| 566 | c1 = Hadrons->length()-1;
|
|---|
| 567 | Hadrons->at(c1)->SetFormationTime(Hadrons->at(c1-1)->GetFormationTime());
|
|---|
| 568 | Hadrons->at(c1)->SetPosition(Hadrons->at(c1-1)->GetPosition());
|
|---|
| 569 | }
|
|---|
| 570 | */
|
|---|
| 571 |
|
|---|
| 572 |
|
|---|
| 573 | //*****************************************************************************
|
|---|
| 574 |
|
|---|
| 575 | void G4VLongitudinalStringDecay::SetSigmaTransverseMomentum(G4double aValue)
|
|---|
| 576 | {
|
|---|
| 577 | if ( PastInitPhase ) {
|
|---|
| 578 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetSigmaTransverseMomentum after FragmentString() not allowed");
|
|---|
| 579 | } else {
|
|---|
| 580 | SigmaQT = aValue;
|
|---|
| 581 | }
|
|---|
| 582 | }
|
|---|
| 583 |
|
|---|
| 584 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 585 |
|
|---|
| 586 | void G4VLongitudinalStringDecay::SetStrangenessSuppression(G4double aValue)
|
|---|
| 587 | {
|
|---|
| 588 | if ( PastInitPhase ) {
|
|---|
| 589 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetStrangenessSuppression after FragmentString() not allowed");
|
|---|
| 590 | } else {
|
|---|
| 591 | StrangeSuppress = aValue;
|
|---|
| 592 | }
|
|---|
| 593 | }
|
|---|
| 594 |
|
|---|
| 595 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 596 |
|
|---|
| 597 | void G4VLongitudinalStringDecay::SetDiquarkSuppression(G4double aValue)
|
|---|
| 598 | {
|
|---|
| 599 | if ( PastInitPhase ) {
|
|---|
| 600 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkSuppression after FragmentString() not allowed");
|
|---|
| 601 | } else {
|
|---|
| 602 | DiquarkSuppress = aValue;
|
|---|
| 603 | }
|
|---|
| 604 | }
|
|---|
| 605 |
|
|---|
| 606 | //----------------------------------------------------------------------------------------
|
|---|
| 607 |
|
|---|
| 608 | void G4VLongitudinalStringDecay::SetDiquarkBreakProbability(G4double aValue)
|
|---|
| 609 | {
|
|---|
| 610 | if ( PastInitPhase ) {
|
|---|
| 611 | throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkBreakProbability after FragmentString() not allowed");
|
|---|
| 612 | } else {
|
|---|
| 613 | DiquarkBreakProb = aValue;
|
|---|
| 614 | }
|
|---|
| 615 | }
|
|---|
| 616 |
|
|---|
| 617 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 618 |
|
|---|
| 619 | void G4VLongitudinalStringDecay::SetVectorMesonProbability(G4double aValue)
|
|---|
| 620 | {
|
|---|
| 621 | if ( PastInitPhase ) {
|
|---|
| 622 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonProbability after FragmentString() not allowed");
|
|---|
| 623 | } else {
|
|---|
| 624 | pspin_meson = aValue;
|
|---|
| 625 | delete hadronizer;
|
|---|
| 626 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 627 | scalarMesonMix,vectorMesonMix);
|
|---|
| 628 | }
|
|---|
| 629 | }
|
|---|
| 630 |
|
|---|
| 631 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 632 |
|
|---|
| 633 | void G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability(G4double aValue)
|
|---|
| 634 | {
|
|---|
| 635 | if ( PastInitPhase ) {
|
|---|
| 636 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability after FragmentString() not allowed");
|
|---|
| 637 | } else {
|
|---|
| 638 | pspin_barion = aValue;
|
|---|
| 639 | delete hadronizer;
|
|---|
| 640 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 641 | scalarMesonMix,vectorMesonMix);
|
|---|
| 642 | }
|
|---|
| 643 | }
|
|---|
| 644 |
|
|---|
| 645 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 646 |
|
|---|
| 647 | void G4VLongitudinalStringDecay::SetScalarMesonMixings(std::vector<G4double> aVector)
|
|---|
| 648 | {
|
|---|
| 649 | if ( PastInitPhase ) {
|
|---|
| 650 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings after FragmentString() not allowed");
|
|---|
| 651 | } else {
|
|---|
| 652 | if ( aVector.size() < 6 )
|
|---|
| 653 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings( argument Vector too small");
|
|---|
| 654 | scalarMesonMix[0] = aVector[0];
|
|---|
| 655 | scalarMesonMix[1] = aVector[1];
|
|---|
| 656 | scalarMesonMix[2] = aVector[2];
|
|---|
| 657 | scalarMesonMix[3] = aVector[3];
|
|---|
| 658 | scalarMesonMix[4] = aVector[4];
|
|---|
| 659 | scalarMesonMix[5] = aVector[5];
|
|---|
| 660 | delete hadronizer;
|
|---|
| 661 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 662 | scalarMesonMix,vectorMesonMix);
|
|---|
| 663 | }
|
|---|
| 664 | }
|
|---|
| 665 |
|
|---|
| 666 | //----------------------------------------------------------------------------------------------------------
|
|---|
| 667 |
|
|---|
| 668 | void G4VLongitudinalStringDecay::SetVectorMesonMixings(std::vector<G4double> aVector)
|
|---|
| 669 | {
|
|---|
| 670 | if ( PastInitPhase ) {
|
|---|
| 671 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings after FragmentString() not allowed");
|
|---|
| 672 | } else {
|
|---|
| 673 | if ( aVector.size() < 6 )
|
|---|
| 674 | throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings( argument Vector too small");
|
|---|
| 675 | vectorMesonMix[0] = aVector[0];
|
|---|
| 676 | vectorMesonMix[1] = aVector[1];
|
|---|
| 677 | vectorMesonMix[2] = aVector[2];
|
|---|
| 678 | vectorMesonMix[3] = aVector[3];
|
|---|
| 679 | vectorMesonMix[4] = aVector[4];
|
|---|
| 680 | vectorMesonMix[5] = aVector[5];
|
|---|
| 681 | delete hadronizer;
|
|---|
| 682 | hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
|
|---|
| 683 | scalarMesonMix,vectorMesonMix);
|
|---|
| 684 |
|
|---|
| 685 | }
|
|---|
| 686 | }
|
|---|
| 687 |
|
|---|
| 688 | //-------------------------------------------------------------------------------------------
|
|---|
| 689 | void G4VLongitudinalStringDecay::SetStringTensionParameter(G4double aValue)// Uzhi 20 June 08
|
|---|
| 690 | {
|
|---|
| 691 | Kappa = aValue * GeV/fermi;
|
|---|
| 692 | }
|
|---|
| 693 | //**************************************************************************************
|
|---|
| 694 |
|
|---|