source: trunk/source/processes/hadronic/models/coherent_elastic/src/G4WHadronElasticProcess.cc @ 1315

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

update geant4-09-04-beta-cand-01 interfaces-V09-03-09 vis-V09-03-08

File size: 7.5 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: G4WHadronElasticProcess.cc,v 1.3 2010/01/13 15:42:06 mkossov Exp $
27// GEANT4 tag $Name: geant4-09-04-beta-cand-01 $
28//
29// Geant4 Hadron Elastic Scattering Process
30//
31// Created 21 April 2006 V.Ivanchenko
32// 
33// Modified:
34// 24.04.06 V.Ivanchenko add neutron scattering on hydrogen from CHIPS
35// 07.06.06 V.Ivanchenko fix problem of rotation of final state
36// 25.07.06 V.Ivanchenko add 19 MeV low energy for CHIPS
37// 26.09.06 V.Ivanchenko add lowestEnergy
38// 20.10.06 V.Ivanchenko initialise lowestEnergy=0 for neitrals, eV for charged
39// 23.01.07 V.Ivanchenko add cross section interfaces with Z and A
40// 02.05.07 V.Ivanchenko add He3
41// 13.01.10: M.Kosov: Commented not used G4QElasticCrossSection & G4QCHIPSWorld
42//
43
44#include "G4WHadronElasticProcess.hh"
45#include "globals.hh"
46#include "G4CrossSectionDataStore.hh"
47#include "G4HadronElasticDataSet.hh"
48#include "G4VQCrossSection.hh"
49//#include "G4QElasticCrossSection.hh"
50//#include "G4QCHIPSWorld.hh"
51#include "G4Element.hh"
52#include "G4ElementVector.hh"
53#include "G4IsotopeVector.hh"
54#include "G4Neutron.hh"
55#include "G4ProductionCutsTable.hh"
56 
57G4WHadronElasticProcess::G4WHadronElasticProcess(const G4String& pName)
58  : G4HadronicProcess(pName) 
59{
60  SetProcessSubType(fHadronElastic);
61  AddDataSet(new G4HadronElasticDataSet);
62  theNeutron  = G4Neutron::Neutron();
63  lowestEnergy = 1.*keV;
64  lowestEnergyNeutron = 1.e-6*eV;
65}
66
67G4WHadronElasticProcess::~G4WHadronElasticProcess()
68{
69}
70
71G4VParticleChange* G4WHadronElasticProcess::PostStepDoIt(
72                                  const G4Track& track, 
73                                  const G4Step& step)
74{
75  aParticleChange.Initialize(track);
76  G4double kineticEnergy = track.GetKineticEnergy();
77  const G4DynamicParticle* dynParticle = track.GetDynamicParticle();
78  const G4ParticleDefinition* part = dynParticle->GetDefinition();
79
80  // protection against numerical problems
81  if(part == theNeutron) {
82    if(kineticEnergy <= lowestEnergyNeutron) 
83      return G4VDiscreteProcess::PostStepDoIt(track,step);
84  } else {
85    if(kineticEnergy <= lowestEnergy)
86      return G4VDiscreteProcess::PostStepDoIt(track,step);
87  }
88
89  G4Material* material = track.GetMaterial();
90  G4CrossSectionDataStore* store = GetCrossSectionDataStore();
91  G4double xsec = store->GetCrossSection(dynParticle,material);
92  if(xsec <= DBL_MIN) return G4VDiscreteProcess::PostStepDoIt(track,step);
93
94  // Select element
95  G4Element* elm = store->SampleZandA(dynParticle,material,targetNucleus);
96
97  G4HadronicInteraction* hadi = 
98    ChooseHadronicInteraction( kineticEnergy, material, elm);
99
100  size_t idx = track.GetMaterialCutsCouple()->GetIndex();
101  G4double tcut = 
102    (*(G4ProductionCutsTable::GetProductionCutsTable()->GetEnergyCutsVector(3)))[idx];
103  hadi->SetRecoilEnergyThreshold(tcut);
104
105  // Initialize the hadronic projectile from the track
106  //  G4cout << "track " << track.GetDynamicParticle()->Get4Momentum()<<G4endl;
107  G4HadProjectile thePro(track);
108  if(verboseLevel>1) {
109    G4cout << "G4WHadronElasticProcess::PostStepDoIt for " 
110           << part->GetParticleName() 
111           << " Target Z= " << targetNucleus.GetZ() 
112           << " A= " << targetNucleus.GetZ() << G4endl; 
113  }
114  G4HadFinalState* result = hadi->ApplyYourself(thePro, targetNucleus);
115  G4ThreeVector indir = track.GetMomentumDirection();
116  G4ThreeVector outdir = (result->GetMomentumChange()).rotateUz(indir);
117 
118  if(verboseLevel>1) {
119    G4cout << "Efin= " << result->GetEnergyChange()
120           << " de= " << result->GetLocalEnergyDeposit()
121           << " nsec= " << result->GetNumberOfSecondaries()
122           << " dir= " << outdir
123           << G4endl;
124  }
125 
126  aParticleChange.ProposeEnergy(result->GetEnergyChange());
127  aParticleChange.ProposeMomentumDirection(outdir);
128  if(result->GetNumberOfSecondaries() > 0) {
129    aParticleChange.SetNumberOfSecondaries(1);
130    G4DynamicParticle* p = result->GetSecondary(0)->GetParticle();
131    G4ThreeVector pdir = p->GetMomentumDirection();
132    // G4cout << "recoil " << pdir << G4endl;
133    pdir = pdir.rotateUz(indir);
134    // G4cout << "recoil rotated " << pdir << G4endl;
135    p->SetMomentumDirection(pdir);
136    aParticleChange.AddSecondary(p);
137  } else {
138    aParticleChange.SetNumberOfSecondaries(0);
139    aParticleChange.ProposeLocalEnergyDeposit(result->GetLocalEnergyDeposit());
140  }
141  result->Clear();
142
143  return G4VDiscreteProcess::PostStepDoIt(track,step);
144}
145
146G4bool G4WHadronElasticProcess::
147IsApplicable(const G4ParticleDefinition& aParticleType)
148{
149   return (&aParticleType == G4PionPlus::PionPlus() ||
150           &aParticleType == G4PionMinus::PionMinus() ||
151           &aParticleType == G4KaonPlus::KaonPlus() ||
152           &aParticleType == G4KaonZeroShort::KaonZeroShort() ||
153           &aParticleType == G4KaonZeroLong::KaonZeroLong() ||
154           &aParticleType == G4KaonMinus::KaonMinus() ||
155           &aParticleType == G4Proton::Proton() ||
156           &aParticleType == G4AntiProton::AntiProton() ||
157           &aParticleType == G4Neutron::Neutron() ||
158           &aParticleType == G4AntiNeutron::AntiNeutron() ||
159           &aParticleType == G4Lambda::Lambda() ||
160           &aParticleType == G4AntiLambda::AntiLambda() ||
161           &aParticleType == G4SigmaPlus::SigmaPlus() ||
162           &aParticleType == G4SigmaZero::SigmaZero() ||
163           &aParticleType == G4SigmaMinus::SigmaMinus() ||
164           &aParticleType == G4AntiSigmaPlus::AntiSigmaPlus() ||
165           &aParticleType == G4AntiSigmaZero::AntiSigmaZero() ||
166           &aParticleType == G4AntiSigmaMinus::AntiSigmaMinus() ||
167           &aParticleType == G4XiZero::XiZero() ||
168           &aParticleType == G4XiMinus::XiMinus() ||
169           &aParticleType == G4AntiXiZero::AntiXiZero() ||
170           &aParticleType == G4AntiXiMinus::AntiXiMinus() ||
171           &aParticleType == G4Deuteron::Deuteron() ||
172           &aParticleType == G4Triton::Triton() ||
173           &aParticleType == G4He3::He3() ||
174           &aParticleType == G4Alpha::Alpha() ||
175           &aParticleType == G4OmegaMinus::OmegaMinus() ||
176           &aParticleType == G4AntiOmegaMinus::AntiOmegaMinus() ||
177           &aParticleType == G4GenericIon::GenericIon());
178}
179
180
Note: See TracBrowser for help on using the repository browser.