| 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 | #include "G4QMDNucleus.hh"
|
|---|
| 27 | #include "G4Proton.hh"
|
|---|
| 28 | #include "G4Neutron.hh"
|
|---|
| 29 | #include "G4NucleiProperties.hh"
|
|---|
| 30 |
|
|---|
| 31 | #include <numeric>
|
|---|
| 32 |
|
|---|
| 33 | G4QMDNucleus::G4QMDNucleus()
|
|---|
| 34 | {
|
|---|
| 35 | G4QMDParameters* parameters = G4QMDParameters::GetInstance();
|
|---|
| 36 | hbc = parameters->Get_hbc();
|
|---|
| 37 | }
|
|---|
| 38 |
|
|---|
| 39 |
|
|---|
| 40 |
|
|---|
| 41 | G4QMDNucleus::~G4QMDNucleus()
|
|---|
| 42 | {
|
|---|
| 43 | ;
|
|---|
| 44 | }
|
|---|
| 45 |
|
|---|
| 46 |
|
|---|
| 47 | G4LorentzVector G4QMDNucleus::Get4Momentum()
|
|---|
| 48 | {
|
|---|
| 49 | G4LorentzVector p( 0 );
|
|---|
| 50 | std::vector< G4QMDParticipant* >::iterator it;
|
|---|
| 51 | for ( it = participants.begin() ; it != participants.end() ; it++ )
|
|---|
| 52 | p += (*it)->Get4Momentum();
|
|---|
| 53 |
|
|---|
| 54 | return p;
|
|---|
| 55 | }
|
|---|
| 56 |
|
|---|
| 57 |
|
|---|
| 58 |
|
|---|
| 59 | G4int G4QMDNucleus::GetMassNumber()
|
|---|
| 60 | {
|
|---|
| 61 |
|
|---|
| 62 | G4int A = 0;
|
|---|
| 63 | std::vector< G4QMDParticipant* >::iterator it;
|
|---|
| 64 | for ( it = participants.begin() ; it != participants.end() ; it++ )
|
|---|
| 65 | {
|
|---|
| 66 | if ( (*it)->GetDefinition() == G4Proton::Proton()
|
|---|
| 67 | || (*it)->GetDefinition() == G4Neutron::Neutron() )
|
|---|
| 68 | A++;
|
|---|
| 69 | }
|
|---|
| 70 |
|
|---|
| 71 | return A;
|
|---|
| 72 | }
|
|---|
| 73 |
|
|---|
| 74 |
|
|---|
| 75 |
|
|---|
| 76 | G4int G4QMDNucleus::GetAtomicNumber()
|
|---|
| 77 | {
|
|---|
| 78 | G4int Z = 0;
|
|---|
| 79 | std::vector< G4QMDParticipant* >::iterator it;
|
|---|
| 80 | for ( it = participants.begin() ; it != participants.end() ; it++ )
|
|---|
| 81 | {
|
|---|
| 82 | if ( (*it)->GetDefinition() == G4Proton::Proton() )
|
|---|
| 83 | Z++;
|
|---|
| 84 | }
|
|---|
| 85 | return Z;
|
|---|
| 86 | }
|
|---|
| 87 |
|
|---|
| 88 |
|
|---|
| 89 |
|
|---|
| 90 | G4double G4QMDNucleus::GetNuclearMass()
|
|---|
| 91 | {
|
|---|
| 92 |
|
|---|
| 93 | G4double mass = G4NucleiPropertiesTable::GetNuclearMass( GetAtomicNumber() , GetMassNumber() );
|
|---|
| 94 |
|
|---|
| 95 | if ( mass == 0.0 )
|
|---|
| 96 | {
|
|---|
| 97 |
|
|---|
| 98 | G4int Z = GetAtomicNumber();
|
|---|
| 99 | G4int A = GetMassNumber();
|
|---|
| 100 | G4int N = A - Z;
|
|---|
| 101 |
|
|---|
| 102 | // Weizsacker-Bethe
|
|---|
| 103 |
|
|---|
| 104 | G4double Av = 16*MeV;
|
|---|
| 105 | G4double As = 17*MeV;
|
|---|
| 106 | G4double Ac = 0.7*MeV;
|
|---|
| 107 | G4double Asym = 23*MeV;
|
|---|
| 108 |
|
|---|
| 109 | G4double BE = Av * A
|
|---|
| 110 | - As * std::pow ( G4double ( A ) , 2.0/3.0 )
|
|---|
| 111 | - Ac * Z*Z/std::pow ( G4double ( A ) , 1.0/3.0 )
|
|---|
| 112 | - Asym * ( N - Z )* ( N - Z ) / A;
|
|---|
| 113 |
|
|---|
| 114 | mass = Z * G4Proton::Proton()->GetPDGMass()
|
|---|
| 115 | + N * G4Neutron::Neutron()->GetPDGMass()
|
|---|
| 116 | - BE;
|
|---|
| 117 |
|
|---|
| 118 | }
|
|---|
| 119 |
|
|---|
| 120 | return mass;
|
|---|
| 121 | }
|
|---|
| 122 |
|
|---|
| 123 |
|
|---|
| 124 |
|
|---|
| 125 | void G4QMDNucleus::CalEnergyAndAngularMomentumInCM()
|
|---|
| 126 | {
|
|---|
| 127 |
|
|---|
| 128 | //G4cout << "CalEnergyAndAngularMomentumInCM " << this->GetAtomicNumber() << " " << GetMassNumber() << G4endl;
|
|---|
| 129 |
|
|---|
| 130 | G4double gamma = Get4Momentum().gamma();
|
|---|
| 131 | G4ThreeVector beta = Get4Momentum().v()/ Get4Momentum().e();
|
|---|
| 132 |
|
|---|
| 133 | G4ThreeVector pcm0( 0.0 ) ;
|
|---|
| 134 |
|
|---|
| 135 | G4int n = GetTotalNumberOfParticipant();
|
|---|
| 136 | pcm.resize( n );
|
|---|
| 137 |
|
|---|
| 138 | for ( G4int i= 0; i < n ; i++ )
|
|---|
| 139 | {
|
|---|
| 140 | G4ThreeVector p_i = GetParticipant( i )->GetMomentum();
|
|---|
| 141 |
|
|---|
| 142 | G4double trans = gamma / ( gamma + 1.0 ) * p_i * beta;
|
|---|
| 143 | pcm[i] = p_i - trans*beta;
|
|---|
| 144 |
|
|---|
| 145 | pcm0 += pcm[i];
|
|---|
| 146 | }
|
|---|
| 147 |
|
|---|
| 148 | pcm0 = pcm0 / double ( n );
|
|---|
| 149 |
|
|---|
| 150 | //G4cout << "pcm0 " << pcm0 << G4endl;
|
|---|
| 151 |
|
|---|
| 152 | for ( G4int i= 0; i < n ; i++ )
|
|---|
| 153 | {
|
|---|
| 154 | pcm[i] += -pcm0;
|
|---|
| 155 | //G4cout << "pcm " << i << " " << pcm[i] << G4endl;
|
|---|
| 156 | }
|
|---|
| 157 |
|
|---|
| 158 |
|
|---|
| 159 | G4double tmass = 0;
|
|---|
| 160 | G4ThreeVector rcm0( 0.0 ) ;
|
|---|
| 161 | rcm.resize( n );
|
|---|
| 162 | es.resize( n );
|
|---|
| 163 |
|
|---|
| 164 | for ( G4int i= 0; i < n ; i++ )
|
|---|
| 165 | {
|
|---|
| 166 | G4ThreeVector ri = GetParticipant( i )->GetPosition();
|
|---|
| 167 | G4double trans = gamma / ( gamma + 1.0 ) * ri * beta;
|
|---|
| 168 |
|
|---|
| 169 | es[i] = std::sqrt ( std::pow ( GetParticipant( i )->GetMass() , 2 ) + pcm[i]*pcm[i] );
|
|---|
| 170 |
|
|---|
| 171 | rcm[i] = ri + trans*beta;
|
|---|
| 172 |
|
|---|
| 173 | rcm0 += rcm[i]*es[i];
|
|---|
| 174 |
|
|---|
| 175 | tmass += es[i];
|
|---|
| 176 | }
|
|---|
| 177 |
|
|---|
| 178 | rcm0 = rcm0/tmass;
|
|---|
| 179 |
|
|---|
| 180 | for ( G4int i= 0; i < n ; i++ )
|
|---|
| 181 | {
|
|---|
| 182 | rcm[i] += -rcm0;
|
|---|
| 183 | //G4cout << "rcm " << i << " " << rcm[i] << G4endl;
|
|---|
| 184 | }
|
|---|
| 185 |
|
|---|
| 186 | // Angluar momentum
|
|---|
| 187 |
|
|---|
| 188 | G4ThreeVector rl ( 0.0 );
|
|---|
| 189 | for ( G4int i= 0; i < n ; i++ )
|
|---|
| 190 | {
|
|---|
| 191 | rl += rcm[i].cross ( pcm[i] );
|
|---|
| 192 | }
|
|---|
| 193 |
|
|---|
| 194 | jj = int ( std::sqrt ( rl*rl / hbc ) + 0.5 );
|
|---|
| 195 |
|
|---|
| 196 |
|
|---|
| 197 | // kinetic energy per nucleon in CM
|
|---|
| 198 |
|
|---|
| 199 | G4double totalMass = 0.0;
|
|---|
| 200 | for ( G4int i= 0; i < n ; i++ )
|
|---|
| 201 | {
|
|---|
| 202 | // following two lines are equivalent
|
|---|
| 203 | //totalMass += GetParticipant( i )->GetDefinition()->GetPDGMass()/GeV;
|
|---|
| 204 | totalMass += GetParticipant( i )->GetMass();
|
|---|
| 205 | }
|
|---|
| 206 |
|
|---|
| 207 | kineticEnergyPerNucleon = ( std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass )/n;
|
|---|
| 208 |
|
|---|
| 209 | // Total (not per nucleion ) Binding Energy
|
|---|
| 210 | bindingEnergy = ( std::accumulate ( es.begin() , es.end() , 0.0 ) -totalMass ) + potentialEnergy;
|
|---|
| 211 |
|
|---|
| 212 | //G4cout << "KineticEnergyPerNucleon in GeV " << kineticEnergyPerNucleon << G4endl;
|
|---|
| 213 | //G4cout << "KineticEnergySum in GeV " << std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass << G4endl;
|
|---|
| 214 | //G4cout << "PotentialEnergy in GeV " << potentialEnergy << G4endl;
|
|---|
| 215 | //G4cout << "BindingEnergy in GeV " << bindingEnergy << G4endl;
|
|---|
| 216 | //G4cout << "G4BindingEnergy in GeV " << G4NucleiPropertiesTable::GetBindingEnergy( GetAtomicNumber() , GetMassNumber() )/GeV << G4endl;
|
|---|
| 217 |
|
|---|
| 218 | excitationEnergy = bindingEnergy + G4NucleiPropertiesTable::GetBindingEnergy( GetAtomicNumber() , GetMassNumber() )/GeV;
|
|---|
| 219 | //G4cout << "excitationEnergy in GeV " << excitationEnergy << G4endl;
|
|---|
| 220 | if ( excitationEnergy < 0 ) excitationEnergy = 0.0;
|
|---|
| 221 |
|
|---|
| 222 | }
|
|---|