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

Last change on this file since 949 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 5.8 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;
63
64class G4BinaryCascade : public G4VIntraNuclearTransportModel
65{
66public:
67
68  G4BinaryCascade();
69  G4BinaryCascade(const G4BinaryCascade & right);
70
71  virtual ~G4BinaryCascade();
72
73  const G4BinaryCascade& operator=(G4BinaryCascade & right);
74  G4int operator==(G4BinaryCascade& right) {return (this == &right);}
75  G4int operator!=(G4BinaryCascade& right) {return (this != &right);}
76
77  G4HadFinalState* ApplyYourself(const G4HadProjectile& aTrack, 
78                                              G4Nucleus& theNucleus);
79  virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *,
80                                              G4V3DNucleus *);
81
82private:
83
84  G4int GetTotalCharge(std::vector<G4KineticTrack *> & aV)
85  {
86    G4int result = 0;
87    std::vector<G4KineticTrack *>::iterator i;
88    for(i = aV.begin(); i != aV.end(); ++i)
89    {
90       result += G4lrint((*i)->GetDefinition()->GetPDGCharge());
91    }
92    return result;
93  }
94  void PrintWelcomeMessage();
95  void BuildTargetList();
96  void FindCollisions(G4KineticTrackVector *);
97  void FindDecayCollision(G4KineticTrack *);
98  void FindLateParticleCollision(G4KineticTrack *);
99  G4bool ApplyCollision(G4CollisionInitialState *);
100  G4bool Capture(G4bool verbose=false);
101  G4bool Absorb();
102  G4bool CheckPauliPrinciple(G4KineticTrackVector *);
103  G4double GetExcitationEnergy();
104  G4bool CheckDecay(G4KineticTrackVector *);
105  void CorrectFinalPandE();
106  void UpdateTracksAndCollisions(G4KineticTrackVector * oldSecondaries,
107                                 G4KineticTrackVector * oldTarget,
108                                 G4KineticTrackVector * newSecondaries);
109  G4bool DoTimeStep(G4double timeStep);
110  G4KineticTrackVector* CorrectBarionsOnBoundary(G4KineticTrackVector *in, 
111                                               G4KineticTrackVector *out);
112  G4Fragment * FindFragments();
113  void StepParticlesOut();
114  G4LorentzVector GetFinal4Momentum();
115  G4LorentzVector GetFinalNucleusMomentum();
116  G4ReactionProductVector * Propagate1H1(G4KineticTrackVector *,
117                                         G4V3DNucleus *);
118  G4double GetIonMass(G4int Z, G4int A);
119 
120// utility methods
121  G4ThreeVector GetSpherePoint(G4double r, const G4LorentzVector & momentumdirection);
122  void ClearAndDestroy(G4KineticTrackVector * ktv);
123  void ClearAndDestroy(G4ReactionProductVector * rpv);
124
125// for debugging purpose
126  void PrintKTVector(G4KineticTrackVector * ktv, std::string comment=std::string(""));
127  void PrintKTVector(G4KineticTrack* kt, std::string comment=std::string(""));
128  void DebugApplyCollision(G4CollisionInitialState * collision, 
129                           G4KineticTrackVector *products);
130private:
131  G4KineticTrackVector theProjectileList;
132  G4KineticTrackVector theTargetList;
133  G4KineticTrackVector theSecondaryList;
134  G4KineticTrackVector theCapturedList;
135  G4KineticTrackVector theFinalState;
136
137  G4ExcitationHandler * theExcitationHandler;
138  G4CollisionManager * theCollisionMgr;
139
140  std::vector<G4BCAction *> theImR;
141  G4BCDecay * theDecay;
142  G4BCLateParticle * theLateParticle;
143  G4VFieldPropagation * thePropagator;
144  G4double theCurrentTime;
145  G4double theBCminP;
146  G4double theCutOnP;
147  G4double theCutOnPAbsorb;
148  G4LorentzVector theInitial4Mom;
149  G4int currentA, currentZ;
150  G4double massInNucleus;
151  G4double currentInitialEnergy;   // for debugging
152  G4LorentzRotation precompoundLorentzboost;
153  G4double theOuterRadius;
154  G4bool thePrimaryEscape;
155  G4ParticleDefinition * thePrimaryType;
156  G4ThreeVector theMomentumTransfer;
157 
158 
159
160};
161
162#endif
163
164
165
166
Note: See TracBrowser for help on using the repository browser.