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

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

update ti head

File size: 8.5 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.8 2010/09/20 12:46:23 vuzhinsk Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
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
55   int operator==(const G4VLongitudinalStringDecay &right) const;
56   int operator!=(const G4VLongitudinalStringDecay &right) const;
57
58public:
59   virtual G4KineticTrackVector* FragmentString(const G4ExcitedString& theString)=0;
60
61protected: 
62
63// For changing Mass Cut used for selection of very small mass strings
64   virtual void SetMassCut(G4double aValue);
65
66// For handling a string with very low mass
67   G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString * const theString);
68
69// To store created quarks or 2 last hadrons
70   typedef std::pair<G4ParticleDefinition*, G4ParticleDefinition*> pDefPair;
71
72// For creation of hadrons from given quark pair
73   typedef G4ParticleDefinition * (G4HadronBuilder::*Pcreate)
74                                        (G4ParticleDefinition*, G4ParticleDefinition*);
75
76//-----------------------------------------------------------------------------
77// Used by LightFragmentationTest for estimation of lowest possible mass of
78// given quark system
79   G4double FragmentationMass(const G4FragmentingString * const string,
80                              Pcreate build=0,
81                              pDefPair * pdefs=0);
82
83   G4ParticleDefinition* FindParticle(G4int Encoding); 
84
85   virtual void Sample4Momentum(G4LorentzVector* Mom,     G4double Mass, 
86                                G4LorentzVector* AntiMom, G4double AntiMass, 
87                                G4double InitialMass)=0; 
88//-----------------------------------------------------------------------------
89// For decision on continue or stop string fragmentation
90   virtual G4bool StopFragmenting(const G4FragmentingString  * const string)=0;
91   virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0;
92
93// If a string can not fragment, make last break into 2 hadrons
94   virtual G4bool SplitLast(G4FragmentingString * string, 
95                    G4KineticTrackVector * LeftVector,
96                    G4KineticTrackVector * RightVector)=0;
97//-----------------------------------------------------------------------------
98
99// If a string fragments, do the following
100
101// For transver of a string to its CMS frame
102   G4ExcitedString *CPExcited(const G4ExcitedString& string);
103
104   G4KineticTrack * Splitup(G4FragmentingString *string, 
105                            G4FragmentingString *&newString);
106
107   G4ParticleDefinition * QuarkSplitup(G4ParticleDefinition* decay,
108                                       G4ParticleDefinition *&created);
109
110   G4ParticleDefinition * DiQuarkSplitup(G4ParticleDefinition* decay,
111                                         G4ParticleDefinition *&created);
112                                       
113   pDefPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true);
114
115public:
116//   used by G4VKinkyStringDecy..
117   G4int SampleQuarkFlavor(void);
118   G4ThreeVector SampleQuarkPt(G4double ptMax=-1.); // -1. no limit on maxpt.
119
120protected:
121
122//-----------------------------------------------------------------------------
123// For determination of kinematical properties of created hadron
124//   virtual G4LorentzVector * SplitEandP(G4ParticleDefinition * pHadron,    // Uzhi
125//                                        G4FragmentingString * string  )=0; // Uzhi
126
127   virtual G4LorentzVector * SplitEandP(G4ParticleDefinition * pHadron,      // Uzhi
128                                        G4FragmentingString * string,        // Uzhi
129                                        G4FragmentingString * newString  )=0;// Uzhi
130
131   virtual G4double GetLightConeZ(G4double zmin, G4double zmax, 
132                                  G4int PartonEncoding, 
133                                  G4ParticleDefinition* pHadron, 
134                                  G4double Px, G4double Py       ) = 0;     
135
136   void CalculateHadronTimePosition(G4double theInitialStringMass, 
137                                    G4KineticTrackVector *);
138
139// Used for some test purposes ------------------------------------------------
140   void ConstructParticle();
141
142   G4ParticleDefinition* CreateHadron(G4int id1, G4int id2, 
143                                      G4bool theGivenSpin, G4int theSpin); 
144   
145//-----------------------------------------------------------------------------
146public:
147
148   G4KineticTrackVector* DecayResonans (G4KineticTrackVector* aHadrons);
149
150   void SetSigmaTransverseMomentum(G4double aQT);
151   void SetStrangenessSuppression(G4double aValue);
152   void SetDiquarkSuppression(G4double aValue);
153   void SetDiquarkBreakProbability(G4double aValue);
154
155   void SetVectorMesonProbability(G4double aValue);
156   void SetSpinThreeHalfBarionProbability(G4double aValue);
157   
158   void SetScalarMesonMixings( std::vector<G4double> aVector);
159   void SetVectorMesonMixings( std::vector<G4double> aVector);
160
161   void SetStringTensionParameter(G4double aValue);            // Uzhi 20 June 08
162
163//private:
164protected: 
165   G4double GetDiquarkSuppress()        {return DiquarkSuppress;};
166   G4double GetDiquarkBreakProb()       {return DiquarkBreakProb;};
167   G4double GetStrangeSuppress()        {return StrangeSuppress;};
168   G4double GetClusterMass()            {return ClusterMass;};
169   G4int    GetClusterLoopInterrupt()   {return ClusterLoopInterrupt;};
170
171   G4double GetStringTensionParameter() {return Kappa;};       // Uzhi 20 June 08
172   
173//private:
174protected: 
175   G4double  MassCut;
176   G4double  ClusterMass;
177   G4double  SigmaQT;          // sigma_q_t is quark transverse momentum distribution parameter
178   G4double  DiquarkSuppress;  // is Diquark suppression parameter 
179   G4double  DiquarkBreakProb; // is Diquark breaking probability
180   G4double  SmoothParam;      // model parameter
181   G4double  StrangeSuppress ;
182   G4int     StringLoopInterrupt;
183   G4int     ClusterLoopInterrupt;
184   G4int     SideOfDecay;
185   G4HadronBuilder *hadronizer;
186
187   G4double pspin_meson;
188   G4double pspin_barion;
189   std::vector<G4double> vectorMesonMix;
190   std::vector<G4double> scalarMesonMix;
191   
192   G4bool    PastInitPhase;
193
194   G4double Kappa; // String tension parameter                 // Uzhi 20 June 08
195
196//   G4double MinFragmentationMass(G4ExcitedString * theString,
197//                              G4ParticleDefinition*& Hadron1,
198//                              G4ParticleDefinition*& Hadron2);
199
200};
201
202//*************************************************************************************
203// Class G4VLongitudinalStringDecay
204#endif
Note: See TracBrowser for help on using the repository browser.