Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4VDNAPTBModel.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 #ifndef G4VDNAPTBMODEL_HH
27 #define G4VDNAPTBMODEL_HH
28 
29 #include "G4VEmModel.hh"
30 
33 #include "G4LogLogInterpolation.hh"
34 #include "G4ParticleTable.hh"
35 
36 class G4VDNAPTBModel : public G4VEmModel
37 {
38 public:
39  G4VDNAPTBModel(const G4String& nam, const G4String& applyToMaterial);
40 
41  virtual ~G4VDNAPTBModel();
42 
43  // ***********************
44  // Initialisation
45  // ***********************
46 
47  virtual void Initialise(const G4ParticleDefinition* particle,
48  const G4DataVector& cuts) =0;
49 
50  G4bool IsMaterialDefine(const G4String& materialName);
51 
52  G4bool IsParticleExistingInModel(const G4String& particleName);
53 
54  G4bool IsMaterialExistingInModelForParticle(const G4String& particleName, const G4String& materialName);
55 
56  void SetHighELimit(const G4String& material, const G4String& particle, G4double lim) {fHighEnergyLimits[particle][material]=lim;}
57  void SetLowELimit(const G4String& material, const G4String& particle, G4double lim) {fLowEnergyLimits[particle][material]=lim;}
58 
59  G4double GetHighELimit(const G4String& material, const G4String& particle) {return fHighEnergyLimits[particle][material];}
60  G4double GetLowELimit(const G4String& material, const G4String& particle) {return fLowEnergyLimits[particle][material];}
61 
62  // ***********************
63  // Runtime
64  // ***********************
65 
66  virtual G4double CrossSectionPerVolume(const G4Material* material,
67  const G4ParticleDefinition* p,
68  G4double ekin,
69  G4double emin,
70  G4double emax) = 0;
71 
72  virtual void SampleSecondaries(std::vector<G4DynamicParticle*>*,
73  const G4MaterialCutsCouple*,
74  const G4DynamicParticle*,
75  G4double tmin = 0,
76  G4double tmax = DBL_MAX) = 0;
77 
78 
79  G4double GetHighELimit(const G4Material* material) {return fHighEnergyLimitsRuntime.at(material->GetIndex() );}
80  G4double GetLowELimit(const G4Material* material) {return fLowEnergyLimitsRuntime.at(material->GetIndex() );}
81 
82  void SetHighELimit(const G4Material* material, G4double lim) {fHighEnergyLimitsRuntime[material->GetIndex()]=lim;}
83  void SetLowELimit(const G4Material* material, G4double lim) {fLowEnergyLimitsRuntime[material->GetIndex()]=lim;}
84 
85 protected:
86 
87  // ***********************
88  // Initialisation variables
89  // ***********************
90 
91  typedef std::map<G4String, std::map<G4String,G4DNACrossSectionDataSet*, std::less<G4String> > > TableMapData;
92 
94 
96 
97  struct MaterialData
98  {
99  MaterialData(const G4String& mat, const G4String& particule, const G4String& CSFile,
100  const G4String& diffCSFile, G4double scaleFactor) :
101  fMaterial(mat),
102  fParticle(particule),
103  fCSFile(CSFile),
104  fDiffCSFile(diffCSFile),
105  fScaleFactor(scaleFactor)
106  {
107 
108  }
109 
110  G4String fMaterial; // materials that can be activated (and will be by default) within the model
111  G4String fParticle; // particles that can be activated within the model
112  G4String fCSFile; // cross section data files
113  G4String fDiffCSFile; // differential corss section data files
114  G4double fScaleFactor; // model scale factors (they could change with material)
115  };
116 
117  std::vector<MaterialData> fModelMaterialData;
118 
119  // Initisation energy limits
120  std::map<G4String, std::map<G4String, G4double> > fLowEnergyLimits; // List the low energy limits
121  std::map<G4String, std::map<G4String, G4double> > fHighEnergyLimits; // List the high energy limits
122 
123  // ***********************
124  // Runtime variables
125  // ***********************
126 
127  // This vector has the same index as G4MaterialTable. If a material is within G4MaterialTable but not declared in the current model, then
128  // this vector registered a nullptr.
129  std::map<G4int, G4DNACrossSectionDataSet*> fTableDataRuntime;
130 
131  // We do not need the particule id since every model instance is associated to one particle
132  std::map<G4int, G4double> fLowEnergyLimitsRuntime;
133  std::map<G4int, G4double> fHighEnergyLimitsRuntime;
134 
135  // ***********************
136  // Methods
137  // ***********************
138 
140  G4DNACrossSectionDataSet* GetSigmaData(const G4Material* material) {return fTableDataRuntime.at(material->GetIndex() );}
141 
142  std::vector<G4String> BuildApplyToMatVect(const G4String& materials);
143 
144  void ReadAndSaveCSFile(const G4String& materialName, const G4String& particleName, const G4String& file, G4double scaleFactor);
145 
146  G4int RandomSelectShell(G4double k, const G4Material* material);
147 
148  void AddCrossSectionData(const G4String& materialName, const G4String& particleName, const G4String& fileCS, const G4String& fileDiffCS, G4double scaleFactor);
149  void AddCrossSectionData(const G4String& materialName, const G4String& particleName, const G4String& fileCS, G4double scaleFactor);
150 
151  void LoadCrossSectionData(const G4String& particleName);
152 
153  virtual void ReadDiffCSFile(const G4String& materialName,
154  const G4String& particleName,
155  const G4String& path,
156  const G4double scaleFactor);
157 
158  void EnableForMaterialAndParticle(const G4String& materialName, const G4String& particleName);
159 };
160 
161 #endif // G4VDNAPTBMODEL_HH
virtual void Initialise(const G4ParticleDefinition *particle, const G4DataVector &cuts)=0
void SetHighELimit(const G4String &material, const G4String &particle, G4double lim)
size_t GetIndex() const
Definition: G4Material.hh:262
void SetLowELimit(const G4Material *material, G4double lim)
const char * p
Definition: xmltok.h:285
std::map< G4String, std::map< G4String, G4double > > fLowEnergyLimits
G4bool IsParticleExistingInModel(const G4String &particleName)
std::map< G4int, G4double > fLowEnergyLimitsRuntime
void LoadCrossSectionData(const G4String &particleName)
G4bool IsMaterialDefine(const G4String &materialName)
int G4int
Definition: G4Types.hh:78
std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > TableMapData
G4bool IsMaterialExistingInModelForParticle(const G4String &particleName, const G4String &materialName)
void SetHighELimit(const G4Material *material, G4double lim)
MaterialData(const G4String &mat, const G4String &particule, const G4String &CSFile, const G4String &diffCSFile, G4double scaleFactor)
void AddCrossSectionData(const G4String &materialName, const G4String &particleName, const G4String &fileCS, const G4String &fileDiffCS, G4double scaleFactor)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0, G4double tmax=DBL_MAX)=0
std::vector< G4String > BuildApplyToMatVect(const G4String &materials)
bool G4bool
Definition: G4Types.hh:79
G4double GetLowELimit(const G4Material *material)
void EnableForMaterialAndParticle(const G4String &materialName, const G4String &particleName)
static const G4double emax
TableMapData * GetTableData()
G4int RandomSelectShell(G4double k, const G4Material *material)
const G4String fStringOfMaterials
G4double GetLowELimit(const G4String &material, const G4String &particle)
std::map< G4String, std::map< G4String, G4double > > fHighEnergyLimits
void ReadAndSaveCSFile(const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor)
std::map< G4int, G4double > fHighEnergyLimitsRuntime
virtual G4double CrossSectionPerVolume(const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)=0
G4double GetHighELimit(const G4Material *material)
TableMapData fTableData
virtual void ReadDiffCSFile(const G4String &materialName, const G4String &particleName, const G4String &path, const G4double scaleFactor)
std::map< G4int, G4DNACrossSectionDataSet * > fTableDataRuntime
void SetLowELimit(const G4String &material, const G4String &particle, G4double lim)
double G4double
Definition: G4Types.hh:76
std::vector< MaterialData > fModelMaterialData
#define DBL_MAX
Definition: templates.hh:83
G4DNACrossSectionDataSet * GetSigmaData(const G4Material *material)
virtual ~G4VDNAPTBModel()
G4VDNAPTBModel(const G4String &nam, const G4String &applyToMaterial)
G4double GetHighELimit(const G4String &material, const G4String &particle)