source: trunk/source/processes/hadronic/models/parton_string/qgsm/include/G4QGSParticipants.hh @ 1198

Last change on this file since 1198 was 1055, checked in by garnier, 15 years ago

maj sur la beta de geant 4.9.3

File size: 4.2 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 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"
34#include "G4QGSDiffractiveExcitation.hh"
35#include "G4SingleDiffractiveExcitation.hh"
36#include "G4PartonPair.hh"
37#include "G4QGSMSplitableHadron.hh"
38
39class G4QGSParticipants : public G4VParticipants
40{
41  public:
42    G4QGSParticipants();
43    G4QGSParticipants(const G4QGSParticipants &right);
44    const G4QGSParticipants & operator=(const G4QGSParticipants &right);
45    virtual ~G4QGSParticipants(); 
46
47    int operator==(const G4QGSParticipants &right) const;
48    int operator!=(const G4QGSParticipants &right) const;
49
50    virtual void DoLorentzBoost(G4ThreeVector aBoost) 
51    {
52      if(theNucleus) theNucleus->DoLorentzBoost(aBoost);
53      theBoost = aBoost;
54    }
55
56    G4PartonPair* GetNextPartonPair();
57    void BuildInteractions(const G4ReactionProduct  &thePrimary);
58    void StartPartonPairLoop();
59
60  protected:
61    virtual G4VSplitableHadron* SelectInteractions(const G4ReactionProduct  &thePrimary);
62        void SplitHadrons(); 
63    void PerformSoftCollisions();
64    void PerformDiffractiveCollisions();
65     
66  protected:
67    struct DeleteInteractionContent {void operator()(G4InteractionContent*aC){delete aC;}};
68    std::vector<G4InteractionContent*> theInteractions;
69    struct DeleteSplitableHadron{void operator()(G4VSplitableHadron*aS){delete aS;}};
70    std::vector<G4VSplitableHadron*>   theTargets; 
71    struct DeletePartonPair{void operator()(G4PartonPair*aP){delete aP;}};
72    std::vector<G4PartonPair*>   thePartonPairs;
73 
74    G4SingleDiffractiveExcitation theSingleDiffExcitation;
75    G4QGSDiffractiveExcitation theDiffExcitaton;
76    G4int ModelMode;
77    G4bool IsSingleDiffractive();
78 
79    G4ThreeVector theBoost;
80
81  protected:
82    // model parameters HPW
83    enum  { SOFT, DIFFRACTIVE };
84    const G4int nCutMax; 
85    const G4double ThresholdParameter; 
86    const G4double QGSMThreshold; 
87    const G4double theNucleonRadius;
88   
89};
90
91
92inline G4bool G4QGSParticipants::IsSingleDiffractive()
93{
94  G4bool result=false;
95  if(G4UniformRand()<1.) result = true;
96  return result;
97}
98
99inline void G4QGSParticipants::StartPartonPairLoop()
100{
101}
102
103inline G4PartonPair* G4QGSParticipants::GetNextPartonPair()
104{
105  if (thePartonPairs.empty()) return 0;
106  G4PartonPair * result = thePartonPairs.back();
107  thePartonPairs.pop_back();
108  return result;
109}
110
111
112inline void G4QGSParticipants::SplitHadrons()
113{
114  unsigned int i;
115  for(i = 0; i < theInteractions.size(); i++) 
116  {
117    theInteractions[i]->SplitHadrons();
118  }
119}
120
121#endif
122
123
Note: See TracBrowser for help on using the repository browser.