Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
G4FTFModel.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//
27// $Id$
28// GEANT4 tag $Name: $
29//
30// Class Description
31// Final state production code for hadron inelastic scattering above 20 GeV
32// based on the modeling ansatz used in FRITIOF.
33// To be used in your physics list in case you need this physics.
34// In this case you want to register an object of this class with an object
35// of G4TheoFSGenerator.
36// Class Description - End
37
38#ifndef G4FTFModel_h
39#define G4FTFModel_h 1
40
41// ------------------------------------------------------------
42// GEANT 4 class header file
43//
44// ---------------- G4FTFModel ----------------
45// by Gunter Folger, May 1998.
46// class implementing the excitation in the FTF Parton String Model
47// ------------------------------------------------------------
48
49
51
53class G4ExcitedString;
54
55#include "G4FTFParameters.hh"
56#include "G4FTFParticipants.hh"
57
61#include "G4FTFAnnihilation.hh"
62
64{
65
66 public:
67 G4FTFModel(const G4String& modelName = "FTF");
68 //G4FTFModel(G4double , G4double , G4double );
69 //G4FTFModel(G4DiffractiveExcitation * anExcitation);
71 private:
72 G4FTFModel(const G4FTFModel &right);
73 const G4FTFModel & operator=(const G4FTFModel &right);
74 int operator==(const G4FTFModel &right) const;
75 int operator!=(const G4FTFModel &right) const;
76 public:
77 void Init(const G4Nucleus & aNucleus, const G4DynamicParticle & aProjectile);
80 virtual void ModelDescription(std::ostream&) const;
81
82 protected:
83
84 private:
85 void StoreInvolvedNucleon();
86 void ReggeonCascade();
87 G4bool PutOnMassShell();
88 G4bool ExciteParticipants();
89 G4ExcitedStringVector * BuildStrings();
90 void GetResidualNucleus();
91 void AjustTargetNucleonForAnnihilation(G4VSplitableHadron *SelectedAntiBaryon,
92 G4VSplitableHadron *SelectedTargetNucleon);
93 G4ThreeVector GaussianPt(G4double AveragePt2, G4double maxPtSquare) const;
94
95 private:
96
97 G4ReactionProduct theProjectile;
98 G4FTFParticipants theParticipants;
99
100 G4Nucleon * TheInvolvedNucleon[250];
101 G4int NumberOfInvolvedNucleon;
102 G4int NumberOfInvolvedTargetNucleon;
103
104 G4Nucleon * TheInvolvedNucleonOfProjectile[250];
105 G4int NumberOfInvolvedNucleonOfProjectile;
106
107 G4FTFParameters *theParameters;
108 G4DiffractiveExcitation * theExcitation;
109 G4ElasticHNScattering * theElastic;
110 G4FTFAnnihilation * theAnnihilation; // Uzhi 17.11.10
111
112 std::vector<G4VSplitableHadron *> theAdditionalString; // Uzhi 17.11.10
113
114 G4LorentzVector Residual4Momentum;
115 G4double ResidualExcitationEnergy;
116
117};
118
119// ------------------------------------------------------------
120inline
122{
123 return theParticipants.GetWoundedNucleus();
124}
125
126#endif
std::vector< G4ExcitedString * > G4ExcitedStringVector
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
bool G4bool
Definition: G4Types.hh:67
G4V3DNucleus * GetWoundedNucleus() const
Definition: G4FTFModel.hh:121
void Init(const G4Nucleus &aNucleus, const G4DynamicParticle &aProjectile)
Definition: G4FTFModel.cc:112
virtual void ModelDescription(std::ostream &) const
Definition: G4FTFModel.cc:1745
G4ExcitedStringVector * GetStrings()
Definition: G4FTFModel.cc:206
virtual G4V3DNucleus * GetWoundedNucleus() const