Geant4 Cross Reference |
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" 43 #include "G4ElasticHNScattering.hh" 44 #include "G4ThreeVector.hh" 45 #include "G4LorentzVector.hh" 46 #include "G4LorentzRotation.hh" 47 48 class G4VSplitableHadron; 49 class G4ExcitedString; 50 51 52 class G4DiffractiveExcitation { 53 public: 54 G4DiffractiveExcitation(); 55 virtual ~G4DiffractiveExcitation(); 56 57 virtual G4bool ExciteParticipants( G4VSplitableHadron* aPartner, 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: 69 G4DiffractiveExcitation( const G4DiffractiveExcitation& right ); 70 const G4DiffractiveExcitation& operator=( const G4DiffractiveExcitation& right ); 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; 77 G4double ChooseP( G4double Pmin, G4double Pmax ) const; 78 G4double GetQuarkFractionOfKink( G4double zmin, G4double zmax ) 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: 84 struct CommonVariables { 85 G4int ProjectilePDGcode = 0, absProjectilePDGcode = 0, TargetPDGcode = 0, 86 absTargetPDGcode = 0; 87 G4double M0projectile = 0.0, M0projectile2 = 0.0, M0target = 0.0, M0target2 = 0.0, 88 ProjMassT = 0.0, ProjMassT2 = 0.0, TargMassT = 0.0, TargMassT2 = 0.0, 89 MminProjectile = 0.0, MminTarget = 0.0, 90 ProjectileDiffStateMinMass = 0.0, ProjectileDiffStateMinMass2 = 0.0, 91 ProjectileNonDiffStateMinMass = 0.0, ProjectileNonDiffStateMinMass2 = 0.0, 92 TargetDiffStateMinMass = 0.0, TargetDiffStateMinMass2 = 0.0, 93 TargetNonDiffStateMinMass = 0.0, TargetNonDiffStateMinMass2 = 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, 99 ProbProjectileDiffraction = 0.0, ProbTargetDiffraction = 0.0, ProbOfDiffraction = 0.0; 100 G4LorentzVector Pprojectile, Ptarget, Qmomentum; 101 G4LorentzRotation toCms, toLab; 102 G4SampleResonance BrW; 103 }; 104 G4int ExciteParticipants_doChargeExchange( G4VSplitableHadron* projectile, 105 G4VSplitableHadron* target, 106 G4FTFParameters* theParameters, 107 G4ElasticHNScattering* theElastic, 108 CommonVariables& common ) const; 109 G4bool ExciteParticipants_doDiffraction( G4VSplitableHadron* projectile, 110 G4VSplitableHadron* target, 111 G4FTFParameters* theParameters, 112 CommonVariables& common ) const; 113 G4bool ExciteParticipants_doNonDiffraction( G4VSplitableHadron* projectile, 114 G4VSplitableHadron* target, 115 G4FTFParameters* theParameters, 116 CommonVariables& common ) const; 117 }; 118 119 #endif 120 121