|
||||
File indexing completed on 2025-01-18 09:59:19
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 0028 #ifndef G4VChannelingFastSimCrystalData_h 0029 #define G4VChannelingFastSimCrystalData_h 1 0030 0031 #include "G4ios.hh" 0032 #include "globals.hh" 0033 #include "G4ThreeVector.hh" 0034 #include "Randomize.hh" 0035 #include "G4LogicalVolume.hh" 0036 #include "G4Material.hh" 0037 #include "G4VSolid.hh" 0038 #include <unordered_map> 0039 0040 #include "G4ChannelingFastSimInterpolation.hh" 0041 0042 /** \file G4VChannelingFastSimCrystalData.hh 0043 * \brief Definition of the G4VChannelingFastSimCrystalData class 0044 * The class contains the data and properties related to the crystal lattice as well as 0045 * functions to simulate of important physical processes, i.e. coulomb scattering on 0046 * screened atomic potential, on single electrons and ionization energy losses; 0047 * functions of electric fields, nuclear and electron densities and minimum energy 0048 * of ionization (the corresponding interpolation coefficients are in 0049 * G4ChannelingFastSimInterpolation). 0050 * The functions related to the crystal geometry (transformation of coordinates and angles 0051 * from the reference system of the bounding box of the local volume to 0052 * the crystal lattice co-rotating reference system and vice versa) and 0053 * initialization function SetMaterialProperties are created as virtual to make 0054 * material data input and geometry functions flexible for modification. 0055 */ 0056 0057 class G4VChannelingFastSimCrystalData{ 0058 public: 0059 0060 G4VChannelingFastSimCrystalData(); 0061 virtual ~G4VChannelingFastSimCrystalData(); 0062 0063 ///electric fields produced by crystal lattice 0064 G4double Ex(G4double x,G4double y) {return (fElectricFieldX->GetIF(x,y))*(-fZ2/fPV);} 0065 G4double Ey(G4double x,G4double y) {return (fElectricFieldY->GetIF(x,y))*(-fZ2/fPV);} 0066 0067 ///electron density function 0068 G4double ElectronDensity(G4double x,G4double y) 0069 { 0070 G4double nel0=fElectronDensity->GetIF(x,y); 0071 if(nel0<0.) {nel0=0.;}//exception, errors of interpolation functions 0072 return nel0; 0073 } 0074 ///minimum energy of ionization function 0075 G4double MinIonizationEnergy(G4double x,G4double y) 0076 {return fMinIonizationEnergy->GetIF(x,y);} 0077 ///nuclear density function (normalized to average nuclear density) 0078 G4double NuclearDensity(G4double x,G4double y, G4int ielement) 0079 {return std::abs(fNucleiDensity[ielement]->GetIF(x,y));} 0080 //abs to describe exception, errors of interpolation functions, 0081 //don't put it =0, otherwise division on 0 in CoulombAtomicScattering 0082 0083 ///Calculate the value of the Lindhard angle (!!! the value for a straight crystal) 0084 G4double GetLindhardAngle(G4double etotal, G4double mass); 0085 ///Calculate the value of the Lindhard angle (!!! the value for a straight crystal) 0086 G4double GetLindhardAngle();//return the Lindhard angle value calculated in 0087 //SetParticleProperties 0088 0089 ///Calculate simulation step (standard value for channeling particles and 0090 ///reduced value for overbarrier particles) 0091 G4double GetSimulationStep(G4double tx,G4double ty); 0092 ///Calculate maximal simulation step (standard value for channeling particles) 0093 G4double GetMaxSimulationStep(G4double etotal, G4double mass); 0094 0095 ///get particle velocity/c 0096 G4double GetBeta(){return fBeta;} 0097 0098 G4int GetNelements() {return fNelements;} 0099 G4int GetModel() {return iModel;}//=1 for planes, =2 for axes 0100 0101 ///get bending angle of the crystal planes/axes 0102 ///(default BendingAngle=0 => straight crystal); 0103 G4double GetBendingAngle(){return fBendingAngle;} 0104 ///fBendingAngle MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME: 0105 ///THE VOLUME OF A BENT CRYSTAL MAY BE G4Box, while the planes/axes inside may be bent 0106 0107 G4double GetMiscutAngle(){return fMiscutAngle;} 0108 0109 ///get crystal curvature 0110 ///for crystalline undulator the curvature is a function, otherwise it's a constant 0111 G4double GetCurv(G4double z){return fCU ? -fCUK2*GetCUx(z) : fCurv;} 0112 0113 ///get crystalline undulator wave function 0114 G4double GetCUx(G4double z){return fCUAmplitude*std::cos(fCUK*z+fCUPhase);} 0115 ///get crystalline undulator wave 1st derivative function 0116 G4double GetCUtetax(G4double z){ 0117 return fCU ? -fCUAmplitudeK*std::sin(fCUK*z+fCUPhase) : 0;} 0118 0119 ///find and upload crystal lattice input files, calculate all the basic values 0120 ///(to do only once) 0121 virtual void SetMaterialProperties(const G4Material* crystal, 0122 const G4String &lattice) = 0; 0123 0124 ///set geometry parameters from current logical volume 0125 void SetGeometryParameters(const G4LogicalVolume *crystallogic); 0126 0127 ///set bending angle of the crystal planes/axes 0128 ///(default fBendingAngle=0 => straight crystal); 0129 ///only non-negative values! crystal is bent in the positive direction of x 0130 void SetBendingAngle(G4double tetab, const G4LogicalVolume *crystallogic); 0131 ///fBendingAngle MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME 0132 ///THE VOLUME OF A BENT CRYSTAL MAY BE G4Box, while the planes/axes inside may be bent 0133 0134 ///set miscut angle (default fMiscutAngle=0), acceptable range +-1 mrad, 0135 ///otherwise geometry routines may be unstable 0136 void SetMiscutAngle(G4double tetam, const G4LogicalVolume *crystallogic); 0137 0138 ///set crystalline undulator parameters: amplitude, period and phase 0139 /// (default: all 3 value = 0) 0140 /// function to use in Detector Construction 0141 void SetCrystallineUndulatorParameters(G4double amplitude, 0142 G4double period, 0143 G4double phase, 0144 const G4LogicalVolume *crystallogic); 0145 0146 ///set crystalline undulator parameters (internal function of the model) 0147 ///for convenience we put amplitude, period and phase in a G4ThreeVector 0148 void SetCUParameters(const G4ThreeVector &litudePeriodPhase, 0149 const G4LogicalVolume *crystallogic); 0150 0151 ///recalculate all the important values 0152 ///(to do both at the trajectory start and after energy loss) 0153 void SetParticleProperties(G4double etotal, 0154 G4double mp, 0155 G4double charge, 0156 G4bool ifhadron); 0157 0158 ///calculate the coordinates in the co-rotating reference system 0159 ///within a channel (periodic cell) 0160 ///(connected with crystal planes/axes either bent or straight) 0161 virtual G4ThreeVector CoordinatesFromBoxToLattice(const G4ThreeVector &pos0) = 0; 0162 0163 ///calculate the coordinates in the Box reference system 0164 ///(connected with the bounding box of the volume) 0165 virtual G4ThreeVector CoordinatesFromLatticeToBox(const G4ThreeVector &pos) = 0; 0166 0167 ///change the channel if necessary, recalculate x o y 0168 virtual G4ThreeVector ChannelChange(G4double& x, G4double& y, G4double& z) = 0; 0169 0170 ///return correction of the longitudinal coordinate 0171 /// (along current plane/axis vs "central plane/axis") 0172 G4double GetCorrectionZ(){return fCorrectionZ;} 0173 0174 ///calculate the horizontal angle in the co-rotating reference system 0175 ///within a channel (periodic cell) 0176 ///(connected with crystal planes/axes either bent or straight) 0177 virtual G4double AngleXFromBoxToLattice(G4double tx, G4double z)=0; 0178 0179 ///calculate the horizontal angle in the Box reference system 0180 ///(connected with the bounding box of the volume) 0181 virtual G4double AngleXFromLatticeToBox(G4double tx, G4double z)=0; 0182 0183 ///auxialiary function to transform the horizontal angle 0184 virtual G4double AngleXShift(G4double z)=0; 0185 0186 ///multiple and single scattering on screened potential 0187 G4ThreeVector CoulombAtomicScattering( 0188 G4double effectiveStep, 0189 G4double step, 0190 G4int ielement); 0191 ///multiple and single scattering on electrons 0192 G4ThreeVector CoulombElectronScattering(G4double eMinIonization, 0193 G4double electronDensity, 0194 G4double step); 0195 ///ionization losses 0196 G4double IonizationLosses(G4double dz, G4int ielement); 0197 0198 void SetVerbosity(G4int ver){fVerbosity = ver;} 0199 0200 protected: 0201 ///classes containing interpolation coefficients 0202 //horizontal electric field data 0203 G4ChannelingFastSimInterpolation* fElectricFieldX{nullptr}; 0204 //vertical electric field data 0205 G4ChannelingFastSimInterpolation* fElectricFieldY{nullptr}; 0206 //electron density data 0207 G4ChannelingFastSimInterpolation* fElectronDensity{nullptr}; 0208 //minimal energy of ionization data 0209 G4ChannelingFastSimInterpolation* fMinIonizationEnergy{nullptr}; 0210 //nuclear density distributions data 0211 std::vector <G4ChannelingFastSimInterpolation*> fNucleiDensity; 0212 0213 ///values related to the crystal geometry 0214 G4ThreeVector fHalfDimBoundingBox;//bounding box half dimensions 0215 G4int fBent=0;//flag of bent crystal, 0216 //=0 for straight and =1 for bent, by default straight crystal 0217 0218 G4double fBendingAngle=0.;// angle of bending of the crystal planes/axes 0219 //inside the crystal volume 0220 //MAY BE NOT THE SAME AS THE BENDING ANGLE OF THE CRYSTAL VOLUME 0221 //THE VOLUME OF A BENT CRYSTAL MAY BE G4Box, 0222 //while the planes/axes inside may be bent 0223 G4double fBendingR = 0.; // bending radius of the crystal planes/axes 0224 G4double fBending2R=0.; // =2*fBendingR 0225 G4double fBendingRsquare=0.; // =fBendingR**2 0226 G4double fCurv=0.; //=1/fBendingR bending curvature of the crystal planes/axes 0227 0228 G4double fMiscutAngle = 0.;// miscut angle, can be of either sign or 0; 0229 //safe values |ThetaMiscut|<0.001 0230 G4double fCosMiscutAngle=1.;// = std::cos(fMiscutAngle), to economy operations 0231 G4double fSinMiscutAngle=0.;// = std::sin(fMiscutAngle), to economy operations 0232 0233 G4double fCorrectionZ = 1.;//correction of the longitudinal coordinate 0234 //(along current plane/axis vs "central plane/axis"), 1 is default value 0235 //(for "central plane/axis" or a straight crystal) 0236 0237 G4bool fCU = false;//flag of crystalline undulator geometry 0238 //(periodically bent crystal) 0239 G4double fCUAmplitude=0.; //Amplitude of a crystalline undulator 0240 G4double fCUK=0.; //2*pi/period of a crystalline undulator 0241 G4double fCUPhase=0.;//Phase of a crystalline undulator 0242 G4double fCUAmplitudeK=0.;//fCUAmplitude*fCUK 0243 G4double fCUK2=0.; //fCUK^2 0244 0245 ///values related to the crystal lattice 0246 G4int fNelements=1;//number of nuclear elements in a crystal 0247 G4int iModel=1;// model type (iModel=1 for interplanar potential, 0248 //iModel=2 for the interaxial one) 0249 0250 G4double fVmax=0; // the height of the potential well 0251 G4double fVmax2=0; // =2*fVmax 0252 G4double fVMinCrystal=0;// non-zero minimal potential inside the crystal, 0253 // necessary for angle recalculation for entrance/exit 0254 //through the crystal lateral surface 0255 0256 G4double fChangeStep=0;// fChannelingStep = fChangeStep/fTetaL 0257 0258 std::vector <G4double> fI0; //Mean excitation energy 0259 0260 std::vector <G4double> fRF;//Thomas-Fermi screening radius 0261 0262 ///angles necessary for multiple and single coulomb scattering 0263 0264 //minimal scattering angle by coulomb scattering on nuclei 0265 //defined by shielding by electrons 0266 std::vector <G4double> fTeta10;//(in the Channeling model 0267 //teta1=fTeta10/fPz*(1.13+fK40/vz**2) 0268 //maximal scattering angle by coulomb scattering on nuclei defined by nucleus radius 0269 std::vector <G4double> fTetamax0;//(in the Channeling model tetamax=fTetamax0/fPz) 0270 std::vector <G4double> fTetamax2;//=tetamax*tetamax 0271 std::vector <G4double> fTetamax12;//=teta1*teta1+tetamax*tetamax 0272 std::vector <G4double> fTeta12; //= teta1*teta1 0273 0274 ///coefficients necessary for multiple and single coulomb scattering 0275 std::vector <G4double> fK20; //a useful coefficient, fK2=fK20/fPV/fPV 0276 std::vector <G4double> fK2; //a useful coefficient, 0277 //fK2=(fZ2*alpha*hdc)**2*4.*pi*fN0*(fZ1/fPV)**2 0278 std::vector <G4double> fK40; //a useful coefficient, fK40=3.76D0*(alpha*fZ1)**2 0279 G4double fK30=0;//a useful coefficient, fK3=fK30/fPV/fPV 0280 G4double fK3=0;//a useful coefficient, fK3=2.*pi*alpha*hdc/electron_mass_c2/(fPV)**2 0281 0282 std::vector <G4double> fKD; //a useful coefficient for dE/dx 0283 0284 ///coefficients for multiple scattering suppression 0285 std::vector <G4double> fPu11;//a useful coefficient for exponent containing u1 0286 std::vector <G4double> fPzu11;//a useful coefficient for exponent containing u1 0287 std::vector <G4double> fBB;//a useful coefficient 0288 std::vector <G4double> fE1XBbb;//a useful coefficient 0289 std::vector <G4double> fBBDEXP;//a useful coefficient 0290 0291 //Variable to control printout 0292 G4int fVerbosity = 1; 0293 0294 0295 private: 0296 0297 //exponential integral 0298 G4double expint(G4double x); 0299 0300 ///private variables 0301 0302 std::unordered_map<G4int, G4double> fMapBendingAngle;//the map fBendingAngle 0303 //for different logical volumes 0304 0305 std::unordered_map<G4int, G4double> fMapMiscutAngle;//the map fMiscutAngle 0306 //for different logical volumes 0307 0308 std::unordered_map<G4int, G4ThreeVector> fMapCUAmplitudePeriodPhase;//the map of 0309 //AmplitudePeriodPhase 0310 //for different logical volumes 0311 0312 G4double fChannelingStep=0;// simulation step under the channeling conditions = 0313 //channeling oscillation length/fNsteps 0314 // channeling oscillation length: Biryukov book Eq. (1.24) 0315 0316 ///energy depended values 0317 G4double fPz=0; // particle momentum absolute value 0318 G4double fPV=0; // pv 0319 G4double fTetaL=0; //Lindhard angle 0320 G4double fBeta=0; //particle (velocity/c) 0321 G4double fV2=0; // particle (velocity/c)^2 0322 G4double fGamma=0; //Lorentz factor 0323 G4double fMe2Gamma=0; // me^2*fGamma 0324 G4double fTmax=0; // max ionization losses 0325 0326 ///particle properties flags 0327 G4bool fHadron=false;//=true (for hadrons); =false (for leptons) 0328 G4double fZ2=0; //particle charge 0329 0330 }; 0331 0332 #endif 0333
[ Source navigation ] | [ Diff markup ] | [ Identifier search ] | [ general search ] |
This page was automatically generated by the 2.3.7 LXR engine. The LXR team |