[819] | 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 | // Hadronic Process: Reaction Dynamics |
---|
| 29 | // original by H.P. Wellisch |
---|
| 30 | // Modified by J.L.Chuma 19-Nov-96 |
---|
| 31 | // Modified by J.L.Chuma 27-Mar-97 |
---|
| 32 | // Modified by J.L.Chuma 30-Apr-97 |
---|
| 33 | // Modified by J.L.Chuma 06-Aug-97 to include the original incident particle |
---|
| 34 | // before Fermi motion and evaporation effects |
---|
| 35 | |
---|
| 36 | #ifndef G4ReactionDynamics_h |
---|
| 37 | #define G4ReactionDynamics_h 1 |
---|
| 38 | |
---|
| 39 | #include "G4ParticleTypes.hh" |
---|
| 40 | #include "G4DynamicParticle.hh" |
---|
| 41 | #include "G4ReactionProduct.hh" |
---|
| 42 | #include "G4Nucleus.hh" |
---|
| 43 | #include "G4FastVector.hh" |
---|
| 44 | #include "G4HadProjectile.hh" |
---|
| 45 | |
---|
| 46 | enum{ GHADLISTSIZE=256}; |
---|
| 47 | |
---|
| 48 | class G4ReactionDynamics |
---|
| 49 | { |
---|
| 50 | public: |
---|
| 51 | |
---|
| 52 | G4ReactionDynamics() {} |
---|
| 53 | |
---|
| 54 | virtual ~G4ReactionDynamics() {} |
---|
| 55 | |
---|
| 56 | virtual G4double FindInelasticity() |
---|
| 57 | { return 0.0; } |
---|
| 58 | |
---|
| 59 | virtual G4double FindTimeDelay() |
---|
| 60 | { return 0.0; } |
---|
| 61 | |
---|
| 62 | G4bool GenerateXandPt( // derived from GENXPT |
---|
| 63 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 64 | G4int &vecLen, |
---|
| 65 | G4ReactionProduct &modifiedOriginal, // Fermi motion & evap. effect included |
---|
| 66 | const G4HadProjectile *originalIncident, |
---|
| 67 | G4ReactionProduct ¤tParticle, |
---|
| 68 | G4ReactionProduct &targetParticle, |
---|
| 69 | const G4DynamicParticle* originalTarget, |
---|
| 70 | const G4Nucleus &targetNucleus, |
---|
| 71 | G4bool &incidentHasChanged, |
---|
| 72 | G4bool &targetHasChanged, |
---|
| 73 | G4bool leadFlag, |
---|
| 74 | G4ReactionProduct &leadingStrangeParticle ); |
---|
| 75 | |
---|
| 76 | void SuppressChargedPions( |
---|
| 77 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 78 | G4int &vecLen, |
---|
| 79 | const G4ReactionProduct &modifiedOriginal, |
---|
| 80 | G4ReactionProduct ¤tParticle, |
---|
| 81 | G4ReactionProduct &targetParticle, |
---|
| 82 | const G4Nucleus &targetNucleus, |
---|
| 83 | G4bool &incidentHasChanged, |
---|
| 84 | G4bool &targetHasChanged ); |
---|
| 85 | |
---|
| 86 | G4bool TwoCluster( // derived from TWOCLU |
---|
| 87 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 88 | G4int &vecLen, |
---|
| 89 | G4ReactionProduct &modifiedOriginal, // Fermi motion & evap. effect included |
---|
| 90 | const G4HadProjectile *originalIncident, |
---|
| 91 | G4ReactionProduct ¤tParticle, |
---|
| 92 | G4ReactionProduct &targetParticle, |
---|
| 93 | const G4DynamicParticle* originalTarget, |
---|
| 94 | const G4Nucleus &targetNucleus, |
---|
| 95 | G4bool &incidentHasChanged, |
---|
| 96 | G4bool &targetHasChanged, |
---|
| 97 | G4bool leadFlag, |
---|
| 98 | G4ReactionProduct &leadingStrangeParticle ); |
---|
| 99 | |
---|
| 100 | void TwoBody( // derived from TWOB |
---|
| 101 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 102 | G4int &vecLen, |
---|
| 103 | G4ReactionProduct &modifiedOriginal, |
---|
| 104 | const G4DynamicParticle *originalTarget, |
---|
| 105 | G4ReactionProduct ¤tParticle, |
---|
| 106 | G4ReactionProduct &targetParticle, |
---|
| 107 | const G4Nucleus &targetNucleus, |
---|
| 108 | G4bool &targetHasChanged ); |
---|
| 109 | |
---|
| 110 | G4int Factorial( G4int n ); |
---|
| 111 | |
---|
| 112 | G4double GenerateNBodyEvent( // derived from PHASP |
---|
| 113 | const G4double totalEnergy, |
---|
| 114 | const G4bool constantCrossSection, |
---|
| 115 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 116 | G4int &vecLen ); |
---|
| 117 | |
---|
| 118 | void ProduceStrangeParticlePairs( |
---|
| 119 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 120 | G4int &vecLen, |
---|
| 121 | const G4ReactionProduct &modifiedOriginal, |
---|
| 122 | const G4DynamicParticle *originalTarget, |
---|
| 123 | G4ReactionProduct ¤tParticle, |
---|
| 124 | G4ReactionProduct &targetParticle, |
---|
| 125 | G4bool &incidentHasChanged, |
---|
| 126 | G4bool &targetHasChanged ); |
---|
| 127 | |
---|
| 128 | void NuclearReaction( // derived from NUCREC |
---|
| 129 | G4FastVector<G4ReactionProduct,4> &vec, |
---|
| 130 | G4int &vecLen, |
---|
| 131 | const G4HadProjectile *originalIncident, |
---|
| 132 | const G4Nucleus &aNucleus, |
---|
| 133 | const G4double theAtomicMass, |
---|
| 134 | const G4double *massVec ); |
---|
| 135 | |
---|
| 136 | private: |
---|
| 137 | |
---|
| 138 | void Rotate( |
---|
| 139 | const G4double numberofFinalStateNucleons, |
---|
| 140 | const G4ThreeVector &temp, |
---|
| 141 | const G4ReactionProduct &modifiedOriginal, // Fermi motion & evap. effect included |
---|
| 142 | const G4HadProjectile *originalIncident, |
---|
| 143 | const G4Nucleus &targetNucleus, |
---|
| 144 | G4ReactionProduct ¤tParticle, |
---|
| 145 | G4ReactionProduct &targetParticle, |
---|
| 146 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 147 | G4int &vecLen ); |
---|
| 148 | |
---|
| 149 | void Defs1( |
---|
| 150 | const G4ReactionProduct &modifiedOriginal, |
---|
| 151 | G4ReactionProduct ¤tParticle, |
---|
| 152 | G4ReactionProduct &targetParticle, |
---|
| 153 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 154 | G4int &vecLen ); |
---|
| 155 | |
---|
| 156 | void AddBlackTrackParticles( |
---|
| 157 | const G4double epnb, |
---|
| 158 | const G4int npnb, |
---|
| 159 | const G4double edta, |
---|
| 160 | const G4int ndta, |
---|
| 161 | const G4double sprob, |
---|
| 162 | const G4double kineticMinimum, |
---|
| 163 | const G4double kineticFactor, |
---|
| 164 | const G4ReactionProduct &modifiedOriginal, |
---|
| 165 | G4int PinNucleus, |
---|
| 166 | G4int NinNucleus, |
---|
| 167 | const G4Nucleus &aNucleus, |
---|
| 168 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 169 | G4int &vecLen ); |
---|
| 170 | |
---|
| 171 | std::pair<G4int, G4int> GetFinalStateNucleons( |
---|
| 172 | const G4DynamicParticle* originalTarget, |
---|
| 173 | const G4FastVector<G4ReactionProduct,GHADLISTSIZE>& vec, |
---|
| 174 | const G4int& vecLen ); |
---|
| 175 | |
---|
| 176 | void MomentumCheck( |
---|
| 177 | const G4ReactionProduct &modifiedOriginal, |
---|
| 178 | G4ReactionProduct ¤tParticle, |
---|
| 179 | G4ReactionProduct &targetParticle, |
---|
| 180 | G4FastVector<G4ReactionProduct,GHADLISTSIZE> &vec, |
---|
| 181 | G4int &vecLen ); |
---|
| 182 | |
---|
| 183 | G4double normal(); |
---|
| 184 | |
---|
| 185 | G4int Poisson( G4double x ); |
---|
| 186 | |
---|
| 187 | }; |
---|
| 188 | |
---|
| 189 | #endif |
---|
| 190 | |
---|