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