Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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"
36#include "G4PartonPair.hh"
38#include "G4V3DNucleus.hh"
39
40
42{
43public:
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);
61
62protected:
64 void SplitHadrons();
67
68protected:
70 std::vector<G4InteractionContent*> theInteractions;
72 std::vector<G4VSplitableHadron*> theTargets;
73 struct DeletePartonPair{void operator()(G4PartonPair*aP){delete aP;}};
74 std::vector<G4PartonPair*> thePartonPairs;
75
80
82
83protected:
84 // model parameters HPW
85 enum { SOFT, DIFFRACTIVE };
90
91};
92
93
95{
96 G4bool result=false;
97 if(G4UniformRand()<1.) result = true;
98 return result;
99}
100
102{
103}
104
106{
107 if (thePartonPairs.empty()) return 0;
108 G4PartonPair * result = thePartonPairs.back();
109 thePartonPairs.pop_back();
110 return result;
111}
112
113
115{
116 unsigned int i;
117 for(i = 0; i < theInteractions.size(); i++)
118 {
119 theInteractions[i]->SplitHadrons();
120 }
121}
122
123#endif
124
125
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
#define G4UniformRand()
Definition: Randomize.hh:53
const G4double ThresholdParameter
G4SingleDiffractiveExcitation theSingleDiffExcitation
virtual void DoLorentzBoost(G4ThreeVector aBoost)
G4ThreeVector theBoost
std::vector< G4InteractionContent * > theInteractions
G4PartonPair * GetNextPartonPair()
virtual G4VSplitableHadron * SelectInteractions(const G4ReactionProduct &thePrimary)
G4QGSDiffractiveExcitation theDiffExcitaton
virtual ~G4QGSParticipants()
const G4double QGSMThreshold
const G4double theNucleonRadius
void BuildInteractions(const G4ReactionProduct &thePrimary)
int operator==(const G4QGSParticipants &right) const
const G4QGSParticipants & operator=(const G4QGSParticipants &right)
int operator!=(const G4QGSParticipants &right) const
std::vector< G4PartonPair * > thePartonPairs
std::vector< G4VSplitableHadron * > theTargets
virtual void DoLorentzBoost(const G4LorentzVector &theBoost)=0
G4V3DNucleus * theNucleus
void operator()(G4InteractionContent *aC)
void operator()(G4VSplitableHadron *aS)