source: trunk/source/processes/hadronic/models/binary_cascade/include/G4BinaryCascade.hh @ 962

Last change on this file since 962 was 962, checked in by garnier, 15 years ago

update processes

File size: 5.9 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//
28// -------------------------------------------------------------------
29//      GEANT4 Class file
30//
31//
32//      File name:    G4BinaryCascade.hh
33//
34//      Author: Alessandro Brunengo (Alessandro.Brunengo@ge.infn.it)       
35//
36//      Creation date: 8 June 2000
37// -----------------------------------------------------------------------------
38
39#ifndef G4BinaryCascade_hh
40#define G4BinaryCascade_hh
41
42#include "G4VIntraNuclearTransportModel.hh"
43#include "G4ReactionProductVector.hh"
44#include "G4KineticTrackVector.hh"
45#include "G4ListOfCollisions.hh"
46#include "G4V3DNucleus.hh"
47#include "G4Fancy3DNucleus.hh"
48#include "G4Fragment.hh"
49#include "G4VFieldPropagation.hh"
50#include "G4VScatterer.hh"
51#include "G4LorentzVector.hh"
52#include "G4LorentzRotation.hh"
53
54#include "G4BCDecay.hh"
55#include "G4BCLateParticle.hh"
56#include "G4BCAction.hh"
57
58class G4CollisionManager;
59
60class G4Track;
61class G4KineticTrack;
62class G43DNucleus;
63class G4Scatterer;
64
65class G4BinaryCascade : public G4VIntraNuclearTransportModel
66{
67public:
68
69  G4BinaryCascade();
70  G4BinaryCascade(const G4BinaryCascade & right);
71
72  virtual ~G4BinaryCascade();
73
74  const G4BinaryCascade& operator=(G4BinaryCascade & right);
75  G4int operator==(G4BinaryCascade& right) {return (this == &right);}
76  G4int operator!=(G4BinaryCascade& right) {return (this != &right);}
77
78  G4HadFinalState* ApplyYourself(const G4HadProjectile& aTrack, 
79                                              G4Nucleus& theNucleus);
80  virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *,
81                                              G4V3DNucleus *);
82
83private:
84
85  G4int GetTotalCharge(std::vector<G4KineticTrack *> & aV)
86  {
87    G4int result = 0;
88    std::vector<G4KineticTrack *>::iterator i;
89    for(i = aV.begin(); i != aV.end(); ++i)
90    {
91       result += G4lrint((*i)->GetDefinition()->GetPDGCharge());
92    }
93    return result;
94  }
95  void PrintWelcomeMessage();
96  void BuildTargetList();
97  void FindCollisions(G4KineticTrackVector *);
98  void FindDecayCollision(G4KineticTrack *);
99  void FindLateParticleCollision(G4KineticTrack *);
100  G4bool ApplyCollision(G4CollisionInitialState *);
101  G4bool Capture(G4bool verbose=false);
102  G4bool Absorb();
103  G4bool CheckPauliPrinciple(G4KineticTrackVector *);
104  G4double GetExcitationEnergy();
105  G4bool CheckDecay(G4KineticTrackVector *);
106  void CorrectFinalPandE();
107  void UpdateTracksAndCollisions(G4KineticTrackVector * oldSecondaries,
108                                 G4KineticTrackVector * oldTarget,
109                                 G4KineticTrackVector * newSecondaries);
110  G4bool DoTimeStep(G4double timeStep);
111  G4KineticTrackVector* CorrectBarionsOnBoundary(G4KineticTrackVector *in, 
112                                               G4KineticTrackVector *out);
113  G4Fragment * FindFragments();
114  void StepParticlesOut();
115  G4LorentzVector GetFinal4Momentum();
116  G4LorentzVector GetFinalNucleusMomentum();
117  G4ReactionProductVector * Propagate1H1(G4KineticTrackVector *,
118                                         G4V3DNucleus *);
119  G4double GetIonMass(G4int Z, G4int A);
120 
121// utility methods
122  G4ThreeVector GetSpherePoint(G4double r, const G4LorentzVector & momentumdirection);
123  void ClearAndDestroy(G4KineticTrackVector * ktv);
124  void ClearAndDestroy(G4ReactionProductVector * rpv);
125
126// for debugging purpose
127  void PrintKTVector(G4KineticTrackVector * ktv, std::string comment=std::string(""));
128  void PrintKTVector(G4KineticTrack* kt, std::string comment=std::string(""));
129  void DebugApplyCollision(G4CollisionInitialState * collision, 
130                           G4KineticTrackVector *products);
131private:
132  G4KineticTrackVector theProjectileList;
133  G4KineticTrackVector theTargetList;
134  G4KineticTrackVector theSecondaryList;
135  G4KineticTrackVector theCapturedList;
136  G4KineticTrackVector theFinalState;
137
138  G4ExcitationHandler * theExcitationHandler;
139  G4CollisionManager * theCollisionMgr;
140 
141  G4Scatterer * theH1Scatterer;
142
143  std::vector<G4BCAction *> theImR;
144  G4BCDecay * theDecay;
145  G4BCLateParticle * theLateParticle;
146  G4VFieldPropagation * thePropagator;
147  G4double theCurrentTime;
148  G4double theBCminP;
149  G4double theCutOnP;
150  G4double theCutOnPAbsorb;
151  G4LorentzVector theInitial4Mom;
152  G4int currentA, currentZ;
153  G4double massInNucleus;
154  G4double currentInitialEnergy;   // for debugging
155  G4LorentzRotation precompoundLorentzboost;
156  G4double theOuterRadius;
157  G4bool thePrimaryEscape;
158  G4ParticleDefinition * thePrimaryType;
159  G4ThreeVector theMomentumTransfer;
160 
161 
162
163};
164
165#endif
166
167
168
169
Note: See TracBrowser for help on using the repository browser.