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 G4FTFAnnihilation_h 29 #define G4FTFAnnihilation_h 1 30 31 // ------------------------------------------------------------ 32 // GEANT 4 class header file 33 // 34 // ---------------- G4FTFAnnihilation -------------- 35 // by Vladimir Uzhinsky, November 2010. 36 // Annihilation used by Fritiof (FTF) model 37 // Takes a projectile and a target 38 // Produces strings (excited hadrons) 39 // ------------------------------------------------------------ 40 41 #include "globals.hh" 42 #include "G4FTFParameters.hh" 43 #include "G4ThreeVector.hh" 44 45 class G4VSplitableHadron; 46 class G4ExcitedString; 47 48 49 class G4FTFAnnihilation { 50 public: 51 G4FTFAnnihilation(); 52 virtual ~G4FTFAnnihilation(); 53 54 virtual G4bool Annihilate( G4VSplitableHadron* aPartner, 55 G4VSplitableHadron* bPartner, 56 G4VSplitableHadron*& AdditionalString, 57 G4FTFParameters* theParameters ) const; 58 59 private: 60 G4FTFAnnihilation( const G4FTFAnnihilation& right ); 61 const G4FTFAnnihilation& operator=( const G4FTFAnnihilation& right ); 62 G4bool operator==( const G4FTFAnnihilation& right ) const; 63 G4bool operator!=( const G4FTFAnnihilation& right ) const; 64 65 // The "Annihilate" method uses the following struct and 4 new utility methods: 66 struct CommonVariables { 67 G4int AQ[3], Q[3]; 68 G4bool RotateStrings = false; 69 G4double S = 0.0, SqrtS = 0.0; 70 G4LorentzVector Pprojectile, Ptarget; 71 G4LorentzRotation toLab, RandomRotation; 72 }; 73 G4bool Create3QuarkAntiQuarkStrings( G4VSplitableHadron* aPartner, 74 G4VSplitableHadron* bPartner, 75 G4VSplitableHadron*& AdditionalString, 76 G4FTFParameters* theParameters, 77 CommonVariables& common ) const; 78 G4int Create1DiquarkAntiDiquarkString( G4VSplitableHadron* aPartner, 79 G4VSplitableHadron* bPartner, 80 CommonVariables& common ) const; 81 G4int Create2QuarkAntiQuarkStrings( G4VSplitableHadron* aPartner, 82 G4VSplitableHadron* bPartner, 83 G4FTFParameters* theParameters, 84 CommonVariables& common ) const; 85 G4bool Create1QuarkAntiQuarkString( G4VSplitableHadron* aPartner, 86 G4VSplitableHadron* bPartner, 87 G4FTFParameters* theParameters, 88 CommonVariables& common ) const; 89 90 G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const; 91 G4double ChooseX( G4double Alpha, G4double Beta ) const; 92 void UnpackBaryon( G4int IdPDG, G4int& Q1, G4int& Q2, G4int& Q3 ) const; 93 }; 94 95 #endif 96 97