![]() |
|
|||
File indexing completed on 2025-02-23 09:22:33
0001 // 0002 // ******************************************************************** 0003 // * License and Disclaimer * 0004 // * * 0005 // * The Geant4 software is copyright of the Copyright Holders of * 0006 // * the Geant4 Collaboration. It is provided under the terms and * 0007 // * conditions of the Geant4 Software License, included in the file * 0008 // * LICENSE and available at http://cern.ch/geant4/license . These * 0009 // * include a list of copyright holders. * 0010 // * * 0011 // * Neither the authors of this software system, nor their employing * 0012 // * institutes,nor the agencies providing financial support for this * 0013 // * work make any representation or warranty, express or implied, * 0014 // * regarding this software system or assume any liability for its * 0015 // * use. Please see the license in the file LICENSE and URL above * 0016 // * for the full disclaimer and the limitation of liability. * 0017 // * * 0018 // * This code implementation is the result of the scientific and * 0019 // * technical work of the GEANT4 collaboration. * 0020 // * By using, copying, modifying or distributing the software (or * 0021 // * any work based on the software) you agree to acknowledge its * 0022 // * use in resulting scientific publications, and indicate your * 0023 // * acceptance of all terms of the Geant4 Software license. * 0024 // ******************************************************************** 0025 // 0026 #ifndef PAR03EMSHOWERMODEL_HH 0027 #define PAR03EMSHOWERMODEL_HH 0028 0029 #include "G4VFastSimulationModel.hh" 0030 0031 class Par03EMShowerMessenger; 0032 class G4FastSimHitMaker; 0033 0034 /** 0035 * @brief Example fast simulation model for EM showers. 0036 * 0037 * Parametrisation of electrons, positrons, and gammas. It is triggered if those 0038 * particles enter the detector so that there is sufficient length for the 0039 * shower development (max depth, controlled by the UI command). 0040 * 0041 * Parametrisation is based on the PDG chapter on the electromagnetic cascades 0042 * (chapter 33.5). Longitudinal profile of the shower is described with Gamma 0043 * distribution, with beta parameter on average equal to 0.5 (default value, 0044 * Fig. 33.21), and alpha parameter calcluated from the incident particle energy 0045 * and material of the detector (critical energy) following Eq.(33.36). 0046 * 0047 * Transverse profile is in this model approximated by the Gaussian 0048 * distribution, with the mean along the shower axis (incident particle momentum 0049 * direction) and the standard deviation calculated from the detector material 0050 * (Moliere radius). This assumes that EM shower is in 90% contained within a 0051 * cylinder of radius equal to Moliere radius, and that area below Gaussian 0052 * distribution from `mean-1.645 sigma` to `mean+1.645 sigma` is also equal to 0053 * 90% of total distribution. 0054 * 0055 * Parameters of both distributions (alpha, beta for Gamma, sigma for Gaussian) 0056 * can be overwritten by UI commands. 0057 * 0058 * Parametrisation creates N hits of same energy (N can be set by UI command), 0059 * using rejection sampling to generate position along shower axis from Gamma 0060 * distribution, and then sampling from uniform and Gaussian distributions to 0061 * sample phi and radius, respectively. Created hits are deposited in the 0062 * detector using its readout geometry, using the helper class G4FastSimHitMaker 0063 * that locates the volume, and calls appropriate sensitive detector class. 0064 * 0065 * PDG Chapter 33: 0066 * https://pdg.lbl.gov/2019/reviews/rpp2018-rev-passage-particles-matter.pdf 0067 * 0068 */ 0069 0070 class Par03EMShowerModel : public G4VFastSimulationModel 0071 { 0072 public: 0073 Par03EMShowerModel(G4String, G4Region*); 0074 Par03EMShowerModel(G4String); 0075 ~Par03EMShowerModel(); 0076 0077 /// There are no kinematics constraints. True is returned. 0078 virtual G4bool ModelTrigger(const G4FastTrack&) final; 0079 /// Model is applicable to electrons, positrons, and photons. 0080 virtual G4bool IsApplicable(const G4ParticleDefinition&) final; 0081 /// Take particle out of the full simulation (kill it at the entrance 0082 /// depositing all the energy). Calculate energy deposited in the detector 0083 /// according to Gamma distribution (along the particle direction) and 0084 /// Gaussian distribution in the transverse direction. Mean of the Gaussian is 0085 /// centred on the shower axis. Create energy deposits on a cylindrical mesh. 0086 /// Parameters of the mesh (size, number of cells) and of the distributions 0087 /// (alpha, beta for Gamma, sigma for Gaussian) can be set with UI commands. 0088 virtual void DoIt(const G4FastTrack&, G4FastStep&) final; 0089 0090 /// Print current settings. 0091 void Print() const; 0092 /// Set standard deviation of a Gaussian distribution that describes the 0093 /// transverse shower profile. 0094 inline void SetSigma(const G4double aSigma) { fSigma = aSigma; }; 0095 /// Get standard deviation of a Gaussian distribution that describes the 0096 /// transverse shower profile. 0097 inline G4double GetSigma() const { return fSigma; }; 0098 /// Set alpha parameter of a Gamma distribution that describes the 0099 /// longitudinal shower profile. 0100 inline void SetAlpha(const G4double aAlpha) { fAlpha = aAlpha; }; 0101 /// Get alpha parameter of a Gamma distribution that describes the 0102 /// longitudinal shower profile. 0103 inline G4double GetAlpha() const { return fAlpha; }; 0104 /// Set beta parameter of a Gamma distribution that describes the longitudinal 0105 /// shower profile. 0106 inline void SetBeta(const G4double aBeta) { fBeta = aBeta; }; 0107 /// Get beta parameter of a Gamma distribution that describes the longitudinal 0108 /// shower profile. 0109 inline G4double GetBeta() const { return fBeta; }; 0110 /// Set number of (same energy) hits created in the parametrisation. 0111 inline void SetNbOfHits(const G4int aNumber) { fNbOfHits = aNumber; }; 0112 /// Get number of (same energy) hits created in the parametrisation.s 0113 inline G4int GetNbOfHits() const { return fNbOfHits; }; 0114 /// Set maximum depth of shower created in fast simulation. It is expressed in 0115 /// units of radiaton length. 0116 inline void SetLongMaxDepth(const G4double aDepth) { fLongMaxDepth = aDepth; }; 0117 /// Get maximum depth of shower created in fast simulation. It is expressed in 0118 /// units of radiaton length. 0119 inline G4double GetLongMaxDepth() const { return fLongMaxDepth; }; 0120 0121 private: 0122 /// Gamma distribution 0123 inline G4double Gamma(G4double x, G4double alpha, G4double beta) 0124 { 0125 return (std::pow(beta, alpha) / std::tgamma(alpha) * std::pow(x, alpha - 1) 0126 * std::exp(-beta * x)); 0127 } 0128 /// Gaussian distribution 0129 inline G4double Gaussian(G4double x, G4double sigma = 1, G4double x0 = 0) 0130 { 0131 G4double tmp = (x - x0) / sigma; 0132 return (1.0 / (std::sqrt(2 * CLHEP::pi) * sigma)) * std::exp(-tmp * tmp / 2); 0133 } 0134 0135 private: 0136 /// Messenger for configuration 0137 Par03EMShowerMessenger* fMessenger; 0138 /// Helper class for creation of hits within the sensitive detector 0139 std::unique_ptr<G4FastSimHitMaker> fHitMaker; 0140 /// Standard deviation of the Gaussian distribution 0141 /// Can be changed with UI command `/Par03/fastSim/transverseProfile/sigma 0142 /// <sigma>` 0143 /// If sigma is smaller than 0, it will be estimated from the detector 0144 /// material (Moliere radius). 0145 G4double fSigma = -1; 0146 /// Alpha parameter of the Gamma distribution 0147 /// Can be changed with UI command `/Par03/fastSim/longitudunalProfile/alpha 0148 /// <alpha>` 0149 /// If alpha is smaller than 0, it will be estimated from particle energy and 0150 /// the detector material. 0151 G4double fAlpha = -1; 0152 /// Beta parameter of the Gamma distribution 0153 /// Can be changed with UI command `/Par03/fastSim/longitudinalProfile/beta 0154 /// <beta>` 0155 G4double fBeta = 0.5; 0156 /// Number of (same energy) hits created by the parametrisation. Can be 0157 /// changed with UI command `/Par03/fastSim/numberOfHits <number>` 0158 G4int fNbOfHits = 100; 0159 /// Maximum depth of a shower created in fast simulation. 0160 /// It is expressed in units of radiation length. Can be changed with UI 0161 /// command `/Par03/fastSim/longitudinalProfile/maxDepth <depth>` 0162 G4double fLongMaxDepth = 30; 0163 }; 0164 #endif /* PAR03EMSHOWERMODEL_HH */
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |
![]() ![]() |