Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LEPTSAttachmentModel.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 // constructor
30  : G4VLEPTSModel( modelName )
31 {
33 
34 } // constructor
35 
36 
38 }
39 
40 
41 
42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
44  const G4DataVector&)
45 {
46  Init();
47  BuildPhysicsTable( *aParticle );
48 
49  fParticleChangeForGamma = GetParticleChangeForGamma();
50 
51 }
52 
53 
54 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56  const G4ParticleDefinition* aParticle,
57  G4double kineticEnergy,
58  G4double,
59  G4double)
60 {
61  return 1./GetMeanFreePath( mate, aParticle, kineticEnergy );
62 
63 }
64 
65 
66 void G4LEPTSAttachmentModel::SampleSecondaries(std::vector<G4DynamicParticle*>* ,
67  const G4MaterialCutsCouple*,
68  const G4DynamicParticle* aDynamicParticle,
69  G4double,
70  G4double)
71 {
72  G4double P0KinEn = aDynamicParticle->GetKineticEnergy();
73 
74  fParticleChangeForGamma->SetProposedKineticEnergy(0.);
75  fParticleChangeForGamma->ProposeLocalEnergyDeposit (P0KinEn);
76  fParticleChangeForGamma->ProposeTrackStatus( fStopAndKill);
77 
78 }
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)
G4double GetKineticEnergy() const
G4LEPTSAttachmentModel(const G4String &processName="G4LEPTSAttachmentModel")
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double GetMeanFreePath(const G4Material *mate, const G4ParticleDefinition *aParticle, G4double kineticEnergy)
void BuildPhysicsTable(const G4ParticleDefinition &aParticleType)
void SetProposedKineticEnergy(G4double proposedKinEnergy)
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132