Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4VLongitudinalStringDecay.hh
이 파일의 문서화 페이지로 가기
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 107869 2017-12-07 14:46:39Z gcosmo $
28 // Maxim Komogorov
29 //
30 // -----------------------------------------------------------------------------
31 // GEANT 4 class implementation file
32 //
33 // History: first implementation, Maxim Komogorov, 1-Jul-1998
34 // -----------------------------------------------------------------------------
35 #ifndef G4VLongitudinalStringDecay_h
36 #define G4VLongitudinalStringDecay_h 1
38 #include "G4DynamicParticle.hh"
39 #include "G4KineticTrack.hh"
40 #include "G4KineticTrackVector.hh"
41 #include "G4HadronBuilder.hh"
42 
44 //**************************************************************************************
45 
47  {
48 public:
51 
52 private:
53  // not implemented to protect/forbid use
56  int operator==(const G4VLongitudinalStringDecay &right) const;
57  int operator!=(const G4VLongitudinalStringDecay &right) const;
58 
59 public:
60  virtual G4KineticTrackVector* FragmentString(const G4ExcitedString& theString)=0;
61 
62 protected:
63 
64 // For changing Mass Cut used for selection of very small mass strings
65  virtual void SetMassCut(G4double aValue);
67 
68 // For handling a string with very low mass
70 
71 // To store created quarks or 2 last hadrons
72  typedef std::pair<G4ParticleDefinition*, G4ParticleDefinition*> pDefPair;
73 
74 // For creation of hadrons from given quark pair
77 
78 //-----------------------------------------------------------------------------
79 // Used by LightFragmentationTest for estimation of lowest possible mass of
80 // given quark system
81  G4double FragmentationMass(const G4FragmentingString * const string,
82  Pcreate build=0,
83  pDefPair * pdefs=0);
84 
86 
87  virtual void Sample4Momentum(G4LorentzVector* Mom, G4double Mass,
88  G4LorentzVector* AntiMom, G4double AntiMass,
89  G4double InitialMass)=0;
90 //-----------------------------------------------------------------------------
91 // For decision on continue or stop string fragmentation
92  virtual G4bool StopFragmenting(const G4FragmentingString * const string)=0;
93  virtual G4bool IsFragmentable(const G4FragmentingString * const string)=0;
94 
95 // If a string can not fragment, make last break into 2 hadrons
96  virtual G4bool SplitLast(G4FragmentingString * string,
97  G4KineticTrackVector * LeftVector,
98  G4KineticTrackVector * RightVector)=0;
99 //-----------------------------------------------------------------------------
100 
101 // If a string fragments, do the following
102 
103 // To make a copy of a string
105 
107  G4ParticleDefinition *&created);
108 
110  G4ParticleDefinition *&created)=0;
111 
112  pDefPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true);
113 
114 public:
115 // used by G4VKinkyStringDecy..
116  G4int SampleQuarkFlavor(void);
117  G4ThreeVector SampleQuarkPt(G4double ptMax=-1.); // -1. no limit on maxpt.
118 
119 protected:
120 
121 //-----------------------------------------------------------------------------
122 // For determination of kinematical properties of created hadron
123 // virtual G4LorentzVector * SplitEandP(G4ParticleDefinition * pHadron,
124 // G4FragmentingString * string )=0;
125  virtual G4KineticTrack * Splitup(G4FragmentingString *string, // Uzhi 28 June 2016
126  G4FragmentingString *&newString)=0;
127 
128  virtual G4LorentzVector * SplitEandP(G4ParticleDefinition * pHadron,
129  G4FragmentingString * string,
130  G4FragmentingString * newString )=0;
131 
132  virtual G4double GetLightConeZ(G4double zmin, G4double zmax,
133  G4int PartonEncoding,
134  G4ParticleDefinition* pHadron,
135  G4double Px, G4double Py ) = 0;
136 
137  void CalculateHadronTimePosition(G4double theInitialStringMass,
139 
140 // Used for some test purposes ------------------------------------------------
141  void ConstructParticle();
142 
144  G4bool theGivenSpin, G4int theSpin);
145 
146 //-----------------------------------------------------------------------------
147 public:
148 
150 
152  void SetStrangenessSuppression(G4double aValue);
153  void SetDiquarkSuppression(G4double aValue);
155 
156  void SetVectorMesonProbability(G4double aValue);
158 
159  void SetScalarMesonMixings( std::vector<G4double> aVector);
160  void SetVectorMesonMixings( std::vector<G4double> aVector);
161 
162  void SetStringTensionParameter(G4double aValue);
163 
164 // G4LorentzVector RestMomentum; // Uzhi June 2016
165 //private:
166 protected:
172 
174 
175 //private:
176 protected:
179  G4double SigmaQT; // sigma_q_t is quark transverse momentum distribution parameter
180  G4double DiquarkSuppress; // is Diquark suppression parameter
181  G4double DiquarkBreakProb; // is Diquark breaking probability
182  G4double SmoothParam; // model parameter
187 
190  std::vector<G4double> vectorMesonMix;
191  std::vector<G4double> scalarMesonMix;
192 
194 
195  G4double Kappa; // String tension parameter
196 
197 // G4double MinFragmentationMass(G4ExcitedString * theString,
198 // G4ParticleDefinition*& Hadron1,
199 // G4ParticleDefinition*& Hadron2);
200 
201 };
202 
203 //*************************************************************************************
204 // Class G4VLongitudinalStringDecay
205 #endif
ParticleList decay(Cluster *const c)
Carries out a cluster decay.
G4KineticTrackVector * DecayResonans(G4KineticTrackVector *aHadrons)
void SetDiquarkSuppression(G4double aValue)
virtual void SetMassCut(G4double aValue)
G4ParticleDefinition * CreateHadron(G4int id1, G4int id2, G4bool theGivenSpin, G4int theSpin)
void SetStringTensionParameter(G4double aValue)
void SetVectorMesonMixings(std::vector< G4double > aVector)
G4double FragmentationMass(const G4FragmentingString *const string, Pcreate build=0, pDefPair *pdefs=0)
void SetScalarMesonMixings(std::vector< G4double > aVector)
int operator==(const G4VLongitudinalStringDecay &right) const
G4ParticleDefinition *(G4HadronBuilder::* Pcreate)(G4ParticleDefinition *, G4ParticleDefinition *)
void SetSpinThreeHalfBarionProbability(G4double aValue)
void SetDiquarkBreakProbability(G4double aValue)
std::vector< G4double > vectorMesonMix
virtual G4bool SplitLast(G4FragmentingString *string, G4KineticTrackVector *LeftVector, G4KineticTrackVector *RightVector)=0
double G4double
Definition: G4Types.hh:76
bool G4bool
Definition: G4Types.hh:79
void SetVectorMesonProbability(G4double aValue)
G4ParticleDefinition * FindParticle(G4int Encoding)
G4ThreeVector SampleQuarkPt(G4double ptMax=-1.)
virtual G4bool StopFragmenting(const G4FragmentingString *const string)=0
void CalculateHadronTimePosition(G4double theInitialStringMass, G4KineticTrackVector *)
virtual G4double GetLightConeZ(G4double zmin, G4double zmax, G4int PartonEncoding, G4ParticleDefinition *pHadron, G4double Px, G4double Py)=0
virtual G4ParticleDefinition * DiQuarkSplitup(G4ParticleDefinition *decay, G4ParticleDefinition *&created)=0
const G4VLongitudinalStringDecay & operator=(const G4VLongitudinalStringDecay &right)
G4KineticTrackVector * LightFragmentationTest(const G4ExcitedString *const theString)
void SetStrangenessSuppression(G4double aValue)
int operator!=(const G4VLongitudinalStringDecay &right) const
G4ParticleDefinition * QuarkSplitup(G4ParticleDefinition *decay, G4ParticleDefinition *&created)
virtual G4bool IsFragmentable(const G4FragmentingString *const string)=0
virtual void Sample4Momentum(G4LorentzVector *Mom, G4double Mass, G4LorentzVector *AntiMom, G4double AntiMass, G4double InitialMass)=0
int G4int
Definition: G4Types.hh:78
virtual G4LorentzVector * SplitEandP(G4ParticleDefinition *pHadron, G4FragmentingString *string, G4FragmentingString *newString)=0
pDefPair CreatePartonPair(G4int NeedParticle, G4bool AllowDiquarks=true)
virtual G4KineticTrack * Splitup(G4FragmentingString *string, G4FragmentingString *&newString)=0
G4ExcitedString * CopyExcited(const G4ExcitedString &string)
std::pair< G4ParticleDefinition *, G4ParticleDefinition * > pDefPair
std::vector< G4double > scalarMesonMix
virtual G4KineticTrackVector * FragmentString(const G4ExcitedString &theString)=0