source: trunk/source/processes/hadronic/models/parton_string/hadronization/src/G4VLongitudinalStringDecay.cc @ 1358

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

update ti head

File size: 23.2 KB
Line 
1//
2// ********************************************************************
3// * License and Disclaimer                                           *
4// *                                                                  *
5// * The  Geant4 software  is  copyright of the Copyright Holders  of *
6// * the Geant4 Collaboration.  It is provided  under  the terms  and *
7// * conditions of the Geant4 Software License,  included in the file *
8// * LICENSE and available at  http://cern.ch/geant4/license .  These *
9// * include a list of copyright holders.                             *
10// *                                                                  *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work  make  any representation or  warranty, express or implied, *
14// * regarding  this  software system or assume any liability for its *
15// * use.  Please see the license in the file  LICENSE  and URL above *
16// * for the full disclaimer and the limitation of liability.         *
17// *                                                                  *
18// * This  code  implementation is the result of  the  scientific and *
19// * technical work of the GEANT4 collaboration.                      *
20// * By using,  copying,  modifying or  distributing the software (or *
21// * any work based  on the software)  you  agree  to acknowledge its *
22// * use  in  resulting  scientific  publications,  and indicate your *
23// * acceptance of all terms of the Geant4 Software license.          *
24// ********************************************************************
25//
26//
27// $Id: G4VLongitudinalStringDecay.cc,v 1.22 2010/09/20 12:46:23 vuzhinsk Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
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
63G4VLongitudinalStringDecay::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;  // 0.5
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
113G4VLongitudinalStringDecay::~G4VLongitudinalStringDecay()
114   {
115   delete hadronizer;
116   }
117
118//=============================================================================
119
120// Operators
121
122//const  & G4VLongitudinalStringDecay::operator=(const G4VLongitudinalStringDecay &)
123//    {
124//    }
125
126//-----------------------------------------------------------------------------
127
128int G4VLongitudinalStringDecay::operator==(const G4VLongitudinalStringDecay &) const
129    {
130        throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::operator== forbidden");
131        return false;
132    }
133
134//-------------------------------------------------------------------------------------
135
136int 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
145void G4VLongitudinalStringDecay::SetMassCut(G4double aValue){MassCut=aValue;}
146
147//-----------------------------------------------------------------------------
148
149// For handling a string with very low mass
150
151G4KineticTrackVector* 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*/
182               G4ThreeVector   Mom3 = string->Get4Momentum().vect();
183               G4LorentzVector Mom(Mom3, 
184                                   std::sqrt(Mom3.mag2() + 
185                                             sqr(hadrons.first->GetPDGMass())));
186               result->push_back(new G4KineticTrack(hadrons.first, 0, 
187                                                  string->GetPosition(),
188                                                          Mom));
189        } else 
190        {
191//... string was qq--qqbar type: Build two stable hadrons,
192
193#ifdef DEBUG_LightFragmentationTest
194               G4cout << "VlongSF Warning replacing qq-qqbar string by TWO hadrons " 
195                      << hadrons.first->GetParticleName() << " / " 
196                      << hadrons.second->GetParticleName()
197                      << "string .. " << string->Get4Momentum() << " " 
198                      << string->Get4Momentum().m() << G4endl;
199#endif               
200
201               G4LorentzVector  Mom1, Mom2;
202               Sample4Momentum(&Mom1, hadrons.first->GetPDGMass(), 
203                               &Mom2,hadrons.second->GetPDGMass(),
204                               string->Get4Momentum().mag());
205
206               result->push_back(new G4KineticTrack(hadrons.first, 0, 
207                                                    string->GetPosition(), 
208                                                            Mom1));
209               result->push_back(new G4KineticTrack(hadrons.second, 0, 
210                                                    string->GetPosition(), 
211                                                    Mom2));
212
213               G4ThreeVector Velocity = string->Get4Momentum().boostVector();
214               result->Boost(Velocity);         
215        }
216
217        return result;
218       
219}
220
221//----------------------------------------------------------------------------------------
222
223G4double G4VLongitudinalStringDecay::FragmentationMass(
224            const G4FragmentingString * const string,
225                Pcreate build, pDefPair * pdefs       )
226{
227       
228        G4double mass;
229        static G4bool NeedInit(true);
230        static std::vector<double> nomix;
231        static G4HadronBuilder * minMassHadronizer;
232        if ( NeedInit ) 
233        {
234           NeedInit = false;
235           nomix.resize(6);
236           for ( G4int i=0; i<6 ; i++ ) nomix[i]=0;
237
238//         minMassHadronizer=new G4HadronBuilder(pspin_meson,pspin_barion,nomix,nomix);
239           minMassHadronizer=hadronizer;
240        }
241
242        if ( build==0 ) build=&G4HadronBuilder::BuildLowSpin;
243
244        G4ParticleDefinition *Hadron1, *Hadron2=0;
245
246        if (!string->FourQuarkString() )
247        {
248           // spin 0 meson or spin 1/2 barion will be built
249
250           Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
251                                                 string->GetRightParton());
252//G4cout<<"Hadron1 "<<Hadron1->GetParticleName()<<G4endl;
253           mass= (Hadron1)->GetPDGMass();
254        } else
255        {
256           //... string is qq--qqbar: Build two stable hadrons,
257           //...    with extra uubar or ddbar quark pair
258           G4int iflc = (G4UniformRand() < 0.5)? 1 : 2;
259           if (string->GetLeftParton()->GetPDGEncoding() < 0) iflc = -iflc;
260
261           //... theSpin = 4; spin 3/2 baryons will be built
262           Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
263                                                 FindParticle(iflc)       );
264           Hadron2 = (minMassHadronizer->*build)(string->GetRightParton(),
265                                                 FindParticle(-iflc)      );
266           mass = (Hadron1)->GetPDGMass() + (Hadron2)->GetPDGMass();
267        }
268       
269        if ( pdefs != 0 ) 
270        { // need to return hadrons as well....
271           pdefs->first  = Hadron1;
272           pdefs->second = Hadron2;
273        }
274           
275        return mass;
276}
277
278//----------------------------------------------------------------------------
279
280G4ParticleDefinition* G4VLongitudinalStringDecay::FindParticle(G4int Encoding) 
281   {
282   G4ParticleDefinition* ptr = G4ParticleTable::GetParticleTable()->FindParticle(Encoding);
283      if (ptr == NULL)
284       {
285       G4cout << "Particle with encoding "<<Encoding<<" does not exist!!!"<<G4endl;
286       throw G4HadronicException(__FILE__, __LINE__, "Check your particle table");
287       }
288   return ptr;   
289   }
290
291//-----------------------------------------------------------------------------
292//   virtual void Sample4Momentum(G4LorentzVector* Mom,     G4double Mass,
293//                                G4LorentzVector* AntiMom, G4double AntiMass,
294//                                G4double InitialMass)=0;
295//-----------------------------------------------------------------------------
296
297//*********************************************************************************
298//   For decision on continue or stop string fragmentation
299//   virtual G4bool StopFragmenting(const G4FragmentingString  * const string)=0;
300//   virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0;
301
302//   If a string can not fragment, make last break into 2 hadrons
303//   virtual G4bool SplitLast(G4FragmentingString * string,
304//                            G4KineticTrackVector * LeftVector,
305//                            G4KineticTrackVector * RightVector)=0;
306//-----------------------------------------------------------------------------
307//
308//   If a string fragments, do the following
309//
310//   For transver of a string to its CMS frame
311//-----------------------------------------------------------------------------
312
313G4ExcitedString *G4VLongitudinalStringDecay::CPExcited(const G4ExcitedString & in)
314{
315        G4Parton *Left=new G4Parton(*in.GetLeftParton());
316        G4Parton *Right=new G4Parton(*in.GetRightParton());
317        return new G4ExcitedString(Left,Right,in.GetDirection());
318}
319
320//-----------------------------------------------------------------------------
321
322G4KineticTrack * G4VLongitudinalStringDecay::Splitup(
323                        G4FragmentingString *string, 
324                        G4FragmentingString *&newString)
325{
326//G4cout<<"Start SplitUP"<<G4endl;
327       //... random choice of string end to use for creating the hadron (decay)   
328       SideOfDecay = (G4UniformRand() < 0.5)? 1: -1;
329       if (SideOfDecay < 0)
330       {
331          string->SetLeftPartonStable();
332       } else
333       {
334          string->SetRightPartonStable();
335       }
336
337       G4ParticleDefinition *newStringEnd;
338       G4ParticleDefinition * HadronDefinition;
339       if (string->DecayIsQuark())
340       {
341           HadronDefinition= QuarkSplitup(string->GetDecayParton(), newStringEnd);
342       } else {
343           HadronDefinition= DiQuarkSplitup(string->GetDecayParton(), newStringEnd);
344       }     
345
346//G4cout<<"New had "<<HadronDefinition->GetParticleName()<<G4endl;
347// create new String from old, ie. keep Left and Right order, but replace decay
348
349       newString=new G4FragmentingString(*string,newStringEnd); // To store possible
350                                                                // quark containt of new string
351//G4cout<<"SplitEandP "<<G4endl;
352       G4LorentzVector* HadronMomentum=SplitEandP(HadronDefinition, string, newString);
353
354       delete newString; newString=0;                         
355       
356       G4KineticTrack * Hadron =0;
357       if ( HadronMomentum != 0 ) {   
358
359           G4ThreeVector   Pos;
360           Hadron = new G4KineticTrack(HadronDefinition, 0,Pos, *HadronMomentum);
361 
362           newString=new G4FragmentingString(*string,newStringEnd,
363                                        HadronMomentum);
364
365           delete HadronMomentum;
366       }     
367//G4cout<<"End SplitUP"<<G4endl;
368       return Hadron;
369}
370
371//--------------------------------------------------------------------------------------
372
373G4ParticleDefinition *
374                G4VLongitudinalStringDecay::QuarkSplitup(G4ParticleDefinition*
375                decay, G4ParticleDefinition *&created)
376{
377    G4int IsParticle=(decay->GetPDGEncoding()>0) ? -1 : +1; // if we have a quark,
378                                                            // we need antiquark
379                                                            // (or diquark)
380    pDefPair QuarkPair = CreatePartonPair(IsParticle);
381    created = QuarkPair.second;
382    return hadronizer->Build(QuarkPair.first, decay);
383   
384}
385
386//-----------------------------------------------------------------------------
387
388G4ParticleDefinition *G4VLongitudinalStringDecay::DiQuarkSplitup(
389                                                        G4ParticleDefinition* decay,
390                                                        G4ParticleDefinition *&created)
391{
392   //... can Diquark break or not?
393   if (G4UniformRand() < DiquarkBreakProb ){
394   //... Diquark break
395
396      G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000;
397      G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10;
398      if (G4UniformRand() < 0.5)
399         {
400         G4int Swap = stableQuarkEncoding;
401         stableQuarkEncoding = decayQuarkEncoding;
402         decayQuarkEncoding = Swap;
403         }
404
405      G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1; 
406                        // if we have a quark, we need antiquark)
407      pDefPair QuarkPair = CreatePartonPair(IsParticle,false);  // no diquarks wanted
408      //... Build new Diquark
409      G4int QuarkEncoding=QuarkPair.second->GetPDGEncoding();
410      G4int i10  = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
411      G4int i20  = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
412      G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3;
413      G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin);
414      created = FindParticle(NewDecayEncoding);
415      G4ParticleDefinition * decayQuark=FindParticle(decayQuarkEncoding);
416      G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decayQuark);
417      return had;
418//      return hadronizer->Build(QuarkPair.first, decayQuark);
419   
420   } else {
421   //... Diquark does not break
422 
423      G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1; 
424                        // if we have a diquark, we need quark)
425      pDefPair QuarkPair = CreatePartonPair(IsParticle,false);  // no diquarks wanted
426      created = QuarkPair.second;
427
428      G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
429      return had;
430//      return G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
431   }
432}
433
434//-----------------------------------------------------------------------------
435
436G4int G4VLongitudinalStringDecay::SampleQuarkFlavor(void)
437   {
438   return (1 + (int)(G4UniformRand()/StrangeSuppress));
439   }
440
441//-----------------------------------------------------------------------------
442
443G4VLongitudinalStringDecay::pDefPair G4VLongitudinalStringDecay::CreatePartonPair(G4int NeedParticle,G4bool AllowDiquarks)
444{
445//  NeedParticle = +1 for Particle, -1 for Antiparticle
446
447    if ( AllowDiquarks && G4UniformRand() < DiquarkSuppress )
448    {
449      // Create a Diquark - AntiDiquark pair , first in pair is anti to IsParticle
450      G4int q1  = SampleQuarkFlavor();
451      G4int q2  = SampleQuarkFlavor();
452      G4int spin = (q1 != q2 && G4UniformRand() <= 0.5)? 1 : 3;
453                                     //   convention: quark with higher PDG number is first
454      G4int PDGcode = (std::max(q1,q2) * 1000 + std::min(q1,q2) * 100 + spin) * NeedParticle;
455      return pDefPair (FindParticle(-PDGcode),FindParticle(PDGcode));
456     
457
458    } else {
459      // Create a Quark - AntiQuark pair, first in pair  IsParticle
460      G4int PDGcode=SampleQuarkFlavor()*NeedParticle;
461      return pDefPair (FindParticle(PDGcode),FindParticle(-PDGcode));
462    }
463
464}
465
466//-----------------------------------------------------------------------------
467G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt(G4double ptMax)
468   {
469   G4double Pt;
470   if ( ptMax < 0 ) {
471      // sample full gaussian
472      Pt = -std::log(G4UniformRand());
473   } else {
474      // sample in limited range
475      Pt = -std::log(CLHEP::RandFlat::shoot(std::exp(-sqr(ptMax)/sqr(SigmaQT)), 1.));
476   }
477   Pt = SigmaQT * std::sqrt(Pt);
478   G4double phi = 2.*pi*G4UniformRand();
479   return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0);
480   }
481
482//******************************************************************************
483
484void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons)
485   {
486
487//   `yo-yo` formation time
488//   const G4double kappa = 1.0 * GeV/fermi/4.;     
489   G4double kappa = GetStringTensionParameter();
490   for(size_t c1 = 0; c1 < Hadrons->size(); c1++)
491      {
492      G4double SumPz = 0; 
493      G4double SumE  = 0;
494      for(size_t c2 = 0; c2 < c1; c2++)
495         {
496         SumPz += Hadrons->operator[](c2)->Get4Momentum().pz();
497         SumE  += Hadrons->operator[](c2)->Get4Momentum().e();   
498         } 
499      G4double HadronE  = Hadrons->operator[](c1)->Get4Momentum().e();
500      G4double HadronPz = Hadrons->operator[](c1)->Get4Momentum().pz();
501      Hadrons->operator[](c1)->SetFormationTime(
502(theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa)/c_light); 
503
504      G4ThreeVector aPosition(0, 0,     
505(theInitialStringMass - 2.*SumE  - HadronE + HadronPz)/(2.*kappa));
506      Hadrons->operator[](c1)->SetPosition(aPosition);
507
508      } 
509   }
510
511//-----------------------------------------------------------------------------
512
513void G4VLongitudinalStringDecay::SetSigmaTransverseMomentum(G4double aValue)
514{
515        if ( PastInitPhase ) {
516                throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetSigmaTransverseMomentum after FragmentString() not allowed");
517        } else {
518                SigmaQT = aValue;
519        }
520}
521
522//----------------------------------------------------------------------------------------------------------
523
524void G4VLongitudinalStringDecay::SetStrangenessSuppression(G4double aValue)
525{
526        if ( PastInitPhase ) {
527                throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetStrangenessSuppression after FragmentString() not allowed");
528        } else {
529                StrangeSuppress = aValue;
530        }
531}
532
533//----------------------------------------------------------------------------------------------------------
534
535void G4VLongitudinalStringDecay::SetDiquarkSuppression(G4double aValue)
536{
537        if ( PastInitPhase ) {
538                throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkSuppression after FragmentString() not allowed");
539        } else {
540                DiquarkSuppress = aValue;
541        }
542}
543
544//----------------------------------------------------------------------------------------
545
546void G4VLongitudinalStringDecay::SetDiquarkBreakProbability(G4double aValue)
547{
548        if ( PastInitPhase ) {
549                throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkBreakProbability after FragmentString() not allowed");
550        } else {
551                DiquarkBreakProb = aValue;
552        }
553}
554
555//----------------------------------------------------------------------------------------------------------
556
557void G4VLongitudinalStringDecay::SetVectorMesonProbability(G4double aValue)
558{
559        if ( PastInitPhase ) {
560                throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonProbability after FragmentString() not allowed");
561        } else {
562                pspin_meson = aValue;
563                delete hadronizer;
564                hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
565                                scalarMesonMix,vectorMesonMix);
566        }
567}
568
569//----------------------------------------------------------------------------------------------------------
570
571void G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability(G4double aValue)
572{
573        if ( PastInitPhase ) {
574                throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability after FragmentString() not allowed");
575        } else {
576                pspin_barion = aValue;
577                delete hadronizer;
578                hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
579                                scalarMesonMix,vectorMesonMix);
580        }
581}
582
583//----------------------------------------------------------------------------------------------------------
584
585void G4VLongitudinalStringDecay::SetScalarMesonMixings(std::vector<G4double> aVector)
586{
587        if ( PastInitPhase ) {
588                throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings after FragmentString() not allowed");
589        } else {
590          if ( aVector.size() < 6 ) 
591              throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings( argument Vector too small");
592          scalarMesonMix[0] = aVector[0];
593          scalarMesonMix[1] = aVector[1];
594          scalarMesonMix[2] = aVector[2];
595          scalarMesonMix[3] = aVector[3];
596          scalarMesonMix[4] = aVector[4];
597          scalarMesonMix[5] = aVector[5];
598          delete hadronizer;
599          hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
600                                scalarMesonMix,vectorMesonMix);
601        }
602}
603
604//----------------------------------------------------------------------------------------------------------
605
606void G4VLongitudinalStringDecay::SetVectorMesonMixings(std::vector<G4double> aVector)
607{
608        if ( PastInitPhase ) {
609                throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings after FragmentString() not allowed");
610        } else {
611          if ( aVector.size() < 6 ) 
612              throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings( argument Vector too small");
613          vectorMesonMix[0] = aVector[0];
614          vectorMesonMix[1] = aVector[1];
615          vectorMesonMix[2] = aVector[2];
616          vectorMesonMix[3] = aVector[3];
617          vectorMesonMix[4] = aVector[4];
618          vectorMesonMix[5] = aVector[5];
619          delete hadronizer;
620          hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
621                                scalarMesonMix,vectorMesonMix);
622 
623        }
624}
625
626//-------------------------------------------------------------------------------------------
627void G4VLongitudinalStringDecay::SetStringTensionParameter(G4double aValue)// Uzhi 20 June 08
628{
629          Kappa = aValue * GeV/fermi;
630}       
631//**************************************************************************************
632
Note: See TracBrowser for help on using the repository browser.