G4KM_OpticalEqRhs.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 //
00027 // -------------------------------------------------------------------
00028 //      GEANT 4 class implementation file 
00029 //
00030 //      CERN, Geneva, Switzerland
00031 //
00032 //      File name:     G4KM_OpticalEqRhs.cc
00033 //
00034 //      Author:        Alessandro Brunengo (Alessandro.Brunengo@ge.infn.it)
00035 // 
00036 //      Creation date: 5 June 2000
00037 // -------------------------------------------------------------------
00038 
00039 #include "G4KM_OpticalEqRhs.hh"
00040 #include "G4PhysicalConstants.hh"
00041 #include "G4NucleiProperties.hh"
00042 #include "G4VNuclearDensity.hh"
00043 #include "G4HadTmpUtil.hh"
00044 
00045 G4KM_OpticalEqRhs::G4KM_OpticalEqRhs(G4KM_DummyField *field,
00046                                      G4V3DNucleus * nucleus) :
00047   G4Mag_EqRhs(field), theNucleus(nucleus)
00048 {
00049   theFactor = 0;
00050   theMass = 0;
00051 }
00052 
00053 
00054 void G4KM_OpticalEqRhs::SetFactor(G4double mass, G4double opticalParameter)
00055 {
00056   G4double A = theNucleus->GetMassNumber();
00057   G4double Z = theNucleus->GetCharge();
00058   G4double bindingEnergy = G4NucleiProperties::GetBindingEnergy(G4lrint(A), G4lrint(Z));
00059   G4double nucleusMass = Z*proton_mass_c2+(A-Z)*neutron_mass_c2+bindingEnergy;
00060   G4double reducedMass = mass*nucleusMass/(mass+nucleusMass);
00061 
00062   G4double nucleonMass = (proton_mass_c2+neutron_mass_c2)/2;
00063 
00064 // _factor in (MeV*fermi)*fermi/MeV = fermi*fermi  -- need to have A as density normalized to 1
00065   theFactor = 2*pi*hbarc*hbarc*(1+mass/nucleonMass)* opticalParameter/reducedMass * A;
00066 
00067   theMass = mass;
00068 }
00069 
00070 
00071 void G4KM_OpticalEqRhs::EvaluateRhsGivenB(const G4double y[], const G4double *,
00072                                           G4double dydx[]) const
00073 {
00074   G4double yMod = std::sqrt(y[0]*y[0]+y[1]*y[1]+y[2]*y[2]);
00075   G4double e = std::sqrt(theMass*theMass+y[3]*y[3]+y[4]*y[4]+y[5]*y[5]);
00076   dydx[0] = c_light*y[3]/e;   //
00077   dydx[1] = c_light*y[4]/e;   //  dq/dt=dH/dp = c*p/e
00078   dydx[2] = c_light*y[5]/e;   // 
00079 
00080 // V=K*rho(r) ==> dydx[3] = -dV/dr*dr/dx = -K*d(rho)/dr*dr/dx.
00081 // Idem for dydx[4] and dydx[5]
00082 
00083   const G4VNuclearDensity * nuclearDensity=theNucleus->GetNuclearDensity();
00084 
00085   G4ThreeVector pos(y[0],y[1],y[2]);
00086   G4double deriv = theFactor*nuclearDensity->GetDeriv(pos);
00087 
00088   dydx[3] = yMod == 0 ? 0 : -deriv*y[0]/yMod*c_light;
00089   dydx[4] = yMod == 0 ? 0 : -deriv*y[1]/yMod*c_light;
00090   dydx[5] = yMod == 0 ? 0 : -deriv*y[2]/yMod*c_light;
00091 }
00092 
00093 
00094 
00095 

Generated on Mon May 27 17:48:43 2013 for Geant4 by  doxygen 1.4.7