| [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: Nuclear De-excitations
|
|---|
| 29 | // by V. Lara (May 1998)
|
|---|
| 30 |
|
|---|
| 31 | #ifndef G4Fragment_h
|
|---|
| 32 | #define G4Fragment_h 1
|
|---|
| 33 |
|
|---|
| 34 | #include "G4ios.hh"
|
|---|
| 35 | #include <iomanip>
|
|---|
| 36 | #include <vector>
|
|---|
| 37 |
|
|---|
| 38 | #include "globals.hh"
|
|---|
| 39 | #include "G4LorentzVector.hh"
|
|---|
| 40 | #include "G4ParticleMomentum.hh"
|
|---|
| 41 | #include "G4ThreeVector.hh"
|
|---|
| 42 | #include "G4NucleiProperties.hh"
|
|---|
| 43 | #include "G4ParticleTable.hh"
|
|---|
| 44 | #include "G4IonTable.hh"
|
|---|
| 45 | #include "Randomize.hh"
|
|---|
| 46 | #include "G4Proton.hh"
|
|---|
| 47 | #include "G4Neutron.hh"
|
|---|
| 48 | #include "G4HadronicException.hh"
|
|---|
| 49 | #include "G4HadTmpUtil.hh"
|
|---|
| 50 |
|
|---|
| 51 |
|
|---|
| 52 | class G4ParticleDefinition;
|
|---|
| 53 |
|
|---|
| 54 | class G4Fragment; // Forward deckaration
|
|---|
| 55 | typedef std::vector<G4Fragment*> G4FragmentVector;
|
|---|
| 56 |
|
|---|
| 57 | class G4Fragment
|
|---|
| 58 | {
|
|---|
| 59 | public:
|
|---|
| 60 |
|
|---|
| 61 | // ============= CONSTRUCTORS ==================
|
|---|
| 62 |
|
|---|
| 63 | // Default constructor
|
|---|
| 64 | G4Fragment();
|
|---|
| 65 |
|
|---|
| 66 | // Destructor
|
|---|
| 67 | ~G4Fragment();
|
|---|
| 68 |
|
|---|
| 69 | // Copy constructor
|
|---|
| 70 | G4Fragment(const G4Fragment &right);
|
|---|
| 71 |
|
|---|
| 72 | // Several constructors
|
|---|
| 73 |
|
|---|
| 74 | // A,Z and 4-momentum
|
|---|
| 75 | G4Fragment(const G4int A, const G4int Z, const G4LorentzVector aMomentum);
|
|---|
| 76 |
|
|---|
| 77 | // 4-momentum and pointer to G4particleDefinition (for gammas)
|
|---|
| 78 | G4Fragment(const G4LorentzVector aMomentum, G4ParticleDefinition * aParticleDefinition);
|
|---|
| 79 |
|
|---|
| 80 | // ============= OPERATORS ==================
|
|---|
| 81 |
|
|---|
| 82 | const G4Fragment & operator=(const G4Fragment &right);
|
|---|
| 83 | G4bool operator==(const G4Fragment &right) const;
|
|---|
| 84 | G4bool operator!=(const G4Fragment &right) const;
|
|---|
| 85 |
|
|---|
| 86 | friend std::ostream& operator<<(std::ostream&, const G4Fragment*);
|
|---|
| 87 | friend std::ostream& operator<<(std::ostream&, const G4Fragment&);
|
|---|
| 88 |
|
|---|
| 89 | // ============= METHODS ==================
|
|---|
| 90 |
|
|---|
| 91 | inline G4double GetA(void) const;
|
|---|
| 92 | void SetA(const G4double value);
|
|---|
| 93 |
|
|---|
| 94 | G4double GetZ(void) const;
|
|---|
| 95 | void SetZ(const G4double value);
|
|---|
| 96 |
|
|---|
| 97 | G4double GetExcitationEnergy(void) const;
|
|---|
| 98 | void SetExcitationEnergy(const G4double value);
|
|---|
| 99 |
|
|---|
| 100 | const G4LorentzVector GetMomentum(void) const;
|
|---|
| 101 | void SetMomentum(const G4LorentzVector value);
|
|---|
| 102 |
|
|---|
| 103 | const G4ThreeVector GetAngularMomentum(void) const;
|
|---|
| 104 | void SetAngularMomentum(const G4ThreeVector value);
|
|---|
| 105 |
|
|---|
| 106 | G4int GetNumberOfExcitons(void) const;
|
|---|
| 107 | // void SetNumberOfExcitons(const G4int value);
|
|---|
| 108 |
|
|---|
| 109 | G4int GetNumberOfHoles(void) const;
|
|---|
| 110 | void SetNumberOfHoles(const G4int value);
|
|---|
| 111 |
|
|---|
| 112 | G4int GetNumberOfCharged(void) const;
|
|---|
| 113 | void SetNumberOfCharged(const G4int value);
|
|---|
| 114 |
|
|---|
| 115 | G4int GetNumberOfParticles(void) const;
|
|---|
| 116 | void SetNumberOfParticles(const G4int value);
|
|---|
| 117 |
|
|---|
| 118 | inline G4ParticleDefinition * GetParticleDefinition(void) const;
|
|---|
| 119 | void SetParticleDefinition(G4ParticleDefinition * aParticleDefinition);
|
|---|
| 120 |
|
|---|
| 121 | G4double GetCreationTime(void) const;
|
|---|
| 122 | void SetCreationTime(const G4double time);
|
|---|
| 123 |
|
|---|
| 124 |
|
|---|
| 125 | // Some utility methods
|
|---|
| 126 |
|
|---|
| 127 | inline G4double GetGroundStateMass(void) const;
|
|---|
| 128 |
|
|---|
| 129 | inline G4double GetBindingEnergy(void) const;
|
|---|
| 130 |
|
|---|
| 131 | #ifdef PRECOMPOUND_TEST
|
|---|
| 132 | G4String GetCreatorModel() const { return theCreatorModel; }
|
|---|
| 133 | void SetCreatorModel(const G4String & aModel)
|
|---|
| 134 | { theCreatorModel = aModel; }
|
|---|
| 135 | #endif
|
|---|
| 136 |
|
|---|
| 137 | private:
|
|---|
| 138 |
|
|---|
| 139 | G4double CalculateExcitationEnergy(const G4LorentzVector value) const;
|
|---|
| 140 |
|
|---|
| 141 | G4ThreeVector IsotropicRandom3Vector(const G4double Magnitude = 1.0) const;
|
|---|
| 142 |
|
|---|
| 143 |
|
|---|
| 144 | // ============= DATA MEMBERS ==================
|
|---|
| 145 |
|
|---|
| 146 | G4double theA;
|
|---|
| 147 |
|
|---|
| 148 | G4double theZ;
|
|---|
| 149 |
|
|---|
| 150 | G4double theExcitationEnergy;
|
|---|
| 151 |
|
|---|
| 152 | G4LorentzVector theMomentum;
|
|---|
| 153 |
|
|---|
| 154 | G4ThreeVector theAngularMomentum;
|
|---|
| 155 |
|
|---|
| 156 | G4int numberOfParticles;
|
|---|
| 157 |
|
|---|
| 158 | G4int numberOfHoles;
|
|---|
| 159 |
|
|---|
| 160 | G4int numberOfCharged;
|
|---|
| 161 |
|
|---|
| 162 |
|
|---|
| 163 | // Gamma evaporation requeriments
|
|---|
| 164 |
|
|---|
| 165 | G4ParticleDefinition * theParticleDefinition;
|
|---|
| 166 |
|
|---|
| 167 | G4double theCreationTime;
|
|---|
| 168 |
|
|---|
| 169 | #ifdef PRECOMPOUND_TEST
|
|---|
| 170 | G4String theCreatorModel;
|
|---|
| 171 | #endif
|
|---|
| 172 | };
|
|---|
| 173 |
|
|---|
| 174 | // Class G4Fragment
|
|---|
| 175 |
|
|---|
| 176 | inline G4double G4Fragment::GetA() const
|
|---|
| 177 | {
|
|---|
| 178 | return theA;
|
|---|
| 179 | }
|
|---|
| 180 |
|
|---|
| 181 | inline void G4Fragment::SetA(const G4double value)
|
|---|
| 182 | {
|
|---|
| 183 | theA = value;
|
|---|
| 184 | }
|
|---|
| 185 |
|
|---|
| 186 | inline G4double G4Fragment::GetZ() const
|
|---|
| 187 | {
|
|---|
| 188 | return theZ;
|
|---|
| 189 | }
|
|---|
| 190 |
|
|---|
| 191 | inline void G4Fragment::SetZ(const G4double value)
|
|---|
| 192 | {
|
|---|
| 193 | theZ = value;
|
|---|
| 194 | }
|
|---|
| 195 |
|
|---|
| 196 | inline G4double G4Fragment::GetExcitationEnergy() const
|
|---|
| 197 | {
|
|---|
| 198 | // temporary fix for what seems to be
|
|---|
| 199 | // a problem with rounding errors for on-shell lorentz-vectors in CLHEP.
|
|---|
| 200 | // HPW Apr 1999 @@@@@@@
|
|---|
| 201 |
|
|---|
| 202 | if(std::abs(theExcitationEnergy)<10*eV) return 0;
|
|---|
| 203 | return theExcitationEnergy;
|
|---|
| 204 | }
|
|---|
| 205 |
|
|---|
| 206 | inline void G4Fragment::SetExcitationEnergy(const G4double )
|
|---|
| 207 | {
|
|---|
| 208 | // theExcitationEnergy = value;
|
|---|
| 209 | G4cout << "Warning: G4Fragment::SetExcitationEnergy() is a dummy method. Please, avoid to use it." << G4endl;
|
|---|
| 210 | }
|
|---|
| 211 |
|
|---|
| 212 | inline const G4LorentzVector G4Fragment::GetMomentum() const
|
|---|
| 213 | {
|
|---|
| 214 | return theMomentum;
|
|---|
| 215 | }
|
|---|
| 216 |
|
|---|
| 217 | inline void G4Fragment::SetMomentum(const G4LorentzVector value)
|
|---|
| 218 | {
|
|---|
| 219 | theMomentum = value;
|
|---|
| 220 | theExcitationEnergy = CalculateExcitationEnergy(value);
|
|---|
| 221 | }
|
|---|
| 222 |
|
|---|
| 223 | inline const G4ThreeVector G4Fragment::GetAngularMomentum() const
|
|---|
| 224 | {
|
|---|
| 225 | return theAngularMomentum;
|
|---|
| 226 | }
|
|---|
| 227 |
|
|---|
| 228 | inline void G4Fragment::SetAngularMomentum(const G4ThreeVector value)
|
|---|
| 229 | {
|
|---|
| 230 | theAngularMomentum = value;
|
|---|
| 231 | }
|
|---|
| 232 |
|
|---|
| 233 | inline G4int G4Fragment::GetNumberOfExcitons() const
|
|---|
| 234 | {
|
|---|
| 235 | return numberOfParticles + numberOfHoles;
|
|---|
| 236 | }
|
|---|
| 237 |
|
|---|
| 238 | inline void G4Fragment::SetNumberOfParticles(const G4int value)
|
|---|
| 239 | {
|
|---|
| 240 | numberOfParticles = value;
|
|---|
| 241 | }
|
|---|
| 242 |
|
|---|
| 243 | inline G4int G4Fragment::GetNumberOfHoles() const
|
|---|
| 244 | {
|
|---|
| 245 | return numberOfHoles;
|
|---|
| 246 | }
|
|---|
| 247 |
|
|---|
| 248 | inline void G4Fragment::SetNumberOfHoles(const G4int value)
|
|---|
| 249 | {
|
|---|
| 250 | numberOfHoles = value;
|
|---|
| 251 | }
|
|---|
| 252 |
|
|---|
| 253 | inline G4int G4Fragment::GetNumberOfCharged() const
|
|---|
| 254 | {
|
|---|
| 255 | return numberOfCharged;
|
|---|
| 256 | }
|
|---|
| 257 |
|
|---|
| 258 | inline void G4Fragment::SetNumberOfCharged(const G4int value)
|
|---|
| 259 | {
|
|---|
| 260 | if (value <= numberOfParticles) numberOfCharged = value;
|
|---|
| 261 | else
|
|---|
| 262 | {
|
|---|
| 263 | G4String text = "G4Fragment::SetNumberOfCharged: Number of charged particles can't be greater than number of particles";
|
|---|
| 264 | throw G4HadronicException(__FILE__, __LINE__, text);
|
|---|
| 265 | }
|
|---|
| 266 | }
|
|---|
| 267 |
|
|---|
| 268 | inline G4int G4Fragment::GetNumberOfParticles() const
|
|---|
| 269 | {
|
|---|
| 270 | return numberOfParticles;
|
|---|
| 271 | }
|
|---|
| 272 |
|
|---|
| 273 |
|
|---|
| 274 | inline G4ParticleDefinition * G4Fragment::GetParticleDefinition(void) const
|
|---|
| 275 | {
|
|---|
| 276 | return theParticleDefinition;
|
|---|
| 277 | }
|
|---|
| 278 |
|
|---|
| 279 | inline void G4Fragment::SetParticleDefinition(G4ParticleDefinition * aParticleDefinition)
|
|---|
| 280 | {
|
|---|
| 281 | theParticleDefinition = aParticleDefinition;
|
|---|
| 282 | }
|
|---|
| 283 |
|
|---|
| 284 |
|
|---|
| 285 | inline G4double G4Fragment::GetCreationTime(void) const
|
|---|
| 286 | {
|
|---|
| 287 | return theCreationTime;
|
|---|
| 288 | }
|
|---|
| 289 |
|
|---|
| 290 | inline void G4Fragment::SetCreationTime(const G4double time)
|
|---|
| 291 | {
|
|---|
| 292 | theCreationTime = time;
|
|---|
| 293 | }
|
|---|
| 294 |
|
|---|
| 295 | inline G4double G4Fragment::GetGroundStateMass(void) const
|
|---|
| 296 | {
|
|---|
| 297 | if (theA == 0) return 0.0; // photon
|
|---|
| 298 | else return G4ParticleTable::GetParticleTable()->
|
|---|
| 299 | GetIonTable()->GetIonMass(G4lrint(theZ),G4lrint(theA));
|
|---|
| 300 | }
|
|---|
| 301 |
|
|---|
| 302 | inline G4double G4Fragment::GetBindingEnergy(void) const
|
|---|
| 303 | {
|
|---|
| 304 | return -GetGroundStateMass()+(theA-theZ)*G4Neutron::Neutron()->GetPDGMass()
|
|---|
| 305 | + theZ*G4Proton::Proton()->GetPDGMass();
|
|---|
| 306 | }
|
|---|
| 307 |
|
|---|
| 308 |
|
|---|
| 309 | #endif
|
|---|
| 310 |
|
|---|
| 311 |
|
|---|