source: trunk/source/processes/hadronic/models/parton_string/hadronization/include/G4VLongitudinalStringDecay.hh@ 880

Last change on this file since 880 was 819, checked in by garnier, 17 years ago

import all except CVS

File size: 6.7 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.hh,v 1.4 2007/04/24 14:55:23 gunter Exp $
28// GEANT4 tag $Name: geant4-09-01-patch-02 $
29// Maxim Komogorov
30//
31// -----------------------------------------------------------------------------
32// GEANT 4 class implementation file
33//
34// History: first implementation, Maxim Komogorov, 1-Jul-1998
35// -----------------------------------------------------------------------------
36#ifndef G4VLongitudinalStringDecay_h
37#define G4VLongitudinalStringDecay_h 1
38#include "G4VStringFragmentation.hh"
39#include "G4DynamicParticle.hh"
40#include "G4KineticTrack.hh"
41#include "G4KineticTrackVector.hh"
42#include "G4HadronBuilder.hh"
43
44class G4FragmentingString;
45//**********************************************************************************************
46
47class G4VLongitudinalStringDecay
48 {
49public:
50 G4VLongitudinalStringDecay();
51 virtual ~G4VLongitudinalStringDecay();
52
53private:
54// G4VLongitudinalStringDecay(const G4VLongitudinalStringDecay &right);
55// const G4VLongitudinalStringDecay & operator=(const G4VLongitudinalStringDecay &right);
56 int operator==(const G4VLongitudinalStringDecay &right) const;
57 int operator!=(const G4VLongitudinalStringDecay &right) const;
58
59public:
60 virtual G4KineticTrackVector* FragmentString(const G4ExcitedString& theString)=0;
61
62 G4KineticTrackVector* DecayResonans (G4KineticTrackVector* aHadrons);
63 void SetSigmaTransverseMomentum(G4double aQT);
64 void SetStrangenessSuppression(G4double aValue);
65 void SetDiquarkSuppression(G4double aValue);
66 void SetDiquarkBreakProbability(G4double aValue);
67
68 void SetVectorMesonProbability(G4double aValue);
69 void SetSpinThreeHalfBarionProbability(G4double aValue);
70
71 void SetScalarMesonMixings( std::vector<G4double> aVector);
72 void SetVectorMesonMixings( std::vector<G4double> aVector);
73
74// used by G4VKinkyStringDecy..
75 G4int SampleQuarkFlavor(void);
76 G4ThreeVector SampleQuarkPt();
77
78//private:
79protected:
80 G4double GetDiquarkSuppress() {return DiquarkSuppress;};
81 G4double GetDiquarkBreakProb() {return DiquarkBreakProb;};
82 G4double GetStrangeSuppress() {return StrangeSuppress;};
83 G4double GetClusterMass() {return ClusterMass;};
84 G4int GetClusterLoopInterrupt() {return ClusterLoopInterrupt;};
85
86 G4ParticleDefinition* CreateHadron(G4int id1, G4int id2, G4bool theGivenSpin, G4int theSpin);
87 virtual void Sample4Momentum(G4LorentzVector* Mom, G4double Mass, G4LorentzVector* AntiMom, G4double AntiMass, G4double InitialMass)=0;
88
89protected:
90 // Additional protected declarations
91 virtual G4double GetLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding, G4ParticleDefinition* pHadron, G4double Px, G4double Py) = 0;
92
93//private:
94protected:
95 G4double MassCut;
96 G4double ClusterMass;
97 G4double SigmaQT; // sigma_q_t is quark transverse momentum distribution parameter
98 G4double DiquarkSuppress; // is Diquark suppression parameter
99 G4double DiquarkBreakProb; // is Diquark breaking probability
100 G4double SmoothParam; // model parameter
101 G4double StrangeSuppress ;
102 G4int StringLoopInterrupt;
103 G4int ClusterLoopInterrupt;
104 G4int SideOfDecay;
105 G4HadronBuilder *hadronizer;
106
107 void ConstructParticle();
108
109 G4double pspin_meson;
110 G4double pspin_barion;
111 std::vector<G4double> vectorMesonMix;
112 std::vector<G4double> scalarMesonMix;
113
114 G4bool PastInitPhase;
115
116
117 G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString * const theString);
118 virtual G4bool StopFragmenting(const G4FragmentingString * const string)=0;
119 virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0;
120// G4double MinFragmentationMass(G4ExcitedString * theString,
121// G4ParticleDefinition*& Hadron1,
122// G4ParticleDefinition*& Hadron2);
123 typedef std::pair<G4ParticleDefinition*, G4ParticleDefinition*> pDefPair;
124 typedef G4ParticleDefinition * (G4HadronBuilder::*Pcreate)
125 (G4ParticleDefinition*, G4ParticleDefinition*);
126 G4double FragmentationMass(
127 const G4FragmentingString * const string,
128 Pcreate build=0,
129 pDefPair * pdefs=0);
130 G4KineticTrack * Splitup(G4FragmentingString *string, G4FragmentingString *&newString);
131 virtual G4LorentzVector * SplitEandP(G4ParticleDefinition * pHadron, G4FragmentingString * string)=0;
132 virtual G4bool SplitLast(G4FragmentingString * string,
133 G4KineticTrackVector * LeftVector,
134 G4KineticTrackVector * RightVector)=0;
135 void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *);
136 G4ExcitedString *CPExcited(const G4ExcitedString& string);
137 G4ParticleDefinition* FindParticle(G4int Encoding);
138
139 // Additional Implementation Declarations
140 G4ParticleDefinition * QuarkSplitup(G4ParticleDefinition* decay,
141 G4ParticleDefinition *&created);
142 G4ParticleDefinition * DiQuarkSplitup(G4ParticleDefinition* decay,
143 G4ParticleDefinition *&created);
144
145 pDefPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true);
146
147
148};
149
150
151//**********************************************************************************************
152// Class G4VLongitudinalStringDecay
153#endif
154
155
Note: See TracBrowser for help on using the repository browser.