source: trunk/source/processes/hadronic/models/parton_string/qgsm/include/G4QGSMSplitableHadron.hh @ 1055

Last change on this file since 1055 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 4.1 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#ifndef G4QGSMSplitableHadron_h
27#define G4QGSMSplitableHadron_h 1
28
29#include "G4VSplitableHadron.hh"
30#include "G4PartonVector.hh"
31#include "G4MesonSplitter.hh"
32#include "G4BaryonSplitter.hh"
33#include "Randomize.hh"
34#include <deque>
35
36// based on prototype by Maxim Komogorov
37// Splitting into methods, and centralizing of model parameters HPW Feb 1999
38// continued clean-up of interfaces and algorithms HPW 1999.
39// Redesign of data structures and algorithms HPW Feb 1999
40
41class G4QGSMSplitableHadron : public G4VSplitableHadron
42{
43
44  public:
45      G4QGSMSplitableHadron();
46      G4QGSMSplitableHadron(const G4ReactionProduct & aPrimary);
47      G4QGSMSplitableHadron(const G4ReactionProduct & aPrimary, G4bool Direction);
48      G4QGSMSplitableHadron(const G4Nucleon & aNucleon); 
49      G4QGSMSplitableHadron(const G4Nucleon & aNucleon, G4bool Direction); 
50
51      virtual ~G4QGSMSplitableHadron();
52     
53      const G4QGSMSplitableHadron & operator=(const G4QGSMSplitableHadron &right);
54
55      virtual void SplitUp();
56      virtual G4Parton * GetNextParton();
57      virtual G4Parton * GetNextAntiParton();
58
59  private:
60      void InitParameters();
61      void DiffractiveSplitUp();
62      void SoftSplitUp();
63      G4ThreeVector GaussianPt(G4double widthSquare, G4double maxPtSquare);
64      void GetValenceQuarkFlavors(const G4ParticleDefinition * aPart, 
65                                  G4Parton *& Parton1, G4Parton *& Parton2);
66      G4Parton * BuildSeaQuark(G4bool isAntiQuark, G4int aPDGCode, G4int nSeaPair);
67      G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta);
68 
69  private:
70  // aggregated data
71    G4bool Direction; // FALSE is target. - candidate for more detailed design. @@@@ HPW
72
73    std::deque<G4Parton *> Color;
74    std::deque<G4Parton *> AntiColor;   
75  private:
76    // associated classes
77    G4MesonSplitter theMesonSplitter;
78    G4BaryonSplitter theBaryonSplitter;
79 
80  private:
81   // model parameters
82    G4double alpha;
83    G4double beta;
84    G4double theMinPz; 
85    G4double StrangeSuppress;
86    G4double sigmaPt;
87    G4double widthOfPtSquare; 
88    G4double minTransverseMass;
89};
90
91inline G4Parton* G4QGSMSplitableHadron::GetNextParton()
92{
93   if(Color.size()==0) return 0;
94   G4Parton * result = Color.back();
95   Color.pop_back();
96   return result;
97}
98
99inline G4Parton* G4QGSMSplitableHadron::GetNextAntiParton()
100{
101   if(AntiColor.size() == 0) return 0;
102   G4Parton * result = AntiColor.front();
103   AntiColor.pop_front();
104   return result;
105}
106
107#endif
108
109
Note: See TracBrowser for help on using the repository browser.