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

Last change on this file since 1338 was 1337, checked in by garnier, 15 years ago

tag geant4.9.4 beta 1 + modifs locales

File size: 23.0 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.20 2010/06/21 17:50:48 vuzhinsk Exp $
28// GEANT4 tag $Name: geant4-09-04-beta-01 $
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;
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 mass= (Hadron1)->GetPDGMass();
253 } else
254 {
255 //... string is qq--qqbar: Build two stable hadrons,
256 //... with extra uubar or ddbar quark pair
257 G4int iflc = (G4UniformRand() < 0.5)? 1 : 2;
258 if (string->GetLeftParton()->GetPDGEncoding() < 0) iflc = -iflc;
259
260 //... theSpin = 4; spin 3/2 baryons will be built
261 Hadron1 = (minMassHadronizer->*build)(string->GetLeftParton(),
262 FindParticle(iflc) );
263 Hadron2 = (minMassHadronizer->*build)(string->GetRightParton(),
264 FindParticle(-iflc) );
265 mass = (Hadron1)->GetPDGMass() + (Hadron2)->GetPDGMass();
266 }
267
268 if ( pdefs != 0 )
269 { // need to return hadrons as well....
270 pdefs->first = Hadron1;
271 pdefs->second = Hadron2;
272 }
273
274 return mass;
275}
276
277//----------------------------------------------------------------------------
278
279G4ParticleDefinition* G4VLongitudinalStringDecay::FindParticle(G4int Encoding)
280 {
281 G4ParticleDefinition* ptr = G4ParticleTable::GetParticleTable()->FindParticle(Encoding);
282 if (ptr == NULL)
283 {
284 G4cout << "Particle with encoding "<<Encoding<<" does not exist!!!"<<G4endl;
285 throw G4HadronicException(__FILE__, __LINE__, "Check your particle table");
286 }
287 return ptr;
288 }
289
290//-----------------------------------------------------------------------------
291// virtual void Sample4Momentum(G4LorentzVector* Mom, G4double Mass,
292// G4LorentzVector* AntiMom, G4double AntiMass,
293// G4double InitialMass)=0;
294//-----------------------------------------------------------------------------
295
296//*********************************************************************************
297// For decision on continue or stop string fragmentation
298// virtual G4bool StopFragmenting(const G4FragmentingString * const string)=0;
299// virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0;
300
301// If a string can not fragment, make last break into 2 hadrons
302// virtual G4bool SplitLast(G4FragmentingString * string,
303// G4KineticTrackVector * LeftVector,
304// G4KineticTrackVector * RightVector)=0;
305//-----------------------------------------------------------------------------
306//
307// If a string fragments, do the following
308//
309// For transver of a string to its CMS frame
310//-----------------------------------------------------------------------------
311
312G4ExcitedString *G4VLongitudinalStringDecay::CPExcited(const G4ExcitedString & in)
313{
314 G4Parton *Left=new G4Parton(*in.GetLeftParton());
315 G4Parton *Right=new G4Parton(*in.GetRightParton());
316 return new G4ExcitedString(Left,Right,in.GetDirection());
317}
318
319//-----------------------------------------------------------------------------
320
321G4KineticTrack * G4VLongitudinalStringDecay::Splitup(
322 G4FragmentingString *string,
323 G4FragmentingString *&newString)
324{
325 //... random choice of string end to use for creating the hadron (decay)
326 SideOfDecay = (G4UniformRand() < 0.5)? 1: -1;
327 if (SideOfDecay < 0)
328 {
329 string->SetLeftPartonStable();
330 } else
331 {
332 string->SetRightPartonStable();
333 }
334
335 G4ParticleDefinition *newStringEnd;
336 G4ParticleDefinition * HadronDefinition;
337 if (string->DecayIsQuark())
338 {
339 HadronDefinition= QuarkSplitup(string->GetDecayParton(), newStringEnd);
340 } else {
341 HadronDefinition= DiQuarkSplitup(string->GetDecayParton(), newStringEnd);
342 }
343// create new String from old, ie. keep Left and Right order, but replace decay
344
345 newString=new G4FragmentingString(*string,newStringEnd); // To store possible
346 // quark containt of new string
347 G4LorentzVector* HadronMomentum=SplitEandP(HadronDefinition, string, newString);
348
349 delete newString; newString=0;
350
351 G4KineticTrack * Hadron =0;
352 if ( HadronMomentum != 0 ) {
353
354 G4ThreeVector Pos;
355 Hadron = new G4KineticTrack(HadronDefinition, 0,Pos, *HadronMomentum);
356
357 newString=new G4FragmentingString(*string,newStringEnd,
358 HadronMomentum);
359
360 delete HadronMomentum;
361 }
362 return Hadron;
363}
364
365//--------------------------------------------------------------------------------------
366
367G4ParticleDefinition *
368 G4VLongitudinalStringDecay::QuarkSplitup(G4ParticleDefinition*
369 decay, G4ParticleDefinition *&created)
370{
371 G4int IsParticle=(decay->GetPDGEncoding()>0) ? -1 : +1; // if we have a quark,
372 // we need antiquark
373 // (or diquark)
374 pDefPair QuarkPair = CreatePartonPair(IsParticle);
375 created = QuarkPair.second;
376 return hadronizer->Build(QuarkPair.first, decay);
377
378}
379
380//-----------------------------------------------------------------------------
381
382G4ParticleDefinition *G4VLongitudinalStringDecay::DiQuarkSplitup(
383 G4ParticleDefinition* decay,
384 G4ParticleDefinition *&created)
385{
386 //... can Diquark break or not?
387 if (G4UniformRand() < DiquarkBreakProb ){
388 //... Diquark break
389
390 G4int stableQuarkEncoding = decay->GetPDGEncoding()/1000;
391 G4int decayQuarkEncoding = (decay->GetPDGEncoding()/100)%10;
392 if (G4UniformRand() < 0.5)
393 {
394 G4int Swap = stableQuarkEncoding;
395 stableQuarkEncoding = decayQuarkEncoding;
396 decayQuarkEncoding = Swap;
397 }
398
399 G4int IsParticle=(decayQuarkEncoding>0) ? -1 : +1;
400 // if we have a quark, we need antiquark)
401 pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
402 //... Build new Diquark
403 G4int QuarkEncoding=QuarkPair.second->GetPDGEncoding();
404 G4int i10 = std::max(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
405 G4int i20 = std::min(std::abs(QuarkEncoding), std::abs(stableQuarkEncoding));
406 G4int spin = (i10 != i20 && G4UniformRand() <= 0.5)? 1 : 3;
407 G4int NewDecayEncoding = -1*IsParticle*(i10 * 1000 + i20 * 100 + spin);
408 created = FindParticle(NewDecayEncoding);
409 G4ParticleDefinition * decayQuark=FindParticle(decayQuarkEncoding);
410 G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decayQuark);
411 return had;
412// return hadronizer->Build(QuarkPair.first, decayQuark);
413
414 } else {
415 //... Diquark does not break
416
417 G4int IsParticle=(decay->GetPDGEncoding()>0) ? +1 : -1;
418 // if we have a diquark, we need quark)
419 pDefPair QuarkPair = CreatePartonPair(IsParticle,false); // no diquarks wanted
420 created = QuarkPair.second;
421
422 G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
423 return had;
424// return G4ParticleDefinition * had=hadronizer->Build(QuarkPair.first, decay);
425 }
426}
427
428//-----------------------------------------------------------------------------
429
430G4int G4VLongitudinalStringDecay::SampleQuarkFlavor(void)
431 {
432 return (1 + (int)(G4UniformRand()/StrangeSuppress));
433 }
434
435//-----------------------------------------------------------------------------
436
437G4VLongitudinalStringDecay::pDefPair G4VLongitudinalStringDecay::CreatePartonPair(G4int NeedParticle,G4bool AllowDiquarks)
438{
439// NeedParticle = +1 for Particle, -1 for Antiparticle
440
441 if ( AllowDiquarks && G4UniformRand() < DiquarkSuppress )
442 {
443 // Create a Diquark - AntiDiquark pair , first in pair is anti to IsParticle
444 G4int q1 = SampleQuarkFlavor();
445 G4int q2 = SampleQuarkFlavor();
446 G4int spin = (q1 != q2 && G4UniformRand() <= 0.5)? 1 : 3;
447 // convention: quark with higher PDG number is first
448 G4int PDGcode = (std::max(q1,q2) * 1000 + std::min(q1,q2) * 100 + spin) * NeedParticle;
449 return pDefPair (FindParticle(-PDGcode),FindParticle(PDGcode));
450
451
452 } else {
453 // Create a Quark - AntiQuark pair, first in pair IsParticle
454 G4int PDGcode=SampleQuarkFlavor()*NeedParticle;
455 return pDefPair (FindParticle(PDGcode),FindParticle(-PDGcode));
456 }
457
458}
459
460//-----------------------------------------------------------------------------
461G4ThreeVector G4VLongitudinalStringDecay::SampleQuarkPt(G4double ptMax)
462 {
463 G4double Pt;
464 if ( ptMax < 0 ) {
465 // sample full gaussian
466 Pt = -std::log(G4UniformRand());
467 } else {
468 // sample in limited range
469 Pt = -std::log(CLHEP::RandFlat::shoot(std::exp(-sqr(ptMax)/sqr(SigmaQT)), 1.));
470 }
471 Pt = SigmaQT * std::sqrt(Pt);
472 G4double phi = 2.*pi*G4UniformRand();
473 return G4ThreeVector(Pt * std::cos(phi),Pt * std::sin(phi),0);
474 }
475
476//******************************************************************************
477
478void G4VLongitudinalStringDecay::CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector* Hadrons)
479 {
480
481// `yo-yo` formation time
482// const G4double kappa = 1.0 * GeV/fermi/4.;
483 G4double kappa = GetStringTensionParameter();
484 for(size_t c1 = 0; c1 < Hadrons->size(); c1++)
485 {
486 G4double SumPz = 0;
487 G4double SumE = 0;
488 for(size_t c2 = 0; c2 < c1; c2++)
489 {
490 SumPz += Hadrons->operator[](c2)->Get4Momentum().pz();
491 SumE += Hadrons->operator[](c2)->Get4Momentum().e();
492 }
493 G4double HadronE = Hadrons->operator[](c1)->Get4Momentum().e();
494 G4double HadronPz = Hadrons->operator[](c1)->Get4Momentum().pz();
495 Hadrons->operator[](c1)->SetFormationTime(
496(theInitialStringMass - 2.*SumPz + HadronE - HadronPz)/(2.*kappa)/c_light);
497
498 G4ThreeVector aPosition(0, 0,
499(theInitialStringMass - 2.*SumE - HadronE + HadronPz)/(2.*kappa));
500 Hadrons->operator[](c1)->SetPosition(aPosition);
501
502 }
503 }
504
505//-----------------------------------------------------------------------------
506
507void G4VLongitudinalStringDecay::SetSigmaTransverseMomentum(G4double aValue)
508{
509 if ( PastInitPhase ) {
510 throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetSigmaTransverseMomentum after FragmentString() not allowed");
511 } else {
512 SigmaQT = aValue;
513 }
514}
515
516//----------------------------------------------------------------------------------------------------------
517
518void G4VLongitudinalStringDecay::SetStrangenessSuppression(G4double aValue)
519{
520 if ( PastInitPhase ) {
521 throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetStrangenessSuppression after FragmentString() not allowed");
522 } else {
523 StrangeSuppress = aValue;
524 }
525}
526
527//----------------------------------------------------------------------------------------------------------
528
529void G4VLongitudinalStringDecay::SetDiquarkSuppression(G4double aValue)
530{
531 if ( PastInitPhase ) {
532 throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkSuppression after FragmentString() not allowed");
533 } else {
534 DiquarkSuppress = aValue;
535 }
536}
537
538//----------------------------------------------------------------------------------------
539
540void G4VLongitudinalStringDecay::SetDiquarkBreakProbability(G4double aValue)
541{
542 if ( PastInitPhase ) {
543 throw G4HadronicException(__FILE__, __LINE__, "4VLongitudinalStringDecay::SetDiquarkBreakProbability after FragmentString() not allowed");
544 } else {
545 DiquarkBreakProb = aValue;
546 }
547}
548
549//----------------------------------------------------------------------------------------------------------
550
551void G4VLongitudinalStringDecay::SetVectorMesonProbability(G4double aValue)
552{
553 if ( PastInitPhase ) {
554 throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonProbability after FragmentString() not allowed");
555 } else {
556 pspin_meson = aValue;
557 delete hadronizer;
558 hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
559 scalarMesonMix,vectorMesonMix);
560 }
561}
562
563//----------------------------------------------------------------------------------------------------------
564
565void G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability(G4double aValue)
566{
567 if ( PastInitPhase ) {
568 throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetSpinThreeHalfBarionProbability after FragmentString() not allowed");
569 } else {
570 pspin_barion = aValue;
571 delete hadronizer;
572 hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
573 scalarMesonMix,vectorMesonMix);
574 }
575}
576
577//----------------------------------------------------------------------------------------------------------
578
579void G4VLongitudinalStringDecay::SetScalarMesonMixings(std::vector<G4double> aVector)
580{
581 if ( PastInitPhase ) {
582 throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings after FragmentString() not allowed");
583 } else {
584 if ( aVector.size() < 6 )
585 throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetScalarMesonMixings( argument Vector too small");
586 scalarMesonMix[0] = aVector[0];
587 scalarMesonMix[1] = aVector[1];
588 scalarMesonMix[2] = aVector[2];
589 scalarMesonMix[3] = aVector[3];
590 scalarMesonMix[4] = aVector[4];
591 scalarMesonMix[5] = aVector[5];
592 delete hadronizer;
593 hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
594 scalarMesonMix,vectorMesonMix);
595 }
596}
597
598//----------------------------------------------------------------------------------------------------------
599
600void G4VLongitudinalStringDecay::SetVectorMesonMixings(std::vector<G4double> aVector)
601{
602 if ( PastInitPhase ) {
603 throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings after FragmentString() not allowed");
604 } else {
605 if ( aVector.size() < 6 )
606 throw G4HadronicException(__FILE__, __LINE__, "G4VLongitudinalStringDecay::SetVectorMesonMixings( argument Vector too small");
607 vectorMesonMix[0] = aVector[0];
608 vectorMesonMix[1] = aVector[1];
609 vectorMesonMix[2] = aVector[2];
610 vectorMesonMix[3] = aVector[3];
611 vectorMesonMix[4] = aVector[4];
612 vectorMesonMix[5] = aVector[5];
613 delete hadronizer;
614 hadronizer = new G4HadronBuilder(pspin_meson,pspin_barion,
615 scalarMesonMix,vectorMesonMix);
616
617 }
618}
619
620//-------------------------------------------------------------------------------------------
621void G4VLongitudinalStringDecay::SetStringTensionParameter(G4double aValue)// Uzhi 20 June 08
622{
623 Kappa = aValue * GeV/fermi;
624}
625//**************************************************************************************
626
Note: See TracBrowser for help on using the repository browser.