Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4ChannelingFastSimModel.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// Author: Alexei Sytov
27// Co-author: Gianfranco Paternò (modifications & testing)
28// On the base of the CRYSTALRAD realization of channeling model:
29// A. I. Sytov, V. V. Tikhomirov, and L. Bandiera PRAB 22, 064601 (2019)
30
31#ifndef G4ChannelingFastSimModel_h
32#define G4ChannelingFastSimModel_h 1
33
35#include "G4Step.hh"
36#include "G4TouchableHandle.hh"
37#include <vector>
40
42#include <unordered_map>
43#include "G4BaierKatkov.hh"
44#include "G4LogicalVolume.hh"
45#include "G4ParticleTable.hh"
46
47/** \file G4ChannelingFastSimModel.hh
48* \brief Definition of the G4ChannelingFastSimModel class
49* FastSimulation Channeling model: calculates charge particle trajectories
50* in oriented crystals in the field of crystal planes/axes either straight or bent.
51* It is also possible to simulate radiation using Baier-Katkov method.
52*/
53
55{
56public:
57 // Constructor, destructor
61
62 /// -- IsApplicable
64 /// -- ModelTrigger
65 G4bool ModelTrigger(const G4FastTrack &) override;
66 /// -- User method DoIt
67 void DoIt(const G4FastTrack&, G4FastStep&) override;
68
69 ///special functions
70 void Input(const G4Material* crystal,
71 const G4String &lattice)
72 {Input(crystal,lattice,"");}
73
74 void Input(const G4Material* crystal,
75 const G4String &lattice,
76 const G4String &filePath);
77
79
81
82 G4BaierKatkov* GetRadiationModel() {return fBaierKatkov;}
83
85
86 ///set cuts
87 void SetLowKineticEnergyLimit(G4double ekinetic, const G4String& particleName)
88 {fLowEnergyLimit[particleTable->FindParticle(particleName)->
89 GetParticleDefinitionID()] = ekinetic;}
90 void SetLindhardAngleNumberHighLimit(G4double angleNumber, const G4String& particleName)
91 {fLindhardAngleNumberHighLimit[particleTable->FindParticle(particleName)->
92 GetParticleDefinitionID()]=angleNumber;}
93 void SetHighAngleLimit(G4double anglemax, const G4String& particleName)
94 {fHighAngleLimit[particleTable->FindParticle(particleName)->
95 GetParticleDefinitionID()] = anglemax;}
96
98 {fDefaultLowEnergyLimit=ekinetic;}
100 {fDefaultLindhardAngleNumberHighLimit=angleNumber;}
102 {fDefaultHighAngleLimit=anglemax;}
103
104
105 /// get the maximal number of photons that can be produced per fastStep
106 /// Caution: is redundant, if the radiation model is not activated
108 {fMaxPhotonsProducedPerStep=nPhotons;}
109
110 ///get cuts
112 {return (fLowEnergyLimit.count(particleDefinitionID) == 1)
113 ? fLowEnergyLimit[particleDefinitionID]
114 : fDefaultLowEnergyLimit;}
116 {return (fLindhardAngleNumberHighLimit.count(particleDefinitionID) == 1)
117 ? fLindhardAngleNumberHighLimit[particleDefinitionID]
118 : fDefaultLindhardAngleNumberHighLimit;}
119 G4double GetHighAngleLimit(G4int particleDefinitionID)
120 {return (fHighAngleLimit.count(particleDefinitionID) == 1)
121 ? fHighAngleLimit[particleDefinitionID]
122 : fDefaultHighAngleLimit;}
123
124 /// get the maximal number of photons that can be produced per fastStep
125 G4int GetMaxPhotonsProducedPerStep(){return fMaxPhotonsProducedPerStep;}
126
127private:
128
129 G4ChannelingFastSimCrystalData* fCrystalData{nullptr};
130 G4BaierKatkov* fBaierKatkov{nullptr};
131
132 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
133
134 ///flag of radiation model
135 G4bool fRad = false;
136
137 /// maps of cuts (angular cuts are chosen as std::max of
138 /// fHighAngleLimit and calculated Lindhard angle)
139 std::unordered_map<G4int, G4double> fLowEnergyLimit;
140 std::unordered_map<G4int, G4double> fLindhardAngleNumberHighLimit;
141 std::unordered_map<G4int, G4double> fHighAngleLimit;
142
143 G4double fDefaultLowEnergyLimit = 200*CLHEP::MeV;
144 G4double fDefaultLindhardAngleNumberHighLimit = 100.;
145 G4double fDefaultHighAngleLimit = 0.;
146
147 /// the maximal number of photons that can be produced per fastStep
148 G4int fMaxPhotonsProducedPerStep=1000.;
149
150};
151#endif
152
153
154
155
Definition of the G4BaierKatkov class This class is designed for the calculation of radiation probabi...
Definition of the G4ChannelingFastSimCrystalData class The class inherits G4VChannelingFastSimCrystal...
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
int G4int
Definition G4Types.hh:85
void SetDefaultLindhardAngleNumberHighLimit(G4double angleNumber)
void SetDefaultHighAngleLimit(G4double anglemax)
G4ChannelingFastSimCrystalData * GetCrystalData()
G4ChannelingFastSimModel(const G4String &, G4Region *)
G4double GetHighAngleLimit(G4int particleDefinitionID)
void SetDefaultLowKineticEnergyLimit(G4double ekinetic)
void DoIt(const G4FastTrack &, G4FastStep &) override
– User method DoIt
void Input(const G4Material *crystal, const G4String &lattice)
special functions
G4bool IsApplicable(const G4ParticleDefinition &) override
– IsApplicable
void SetLindhardAngleNumberHighLimit(G4double angleNumber, const G4String &particleName)
G4double GetLindhardAngleNumberHighLimit(G4int particleDefinitionID)
void SetMaxPhotonsProducedPerStep(G4double nPhotons)
G4bool ModelTrigger(const G4FastTrack &) override
– ModelTrigger
void SetLowKineticEnergyLimit(G4double ekinetic, const G4String &particleName)
set cuts
void SetHighAngleLimit(G4double anglemax, const G4String &particleName)
G4int GetMaxPhotonsProducedPerStep()
get the maximal number of photons that can be produced per fastStep
G4double GetLowKineticEnergyLimit(G4int particleDefinitionID)
get cuts
static G4ParticleTable * GetParticleTable()
G4VFastSimulationModel(const G4String &aName)