Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNADummyModel.hh
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// Contact authors: S. Meylan, C. Villagrasa
28//
30
31#ifndef G4DNADUMMYMODEL_HH
32#define G4DNADUMMYMODEL_HH
33
34#include "G4VDNAModel.hh"
35#include "G4VEmModel.hh"
36#include "G4Electron.hh"
37#include "G4Proton.hh"
39
41{
42public:
43 G4DNADummyModel(const G4String& applyToMaterial,
44 const G4ParticleDefinition* p,
45 const G4String& nam,
46 G4VEmModel* emModel);
48
49 virtual void Initialise(const G4ParticleDefinition* particle, const G4DataVector& = *(new G4DataVector()), G4ParticleChangeForGamma* changeForGamme=nullptr);
50
51 virtual G4double CrossSectionPerVolume(const G4Material* material,
52 const G4String& materialName,
53 const G4ParticleDefinition* p,
54 G4double ekin,
55 G4double emin,
56 G4double emax);
57
58 virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
60 const G4String& materialName,
61 const G4DynamicParticle*,
62 G4ParticleChangeForGamma *particleChangeForGamma,
63 G4double tmin,
64 G4double tmax);
65
66 const G4VEmModel* GetEmModel() const {return fpEmModel;}
67 G4VEmModel* GetEmModel() {return fpEmModel;}
68
69private:
70 G4VEmModel* fpEmModel;
71 const G4ParticleDefinition* fpParticleDef;
72 const std::vector<double>* fMaterialMolPerVol;
73
74 G4double GetNumMoleculePerVolumeUnitForMaterial(const G4Material *mat);
75};
76
77#endif // G4DNADUMMYMODEL_HH
double G4double
Definition: G4Types.hh:83
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &= *(new G4DataVector()), G4ParticleChangeForGamma *changeForGamme=nullptr)
Initialise Each model must implement an Initialize method.
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method....
const G4VEmModel * GetEmModel() const
G4VEmModel * GetEmModel()
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax)
SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created...
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50