Geant4-11
G4DiffractiveExcitation.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
28#ifndef G4DiffractiveExcitation_h
29#define G4DiffractiveExcitation_h 1
30
31// ------------------------------------------------------------
32// GEANT 4 class header file
33//
34// ---------------- G4DiffractiveExcitation --------------
35// by Gunter Folger, October 1998.
36// diffractive Excitation used by strings models
37// Take a projectile and a target
38// excite the projectile and target
39// ------------------------------------------------------------
40
41#include "globals.hh"
42#include "G4FTFParameters.hh"
44#include "G4ThreeVector.hh"
45#include "G4LorentzVector.hh"
46#include "G4LorentzRotation.hh"
47
49class G4ExcitedString;
50
51
53 public:
56
58 G4VSplitableHadron* bPartner,
59 G4FTFParameters* theParameters,
60 G4ElasticHNScattering* theElastic ) const;
61
62 virtual void CreateStrings( G4VSplitableHadron* aHadron,
63 G4bool isProjectile,
64 G4ExcitedString*& FirstString,
65 G4ExcitedString*& SecondString,
66 G4FTFParameters* theParameters ) const;
67
68 private:
71 G4bool operator==( const G4DiffractiveExcitation& right ) const;
72 G4bool operator!=( const G4DiffractiveExcitation& right ) const;
73
74 G4double LambdaF(G4double sqrM, G4double sqrM1, G4double sqrM2) const;
75
76 G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
79 void UnpackMeson( G4int IdPDG, G4int& Q1, G4int& Q2 ) const;
80 void UnpackBaryon( G4int IdPDG, G4int& Q1, G4int& Q2, G4int& Q3 ) const;
81 G4int NewNucleonId( G4int Q1, G4int Q2, G4int Q3 ) const;
82
83 // The "ExciteParticipants" method uses the following struct and 3 new utility methods:
88 ProjMassT = 0.0, ProjMassT2 = 0.0, TargMassT = 0.0, TargMassT2 = 0.0,
94 S = 0.0, SqrtS = 0.0, Pt2 = 0.0, PZcms = 0.0, PZcms2 = 0.0,
95 AveragePt2 = 0.0, maxPtSquare = 0.0,
96 ProbExc = 0.0, Qminus = 0.0, Qplus = 0.0,
97 PMinusNew = 0.0, PPlusNew = 0.0, TMinusNew = 0.0, TPlusNew = 0.0,
98 PMinusMin = 0.0, PMinusMax = 0.0, TPlusMin = 0.0, TPlusMax = 0.0,
103 };
105 G4VSplitableHadron* target,
106 G4FTFParameters* theParameters,
107 G4ElasticHNScattering* theElastic,
108 CommonVariables& common ) const;
110 G4VSplitableHadron* target,
111 G4FTFParameters* theParameters,
112 CommonVariables& common ) const;
114 G4VSplitableHadron* target,
115 G4FTFParameters* theParameters,
116 CommonVariables& common ) const;
117};
118
119#endif
120
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4int NewNucleonId(G4int Q1, G4int Q2, G4int Q3) const
void UnpackMeson(G4int IdPDG, G4int &Q1, G4int &Q2) const
G4double ChooseP(G4double Pmin, G4double Pmax) const
G4double LambdaF(G4double sqrM, G4double sqrM1, G4double sqrM2) const
G4ThreeVector GaussianPt(G4double AveragePt2, G4double maxPtSquare) const
G4bool operator!=(const G4DiffractiveExcitation &right) const
G4double GetQuarkFractionOfKink(G4double zmin, G4double zmax) const
virtual G4bool ExciteParticipants(G4VSplitableHadron *aPartner, G4VSplitableHadron *bPartner, G4FTFParameters *theParameters, G4ElasticHNScattering *theElastic) const
void UnpackBaryon(G4int IdPDG, G4int &Q1, G4int &Q2, G4int &Q3) const
G4bool ExciteParticipants_doDiffraction(G4VSplitableHadron *projectile, G4VSplitableHadron *target, G4FTFParameters *theParameters, CommonVariables &common) const
G4int ExciteParticipants_doChargeExchange(G4VSplitableHadron *projectile, G4VSplitableHadron *target, G4FTFParameters *theParameters, G4ElasticHNScattering *theElastic, CommonVariables &common) const
G4bool ExciteParticipants_doNonDiffraction(G4VSplitableHadron *projectile, G4VSplitableHadron *target, G4FTFParameters *theParameters, CommonVariables &common) const
virtual void CreateStrings(G4VSplitableHadron *aHadron, G4bool isProjectile, G4ExcitedString *&FirstString, G4ExcitedString *&SecondString, G4FTFParameters *theParameters) const
const G4DiffractiveExcitation & operator=(const G4DiffractiveExcitation &right)
G4bool operator==(const G4DiffractiveExcitation &right) const
static int FASTCALL common(PROLOG_STATE *state, int tok)
Definition: xmlrole.cc:1309