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

Last change on this file was 1340, checked in by garnier, 14 years ago

update ti head

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