source: trunk/source/processes/hadronic/models/chiral_inv_phase_space/interface/include/G4ElectroNuclearReaction.hh @ 836

Last change on this file since 836 was 819, checked in by garnier, 16 years ago

import all except CVS

File size: 10.1 KB
Line 
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// $Id: G4ElectroNuclearReaction.hh,v 1.23 2006/06/29 20:07:46 gunter Exp $
29// GEANT4 tag $Name:  $
30//
31//
32// GEANT4 physics class: G4ElectroNuclearReaction -- header file
33// Created: J.P. Wellisch, 12/11/2001
34// The last update: J.P. Wellisch, 06-June-02
35//
36#ifndef G4ElectroNuclearReaction_h
37#define G4ElectroNuclearReaction_h
38
39#include "globals.hh"
40#include "G4HadronicInteraction.hh"
41#include "G4ChiralInvariantPhaseSpace.hh"
42#include "G4ElectroNuclearCrossSection.hh"
43#include "G4PhotoNuclearCrossSection.hh"
44#include "G4Electron.hh"
45#include "G4Positron.hh"
46#include "G4Gamma.hh"
47#include "G4GammaParticipants.hh"
48#include "G4QGSModel.hh"
49#include "G4TheoFSGenerator.hh"
50#include "G4GeneratorPrecompoundInterface.hh"
51#include "G4QGSMFragmentation.hh"
52#include "G4ExcitedStringDecay.hh"
53
54class G4ElectroNuclearReaction : public G4HadronicInteraction
55{
56  public: 
57    virtual ~G4ElectroNuclearReaction(){}
58    G4ElectroNuclearReaction()
59    {
60      SetMinEnergy(0*GeV);
61      SetMaxEnergy(30*TeV);
62     
63      theHEModel = new G4TheoFSGenerator;
64      theCascade = new G4GeneratorPrecompoundInterface;
65    }
66   
67    virtual G4HadFinalState* ApplyYourself(const G4HadProjectile& aTrack, 
68    G4Nucleus& aTargetNucleus);
69
70  private:
71    G4ChiralInvariantPhaseSpace theLEModel;
72    G4TheoFSGenerator * theHEModel;
73    G4GeneratorPrecompoundInterface * theCascade;
74    G4QGSModel< G4GammaParticipants > theStringModel;
75    G4QGSMFragmentation theFragmentation;
76    G4ExcitedStringDecay * theStringDecay;
77    G4ElectroNuclearCrossSection theElectronData;
78    G4PhotoNuclearCrossSection thePhotonData;
79    G4HadFinalState theResult;
80};
81
82inline G4HadFinalState* G4ElectroNuclearReaction::
83ApplyYourself(const G4HadProjectile& aTrack, G4Nucleus& aTargetNucleus)
84{
85  theResult.Clear();
86  static const G4double dM=G4Proton::Proton()->GetPDGMass()+G4Neutron::Neutron()->GetPDGMass(); // Mean double nucleon mass = m_n+m_p (@@ no binding)
87  static const G4double me=G4Electron::Electron()->GetPDGMass();      // electron mass
88  static const G4double me2=me*me;        // squared electron mass
89  static const G4double dpi=twopi; // 2*pi
90  G4DynamicParticle theTempEl(const_cast<G4ParticleDefinition *>(aTrack.GetDefinition()), 
91                              aTrack.Get4Momentum().vect());
92  const G4DynamicParticle* theElectron=&theTempEl;
93  const G4ParticleDefinition* aD = theElectron->GetDefinition();
94  if((aD != G4Electron::ElectronDefinition()) && (aD != G4Positron::PositronDefinition()))
95    throw G4HadronicException(__FILE__, __LINE__, "G4ElectroNuclearReaction::ApplyYourself called for neither electron or positron");
96 
97  const G4ElementTable* aTab = G4Element::GetElementTable();
98  G4Element * anElement = 0;
99  G4int aZ = static_cast<G4int>(aTargetNucleus.GetZ()+.1);
100  for(size_t ii=0; ii<aTab->size(); ii++)
101  {
102    if ( std::abs((*aTab)[ii]->GetZ()-aZ) < .1)
103    {
104      anElement = (*aTab)[ii];
105      break;
106    }
107  }
108  if(0==anElement) 
109  {
110    G4cerr<<"***G4ElectroNuclearReaction::ApplyYourself: element with Z="<<aTargetNucleus.GetZ()<<
111          " is not in the element table"<<G4endl; // @@ how to retrieve A or N for the isotop?
112    throw G4HadronicException(__FILE__, __LINE__, "Anomalous element error.");
113  }
114
115  // Note: high energy gamma nuclear now implemented.
116 
117  G4double xSec = theElectronData.GetCrossSection(theElectron, anElement); // Check cross section
118  if(xSec<=0.) 
119  {
120    theResult.SetStatusChange(isAlive);
121    theResult.SetEnergyChange(theElectron->GetKineticEnergy());
122    theResult.SetMomentumChange(theElectron->GetMomentumDirection()); // new direction for the electron
123    return &theResult;        // DO-NOTHING condition
124  }
125  G4double photonEnergy = theElectronData.GetEquivalentPhotonEnergy();
126  G4double theElectronKinEnergy=theElectron->GetKineticEnergy();
127  if( theElectronKinEnergy < photonEnergy )
128  {
129    G4cout << "G4ElectroNuclearReaction::ApplyYourself: photonEnergy is very high"<<G4endl;
130    G4cout << "If this condition appears frequently, please contact Hans-Peter.Wellisch@cern.ch"<<G4endl;
131    theResult.SetStatusChange(isAlive);
132    theResult.SetEnergyChange(theElectron->GetKineticEnergy());
133    theResult.SetMomentumChange(theElectron->GetMomentumDirection()); // new direction for the electron
134    return &theResult;        // DO-NOTHING condition
135  }
136  G4double photonQ2 = theElectronData.GetEquivalentPhotonQ2(photonEnergy);
137  G4double W=photonEnergy-photonQ2/dM;   // Hadronic energy flow (W-energy) from the virtual photon
138  if(getenv("debug_G4ElectroNuclearReaction") )
139  {
140    G4cout << "G4ElectroNuclearReaction: Equivalent Energy = "<<W<<G4endl;
141  }
142  if(W<0.) 
143  {
144    theResult.SetStatusChange(isAlive);
145    theResult.SetEnergyChange(theElectron->GetKineticEnergy());
146    theResult.SetMomentumChange(theElectron->GetMomentumDirection()); // new direction for the electron
147    return &theResult;        // DO-NOTHING condition
148    throw G4HadronicException(__FILE__, __LINE__, "G4ElectroNuclearReaction::ApplyYourself: negative equivalent energy");
149  }
150  G4DynamicParticle* theDynamicPhoton = new 
151                     G4DynamicParticle(G4Gamma::GammaDefinition(), 
152                     G4ParticleMomentum(1.,0.,0.), photonEnergy*MeV);            //->-*
153  G4double sigNu=thePhotonData.GetCrossSection(theDynamicPhoton, anElement); //       |
154  theDynamicPhoton->SetKineticEnergy(W); // Redefine photon with equivalent energy    |
155  G4double sigK =thePhotonData.GetCrossSection(theDynamicPhoton, anElement); //       |
156  delete theDynamicPhoton; // <-------------------------------------------------------*
157  G4double rndFraction = theElectronData.GetVirtualFactor(photonEnergy, photonQ2);
158  if(sigNu*G4UniformRand()>sigK*rndFraction) 
159  {
160    theResult.SetStatusChange(isAlive);
161    theResult.SetEnergyChange(theElectron->GetKineticEnergy());
162    theResult.SetMomentumChange(theElectron->GetMomentumDirection()); // new direction for the electron
163    return &theResult; // DO-NOTHING condition
164  }
165  theResult.SetStatusChange(isAlive);
166  // Scatter an electron and make gamma+A reaction
167  G4double iniE=theElectronKinEnergy+me; // Initial total energy of electron
168  G4double finE=iniE-photonEnergy;       // Final total energy of electron
169  theResult.SetEnergyChange(std::max(0.,finE-me));    // Modifies the KINETIC ENERGY (Why not in the name?)
170  G4double EEm=iniE*finE-me2;            // Just an intermediate value to avoid "2*"
171  G4double iniP=std::sqrt(iniE*iniE-me2);     // Initial momentum of the electron
172  G4double finP=std::sqrt(finE*finE-me2);     // Final momentum of the electron
173  G4double cost=(EEm+EEm-photonQ2)/iniP/finP; // std::cos(theta) for the electron scattering
174  if(cost>1.) cost=1.;
175  if(cost<-1.) cost=-1.;
176  G4ThreeVector dir=theElectron->GetMomentumDirection(); // Direction of primary electron
177  G4ThreeVector ort=dir.orthogonal();    // Not normed orthogonal vector (!) (to dir)
178  G4ThreeVector ortx = ort.unit();       // First unit vector orthogonal to the direction
179  G4ThreeVector orty = dir.cross(ortx);  // Second unit vector orthoganal to the direction
180  G4double sint=std::sqrt(1.-cost*cost);      // Perpendicular component
181  G4double phi=dpi*G4UniformRand();      // phi of scattered electron
182  G4double sinx=sint*std::sin(phi);           // x-component
183  G4double siny=sint*std::cos(phi);           // y-component
184  G4ThreeVector findir=cost*dir+sinx*ortx+siny*orty;
185  theResult.SetMomentumChange(findir); // new direction for the electron
186  G4ThreeVector photonMomentum=iniP*dir-finP*findir;
187  G4DynamicParticle localGamma(G4Gamma::GammaDefinition(), photonEnergy, photonMomentum);
188  //G4DynamicParticle localGamma(G4Gamma::GammaDefinition(), photonDirection, photonEnergy);
189  //G4DynamicParticle localGamma(G4Gamma::GammaDefinition(), photonLorentzVector);
190  G4ThreeVector position(0,0,0);
191  G4HadProjectile localTrack(localGamma);
192  G4HadFinalState * result;
193  if(photonEnergy < 3*GeV) 
194  {
195    result = theLEModel.ApplyYourself(localTrack, aTargetNucleus, &theResult);
196  } 
197  else 
198  {
199    // G4cout << "0) Getting a high energy electro-nuclear reaction"<<G4endl;
200    theHEModel->SetTransport(theCascade);
201    theHEModel->SetHighEnergyGenerator(&theStringModel);
202    theStringDecay = new G4ExcitedStringDecay(&theFragmentation);
203    theStringModel.SetFragmentationModel(theStringDecay);
204    theHEModel->SetMinEnergy(2.5*GeV);
205    theHEModel->SetMaxEnergy(100*TeV);
206
207    G4HadFinalState * aResult = theHEModel->ApplyYourself(localTrack, aTargetNucleus);
208    for(G4int all = 0; all < aResult->GetNumberOfSecondaries(); all++)
209    {
210      theResult.AddSecondary(aResult->GetSecondary(all));
211    }
212    aResult->SecondariesAreStale();
213    result = &theResult;
214  }
215  return result;
216}
217
218#endif
Note: See TracBrowser for help on using the repository browser.