Geant4  10.03
CexmcHadronicProcess.cc
Go to the documentation of this file.
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  *
29  * Filename: CexmcHadronicProcess.cc
30  *
31  * Description: hadronic process with production model
32  *
33  * Version: 1.0
34  * Created: 31.10.2009 23:54:38
35  * Revision: none
36  * Compiler: gcc
37  *
38  * Author: Alexey Radkov (),
39  * Company: PNPI
40  *
41  * ============================================================================
42  */
43 
44 #include <G4ParticleChange.hh>
45 #include <G4ParticleDefinition.hh>
46 #include <G4HadronicInteraction.hh>
47 #include <G4Track.hh>
48 #include <G4Step.hh>
49 #include <G4Element.hh>
50 #include <G4StableIsotopes.hh>
51 #include <G4TrackStatus.hh>
52 #include "CexmcHadronicProcess.hh"
53 #include "CexmcProductionModel.hh"
55 #include "CexmcException.hh"
56 
57 
59  G4HadronicProcess( name ), productionModel( NULL ), interaction( NULL ),
60  theTotalResult( NULL ), isInitialized( false )
61 {
63 }
64 
65 
67 {
68  delete theTotalResult;
69 }
70 
71 
73  CexmcProductionModel * model )
74 {
75  productionModel = model;
76 
77  interaction = dynamic_cast< G4HadronicInteraction * >( productionModel );
78 
79  if ( ! interaction )
81 
83 }
84 
85 
87  const G4Material * material )
88 {
89  G4int numberOfElements( material->GetNumberOfElements() );
90  if ( numberOfElements > 1 )
91  {
92  G4cout << CEXMC_LINE_START "WARNING: Number of elements in target "
93  "material is more than 1.\n Only the first "
94  "element will be chosen for target nucleus" << G4endl;
95  }
96 
97  const G4Element * element( material->GetElement( 0 ) );
98  G4double ZZ( element->GetZ() );
99  G4int Z( G4int( ZZ + 0.5 ) );
100 
101  G4StableIsotopes stableIsotopes;
102  G4int index( stableIsotopes.GetFirstIsotope( Z ) );
103  G4double AA( stableIsotopes.GetIsotopeNucleonCount( index ) );
104 
105  targetNucleus.SetParameters( AA, ZZ );
106 }
107 
108 
110  const G4Track & track )
111 {
112  G4int numberOfSecondaries( hadFinalState->GetNumberOfSecondaries() );
113 
115  theTotalResult->Initialize( track );
118  hadFinalState->GetLocalEnergyDeposit() );
119  theTotalResult->SetNumberOfSecondaries( numberOfSecondaries );
120  theTotalResult->ProposeEnergy( hadFinalState->GetEnergyChange() );
122  if ( hadFinalState->GetStatusChange() == stopAndKill )
124 
125  for ( G4int i( 0 ); i < numberOfSecondaries; ++i )
126  {
127  G4double time( hadFinalState->GetSecondary( i )->GetTime() );
128  if ( time < 0 )
129  time = track.GetGlobalTime();
130 
131  G4Track * newTrack( new G4Track(
132  hadFinalState->GetSecondary( i )->GetParticle(),
133  time, track.GetPosition() ) );
134 
135  G4double newWeight( track.GetWeight() *
136  hadFinalState->GetSecondary( i )->GetWeight() );
137  newTrack->SetWeight( newWeight );
138  newTrack->SetTouchableHandle( track.GetTouchableHandle() );
139  theTotalResult->AddSecondary( newTrack );
140  }
141 
142  hadFinalState->Clear();
143 }
144 
145 
147  const G4Step & )
148 {
149  G4TrackStatus trackStatus( track.GetTrackStatus() );
150 
151  if ( trackStatus != fAlive && trackStatus != fSuspend )
152  {
154  theTotalResult->Initialize( track );
155 
156  return theTotalResult;
157  }
158 
159  /* NB: the target nucleus is chosen only once, it means that it will always
160  * have same Z and A, practically the first stable isotope of the first
161  * element in elements vector will be chosen. This simplification prompts
162  * the user to choose simple single-element material for the target, for
163  * example liquid hydrogen. On the other hand target nucleus is supposedly
164  * only needed if user decides to turn Fermi motion on, so this
165  * simplification should not be very harmful */
166  if ( ! isInitialized )
167  {
169  isInitialized = true;
170  }
171 
172  G4HadProjectile projectile( track );
173  G4HadFinalState * result( interaction->ApplyYourself( projectile,
174  targetNucleus ) );
175  FillTotalResult( result, track );
176 
178  {
179  CexmcTrackInfo * trackInfo( static_cast< CexmcTrackInfo * >(
180  track.GetUserInformation() ) );
181 
182  if ( trackInfo &&
184  {
185  CexmcIncidentParticleTrackInfo * theTrackInfo(
186  static_cast< CexmcIncidentParticleTrackInfo * >( trackInfo ) );
187  theTrackInfo->SetNeedsTrackLengthResampling();
188  }
189  }
190 
191  return theTotalResult;
192 }
193 
194 
196  const G4ParticleDefinition & particle )
197 {
198  if ( ! productionModel )
199  return false;
200 
201  G4ParticleDefinition * incidentParticle(
203 
204  if ( ! incidentParticle )
205  return false;
206 
207  return particle == *incidentParticle;
208 }
209 
G4VParticleChange * PostStepDoIt(const G4Track &track, const G4Step &step)
G4HadronicInteraction * interaction
G4HadSecondary * GetSecondary(size_t i)
G4int GetFirstIsotope(G4int Z)
CexmcHadronicProcess(const G4String &name=CexmcStudiedProcessLastName)
const G4ThreeVector & GetPosition() const
G4TrackStatus GetTrackStatus() const
G4ParticleDefinition * GetIncidentParticle(void) const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetEnergyChange() const
const char * name(G4int ptype)
G4double GetTime() const
const G4Element * GetElement(G4int iel) const
Definition: G4Material.hh:202
void FillTotalResult(G4HadFinalState *hadFinalState, const G4Track &track)
int G4int
Definition: G4Types.hh:78
void SetWeight(G4double aValue)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void RegisterProductionModel(CexmcProductionModel *model)
G4VUserTrackInformation * GetUserInformation() const
void RegisterMe(G4HadronicInteraction *a)
void SetSecondaryWeightByProcess(G4bool)
G4bool IsApplicable(const G4ParticleDefinition &particle)
G4GLOB_DLL std::ostream G4cout
#define CEXMC_LINE_START
Definition: CexmcCommon.hh:52
bool G4bool
Definition: G4Types.hh:79
CexmcProductionModel * productionModel
void CalculateTargetNucleus(const G4Material *material)
Definition: G4Step.hh:76
G4double GetGlobalTime() const
virtual G4int GetTypeInfo(void) const
const G4TouchableHandle & GetTouchableHandle() const
G4Material * GetMaterial() const
virtual void Initialize(const G4Track &)
G4int GetIsotopeNucleonCount(G4int number)
void SetNumberOfSecondaries(G4int totSecondaries)
G4DynamicParticle * GetParticle()
void ProposeEnergy(G4double finalEnergy)
void AddSecondary(G4Track *aSecondary)
G4double GetWeight() const
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4TrackStatus GetTrackStatus() const
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4ParticleChange * theTotalResult
void SetParameters(const G4double A, const G4double Z)
Definition: G4Nucleus.cc:212
G4TrackStatus
G4bool isInitialized()
Check if the generator is initialized.
G4int GetNumberOfSecondaries() const
virtual G4HadFinalState * ApplyYourself(const G4HadProjectile &aTrack, G4Nucleus &targetNucleus)=0
G4double GetLocalEnergyDeposit() const
G4HadFinalStateStatus GetStatusChange() const
G4double GetWeight() const