Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-09-16 08:55:51

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 // Author:      Alexei Sytov
0027 // Co-author:   Gianfranco PaternĂ² (modifications & testing)
0028 // On the base of the CRYSTALRAD realization of channeling model:
0029 // A. I. Sytov, V. V. Tikhomirov, and L. Bandiera PRAB 22, 064601 (2019)
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 /** \file G4ChannelingFastSimModel.hh
0048 * \brief Definition of the G4ChannelingFastSimModel class
0049 * FastSimulation Channeling model: calculates charge particle trajectories
0050 * in oriented crystals in the field of crystal planes/axes either straight or bent.
0051 * It is also possible to simulate radiation using Baier-Katkov method.
0052 */
0053 
0054 class G4ChannelingFastSimModel : public G4VFastSimulationModel
0055 {
0056 public:
0057   // Constructor, destructor
0058   G4ChannelingFastSimModel (const G4String&, G4Region*);
0059   G4ChannelingFastSimModel (const G4String&);
0060   ~G4ChannelingFastSimModel ();
0061 
0062   /// -- IsApplicable
0063   G4bool IsApplicable(const G4ParticleDefinition&) override;
0064   /// -- ModelTrigger
0065   G4bool ModelTrigger(const G4FastTrack &) override;
0066   /// -- User method DoIt
0067   void DoIt(const G4FastTrack&, G4FastStep&) override;
0068 
0069   ///special functions
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   ///set cuts
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   /// get the maximal number of photons that can be produced per fastStep
0106   /// Caution: is redundant, if the radiation model is not activated
0107   void SetMaxPhotonsProducedPerStep(G4double nPhotons)
0108            {fMaxPhotonsProducedPerStep=nPhotons;}
0109 
0110   ///get cuts
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   /// get the maximal number of photons that can be produced per fastStep
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   ///flag of radiation model
0135   G4bool fRad = false;
0136 
0137   /// maps of cuts (angular cuts are chosen as std::max of
0138   /// fHighAngleLimit and calculated Lindhard angle)
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   /// the maximal number of photons that can be produced per fastStep
0148   G4int fMaxPhotonsProducedPerStep=1000.;
0149 
0150 };
0151 #endif
0152 
0153 
0154 
0155