00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036 #ifndef G4ElectroNuclearReaction_h
00037 #define G4ElectroNuclearReaction_h 1
00038
00039 #include <iostream>
00040 #include <CLHEP/Units/PhysicalConstants.h>
00041
00042 #include "globals.hh"
00043 #include "G4HadronicInteraction.hh"
00044 #include "G4ChiralInvariantPhaseSpace.hh"
00045 #include "G4ElectroNuclearCrossSection.hh"
00046 #include "G4PhotoNuclearCrossSection.hh"
00047 #include "G4GammaParticipants.hh"
00048 #include "G4QGSModel.hh"
00049 #include "G4QGSMFragmentation.hh"
00050 #include "G4Nucleus.hh"
00051 #include "G4HadFinalState.hh"
00052 #include "G4HadProjectile.hh"
00053 #include "G4Electron.hh"
00054 #include "G4Positron.hh"
00055 #include "G4Gamma.hh"
00056 #include "G4TheoFSGenerator.hh"
00057 #include "G4GeneratorPrecompoundInterface.hh"
00058 #include "G4ExcitedStringDecay.hh"
00059
00060
00061 class G4ElectroNuclearReaction : public G4HadronicInteraction
00062 {
00063 public:
00064 G4ElectroNuclearReaction():G4HadronicInteraction("CHIPSElectroNuclear")
00065 {
00066 SetMinEnergy(0*CLHEP::GeV);
00067 SetMaxEnergy(100*CLHEP::TeV);
00068
00069 theHEModel = new G4TheoFSGenerator;
00070 theCascade = new G4GeneratorPrecompoundInterface;
00071 theHEModel->SetTransport(theCascade);
00072 theHEModel->SetHighEnergyGenerator(&theStringModel);
00073 theStringDecay = new G4ExcitedStringDecay(&theFragmentation);
00074 theStringModel.SetFragmentationModel(theStringDecay);
00075 theHEModel->SetMinEnergy(2.5*CLHEP::GeV);
00076 theHEModel->SetMaxEnergy(100*CLHEP::TeV);
00077 theElectronData = new G4ElectroNuclearCrossSection;
00078 thePhotonData = new G4PhotoNuclearCrossSection;
00079 Description();
00080 }
00081
00082 ~G4ElectroNuclearReaction()
00083 {
00084 delete theHEModel;
00085 delete theCascade;
00086 delete theStringDecay;
00087 delete theElectronData;
00088 delete thePhotonData;
00089 }
00090
00091 G4HadFinalState*
00092 ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aTargetNucleus);
00093
00094 void Description() const
00095 {
00096 char* dirName = getenv("G4PhysListDocDir");
00097 if (dirName) {
00098 std::ofstream outFile;
00099 G4String outFileName = GetModelName() + ".html";
00100 G4String pathName = G4String(dirName) + "/" + outFileName;
00101
00102 outFile.open(pathName);
00103 outFile << "<html>\n";
00104 outFile << "<head>\n";
00105
00106 outFile << "<title>Description of CHIPS ElectroNuclear Model</title>\n";
00107 outFile << "</head>\n";
00108 outFile << "<body>\n";
00109
00110 outFile << "G4ElectroNuclearReaction handles the inelastic scattering\n"
00111 << "of e- and e+ from nuclei using the Chiral Invariant Phase\n"
00112 << "Space (CHIPS) model of M. Kossov. This model uses the\n"
00113 << "Equivalent Photon Approximation in which the incoming\n"
00114 << "electron generates a virtual photon at the electromagnetic\n"
00115 << "vertex, and the virtual photon is converted to a real photon\n"
00116 << "before it interacts with the nucleus. The real photon\n"
00117 << "interacts with the hadrons in the target by producing\n"
00118 << "quasmons (or generalized excited hadrons) which then decay\n"
00119 << "into final state hadrons. This model is valid for e- and\n"
00120 << "e+ of all incident energies.\n";
00121
00122 outFile << "</body>\n";
00123 outFile << "</html>\n";
00124 outFile.close();
00125 }
00126 }
00127
00128 private:
00129
00130 G4ChiralInvariantPhaseSpace theLEModel;
00131 G4TheoFSGenerator* theHEModel;
00132 G4GeneratorPrecompoundInterface* theCascade;
00133 G4QGSModel<G4GammaParticipants> theStringModel;
00134 G4QGSMFragmentation theFragmentation;
00135 G4ExcitedStringDecay* theStringDecay;
00136 G4ElectroNuclearCrossSection* theElectronData;
00137 G4PhotoNuclearCrossSection* thePhotonData;
00138 G4HadFinalState theResult;
00139 };
00140
00141 inline
00142 G4HadFinalState* G4ElectroNuclearReaction::ApplyYourself(const G4HadProjectile& aTrack,
00143 G4Nucleus& aTargetNucleus)
00144 {
00145 theResult.Clear();
00146 static const G4double dM=G4Proton::Proton()->GetPDGMass() +
00147 G4Neutron::Neutron()->GetPDGMass();
00148 static const G4double me=G4Electron::Electron()->GetPDGMass();
00149 static const G4double me2=me*me;
00150 G4DynamicParticle theTempEl(const_cast<G4ParticleDefinition *>(aTrack.GetDefinition()),
00151 aTrack.Get4Momentum().vect());
00152 const G4DynamicParticle* theElectron=&theTempEl;
00153 const G4ParticleDefinition* aD = theElectron->GetDefinition();
00154 if((aD != G4Electron::ElectronDefinition()) && (aD != G4Positron::PositronDefinition()))
00155 throw G4HadronicException(__FILE__, __LINE__,
00156 "G4ElectroNuclearReaction::ApplyYourself called for neither electron or positron");
00157 const G4ElementTable* aTab = G4Element::GetElementTable();
00158 G4Element * anElement = 0;
00159 G4int aZ = static_cast<G4int>(aTargetNucleus.GetZ_asInt()+.1);
00160 for(size_t ii=0; ii<aTab->size(); ++ii) if ( std::abs((*aTab)[ii]->GetZ()-aZ) < .1)
00161 {
00162 anElement = (*aTab)[ii];
00163 break;
00164 }
00165 if(0==anElement)
00166 {
00167 G4cerr<<"***G4ElectroNuclearReaction::ApplyYourself: element with Z="
00168 <<aTargetNucleus.GetZ_asInt()<<" is not in the element table"<<G4endl;
00169 throw G4HadronicException(__FILE__, __LINE__, "Anomalous element error.");
00170 }
00171
00172
00173 G4double xSec = theElectronData->GetCrossSection(theElectron, anElement);
00174 if(xSec<=0.)
00175 {
00176 theResult.SetStatusChange(isAlive);
00177 theResult.SetEnergyChange(theElectron->GetKineticEnergy());
00178
00179 theResult.SetMomentumChange(theElectron->GetMomentumDirection());
00180 return &theResult;
00181 }
00182 G4double photonEnergy = theElectronData->GetEquivalentPhotonEnergy();
00183 G4double theElectronKinEnergy=theElectron->GetKineticEnergy();
00184 if( theElectronKinEnergy < photonEnergy )
00185 {
00186 G4cout<<"G4ElectroNuclearReaction::ApplyYourself: photonEnergy is very high"<<G4endl;
00187 G4cout<<">>> If this condition persists, please contact Geant4 group"<<G4endl;
00188 theResult.SetStatusChange(isAlive);
00189 theResult.SetEnergyChange(theElectron->GetKineticEnergy());
00190
00191 theResult.SetMomentumChange(theElectron->GetMomentumDirection());
00192 return &theResult;
00193 }
00194 G4double photonQ2 = theElectronData->GetEquivalentPhotonQ2(photonEnergy);
00195 G4double W=photonEnergy-photonQ2/dM;
00196 if(getenv("debug_G4ElectroNuclearReaction") )
00197 {
00198 G4cout << "G4ElectroNuclearReaction: Equivalent Energy = "<<W<<G4endl;
00199 }
00200 if(W<0.)
00201 {
00202 theResult.SetStatusChange(isAlive);
00203 theResult.SetEnergyChange(theElectron->GetKineticEnergy());
00204
00205 theResult.SetMomentumChange(theElectron->GetMomentumDirection());
00206 return &theResult;
00207
00208
00209 }
00210 G4DynamicParticle* theDynamicPhoton = new
00211 G4DynamicParticle(G4Gamma::GammaDefinition(),
00212 G4ParticleMomentum(1.,0.,0.), photonEnergy*CLHEP::MeV);
00213 G4double sigNu=thePhotonData->GetCrossSection(theDynamicPhoton, anElement);
00214 theDynamicPhoton->SetKineticEnergy(W);
00215 G4double sigK =thePhotonData->GetCrossSection(theDynamicPhoton, anElement);
00216 delete theDynamicPhoton;
00217 G4double rndFraction = theElectronData->GetVirtualFactor(photonEnergy, photonQ2);
00218 if(sigNu*G4UniformRand()>sigK*rndFraction)
00219 {
00220 theResult.SetStatusChange(isAlive);
00221 theResult.SetEnergyChange(theElectron->GetKineticEnergy());
00222
00223 theResult.SetMomentumChange(theElectron->GetMomentumDirection());
00224 return &theResult;
00225 }
00226 theResult.SetStatusChange(isAlive);
00227
00228 G4double iniE=theElectronKinEnergy+me;
00229 G4double finE=iniE-photonEnergy;
00230 theResult.SetEnergyChange(std::max(0.,finE-me));
00231 G4double EEm=iniE*finE-me2;
00232 G4double iniP=std::sqrt(iniE*iniE-me2);
00233 G4double finP=std::sqrt(finE*finE-me2);
00234 G4double cost=(EEm+EEm-photonQ2)/iniP/finP;
00235 if(cost> 1.) cost= 1.;
00236 if(cost<-1.) cost=-1.;
00237 G4ThreeVector dir=theElectron->GetMomentumDirection();
00238 G4ThreeVector ort=dir.orthogonal();
00239 G4ThreeVector ortx = ort.unit();
00240 G4ThreeVector orty = dir.cross(ortx);
00241 G4double sint=std::sqrt(1.-cost*cost);
00242 G4double phi=CLHEP::twopi*G4UniformRand();
00243 G4double sinx=sint*std::sin(phi);
00244 G4double siny=sint*std::cos(phi);
00245 G4ThreeVector findir=cost*dir+sinx*ortx+siny*orty;
00246 theResult.SetMomentumChange(findir);
00247 G4ThreeVector photonMomentum=iniP*dir-finP*findir;
00248 G4DynamicParticle localGamma(G4Gamma::GammaDefinition(), photonEnergy, photonMomentum);
00249
00250
00251 G4ThreeVector position(0,0,0);
00252 G4HadProjectile localTrack(localGamma);
00253 G4HadFinalState * result;
00254 if (photonEnergy < 3*CLHEP::GeV)
00255 result = theLEModel.ApplyYourself(localTrack, aTargetNucleus, &theResult);
00256 else {
00257
00258 G4HadFinalState * aResult = theHEModel->ApplyYourself(localTrack, aTargetNucleus);
00259 theResult.AddSecondaries(aResult);
00260 aResult->Clear();
00261 result = &theResult;
00262 }
00263 return result;
00264 }
00265
00266 #endif