Geant4  v4-10.4-release
 모두 클래스 네임스페이스들 파일들 함수 변수 타입정의 열거형 타입 열거형 멤버 Friends 매크로 그룹들 페이지들
G4DNADummyModel.cc
이 파일의 문서화 페이지로 가기
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 //
29 // email: sylvain.meylan@symalgo-tech.com, carmen.villagrasa@irsn.fr
30 
31 #include "../include/G4DNADummyModel.hh"
32 
33 #include "G4SystemOfUnits.hh"
34 
35 G4DNADummyModel::G4DNADummyModel(const G4String& applyToMaterial, const G4ParticleDefinition* p, const G4String& nam, G4VEmModel* emModel)
36  : G4VDNAModel(nam, applyToMaterial)
37 {
38  fpEmModel = emModel;
39  fpParticleDef = p;
40 }
41 
43 {
44  // There is no need to delete the model because it will be done in some G4 class.
45  //if(fpEmModel) delete fpEmModel;
46 }
47 
49 {
51 
52  fpEmModel->SetParticleChange(changeForGamme, nullptr);
53  fpEmModel->Initialise(particle, v);
54 
55  // MatManagSys
59 }
60 
62 {
63  G4double crossSectionTimesDensity = fpEmModel->CrossSectionPerVolume(material, p, ekin, emin, emax);
64  G4double crossSection = crossSectionTimesDensity / GetNumMoleculePerVolumeUnitForMaterial(G4Material::GetMaterial("G4_WATER") );
65 
66  return crossSection;
67 }
68 
69 void G4DNADummyModel::SampleSecondaries(std::vector<G4DynamicParticle*>* a, const G4MaterialCutsCouple* b, const G4String& /*materialName*/, const G4DynamicParticle* c, G4ParticleChangeForGamma* /*particleChangeForGamma*/, G4double tmin, G4double tmax)
70 {
71  fpEmModel->SampleSecondaries(a, b, c, tmin, tmax);
72 }
73 
75 {
76  return fMaterialMolPerVol->at(mat->GetIndex() );
77 }
78 
79 
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
void SetHighELimit(const G4String &material, const G4String &particle, G4double lim)
SetHighEnergyLimit.
Definition: G4VDNAModel.hh:169
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...
G4double GetNumMoleculePerVolumeUnitForMaterial(const G4Material *mat)
size_t GetIndex() const
Definition: G4Material.hh:262
const char * p
Definition: xmltok.h:285
const G4String & GetParticleName() const
static const G4double emax
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:608
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:616
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
G4VEmModel * fpEmModel
double G4double
Definition: G4Types.hh:76
void SetLowELimit(const G4String &material, const G4String &particle, G4double lim)
SetLowEnergyLimit.
Definition: G4VDNAModel.hh:177
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &=*(new G4DataVector()), G4ParticleChangeForGamma *changeForGamme=nullptr)
Initialise Each model must implement an Initialize method.
static G4DNAMolecularMaterial * Instance()
Float_t mat
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0
G4DNADummyModel(const G4String &applyToMaterial, const G4ParticleDefinition *p, const G4String &nam, G4VEmModel *emModel)
const G4ParticleDefinition * fpParticleDef
const std::vector< double > * fMaterialMolPerVol
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:252
void EnableForMaterialAndParticle(const G4String &materialName, const G4String &particleName)
EnableMaterialAndParticle.
Definition: G4VDNAModel.cc:134
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:454
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. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.
The G4VDNAModel class.
Definition: G4VDNAModel.hh:49
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:609