File indexing completed on 2025-01-18 09:58:01
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 #ifndef G4ChannelingFastSimModel_h
0028 #define G4ChannelingFastSimModel_h 1
0029
0030 #include "G4VFastSimulationModel.hh"
0031 #include "G4Step.hh"
0032 #include "G4TouchableHandle.hh"
0033 #include <vector>
0034 #include <CLHEP/Units/SystemOfUnits.h>
0035 #include <CLHEP/Units/PhysicalConstants.h>
0036
0037 #include "G4ChannelingFastSimCrystalData.hh"
0038 #include <unordered_map>
0039 #include "G4BaierKatkov.hh"
0040 #include "G4LogicalVolume.hh"
0041 #include "G4ParticleTable.hh"
0042
0043
0044
0045
0046
0047
0048
0049
0050 class G4ChannelingFastSimModel : public G4VFastSimulationModel
0051 {
0052 public:
0053
0054 G4ChannelingFastSimModel (const G4String&, G4Region*);
0055 G4ChannelingFastSimModel (const G4String&);
0056 ~G4ChannelingFastSimModel ();
0057
0058
0059 G4bool IsApplicable(const G4ParticleDefinition&) override;
0060
0061 G4bool ModelTrigger(const G4FastTrack &) override;
0062
0063 void DoIt(const G4FastTrack&, G4FastStep&) override;
0064
0065
0066 void Input(const G4Material* crystal, const G4String &lattice);
0067
0068 void RadiationModelActivate();
0069
0070 G4ChannelingFastSimCrystalData* GetCrystalData() {return fCrystalData;}
0071
0072 G4BaierKatkov* GetRadiationModel() {return fBaierKatkov;}
0073
0074 G4bool GetIfRadiationModelActive(){return fRad;}
0075
0076
0077 void SetLowKineticEnergyLimit(G4double ekinetic, const G4String& particleName)
0078 {fLowEnergyLimit[particleTable->FindParticle(particleName)->
0079 GetParticleDefinitionID()] = ekinetic;}
0080 void SetLindhardAngleNumberHighLimit(G4double angleNumber, const G4String& particleName)
0081 {fLindhardAngleNumberHighLimit[particleTable->FindParticle(particleName)->
0082 GetParticleDefinitionID()]=angleNumber;}
0083
0084 void SetDefaultLowKineticEnergyLimit(G4double ekinetic)
0085 {fDefaultLowEnergyLimit=ekinetic;}
0086 void SetDefaultLindhardAngleNumberHighLimit(G4double angleNumber)
0087 {fDefaultLindhardAngleNumberHighLimit=angleNumber;}
0088
0089
0090
0091
0092 void SetMaxPhotonsProducedPerStep(G4double nPhotons)
0093 {fMaxPhotonsProducedPerStep=nPhotons;}
0094
0095
0096 G4double GetLowKineticEnergyLimit(const G4String& particleName)
0097 {return GetLowKineticEnergyLimit(particleTable->
0098 FindParticle(particleName)->
0099 GetParticleDefinitionID());}
0100 G4double GetLindhardAngleNumberHighLimit(const G4String& particleName)
0101 {return GetLindhardAngleNumberHighLimit(particleTable->
0102 FindParticle(particleName)->
0103 GetParticleDefinitionID());}
0104
0105 G4double GetLowKineticEnergyLimit(G4int particleDefinitionID)
0106 {return (fLowEnergyLimit.count(particleDefinitionID) == 1)
0107 ? fLowEnergyLimit[particleDefinitionID]
0108 : fDefaultLowEnergyLimit;}
0109 G4double GetLindhardAngleNumberHighLimit(G4int particleDefinitionID)
0110 {return (fLindhardAngleNumberHighLimit.count(particleDefinitionID) == 1)
0111 ? fLindhardAngleNumberHighLimit[particleDefinitionID]
0112 : fDefaultLindhardAngleNumberHighLimit;}
0113
0114
0115 G4int GetMaxPhotonsProducedPerStep(){return fMaxPhotonsProducedPerStep;}
0116
0117 private:
0118
0119 G4ChannelingFastSimCrystalData* fCrystalData{nullptr};
0120 G4BaierKatkov* fBaierKatkov{nullptr};
0121
0122 G4ParticleTable* particleTable = G4ParticleTable::GetParticleTable();
0123
0124
0125 G4bool fRad = false;
0126
0127
0128 std::unordered_map<G4int, G4double> fLowEnergyLimit;
0129 std::unordered_map<G4int, G4double> fLindhardAngleNumberHighLimit;
0130
0131 G4double fDefaultLowEnergyLimit = 200*CLHEP::MeV;
0132 G4double fDefaultLindhardAngleNumberHighLimit = 100.;
0133
0134
0135 G4int fMaxPhotonsProducedPerStep=1000.;
0136
0137 };
0138 #endif
0139
0140
0141
0142