Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LEPTSElasticModel.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//
27
28//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
30 : G4VLEPTSModel( modelName )
31{
33} // constructor
34
35
36//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
38}
39
40
41//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43 const G4DataVector&)
44{
45 Init();
46 BuildPhysicsTable( *aParticle );
47
48 fParticleChangeForGamma = GetParticleChangeForGamma();
49
50 // static const G4double proton_mass_c2 = 938.272013 * MeV;
51 // static const G4double neutron_mass_c2 = 939.56536 * MeV;
52 // static const G4double h2o_mass_c2 = 8*neutron_mass_c2 + 10*(proton_mass_c2 + electron_mass_c2);
53 // G4cout << "mme " << h2o_mass_c2/MeV << " " << H2o_mass_c2/MeV << G4endl;
54
55 const G4MaterialTable * materialTable = G4Material::GetMaterialTable() ;
56 std::vector<G4Material*>::const_iterator matite;
57 for( matite = materialTable->begin(); matite != materialTable->end(); matite++ ) {
58 const G4Material * aMaterial = (*matite);
59 theMassTarget[aMaterial] = theMolecularMass[aMaterial] / (6.02214179e+23/CLHEP::mole) *CLHEP::c_light * CLHEP::c_light;
60 theMassProjectile[aMaterial] = CLHEP::electron_mass_c2;
61
62 if( verboseLevel >= 1) G4cout << "Material: " << aMaterial->GetName() << " MolecularMass: " << theMolecularMass[aMaterial]/(CLHEP::g/CLHEP::mole) << " g/mole "
63 << " MTarget: " << theMassTarget[aMaterial]/CLHEP::MeV << " MeV" << G4endl;
64 }
65
66
67}
68
69//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
71 const G4ParticleDefinition* aParticle,
72 G4double kineticEnergy,
75{
76 if( kineticEnergy < theLowestEnergyLimit ) return DBL_MAX;
77 return 1./GetMeanFreePath( mate, aParticle, kineticEnergy );
78
79}
80
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83void G4LEPTSElasticModel::SampleSecondaries(std::vector<G4DynamicParticle*>*,
84 const G4MaterialCutsCouple* mateCuts,
85 const G4DynamicParticle* aDynamicParticle,
88{
89 G4double P0KinEn = aDynamicParticle->GetKineticEnergy();
90 G4ThreeVector P0Dir = aDynamicParticle->GetMomentumDirection();
91
92 if( P0KinEn < theLowestEnergyLimit ) {
93 fParticleChangeForGamma->ProposeMomentumDirection( P0Dir );
94 fParticleChangeForGamma->SetProposedKineticEnergy( 0.);
95 fParticleChangeForGamma->ProposeLocalEnergyDeposit( P0KinEn);
96 fParticleChangeForGamma->ProposeTrackStatus(fStopAndKill);
97 if( verboseLevel > 2 ) G4cout << " ENERGY LOW " << P0KinEn - theLowestEnergyLimit << G4endl;
98 return;
99 }
100
101 //- G4ParticleDefinition * particleDefDef = aTrack.GetDefinition();
102 //- G4String partName = particleDefDef->GetParticleName();
103
104 // G4ThreeVector pos, pos0, dpos;
105
106 //- G4StepPoG4int * PostPoG4int = aStep.GetPostStepPoG4int();
107 //- G4ThreeVector r = PostPoG4int->GetPosition();
108
109 //TypeOfInteraction=-10;
110
111 const G4Material* aMaterial = mateCuts->GetMaterial();
112 G4double ang = SampleAngle(aMaterial, P0KinEn/CLHEP::eV, 0.0);
113 G4ThreeVector P1Dir = SampleNewDirection(P0Dir, ang);
114#ifdef DEBUG_LEPTS
115 if( verboseLevel >= 2 ) G4cout << " G4LEPTSElasticModel::SampleSecondaries( P1Dir " << P1Dir << " P0Dir " << P0Dir << " ang " << ang << G4endl;
116#endif
117
118 //G4ThreeVector P1Dir = SampleNewDirection(P0Dir, P0KinEn/eV, 0.0);
119 //G4double Energylost1= ElasticEnergyTransferWater2(P0KinEn, ang);
120 G4double Energylost = EnergyTransfer(P0KinEn, ang, theMassTarget[aMaterial], theMassProjectile[aMaterial]);
121 if( verboseLevel >= 3 ) G4cout << " ELASTIC Energylost "<< Energylost << " = " << P0KinEn << " " <<ang << " " << theMassTarget[aMaterial] << " " << theMassProjectile[aMaterial] << G4endl;
122
123 G4double P1KinEn = P0KinEn - Energylost;
124 if( verboseLevel >= 3 ) G4cout << " ELASTIC " << P1KinEn << " = " << P0KinEn << " - " << Energylost << G4endl;
125#ifdef DEBUG_LEPTS
126 if( verboseLevel >= 2 ) G4cout << " G4LEPTSElasticModel::SampleSecondaries( SetProposedKineticEnergy " << P1KinEn << " " << P0KinEn << " - " << Energylost << G4endl;
127#endif
128 fParticleChangeForGamma->ProposeMomentumDirection( P1Dir );
129 fParticleChangeForGamma->SetProposedKineticEnergy( P1KinEn);
130 fParticleChangeForGamma->ProposeLocalEnergyDeposit( Energylost);
131 //G4cout << "elasticEnergyLost: " << Energylost << G4endl;
132
133#ifdef DEBUG_LEPTS
134 if( verboseLevel >= 2 ) G4cout << " G4LEPTSElasticModel::SampleSecondaries( ProposeMomentumDirection " << fParticleChangeForGamma->GetProposedMomentumDirection() << G4endl;
135#endif
136}
137
138//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
139G4double G4LEPTSElasticModel::EnergyTransfer(G4double E, G4double ang, G4double MT, G4double MP)
140{
141 G4double co = std::cos(ang);
142 G4double si = std::sin(ang);
143
144 G4double W = ( (E+MP)*si*si + MT - co*std::sqrt(MT*MT-MP*MP*si*si) ) * E*(E+2*MP)
145 / ( std::pow((E+MP+MT),2) - E*co*co*(E+2*MP) );
146
147 //G4double W2 = 2*MP/MT*(1-co)*E;
148 //G4cout << "WWWWWWWWW: " << W/E << " " << E/W << " " << W2/W << G4endl;
149 //G4cout << "Mm " << MT/MeV << " " << MP/MeV << G4endl;
150
151 return W;
152}
153
std::vector< G4Material * > G4MaterialTable
@ fStopAndKill
double G4double
Definition: G4Types.hh:83
@ XSElastic
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4LEPTSElasticModel(const G4String &modelName="G4LEPTSElasticModel")
const G4Material * GetMaterial() const
static G4MaterialTable * GetMaterialTable()
Definition: G4Material.cc:677
const G4String & GetName() const
Definition: G4Material.hh:172
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
const G4ThreeVector & GetProposedMomentumDirection() const
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4ThreeVector SampleNewDirection(const G4Material *aMaterial, G4ThreeVector Dir, G4double e, G4double el)
std::map< const G4Material *, G4double > theMolecularMass
G4double GetMeanFreePath(const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
G4double theLowestEnergyLimit
G4double SampleAngle(const G4Material *aMaterial, G4double e, G4double el)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define W
Definition: crc32.c:84
#define DBL_MAX
Definition: templates.hh:62