File indexing completed on 2025-01-18 09:57:58
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039
0040
0041
0042 #ifndef G4BinaryCascade_hh
0043 #define G4BinaryCascade_hh
0044
0045 #include "G4VIntraNuclearTransportModel.hh"
0046 #include "G4ReactionProductVector.hh"
0047 #include "G4KineticTrackVector.hh"
0048 #include "G4ListOfCollisions.hh"
0049 #include "G4V3DNucleus.hh"
0050 #include "G4Fancy3DNucleus.hh"
0051 #include "G4Fragment.hh"
0052 #include "G4VFieldPropagation.hh"
0053 #include "G4VScatterer.hh"
0054 #include "G4LorentzVector.hh"
0055 #include "G4LorentzRotation.hh"
0056
0057 #include "G4BCDecay.hh"
0058 #include "G4BCLateParticle.hh"
0059 #include "G4BCAction.hh"
0060
0061 #include "G4DecayKineticTracks.hh"
0062
0063 #include "G4Threading.hh"
0064
0065 class G4CollisionManager;
0066
0067 class G4Track;
0068 class G4KineticTrack;
0069 class G43DNucleus;
0070 class G4Scatterer;
0071
0072 class G4BinaryCascade : public G4VIntraNuclearTransportModel
0073 {
0074 public:
0075
0076 G4BinaryCascade(G4VPreCompoundModel* ptr = 0);
0077
0078 virtual ~G4BinaryCascade();
0079
0080 G4HadFinalState* ApplyYourself(const G4HadProjectile& aTrack,
0081 G4Nucleus& theNucleus);
0082 virtual G4ReactionProductVector * Propagate(G4KineticTrackVector *,
0083 G4V3DNucleus *);
0084
0085 virtual void ModelDescription(std::ostream&) const;
0086 virtual void PropagateModelDescription(std::ostream&) const;
0087
0088 private:
0089
0090 G4BinaryCascade(const G4BinaryCascade & right);
0091 const G4BinaryCascade& operator=(G4BinaryCascade & right);
0092 G4bool operator==(G4BinaryCascade& right) {return (this == &right);}
0093 G4bool operator!=(G4BinaryCascade& right) {return (this != &right);}
0094
0095
0096 void PrintWelcomeMessage();
0097 void BuildTargetList();
0098 G4bool BuildLateParticleCollisions(G4KineticTrackVector * secondaries);
0099 void FindCollisions(G4KineticTrackVector *);
0100 void FindDecayCollision(G4KineticTrack *);
0101 void FindLateParticleCollision(G4KineticTrack *);
0102
0103
0104
0105 G4ReactionProductVector * DeExcite();
0106 G4ReactionProductVector * DecayVoidNucleus();
0107 G4ReactionProductVector * ProductsAddFinalState (G4ReactionProductVector * products, G4KineticTrackVector & finalState );
0108 G4ReactionProductVector * ProductsAddPrecompound(G4ReactionProductVector * products ,G4ReactionProductVector * preco );
0109
0110 G4bool ApplyCollision(G4CollisionInitialState *);
0111 G4bool Capture(G4bool verbose=false);
0112 G4bool Absorb();
0113 G4bool CheckPauliPrinciple(G4KineticTrackVector *);
0114 G4double GetExcitationEnergy();
0115 void CorrectFinalPandE();
0116
0117 G4double CorrectShortlivedPrimaryForFermi(
0118 G4KineticTrack* primary,G4KineticTrackVector target_collection);
0119 G4bool CorrectShortlivedFinalsForFermi(
0120 G4KineticTrackVector * products, G4double initial_Efermi);
0121
0122 void UpdateTracksAndCollisions(G4KineticTrackVector * oldSecondaries,
0123 G4KineticTrackVector * oldTarget,
0124 G4KineticTrackVector * newSecondaries);
0125 G4bool DoTimeStep(G4double timeStep);
0126 G4KineticTrackVector* CorrectBarionsOnBoundary(G4KineticTrackVector *in,
0127 G4KineticTrackVector *out);
0128 G4Fragment * FindFragments();
0129 void StepParticlesOut();
0130
0131 G4LorentzVector GetFinal4Momentum();
0132 G4LorentzVector GetFinalNucleusMomentum();
0133
0134 G4ReactionProductVector * Propagate1H1(G4KineticTrackVector *,
0135 G4V3DNucleus *);
0136 G4double GetIonMass(G4int Z, G4int A);
0137
0138 G4int GetTotalCharge(std::vector<G4KineticTrack *> & aV)
0139 {
0140 G4int result = 0;
0141 std::vector<G4KineticTrack *>::iterator i;
0142 for(i = aV.begin(); i != aV.end(); ++i)
0143 {
0144 result += G4lrint((*i)->GetDefinition()->GetPDGCharge());
0145 }
0146 return result;
0147 }
0148 G4int GetTotalBaryonCharge(std::vector<G4KineticTrack *> & aV)
0149 {
0150 G4int result = 0;
0151 std::vector<G4KineticTrack *>::iterator i;
0152 for(i = aV.begin(); i != aV.end(); ++i)
0153 {
0154 if ( (*i)->GetDefinition()->GetBaryonNumber() != 0 ){
0155 result += G4lrint((*i)->GetDefinition()->GetPDGCharge());
0156 }
0157 }
0158 return result;
0159 }
0160
0161 G4ReactionProductVector * HighEnergyModelFSProducts(G4ReactionProductVector *,
0162 G4KineticTrackVector * secondaries);
0163 G4ReactionProductVector * FillVoidNucleusProducts(G4ReactionProductVector * );
0164
0165
0166 G4ThreeVector GetSpherePoint(G4double r, const G4LorentzVector & momentumdirection);
0167
0168 void ClearAndDestroy(G4KineticTrackVector * ktv);
0169 void ClearAndDestroy(G4ReactionProductVector * rpv);
0170
0171
0172
0173 G4ReactionProductVector * ProductsAddFakeGamma(G4ReactionProductVector * products );
0174
0175 void PrintKTVector(G4KineticTrackVector * ktv, std::string comment=std::string(""));
0176 void PrintKTVector(G4KineticTrack* kt, std::string comment=std::string(""));
0177 void DebugApplyCollisionFail(G4CollisionInitialState * collision,
0178 G4KineticTrackVector * products);
0179 void DebugApplyCollision(G4CollisionInitialState * collision,
0180 G4KineticTrackVector *products);
0181 G4bool DebugFinalEpConservation(const G4HadProjectile & aTrack, G4ReactionProductVector* products);
0182 G4bool DebugEpConservation(const G4String where);
0183 G4bool CheckChargeAndBaryonNumber(G4String where);
0184
0185 private:
0186 G4KineticTrackVector theProjectileList;
0187 G4KineticTrackVector theTargetList;
0188 G4KineticTrackVector theSecondaryList;
0189 G4KineticTrackVector theCapturedList;
0190 G4KineticTrackVector theFinalState;
0191
0192
0193 G4ExcitationHandler * theExcitationHandler;
0194 G4CollisionManager * theCollisionMgr;
0195
0196 G4Scatterer * theH1Scatterer;
0197
0198 std::vector<G4BCAction *> theImR;
0199 G4BCDecay * theDecay;
0200 G4BCLateParticle * theLateParticle;
0201 G4VFieldPropagation * thePropagator;
0202 G4DecayKineticTracks decayKTV;
0203
0204 G4double theCurrentTime;
0205 G4double theBCminP;
0206 G4double theCutOnP;
0207 G4double theCutOnPAbsorb;
0208 G4LorentzVector theInitial4Mom;
0209 G4LorentzVector theProjectile4Momentum;
0210 G4int currentA, currentZ, lateA, lateZ;
0211 G4int initialZ, initialA,projectileA,projectileZ;
0212 G4double massInNucleus, initial_nuclear_mass;
0213 G4double currentInitialEnergy;
0214 G4LorentzRotation precompoundLorentzboost;
0215 G4double theOuterRadius;
0216 G4bool thePrimaryEscape;
0217 const G4ParticleDefinition * thePrimaryType;
0218 G4ThreeVector theMomentumTransfer;
0219 static G4int theBIC_ID;
0220 G4bool fBCDEBUG;
0221 };
0222
0223 #endif
0224
0225
0226
0227