Geant4 Cross Reference

Cross-Referencing   Geant4
Geant4/processes/hadronic/models/parton_string/qgsm/include/G4QGSParticipants.hh

Version: [ ReleaseNotes ] [ 1.0 ] [ 1.1 ] [ 2.0 ] [ 3.0 ] [ 3.1 ] [ 3.2 ] [ 4.0 ] [ 4.0.p1 ] [ 4.0.p2 ] [ 4.1 ] [ 4.1.p1 ] [ 5.0 ] [ 5.0.p1 ] [ 5.1 ] [ 5.1.p1 ] [ 5.2 ] [ 5.2.p1 ] [ 5.2.p2 ] [ 6.0 ] [ 6.0.p1 ] [ 6.1 ] [ 6.2 ] [ 6.2.p1 ] [ 6.2.p2 ] [ 7.0 ] [ 7.0.p1 ] [ 7.1 ] [ 7.1.p1 ] [ 8.0 ] [ 8.0.p1 ] [ 8.1 ] [ 8.1.p1 ] [ 8.1.p2 ] [ 8.2 ] [ 8.2.p1 ] [ 8.3 ] [ 8.3.p1 ] [ 8.3.p2 ] [ 9.0 ] [ 9.0.p1 ] [ 9.0.p2 ] [ 9.1 ] [ 9.1.p1 ] [ 9.1.p2 ] [ 9.1.p3 ] [ 9.2 ] [ 9.2.p1 ] [ 9.2.p2 ] [ 9.2.p3 ] [ 9.2.p4 ] [ 9.3 ] [ 9.3.p1 ] [ 9.3.p2 ] [ 9.4 ] [ 9.4.p1 ] [ 9.4.p2 ] [ 9.4.p3 ] [ 9.4.p4 ] [ 9.5 ] [ 9.5.p1 ] [ 9.5.p2 ] [ 9.6 ] [ 9.6.p1 ] [ 9.6.p2 ] [ 9.6.p3 ] [ 9.6.p4 ] [ 10.0 ] [ 10.0.p1 ] [ 10.0.p2 ] [ 10.0.p3 ] [ 10.0.p4 ] [ 10.1 ] [ 10.1.p1 ] [ 10.1.p2 ] [ 10.1.p3 ] [ 10.2 ] [ 10.2.p1 ] [ 10.2.p2 ] [ 10.2.p3 ] [ 10.3 ] [ 10.3.p1 ] [ 10.3.p2 ] [ 10.3.p3 ] [ 10.4 ] [ 10.4.p1 ] [ 10.4.p2 ] [ 10.4.p3 ] [ 10.5 ] [ 10.5.p1 ] [ 10.6 ] [ 10.6.p1 ] [ 10.6.p2 ] [ 10.6.p3 ] [ 10.7 ] [ 10.7.p1 ] [ 10.7.p2 ] [ 10.7.p3 ] [ 10.7.p4 ] [ 11.0 ] [ 11.0.p1 ] [ 11.0.p2 ] [ 11.0.p3, ] [ 11.0.p4 ] [ 11.1 ] [ 11.1.1 ] [ 11.1.2 ] [ 11.1.3 ] [ 11.2 ] [ 11.2.1 ] [ 11.2.2 ] [ 11.3.0 ]

  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 #ifndef G4QGSParticipants_h
 27 #define G4QGSParticipants_h 1
 28 
 29 #include "Randomize.hh"
 30 #include "G4VParticipants.hh"
 31 #include "G4Nucleon.hh"
 32 #include "G4InteractionContent.hh"
 33 #include "G4QGSDiffractiveExcitation.hh"
 34 #include "G4SingleDiffractiveExcitation.hh"
 35 #include "G4PartonPair.hh" 
 36 #include "G4QGSMSplitableHadron.hh" 
 37 #include "G4V3DNucleus.hh"
 38 
 39 #include "G4VSplitableHadron.hh"
 40 
 41 #include "G4Reggeons.hh"
 42 #include "G4QuarkExchange.hh"
 43 
 44 class G4QGSParticipants : public G4VParticipants
 45 {
 46   public:
 47   G4QGSParticipants();
 48   G4QGSParticipants(const G4QGSParticipants &right);
 49   const G4QGSParticipants & operator=(const G4QGSParticipants &right);
 50   virtual ~G4QGSParticipants();
 51 
 52   G4bool operator==(const G4QGSParticipants &right) const;
 53   G4bool operator!=(const G4QGSParticipants &right) const;
 54 
 55   virtual void DoLorentzBoost(G4ThreeVector aBoost)
 56   {
 57                 theCurrentVelocity = -aBoost;
 58     if(theNucleus) theNucleus->DoLorentzBoost(aBoost);
 59     theBoost = aBoost;
 60   }
 61 
 62   G4PartonPair* GetNextPartonPair();
 63   void BuildInteractions(const G4ReactionProduct  &thePrimary);
 64   void StartPartonPairLoop();
 65 
 66   private:
 67     G4V3DNucleus* GetTargetNucleus() const;
 68     G4V3DNucleus* GetProjectileNucleus() const;
 69 
 70     void PrepareInitialState( const G4ReactionProduct& thePrimary );
 71     void GetList( const G4ReactionProduct& thePrimary );
 72 
 73     void StoreInvolvedNucleon();              
 74     void ReggeonCascade();
 75     G4bool PutOnMassShell();
 76     void GetResiduals();                  
 77 
 78     G4ThreeVector GaussianPt( G4double AveragePt2, G4double maxPtSquare ) const;
 79 
 80     G4bool ComputeNucleusProperties( G4V3DNucleus* nucleus, G4LorentzVector& nucleusMomentum, 
 81                                      G4LorentzVector& residualMomentum, G4double& sumMasses,   
 82                                      G4double& residualExcitationEnergy, G4double& residualMass,
 83                                      G4int& residualMassNumber, G4int& residualCharge );
 84     // Utility methods used by PutOnMassShell.
 85 
 86     G4bool GenerateDeltaIsobar( const G4double sqrtS, const G4int numberOfInvolvedNucleons,
 87                                 G4Nucleon* involvedNucleons[], G4double& sumMasses );
 88 
 89     G4bool SamplingNucleonKinematics( G4double averagePt2, const G4double maxPt2,
 90                                       G4double dCor, G4V3DNucleus* nucleus, 
 91                                       const G4LorentzVector& pResidual, 
 92                                       const G4double residualMass, const G4int residualMassNumber,
 93                                       const G4int numberOfInvolvedNucleons,
 94                                       G4Nucleon* involvedNucleons[], G4double& mass2 );
 95 
 96     G4bool CheckKinematics( const G4double sValue, const G4double sqrtS, 
 97                             const G4double projectileMass2, const G4double targetMass2,
 98                             const G4double nucleusY, const G4bool isProjectileNucleus,
 99                             const G4int numberOfInvolvedNucleons, G4Nucleon* involvedNucleons[],
100                             G4double& targetWminus, G4double& projectileWplus, G4bool& success );
101 
102     G4bool FinalizeKinematics( const G4double w, const G4bool isProjectileNucleus, 
103                                const G4LorentzRotation& boostFromCmsToLab,
104                                const G4double residualMass, const G4int residualMassNumber,
105                                const G4int numberOfInvolvedNucleons, 
106                                G4Nucleon* involvedNucleons[],
107                          G4LorentzVector& residual4Momentum );
108 
109     void CreateStrings();
110 
111   private:
112     // Set parameters of nuclear destruction
113     void SetCofNuclearDestruction( const G4double aValue );
114     void SetR2ofNuclearDestruction( const G4double aValue );
115 
116     void SetExcitationEnergyPerWoundedNucleon( const G4double aValue );
117 
118     void SetDofNuclearDestruction( const G4double aValue );
119     void SetPt2ofNuclearDestruction( const G4double aValue );
120     void SetMaxPt2ofNuclearDestruction( const G4double aValue );
121 
122     // Get parameters of nuclear destruction
123     G4double GetCofNuclearDestruction();
124     G4double GetR2ofNuclearDestruction();
125 
126     G4double GetExcitationEnergyPerWoundedNucleon();
127 
128     G4double GetDofNuclearDestruction();
129     G4double GetPt2ofNuclearDestruction();
130     G4double GetMaxPt2ofNuclearDestruction();
131 
132   protected:
133   virtual G4VSplitableHadron* SelectInteractions(const G4ReactionProduct  &thePrimary);
134   void SplitHadrons(); 
135   void PerformSoftCollisions();
136   void PerformDiffractiveCollisions();
137         G4bool DeterminePartonMomenta();
138 
139   protected:
140   struct DeleteInteractionContent {void operator()(G4InteractionContent*aC){delete aC;}};
141   std::vector<G4InteractionContent*> theInteractions;
142   struct DeleteSplitableHadron{void operator()(G4VSplitableHadron*aS){delete aS;}};
143   std::vector<G4VSplitableHadron*>   theTargets;
144   struct DeletePartonPair{void operator()(G4PartonPair*aP){delete aP;}};
145   std::vector<G4PartonPair*>   thePartonPairs;
146 
147   G4QuarkExchange         theQuarkExchange;
148   G4SingleDiffractiveExcitation theSingleDiffExcitation;
149   G4QGSDiffractiveExcitation    theDiffExcitaton;
150   G4int ModelMode;
151 
152   G4ThreeVector theBoost;
153   G4double SampleX(G4double anXmin, G4int nSea, G4int theTotalSea, G4double aBeta);
154 
155   protected:
156   // model parameters HPW
157   enum  { SOFT, DIFFRACTIVE };
158   enum  { ALL, WITHOUT_R, NON_DIFF };   // Interaction modes
159   enum  { PrD, TrD, DD, NonD, Qexc };   // Interaction types
160 
161   const G4int nCutMax;
162   const G4double ThresholdParameter;
163   const G4double QGSMThreshold;
164   const G4double theNucleonRadius;
165 
166         G4ThreeVector theCurrentVelocity;
167         G4QGSMSplitableHadron* theProjectileSplitable;
168 
169   private:
170     G4ReactionProduct theProjectile;       
171 
172     G4Reggeons* Regge;
173     G4int InteractionMode;
174 
175     G4double alpha;
176     G4double beta;
177 
178     G4double sigmaPt;   
179 
180     G4Nucleon* TheInvolvedNucleonsOfTarget[250];
181     G4int NumberOfInvolvedNucleonsOfTarget;
182 
183     G4Nucleon* TheInvolvedNucleonsOfProjectile[250];
184     G4int NumberOfInvolvedNucleonsOfProjectile;
185 
186     G4LorentzVector ProjectileResidual4Momentum;
187     G4int           ProjectileResidualMassNumber;
188     G4int           ProjectileResidualCharge;
189     G4double        ProjectileResidualExcitationEnergy;
190 
191     G4LorentzVector TargetResidual4Momentum;
192     G4int           TargetResidualMassNumber;
193     G4int           TargetResidualCharge;
194     G4double        TargetResidualExcitationEnergy;
195 
196   private:
197     // Parameters of nuclear destruction
198     G4double CofNuclearDestruction;   // Cnd of nuclear destruction
199     G4double R2ofNuclearDestruction;  // R2nd
200 
201     G4double ExcitationEnergyPerWoundedNucleon;
202 
203     G4double DofNuclearDestruction;       // D for momentum sampling
204     G4double Pt2ofNuclearDestruction;     // Pt2
205     G4double MaxPt2ofNuclearDestruction;  // Max Pt2
206 };
207 
208 inline void G4QGSParticipants::StartPartonPairLoop()
209 {
210 }
211 
212 inline G4PartonPair* G4QGSParticipants::GetNextPartonPair()
213 {
214   if (thePartonPairs.empty()) return 0;
215   G4PartonPair * result = thePartonPairs.back();
216   thePartonPairs.pop_back();
217   return result;
218 }
219 
220 inline void G4QGSParticipants::SplitHadrons()
221 {
222   unsigned int i;
223   for(i = 0; i < theInteractions.size(); i++)
224   {
225     theInteractions[i]->SplitHadrons();
226   }
227 }
228 
229 inline G4V3DNucleus* G4QGSParticipants::GetTargetNucleus() const {
230   return theNucleus;
231 }
232 
233 inline G4V3DNucleus* G4QGSParticipants::GetProjectileNucleus() const {
234   return 0;
235 }
236 
237 // Uzhi Start copy from FTFparameters
238 // Set parameters of nuclear destruction
239 
240 inline void G4QGSParticipants::SetCofNuclearDestruction( const G4double aValue ) {
241   CofNuclearDestruction = aValue;
242 }
243 
244 inline void G4QGSParticipants::SetR2ofNuclearDestruction( const G4double aValue ) {
245   R2ofNuclearDestruction = aValue;
246 }
247 
248 inline void G4QGSParticipants::SetExcitationEnergyPerWoundedNucleon( const G4double aValue ) {
249   ExcitationEnergyPerWoundedNucleon = aValue;
250 }
251 
252 inline void G4QGSParticipants::SetDofNuclearDestruction( const G4double aValue ) {
253   DofNuclearDestruction = aValue;
254 }
255 
256 inline void G4QGSParticipants::SetPt2ofNuclearDestruction( const G4double aValue ) {
257   Pt2ofNuclearDestruction = aValue;
258 }
259 
260 inline void G4QGSParticipants::SetMaxPt2ofNuclearDestruction( const G4double aValue ) {
261   MaxPt2ofNuclearDestruction = aValue;
262 }
263 
264 // Get parameters of nuclear destruction
265 inline G4double G4QGSParticipants::GetCofNuclearDestruction() {
266   return CofNuclearDestruction;
267 }
268 
269 inline G4double G4QGSParticipants::GetR2ofNuclearDestruction() {
270   return R2ofNuclearDestruction;
271 }
272 
273 inline G4double G4QGSParticipants::GetExcitationEnergyPerWoundedNucleon() {
274   return ExcitationEnergyPerWoundedNucleon;
275 }
276 
277 inline G4double G4QGSParticipants::GetDofNuclearDestruction() {
278   return DofNuclearDestruction;
279 }
280 
281 inline G4double G4QGSParticipants::GetPt2ofNuclearDestruction() {
282   return Pt2ofNuclearDestruction;
283 }
284 
285 inline G4double G4QGSParticipants::GetMaxPt2ofNuclearDestruction() {
286   return MaxPt2ofNuclearDestruction;
287 }
288 //Uzhi End copy from FTFparameters
289 #endif
290 
291