File indexing completed on 2025-09-16 08:55:51
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029
0030
0031 #ifndef G4ChannelingFastSimModel_h
0032 #define G4ChannelingFastSimModel_h 1
0033
0034 #include "G4VFastSimulationModel.hh"
0035 #include "G4Step.hh"
0036 #include "G4TouchableHandle.hh"
0037 #include <vector>
0038 #include <CLHEP/Units/SystemOfUnits.h>
0039 #include <CLHEP/Units/PhysicalConstants.h>
0040
0041 #include "G4ChannelingFastSimCrystalData.hh"
0042 #include <unordered_map>
0043 #include "G4BaierKatkov.hh"
0044 #include "G4LogicalVolume.hh"
0045 #include "G4ParticleTable.hh"
0046
0047
0048
0049
0050
0051
0052
0053
0054 class G4ChannelingFastSimModel : public G4VFastSimulationModel
0055 {
0056 public:
0057
0058 G4ChannelingFastSimModel (const G4String&, G4Region*);
0059 G4ChannelingFastSimModel (const G4String&);
0060 ~G4ChannelingFastSimModel ();
0061
0062
0063 G4bool IsApplicable(const G4ParticleDefinition&) override;
0064
0065 G4bool ModelTrigger(const G4FastTrack &) override;
0066
0067 void DoIt(const G4FastTrack&, G4FastStep&) override;
0068
0069
0070 void Input(const G4Material* crystal,
0071 const G4String &lattice)
0072 {Input(crystal,lattice,"");}
0073
0074 void Input(const G4Material* crystal,
0075 const G4String &lattice,
0076 const G4String &filePath);
0077
0078 void RadiationModelActivate();
0079
0080 G4ChannelingFastSimCrystalData* GetCrystalData() {return fCrystalData;}
0081
0082 G4BaierKatkov* GetRadiationModel() {return fBaierKatkov;}
0083
0084 G4bool GetIfRadiationModelActive(){return fRad;}
0085
0086
0087 void SetLowKineticEnergyLimit(G4double ekinetic, const G4String& particleName)
0088 {fLowEnergyLimit[particleTable->FindParticle(particleName)->
0089 GetParticleDefinitionID()] = ekinetic;}
0090 void SetLindhardAngleNumberHighLimit(G4double angleNumber, const G4String& particleName)
0091 {fLindhardAngleNumberHighLimit[particleTable->FindParticle(particleName)->
0092 GetParticleDefinitionID()]=angleNumber;}
0093 void SetHighAngleLimit(G4double anglemax, const G4String& particleName)
0094 {fHighAngleLimit[particleTable->FindParticle(particleName)->
0095 GetParticleDefinitionID()] = anglemax;}
0096
0097 void SetDefaultLowKineticEnergyLimit(G4double ekinetic)
0098 {fDefaultLowEnergyLimit=ekinetic;}
0099 void SetDefaultLindhardAngleNumberHighLimit(G4double angleNumber)
0100 {fDefaultLindhardAngleNumberHighLimit=angleNumber;}
0101 void SetDefaultHighAngleLimit(G4double anglemax)
0102 {fDefaultHighAngleLimit=anglemax;}
0103
0104
0105
0106
0107 void SetMaxPhotonsProducedPerStep(G4double nPhotons)
0108 {fMaxPhotonsProducedPerStep=nPhotons;}
0109
0110
0111 G4double GetLowKineticEnergyLimit(G4int particleDefinitionID)
0112 {return (fLowEnergyLimit.count(particleDefinitionID) == 1)
0113 ? fLowEnergyLimit[particleDefinitionID]
0114 : fDefaultLowEnergyLimit;}
0115 G4double GetLindhardAngleNumberHighLimit(G4int particleDefinitionID)
0116 {return (fLindhardAngleNumberHighLimit.count(particleDefinitionID) == 1)
0117 ? fLindhardAngleNumberHighLimit[particleDefinitionID]
0118 : fDefaultLindhardAngleNumberHighLimit;}
0119 G4double GetHighAngleLimit(G4int particleDefinitionID)
0120 {return (fHighAngleLimit.count(particleDefinitionID) == 1)
0121 ? fHighAngleLimit[particleDefinitionID]
0122 : fDefaultHighAngleLimit;}
0123
0124
0125 G4int GetMaxPhotonsProducedPerStep(){return fMaxPhotonsProducedPerStep;}
0126
0127 private:
0128
0129 G4ChannelingFastSimCrystalData* fCrystalData{nullptr};
0130 G4BaierKatkov* fBaierKatkov{nullptr};
0131
0132 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
0133
0134
0135 G4bool fRad = false;
0136
0137
0138
0139 std::unordered_map<G4int, G4double> fLowEnergyLimit;
0140 std::unordered_map<G4int, G4double> fLindhardAngleNumberHighLimit;
0141 std::unordered_map<G4int, G4double> fHighAngleLimit;
0142
0143 G4double fDefaultLowEnergyLimit = 200*CLHEP::MeV;
0144 G4double fDefaultLindhardAngleNumberHighLimit = 100.;
0145 G4double fDefaultHighAngleLimit = 0.;
0146
0147
0148 G4int fMaxPhotonsProducedPerStep=1000.;
0149
0150 };
0151 #endif
0152
0153
0154
0155