Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNAVacuumModel.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#include "G4DNAVacuumModel.hh"
30
32 const G4String& applyToMaterial, const G4ParticleDefinition*, const G4String& nam)
33 : G4VDNAModel(nam, applyToMaterial)
34{
35 verboseLevel = 0;
36
37 if (verboseLevel > 0) {
38 G4cout << "G4DNAVacuumModel is constructed " << G4endl;
39 }
40}
41
42//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
43
45{
46 if (verboseLevel > 3) G4cout << "Calling G4DNAVacuumModel::Initialise()" << G4endl;
47}
48
49//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50
52 const G4ParticleDefinition* particle, const G4DataVector& /*cuts*/)
53{
54 if (verboseLevel > 3) {
55 G4cout << "Calling G4DNAVacuumModel::Initialise()" << G4endl;
56 }
57 if(G4Material::GetMaterial("G4_Galactic",false) != nullptr)
58 {
59 auto index = (G4int)G4Material::GetMaterial("G4_Galactic")->GetIndex();
60 EnableForMaterialAndParticle(index, particle);
61 }
62}
63
64//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
65
66G4double G4DNAVacuumModel::CrossSectionPerVolume(const G4Material* /*material*/, const G4ParticleDefinition* /*particleDefinition*/,
67 G4double /*ekin*/, G4double /*emin*/, G4double /*emax*/)
68{
69 if (verboseLevel > 3) {
70 G4cout << "Calling CrossSectionPerVolume() of G4DNAVacuumModel" << G4endl;
71 }
72
73 return 0;
74}
75
76//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
77
78void G4DNAVacuumModel::SampleSecondaries(std::vector<G4DynamicParticle*>* /*fvect*/,
79 const G4MaterialCutsCouple* /*couple*/,
80 const G4DynamicParticle* /*aDynamicParticle*/, G4double /*tmin*/, G4double /*tmax*/)
81{
82 if (verboseLevel > 3) {
83 G4cout << "Calling SampleSecondaries() of G4DNAVacuumModel" << G4endl;
84 }
85}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double tmax) override
SampleSecondaries.
G4DNAVacuumModel(const G4String &applyToMaterial="all", const G4ParticleDefinition *p=nullptr, const G4String &nam="DNAPTBVacuumModel")
G4DNAVacuumModel Constructor.
~G4DNAVacuumModel() override
~G4DNAVacuumModel Destructor
void Initialise(const G4ParticleDefinition *, const G4DataVector &) override
Initialise Registers the G4_Galactic material as "void material" for every particle.
G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
CrossSectionPerVolume.
std::size_t GetIndex() const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
The G4VDNAModel class.
void EnableForMaterialAndParticle(const size_t &materialID, const G4ParticleDefinition *p)
EnableMaterialAndParticle.