Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:58:01

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 
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 /** \file G4ChannelingFastSimModel.hh
0044 * \brief Definition of the G4ChannelingFastSimModel class
0045 * FastSimulation Channeling model: calculates charge particle trajectories
0046 * in oriented crystals in the field of crystal planes/axes either straight or bent.
0047 * It is also possible to simulate radiation using Baier-Katkov method.
0048 */
0049 
0050 class G4ChannelingFastSimModel : public G4VFastSimulationModel
0051 {
0052 public:
0053   // Constructor, destructor
0054   G4ChannelingFastSimModel (const G4String&, G4Region*);
0055   G4ChannelingFastSimModel (const G4String&);
0056   ~G4ChannelingFastSimModel ();
0057 
0058   /// -- IsApplicable
0059   G4bool IsApplicable(const G4ParticleDefinition&) override;
0060   /// -- ModelTrigger
0061   G4bool ModelTrigger(const G4FastTrack &) override;
0062   /// -- User method DoIt
0063   void DoIt(const G4FastTrack&, G4FastStep&) override;
0064 
0065   ///special functions
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   ///set cuts
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   /// get the maximal number of photons that can be produced per fastStep
0091   /// Caution: is redundant, if the radiation model is not activated
0092   void SetMaxPhotonsProducedPerStep(G4double nPhotons)
0093            {fMaxPhotonsProducedPerStep=nPhotons;}
0094 
0095   ///get cuts
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   //the same functions but using particleDefinitionID (needed for faster model execution)
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   /// get the maximal number of photons that can be produced per fastStep
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   ///flag of radiation model
0125   G4bool fRad = false;
0126 
0127   /// maps of cuts
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   /// the maximal number of photons that can be produced per fastStep
0135   G4int fMaxPhotonsProducedPerStep=1000.;
0136 
0137 };
0138 #endif
0139 
0140 
0141 
0142