G4QMDNucleus.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // 081024 G4NucleiPropertiesTable:: to G4NucleiProperties::
00027 //
00028 #include <numeric>
00029 
00030 #include "G4QMDNucleus.hh"
00031 #include "G4SystemOfUnits.hh"
00032 #include "G4Proton.hh"
00033 #include "G4Neutron.hh"
00034 #include "G4NucleiProperties.hh"
00035 
00036 G4QMDNucleus::G4QMDNucleus()
00037 {
00038    G4QMDParameters* parameters = G4QMDParameters::GetInstance();
00039    hbc = parameters->Get_hbc();
00040 }
00041 
00042 
00043 
00044 G4QMDNucleus::~G4QMDNucleus()
00045 {
00046    ;
00047 }
00048 
00049 
00050 G4LorentzVector G4QMDNucleus::Get4Momentum()
00051 {
00052    G4LorentzVector p( 0 );
00053    std::vector< G4QMDParticipant* >::iterator it; 
00054    for ( it = participants.begin() ; it != participants.end() ; it++ ) 
00055       p += (*it)->Get4Momentum();   
00056 
00057    return p;
00058 }
00059 
00060 
00061 
00062 G4int G4QMDNucleus::GetMassNumber()
00063 {
00064 
00065    G4int A = 0; 
00066    std::vector< G4QMDParticipant* >::iterator it; 
00067    for ( it = participants.begin() ; it != participants.end() ; it++ ) 
00068    {
00069       if ( (*it)->GetDefinition() == G4Proton::Proton() 
00070         || (*it)->GetDefinition() == G4Neutron::Neutron() ) 
00071          A++; 
00072    }
00073 
00074    return A;
00075 }
00076 
00077 
00078 
00079 G4int G4QMDNucleus::GetAtomicNumber()
00080 {
00081    G4int Z = 0; 
00082    std::vector< G4QMDParticipant* >::iterator it; 
00083    for ( it = participants.begin() ; it != participants.end() ; it++ ) 
00084    {
00085       if ( (*it)->GetDefinition() == G4Proton::Proton() ) 
00086          Z++; 
00087    }
00088    return Z;
00089 }
00090 
00091 
00092 
00093 G4double G4QMDNucleus::GetNuclearMass()
00094 {
00095 
00096    G4double mass = G4NucleiProperties::GetNuclearMass( GetMassNumber() , GetAtomicNumber() );
00097    
00098    if ( mass == 0.0 )
00099    {
00100 
00101       G4int Z = GetAtomicNumber();
00102       G4int A = GetMassNumber();
00103       G4int N = A - Z;
00104 
00105 // Weizsacker-Bethe 
00106 
00107       G4double Av = 16*MeV; 
00108       G4double As = 17*MeV; 
00109       G4double Ac = 0.7*MeV; 
00110       G4double Asym = 23*MeV; 
00111 
00112       G4double BE = Av * A 
00113                   - As * std::pow ( G4double ( A ) , 2.0/3.0 ) 
00114                   - Ac * Z*Z/std::pow ( G4double ( A ) , 1.0/3.0 )
00115                   - Asym * ( N - Z )* ( N - Z ) / A; 
00116 
00117       mass = Z * G4Proton::Proton()->GetPDGMass() 
00118            + N * G4Neutron::Neutron()->GetPDGMass()
00119            - BE;
00120 
00121    }
00122 
00123    return mass;
00124 }
00125 
00126 
00127 
00128 void G4QMDNucleus::CalEnergyAndAngularMomentumInCM()
00129 {
00130 
00131    //G4cout << "CalEnergyAndAngularMomentumInCM " << this->GetAtomicNumber() << " " << GetMassNumber() << G4endl;
00132 
00133    G4double gamma = Get4Momentum().gamma();
00134    G4ThreeVector beta = Get4Momentum().v()/ Get4Momentum().e();
00135 
00136    G4ThreeVector pcm0( 0.0 ) ;
00137 
00138    G4int n = GetTotalNumberOfParticipant();
00139    pcm.resize( n );
00140 
00141    for ( G4int i= 0; i < n ; i++ ) 
00142    {
00143       G4ThreeVector p_i = GetParticipant( i )->GetMomentum();
00144 
00145       G4double trans = gamma / ( gamma + 1.0 ) * p_i * beta; 
00146       pcm[i] = p_i - trans*beta;
00147 
00148       pcm0 += pcm[i];
00149    }
00150 
00151    pcm0 = pcm0 / double ( n );
00152 
00153    //G4cout << "pcm0 " << pcm0 << G4endl;
00154 
00155    for ( G4int i= 0; i < n ; i++ ) 
00156    {
00157       pcm[i] += -pcm0;
00158       //G4cout << "pcm " << i << " " << pcm[i] << G4endl;
00159    }
00160 
00161 
00162    G4double tmass = 0;
00163    G4ThreeVector rcm0( 0.0 ) ;
00164    rcm.resize( n );
00165    es.resize( n );
00166 
00167    for ( G4int i= 0; i < n ; i++ ) 
00168    {
00169       G4ThreeVector ri = GetParticipant( i )->GetPosition();
00170       G4double trans = gamma / ( gamma + 1.0 ) * ri * beta; 
00171 
00172       es[i] = std::sqrt ( std::pow ( GetParticipant( i )->GetMass() , 2 ) + pcm[i]*pcm[i] );
00173 
00174       rcm[i] = ri + trans*beta;
00175 
00176       rcm0 += rcm[i]*es[i];
00177 
00178       tmass += es[i];
00179    }
00180 
00181    rcm0 = rcm0/tmass;
00182 
00183    for ( G4int i= 0; i < n ; i++ ) 
00184    {
00185       rcm[i] += -rcm0;
00186       //G4cout << "rcm " << i << " " << rcm[i] << G4endl;
00187    }
00188 
00189 // Angluar momentum
00190 
00191    G4ThreeVector rl ( 0.0 ); 
00192    for ( G4int i= 0; i < n ; i++ ) 
00193    {
00194       rl += rcm[i].cross ( pcm[i] );
00195    }
00196 
00197    jj = int ( std::sqrt ( rl*rl / hbc ) + 0.5 );
00198 
00199 
00200 // kinetic energy per nucleon in CM
00201 
00202    G4double totalMass = 0.0;
00203    for ( G4int i= 0; i < n ; i++ ) 
00204    {
00205       // following two lines are equivalent
00206       //totalMass += GetParticipant( i )->GetDefinition()->GetPDGMass()/GeV;
00207       totalMass += GetParticipant( i )->GetMass();
00208    }
00209 
00210    kineticEnergyPerNucleon = ( std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass )/n;
00211 
00212 // Total (not per nucleion ) Binding Energy 
00213    bindingEnergy =  ( std::accumulate ( es.begin() , es.end() , 0.0 ) -totalMass ) + potentialEnergy;
00214 
00215    //G4cout << "KineticEnergyPerNucleon in GeV " << kineticEnergyPerNucleon << G4endl;
00216    //G4cout << "KineticEnergySum in GeV " << std::accumulate ( es.begin() , es.end() , 0.0 ) - totalMass << G4endl;
00217    //G4cout << "PotentialEnergy in GeV " << potentialEnergy << G4endl;
00218    //G4cout << "BindingEnergy in GeV " << bindingEnergy << G4endl;
00219    //G4cout << "G4BindingEnergy in GeV " << G4NucleiProperties::GetBindingEnergy( GetAtomicNumber() , GetMassNumber() )/GeV << G4endl;
00220 
00221    excitationEnergy = bindingEnergy + G4NucleiProperties::GetBindingEnergy( GetMassNumber() , GetAtomicNumber() )/GeV;
00222    //G4cout << "excitationEnergy in GeV " << excitationEnergy << G4endl;
00223    if ( excitationEnergy < 0 ) excitationEnergy = 0.0; 
00224 
00225 }

Generated on Mon May 27 17:49:37 2013 for Geant4 by  doxygen 1.4.7