source: trunk/source/processes/hadronic/models/theo_high_energy/src/G4QuasiElasticChannel.cc @ 1340

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

update ti head

File size: 7.8 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// $Id: G4QuasiElasticChannel.cc,v 1.9 2010/09/17 11:34:29 gunter Exp $
28// GEANT4 tag $Name: geant4-09-03-ref-09 $
29//
30
31// Author : Gunter Folger March 2007
32// Modified by Mikhail Kossov. Apr2009, E/M conservation: ResidualNucleus is added (ResNuc)
33// Class Description
34// Final state production model for theoretical models of hadron inelastic
35// quasi elastic scattering in geant4;
36// Class Description - End
37//
38// Modified:
39
40#include "G4QuasiElasticChannel.hh"
41
42#include "G4Fancy3DNucleus.hh"
43
44
45#include "G4HadTmpUtil.hh"  //lrint
46
47//#define debug_scatter
48
49G4QuasiElasticChannel::G4QuasiElasticChannel()
50{
51 theQuasiElastic=G4QuasiFreeRatios::GetPointer();
52}
53
54G4QuasiElasticChannel::~G4QuasiElasticChannel()
55{}
56
57G4double G4QuasiElasticChannel::GetFraction(G4Nucleus &theNucleus,
58    const G4DynamicParticle & thePrimary)
59{
60    #ifdef debug_scatter   
61      G4cout << "G4QuasiElasticChannel:: P=" << thePrimary.GetTotalMomentum()
62             << ", pPDG=" << thePrimary.GetDefinition()->GetPDGEncoding()
63             << ", Z = "  << theNucleus.GetZ_asInt())
64             << ", N = "  << theNucleus.GetN_asInt())
65             << ", A = "  << theNucleus.GetA_asInt() << G4endl;
66    #endif
67
68  std::pair<G4double,G4double> ratios;
69  ratios=theQuasiElastic->GetRatios(thePrimary.GetTotalMomentum(),
70                                    thePrimary.GetDefinition()->GetPDGEncoding(),
71                                    theNucleus.GetZ_asInt(),
72                                    theNucleus.GetN_asInt());
73    #ifdef debug_scatter   
74      G4cout << "G4QuasiElasticChannel::ratios " << ratios.first << " x " <<ratios.second
75             << "  = " << ratios.first*ratios.second << G4endl;
76    #endif
77       
78  return ratios.first*ratios.second;
79}
80
81G4KineticTrackVector * G4QuasiElasticChannel::Scatter(G4Nucleus &theNucleus,
82                                                      const G4DynamicParticle & thePrimary)
83{
84  G4int A=theNucleus.GetA_asInt();
85  G4int Z=theNucleus.GetZ_asInt();
86  //   build Nucleus and choose random nucleon to scatter with
87  the3DNucleus.Init(theNucleus.GetA_asInt(),theNucleus.GetZ_asInt());
88  const std::vector<G4Nucleon *> nucleons=the3DNucleus.GetNucleons();
89  G4double targetNucleusMass=the3DNucleus.GetMass();                  // M.K. ResNuc
90  G4LorentzVector targetNucleus4Mom(0.,0.,0.,targetNucleusMass);      // M.K. ResNuc
91  G4int index;
92  do
93  {
94    index=G4lrint((A-1)*G4UniformRand());
95  } while (index < 0 || index >= static_cast<G4int>(nucleons.size()));
96  G4ParticleDefinition * pDef= nucleons[index]->GetDefinition();
97
98  G4int resA=A-1;                                                     // M.K. ResNuc
99  G4int resZ=Z-static_cast<int>(pDef->GetPDGCharge());                // M.K. ResNuc
100  G4ParticleDefinition* resDef=G4Neutron::Neutron(); // Resolve t-p=nn problem M.K. ResNuc
101  G4double residualNucleusMass=resDef->GetPDGMass();                  // M.K. ResNuc
102  if(resZ)
103  {
104    resDef=G4ParticleTable::GetParticleTable()->FindIon(resZ,resA,0,resZ);// M.K. ResNuc
105    residualNucleusMass=resDef->GetPDGMass();                         // M.K. ResNuc
106  }
107  else {
108    residualNucleusMass*=resA;
109  }
110   #ifdef debug_scatter
111     G4cout<<"G4QElChan::Scatter: neutron - proton? A ="<<A<<", Z="<<Z<<", projName="
112           <<pDef->GetParticleName()<<G4endl;
113   #endif
114
115  G4LorentzVector pNucleon=nucleons[index]->Get4Momentum();
116  G4double residualNucleusEnergy=std::sqrt(residualNucleusMass*residualNucleusMass+
117                                           pNucleon.vect().mag2());   // M.K. ResNuc
118  pNucleon.setE(targetNucleusMass-residualNucleusEnergy);             // M.K. ResNuc
119  G4LorentzVector residualNucleus4Mom=targetNucleus4Mom-pNucleon;     // M.K. ResNuc
120 
121  std::pair<G4LorentzVector,G4LorentzVector> result;
122 
123  result=theQuasiElastic->Scatter(pDef->GetPDGEncoding(),pNucleon,
124                                  thePrimary.GetDefinition()->GetPDGEncoding(),
125  thePrimary.Get4Momentum());
126  G4LorentzVector scatteredHadron4Mom=result.second;                  // M.K. ResNuc
127  if (result.first.e() <= 0.)
128  {
129    //G4cout << "Warning - G4QuasiElasticChannel::Scatter no scattering" << G4endl;
130    //return 0;       //no scatter
131    G4LorentzVector scatteredHadron4Mom=thePrimary.Get4Momentum();   // M.K. ResNuc
132    residualNucleus4Mom=G4LorentzVector(0.,0.,0.,targetNucleusMass); // M.K. ResNuc
133    resDef=G4ParticleTable::GetParticleTable()->FindIon(Z,A,0,Z);    // M.K. ResNuc
134  }
135
136#ifdef debug_scatter
137  G4LorentzVector EpConservation=pNucleon+thePrimary.Get4Momentum() 
138                                 - result.first - result.second;
139  if (   (EpConservation.vect().mag2() > .01*MeV*MeV )
140      || (std::abs(EpConservation.e()) > 0.1 * MeV ) ) 
141  {
142    G4cout << "Warning - G4QuasiElasticChannel::Scatter E-p non conservation : "
143           << EpConservation << G4endl;
144  }   
145#endif
146
147  G4KineticTrackVector * ktv;
148  ktv=new G4KineticTrackVector();
149  G4KineticTrack * sPrim=new G4KineticTrack(thePrimary.GetDefinition(),
150                                            0.,G4ThreeVector(0), scatteredHadron4Mom);
151  ktv->push_back(sPrim);
152  if (result.first.e() > 0.)
153  {
154    G4KineticTrack * sNuc=new G4KineticTrack(pDef, 0.,G4ThreeVector(0), result.first);
155    ktv->push_back(sNuc);
156  }
157  if(resZ || resA==1) // For the only neutron or for tnuclei with Z>0    M.K. ResNuc
158  {
159    G4KineticTrack * rNuc=new G4KineticTrack(resDef,
160                           0.,G4ThreeVector(0), residualNucleus4Mom); // M.K. ResNuc
161    ktv->push_back(rNuc);                                             // M.K. ResNuc
162  }
163  else // The residual nucleus consists of only neutrons                 M.K. ResNuc
164  {
165    residualNucleus4Mom/=resA;     // Split 4-mom of A*n system equally  M.K. ResNuc
166    for(G4int in=0; in<resA; in++) // Loop over neutrons in A*n system.  M.K. ResNuc
167    {
168      G4KineticTrack* rNuc=new G4KineticTrack(resDef,
169                           0.,G4ThreeVector(0), residualNucleus4Mom); // M.K. ResNuc
170      ktv->push_back(rNuc);                                           // M.K. ResNuc
171    }
172  }
173#ifdef debug_scatter
174  G4cout<<"G4QElC::Scat: Nucleon: "<<result.first <<" mass "<<result.first.mag() << G4endl;
175  G4cout<<"G4QElC::Scat: Project: "<<result.second<<" mass "<<result.second.mag()<< G4endl;
176#endif
177  return ktv;
178}
Note: See TracBrowser for help on using the repository browser.