Geant4 11.4.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 "globals.hh"
36#include "G4ios.hh"
37
39#include <unordered_map>
40#include "G4BaierKatkov.hh"
41#include "G4LogicalVolume.hh"
42#include "G4ParticleTable.hh"
43
44/** \file G4ChannelingFastSimModel.hh
45* \brief Definition of the G4ChannelingFastSimModel class
46* FastSimulation Channeling model: calculates charge particle trajectories
47* in oriented crystals in the field of crystal planes/axes either straight or bent.
48* It is also possible to simulate radiation using Baier-Katkov method.
49*/
50
52{
53public:
54 // Constructor, destructor
58
59 /// -- IsApplicable
61 /// -- ModelTrigger
62 G4bool ModelTrigger(const G4FastTrack &) override;
63 /// -- User method DoIt
64 void DoIt(const G4FastTrack&, G4FastStep&) override;
65
66 ///special functions
67 void Input(const G4Material* crystal,
68 const G4String &lattice)
69 {Input(crystal,lattice,"");}
70
71 void Input(const G4Material* crystal,
72 const G4String &lattice,
73 const G4String &filePath);
74
76
78
79 G4BaierKatkov* GetRadiationModel() {return fBaierKatkov;}
80
82
83 ///set cuts
84 void SetLowKineticEnergyLimit(G4double ekinetic, const G4String& particleName)
85 {fLowEnergyLimit[particleTable->FindParticle(particleName)->
86 GetParticleDefinitionID()] = ekinetic;}
87 void SetLindhardAngleNumberHighLimit(G4double angleNumber, const G4String& particleName)
88 {fLindhardAngleNumberHighLimit[particleTable->FindParticle(particleName)->
89 GetParticleDefinitionID()]=angleNumber;}
90 void SetHighAngleLimit(G4double anglemax, const G4String& particleName)
91 {fHighAngleLimit[particleTable->FindParticle(particleName)->
92 GetParticleDefinitionID()] = anglemax;}
93
95 {fDefaultLowEnergyLimit=ekinetic;}
97 {fDefaultLindhardAngleNumberHighLimit=angleNumber;}
99 {fDefaultHighAngleLimit=anglemax;}
100
101
102 /// get the maximal number of photons that can be produced per fastStep
103 /// Caution: is redundant, if the radiation model is not activated
105 {fMaxPhotonsProducedPerStep=nPhotons;}
106
107 ///get cuts
109 {return (fLowEnergyLimit.count(particleDefinitionID) == 1)
110 ? fLowEnergyLimit[particleDefinitionID]
111 : fDefaultLowEnergyLimit;}
113 {return (fLindhardAngleNumberHighLimit.count(particleDefinitionID) == 1)
114 ? fLindhardAngleNumberHighLimit[particleDefinitionID]
115 : fDefaultLindhardAngleNumberHighLimit;}
116 G4double GetHighAngleLimit(G4int particleDefinitionID)
117 {return (fHighAngleLimit.count(particleDefinitionID) == 1)
118 ? fHighAngleLimit[particleDefinitionID]
119 : fDefaultHighAngleLimit;}
120
121 /// get the maximal number of photons that can be produced per fastStep
122 G4int GetMaxPhotonsProducedPerStep(){return fMaxPhotonsProducedPerStep;}
123
124private:
125
126 G4ChannelingFastSimCrystalData* fCrystalData{nullptr};
127 G4BaierKatkov* fBaierKatkov{nullptr};
128
129 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
130
131 ///flag of radiation model
132 G4bool fRad = false;
133
134 /// maps of cuts (angular cuts are chosen as std::max of
135 /// fHighAngleLimit and calculated Lindhard angle)
136 std::unordered_map<G4int, G4double> fLowEnergyLimit;
137 std::unordered_map<G4int, G4double> fLindhardAngleNumberHighLimit;
138 std::unordered_map<G4int, G4double> fHighAngleLimit;
139
140 G4double fDefaultLowEnergyLimit = 200*CLHEP::MeV;
141 G4double fDefaultLindhardAngleNumberHighLimit = 100.;
142 G4double fDefaultHighAngleLimit = 0.;
143
144 /// the maximal number of photons that can be produced per fastStep
145 G4int fMaxPhotonsProducedPerStep=1000.;
146
147};
148#endif
149
150
151
152
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)
~G4ChannelingFastSimModel()=default
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()
G4Region defines a region or a group of regions in the detector geometry setup, sharing properties as...
Definition G4Region.hh:90
G4VFastSimulationModel(const G4String &aName)