source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/fragmentation/include/G4QIonIonCollision.hh @ 1340

Last change on this file since 1340 was 1340, checked in by garnier, 14 years ago

update ti head

File size: 5.1 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//
27// $Id: G4QIonIonCollision.hh,v 1.1 2009/11/16 18:16:04 mkossov Exp $
28// GEANT4 tag $Name: hadr-chips-V09-03-08 $
29//
30// -----------------------------------------------------------------------------
31//      GEANT 4 class header file
32//
33//                 History:
34//     Created by Mikhail Kossov, October 2006
35//     CHIPS QGS fragmentation class
36//     For comparison mirror member functions are taken from G4 classes:
37//     G4QGSParticipants
38//     G4QGSModels
39//     G4ExcitedStringDecay
40// -----------------------------------------------------------------------------
41// Short description: CHIPS string fragmentation class for ion-ion inelastic
42// -----------------------------------------------------------------------------
43//
44#ifndef G4QIonIonCollision_h
45#define G4QIonIonCollision_h 1
46
47#include "globals.hh"
48#include "G4LorentzVector.hh"
49#include "Randomize.hh"
50#include "G4QNucleus.hh"
51#include "G4Quasmon.hh"
52#include "G4QHadronVector.hh"
53#include "G4QEnvironment.hh"
54#include "G4QInteractionVector.hh"
55#include "G4QProbability.hh"
56#include "G4QPartonPairVector.hh"
57#include "G4QuasiFreeRatios.hh"
58#include "G4QStringVector.hh"
59
60class G4QIonIonCollision
61{
62 public:
63  G4QIonIonCollision(G4QNucleus& pNucleus, const G4QNucleus& tNucleus);
64  ~G4QIonIonCollision(); 
65
66  G4QHadronVector* Fragment(); // Calls Breeder & fragments Quasmons inside ResidualNucleus
67
68  // Static functions
69  static void SetParameters(G4int nC, G4double strTens, G4double tubeDens, G4double SigPt);
70
71 protected:
72  G4bool ExciteDiffParticipants(G4QHadron* aPartner, G4QHadron* bPartner) const; //@@Once
73  G4bool ExciteSingDiffParticipants(G4QHadron* aPartner, G4QHadron* bPartner) const;//@@Onc
74  std::pair<G4int,G4int> ReducePair(G4int P1, G4int P2) const; // Reduce Q-pairs to singles
75  void Breeder(); // String fragmentation algoritm, which makes Hadrons & Quasmons from Str
76  // @@ At present always SingleDiffractive, so this is a fake member function
77  G4bool IsSingleDiffractive() {G4bool r=false; if(G4UniformRand()<1.) r=true; return r;}
78  G4int SumPartonPDG(G4int PDG1, G4int PFG2) const;
79  G4double ChooseX(G4double Xmin, G4double Xmax) const;
80  G4ThreeVector GaussianPt(G4double widthSquare, G4double maxPtSquare) const;
81  G4int AnnihilationOrder(G4int LS, G4int MS, G4int uP, G4int mP, G4int sP, G4int nP);
82  void SwapPartons(); // Try to swap partons of strings, if one of strings has negative M2
83
84 private:
85  enum {SOFT, DIFFRACTIVE};
86  // static model parameters
87  static G4int    nCutMax;                               // Maximum number of Soft Cuts
88  static G4double stringTension;                         // String Tension to absorb energy
89  static G4double tubeDensity;                           // Nucleon density in the FluxTube
90  static G4double widthOfPtSquare;                       // width^2 of pt(StringExcitation)
91
92  // Body
93  G4QNucleus      theProjNucleus;                        // ProjectileNucleus moving LS->CM
94  G4QNucleus      theTargNucleus;                        // TargetNucleus moving fromLStoCM
95  G4QStringVector strings;                               // Vector of created strings
96  G4QuasmonVector theQuasmons;                           // Strings converter to Quasmons
97  G4QHadronVector* theResult;                            // Pointer to the OUTPUT Result
98  G4double        maxEn;                                 // Energy absorbed by the nucleus
99  G4double        maxNuc;                                // #0fNucleons in the Flux Tube
100  G4QCHIPSWorld*  theWorld;                              // Pointer to the CHIPS World
101};
102
103#endif
Note: See TracBrowser for help on using the repository browser.