Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrossSectionInelastic.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// GEANT4 Class file
30//
31//
32// File name: G4CrossSectionInelastic
33//
34// Author: Vladimir Ivanchenko
35//
36// Creation date: 19.11.2010
37// Modifications:
38//
39// Class Description:
40//
41// Wrapper for inelastic cross section build from a component
42//
43// -------------------------------------------------------------------
44//
45
49#include "G4DynamicParticle.hh"
50#include "G4Element.hh"
51#include "G4NistManager.hh"
53
55 G4int zmin, G4int zmax,
56 G4double Emin, G4double Emax)
57 : G4VCrossSectionDataSet(c->GetName()), component(c),
58 Zmin(zmin),Zmax(zmax)
59{
61 SetMinKinEnergy(Emin);
62 SetMaxKinEnergy(Emax);
63}
64
66{}
67
69 G4int Z, const G4Material*)
70{
72 return
73 (Z >= Zmin && Z <= Zmax && e >= GetMinKinEnergy() && e <= GetMaxKinEnergy());
74}
75
78 G4int Z,
79 const G4Material*)
80{
81 return component->GetInelasticElementCrossSection(p->GetDefinition(),
82 p->GetKineticEnergy(),
83 Z, nist->GetAtomicMassAmu(Z));
84}
85
87{
88 component->BuildPhysicsTable(p);
89 // For ions, the max energy of applicability of the cross sections must scale
90 // with the absolute baryonic number; however, the cross sections objects are
91 // often shared between the different types of ions (d, t, He3, alpha, and
92 // genericIon) therefore we scale by Zmax - which is safely larger than the
93 // number of nucleons of the heaviest nuclides.
95 ( std::abs( p.GetBaryonNumber() ) > 1 ? Zmax : 1 ) );
96}
97
99{
100 component->DumpPhysicsTable(p);
101}
102
103void
105{
106 component->Description(outFile);
107}
108
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
G4CrossSectionInelastic(G4VComponentCrossSection *, G4int zmin=1, G4int zmax=256, G4double Emin=0.0, G4double Emax=DBL_MAX)
virtual G4double GetElementCrossSection(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void CrossSectionDescription(std::ostream &) const
virtual G4bool IsElementApplicable(const G4DynamicParticle *, G4int Z, const G4Material *mat=0)
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4HadronicParameters * Instance()
static G4NistManager * Instance()
G4double GetAtomicMassAmu(const G4String &symb) const
virtual void DumpPhysicsTable(const G4ParticleDefinition &)
virtual void Description(std::ostream &) const
G4double GetInelasticElementCrossSection(const G4ParticleDefinition *, G4double kinEnergy, const G4Element *)
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
void SetMaxKinEnergy(G4double value)
void SetMinKinEnergy(G4double value)