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

Last change on this file since 1200 was 1197, checked in by garnier, 16 years ago

nvx fichiers dans CVS

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