Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4QGSParticipants.hh
Go to the documentation of this file.
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 G4QGSParticipants_h
27 #define G4QGSParticipants_h 1
28 
29 #include "Randomize.hh"
30 #include "G4VParticipants.hh"
31 #include "G4Nucleon.hh"
32 #include "G4InteractionContent.hh"
33 #include "G4PomeronCrossSection.hh"
36 #include "G4PartonPair.hh"
37 #include "G4QGSMSplitableHadron.hh"
38 #include "G4V3DNucleus.hh"
39 
40 
42 {
43  public:
46  const G4QGSParticipants & operator=(const G4QGSParticipants &right);
47  virtual ~G4QGSParticipants();
48 
49  int operator==(const G4QGSParticipants &right) const;
50  int operator!=(const G4QGSParticipants &right) const;
51 
52  virtual void DoLorentzBoost(G4ThreeVector aBoost)
53  {
55  theBoost = aBoost;
56  }
57 
59  void BuildInteractions(const G4ReactionProduct &thePrimary);
60  void StartPartonPairLoop();
61 
62  protected:
63  virtual G4VSplitableHadron* SelectInteractions(const G4ReactionProduct &thePrimary);
64  void SplitHadrons();
65  void PerformSoftCollisions();
67 
68  protected:
70  std::vector<G4InteractionContent*> theInteractions;
71  struct DeleteSplitableHadron{void operator()(G4VSplitableHadron*aS){delete aS;}};
72  std::vector<G4VSplitableHadron*> theTargets;
73  struct DeletePartonPair{void operator()(G4PartonPair*aP){delete aP;}};
74  std::vector<G4PartonPair*> thePartonPairs;
75 
80 
82 
83  protected:
84  // model parameters HPW
85  enum { SOFT, DIFFRACTIVE };
86  const G4int nCutMax;
90 };
91 
92 
94 {
95  G4bool result=false;
96  if(G4UniformRand()<1.) result = true;
97  return result;
98 }
99 
101 {
102 }
103 
105 {
106  if (thePartonPairs.empty()) return 0;
108  thePartonPairs.pop_back();
109  return result;
110 }
111 
113 {
114  unsigned int i;
115  for(i = 0; i < theInteractions.size(); i++)
116  {
117  theInteractions[i]->SplitHadrons();
118  }
119 }
120 
121 #endif
122 
G4double G4ParticleHPJENDLHEData::G4double result
int operator==(const G4QGSParticipants &right) const
int operator!=(const G4QGSParticipants &right) const
const G4double QGSMThreshold
int G4int
Definition: G4Types.hh:78
void operator()(G4VSplitableHadron *aS)
virtual ~G4QGSParticipants()
#define G4UniformRand()
Definition: Randomize.hh:97
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
bool G4bool
Definition: G4Types.hh:79
G4QGSDiffractiveExcitation theDiffExcitaton
G4V3DNucleus * theNucleus
const G4QGSParticipants & operator=(const G4QGSParticipants &right)
G4ThreeVector theBoost
std::vector< G4PartonPair * > thePartonPairs
virtual void DoLorentzBoost(const G4LorentzVector &theBoost)=0
std::vector< G4VSplitableHadron * > theTargets
void BuildInteractions(const G4ReactionProduct &thePrimary)
void operator()(G4InteractionContent *aC)
double G4double
Definition: G4Types.hh:76
G4SingleDiffractiveExcitation theSingleDiffExcitation
std::vector< G4InteractionContent * > theInteractions
const G4double theNucleonRadius
const G4double ThresholdParameter
virtual void DoLorentzBoost(G4ThreeVector aBoost)
G4PartonPair * GetNextPartonPair()